MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex10.cpp
Go to the documentation of this file.
1 // MFEM Example 10
2 //
3 // Compile with: make ex10
4 //
5 // Sample runs:
6 // ex10 -m ../data/beam-quad.mesh -s 3 -r 2 -o 2 -dt 3
7 // ex10 -m ../data/beam-tri.mesh -s 3 -r 2 -o 2 -dt 3
8 // ex10 -m ../data/beam-hex.mesh -s 2 -r 1 -o 2 -dt 3
9 // ex10 -m ../data/beam-tet.mesh -s 2 -r 1 -o 2 -dt 3
10 // ex10 -m ../data/beam-wedge.mesh -s 2 -r 1 -o 2 -dt 3
11 // ex10 -m ../data/beam-quad.mesh -s 14 -r 2 -o 2 -dt 0.03 -vs 20
12 // ex10 -m ../data/beam-hex.mesh -s 14 -r 1 -o 2 -dt 0.05 -vs 20
13 // ex10 -m ../data/beam-quad-amr.mesh -s 3 -r 2 -o 2 -dt 3
14 //
15 // Description: This examples solves a time dependent nonlinear elasticity
16 // problem of the form dv/dt = H(x) + S v, dx/dt = v, where H is a
17 // hyperelastic model and S is a viscosity operator of Laplacian
18 // type. The geometry of the domain is assumed to be as follows:
19 //
20 // +---------------------+
21 // boundary --->| |
22 // attribute 1 | |
23 // (fixed) +---------------------+
24 //
25 // The example demonstrates the use of nonlinear operators (the
26 // class HyperelasticOperator defining H(x)), as well as their
27 // implicit time integration using a Newton method for solving an
28 // associated reduced backward-Euler type nonlinear equation
29 // (class ReducedSystemOperator). Each Newton step requires the
30 // inversion of a Jacobian matrix, which is done through a
31 // (preconditioned) inner solver. Note that implementing the
32 // method HyperelasticOperator::ImplicitSolve is the only
33 // requirement for high-order implicit (SDIRK) time integration.
34 //
35 // We recommend viewing examples 2 and 9 before viewing this
36 // example.
37 
38 #include "mfem.hpp"
39 #include <memory>
40 #include <iostream>
41 #include <fstream>
42 
43 using namespace std;
44 using namespace mfem;
45 
46 class ReducedSystemOperator;
47 
48 /** After spatial discretization, the hyperelastic model can be written as a
49  * system of ODEs:
50  * dv/dt = -M^{-1}*(H(x) + S*v)
51  * dx/dt = v,
52  * where x is the vector representing the deformation, v is the velocity field,
53  * M is the mass matrix, S is the viscosity matrix, and H(x) is the nonlinear
54  * hyperelastic operator.
55  *
56  * Class HyperelasticOperator represents the right-hand side of the above
57  * system of ODEs. */
58 class HyperelasticOperator : public TimeDependentOperator
59 {
60 protected:
61  FiniteElementSpace &fespace;
62 
63  BilinearForm M, S;
64  NonlinearForm H;
65  double viscosity;
66  HyperelasticModel *model;
67 
68  CGSolver M_solver; // Krylov solver for inverting the mass matrix M
69  DSmoother M_prec; // Preconditioner for the mass matrix M
70 
71  /** Nonlinear operator defining the reduced backward Euler equation for the
72  velocity. Used in the implementation of method ImplicitSolve. */
73  ReducedSystemOperator *reduced_oper;
74 
75  /// Newton solver for the reduced backward Euler equation
76  NewtonSolver newton_solver;
77 
78  /// Solver for the Jacobian solve in the Newton method
79  Solver *J_solver;
80  /// Preconditioner for the Jacobian solve in the Newton method
81  Solver *J_prec;
82 
83  mutable Vector z; // auxiliary vector
84 
85 public:
86  HyperelasticOperator(FiniteElementSpace &f, Array<int> &ess_bdr,
87  double visc, double mu, double K);
88 
89  /// Compute the right-hand side of the ODE system.
90  virtual void Mult(const Vector &vx, Vector &dvx_dt) const;
91  /** Solve the Backward-Euler equation: k = f(x + dt*k, t), for the unknown k.
92  This is the only requirement for high-order SDIRK implicit integration.*/
93  virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k);
94 
95  double ElasticEnergy(const Vector &x) const;
96  double KineticEnergy(const Vector &v) const;
97  void GetElasticEnergyDensity(const GridFunction &x, GridFunction &w) const;
98 
99  virtual ~HyperelasticOperator();
100 };
101 
102 /** Nonlinear operator of the form:
103  k --> (M + dt*S)*k + H(x + dt*v + dt^2*k) + S*v,
104  where M and S are given BilinearForms, H is a given NonlinearForm, v and x
105  are given vectors, and dt is a scalar. */
106 class ReducedSystemOperator : public Operator
107 {
108 private:
109  BilinearForm *M, *S;
110  NonlinearForm *H;
111  mutable SparseMatrix *Jacobian;
112  double dt;
113  const Vector *v, *x;
114  mutable Vector w, z;
115 
116 public:
117  ReducedSystemOperator(BilinearForm *M_, BilinearForm *S_, NonlinearForm *H_);
118 
119  /// Set current dt, v, x values - needed to compute action and Jacobian.
120  void SetParameters(double dt_, const Vector *v_, const Vector *x_);
121 
122  /// Compute y = H(x + dt (v + dt k)) + M k + S (v + dt k).
123  virtual void Mult(const Vector &k, Vector &y) const;
124 
125  /// Compute J = M + dt S + dt^2 grad_H(x + dt (v + dt k)).
126  virtual Operator &GetGradient(const Vector &k) const;
127 
128  virtual ~ReducedSystemOperator();
129 };
130 
131 
132 /** Function representing the elastic energy density for the given hyperelastic
133  model+deformation. Used in HyperelasticOperator::GetElasticEnergyDensity. */
134 class ElasticEnergyCoefficient : public Coefficient
135 {
136 private:
137  HyperelasticModel &model;
138  const GridFunction &x;
139  DenseMatrix J;
140 
141 public:
142  ElasticEnergyCoefficient(HyperelasticModel &m, const GridFunction &x_)
143  : model(m), x(x_) { }
144  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
145  virtual ~ElasticEnergyCoefficient() { }
146 };
147 
148 void InitialDeformation(const Vector &x, Vector &y);
149 
150 void InitialVelocity(const Vector &x, Vector &v);
151 
152 void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes,
153  GridFunction *field, const char *field_name = NULL,
154  bool init_vis = false);
155 
156 
157 int main(int argc, char *argv[])
158 {
159  // 1. Parse command-line options.
160  const char *mesh_file = "../data/beam-quad.mesh";
161  int ref_levels = 2;
162  int order = 2;
163  int ode_solver_type = 3;
164  double t_final = 300.0;
165  double dt = 3.0;
166  double visc = 1e-2;
167  double mu = 0.25;
168  double K = 5.0;
169  bool visualization = true;
170  int vis_steps = 1;
171 
172  OptionsParser args(argc, argv);
173  args.AddOption(&mesh_file, "-m", "--mesh",
174  "Mesh file to use.");
175  args.AddOption(&ref_levels, "-r", "--refine",
176  "Number of times to refine the mesh uniformly.");
177  args.AddOption(&order, "-o", "--order",
178  "Order (degree) of the finite elements.");
179  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
180  "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
181  " 11 - Forward Euler, 12 - RK2,\n\t"
182  " 13 - RK3 SSP, 14 - RK4.");
183  args.AddOption(&t_final, "-tf", "--t-final",
184  "Final time; start time is 0.");
185  args.AddOption(&dt, "-dt", "--time-step",
186  "Time step.");
187  args.AddOption(&visc, "-v", "--viscosity",
188  "Viscosity coefficient.");
189  args.AddOption(&mu, "-mu", "--shear-modulus",
190  "Shear modulus in the Neo-Hookean hyperelastic model.");
191  args.AddOption(&K, "-K", "--bulk-modulus",
192  "Bulk modulus in the Neo-Hookean hyperelastic model.");
193  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
194  "--no-visualization",
195  "Enable or disable GLVis visualization.");
196  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
197  "Visualize every n-th timestep.");
198  args.Parse();
199  if (!args.Good())
200  {
201  args.PrintUsage(cout);
202  return 1;
203  }
204  args.PrintOptions(cout);
205 
206  // 2. Read the mesh from the given mesh file. We can handle triangular,
207  // quadrilateral, tetrahedral and hexahedral meshes with the same code.
208  Mesh *mesh = new Mesh(mesh_file, 1, 1);
209  int dim = mesh->Dimension();
210 
211  // 3. Define the ODE solver used for time integration. Several implicit
212  // singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
213  // explicit Runge-Kutta methods are available.
214  ODESolver *ode_solver;
215  switch (ode_solver_type)
216  {
217  // Implicit L-stable methods
218  case 1: ode_solver = new BackwardEulerSolver; break;
219  case 2: ode_solver = new SDIRK23Solver(2); break;
220  case 3: ode_solver = new SDIRK33Solver; break;
221  // Explicit methods
222  case 11: ode_solver = new ForwardEulerSolver; break;
223  case 12: ode_solver = new RK2Solver(0.5); break; // midpoint method
224  case 13: ode_solver = new RK3SSPSolver; break;
225  case 14: ode_solver = new RK4Solver; break;
226  case 15: ode_solver = new GeneralizedAlphaSolver(0.5); break;
227  // Implicit A-stable methods (not L-stable)
228  case 22: ode_solver = new ImplicitMidpointSolver; break;
229  case 23: ode_solver = new SDIRK23Solver; break;
230  case 24: ode_solver = new SDIRK34Solver; break;
231  default:
232  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
233  delete mesh;
234  return 3;
235  }
236 
237  // 4. Refine the mesh to increase the resolution. In this example we do
238  // 'ref_levels' of uniform refinement, where 'ref_levels' is a
239  // command-line parameter.
240  for (int lev = 0; lev < ref_levels; lev++)
241  {
242  mesh->UniformRefinement();
243  }
244 
245  // 5. Define the vector finite element spaces representing the mesh
246  // deformation x, the velocity v, and the initial configuration, x_ref.
247  // Define also the elastic energy density, w, which is in a discontinuous
248  // higher-order space. Since x and v are integrated in time as a system,
249  // we group them together in block vector vx, with offsets given by the
250  // fe_offset array.
251  H1_FECollection fe_coll(order, dim);
252  FiniteElementSpace fespace(mesh, &fe_coll, dim);
253 
254  int fe_size = fespace.GetTrueVSize();
255  cout << "Number of velocity/deformation unknowns: " << fe_size << endl;
256  Array<int> fe_offset(3);
257  fe_offset[0] = 0;
258  fe_offset[1] = fe_size;
259  fe_offset[2] = 2*fe_size;
260 
261  BlockVector vx(fe_offset);
262  GridFunction v, x;
263  v.MakeTRef(&fespace, vx.GetBlock(0), 0);
264  x.MakeTRef(&fespace, vx.GetBlock(1), 0);
265 
266  GridFunction x_ref(&fespace);
267  mesh->GetNodes(x_ref);
268 
269  L2_FECollection w_fec(order + 1, dim);
270  FiniteElementSpace w_fespace(mesh, &w_fec);
271  GridFunction w(&w_fespace);
272 
273  // 6. Set the initial conditions for v and x, and the boundary conditions on
274  // a beam-like mesh (see description above).
276  v.ProjectCoefficient(velo);
277  v.SetTrueVector();
279  x.ProjectCoefficient(deform);
280  x.SetTrueVector();
281 
282  Array<int> ess_bdr(fespace.GetMesh()->bdr_attributes.Max());
283  ess_bdr = 0;
284  ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed
285 
286  // 7. Initialize the hyperelastic operator, the GLVis visualization and print
287  // the initial energies.
288  HyperelasticOperator oper(fespace, ess_bdr, visc, mu, K);
289 
290  socketstream vis_v, vis_w;
291  if (visualization)
292  {
293  char vishost[] = "localhost";
294  int visport = 19916;
295  vis_v.open(vishost, visport);
296  vis_v.precision(8);
298  visualize(vis_v, mesh, &x, &v, "Velocity", true);
299  vis_w.open(vishost, visport);
300  if (vis_w)
301  {
302  oper.GetElasticEnergyDensity(x, w);
303  vis_w.precision(8);
304  visualize(vis_w, mesh, &x, &w, "Elastic energy density", true);
305  }
306  }
307 
308  double ee0 = oper.ElasticEnergy(x.GetTrueVector());
309  double ke0 = oper.KineticEnergy(v.GetTrueVector());
310  cout << "initial elastic energy (EE) = " << ee0 << endl;
311  cout << "initial kinetic energy (KE) = " << ke0 << endl;
312  cout << "initial total energy (TE) = " << (ee0 + ke0) << endl;
313 
314  double t = 0.0;
315  oper.SetTime(t);
316  ode_solver->Init(oper);
317 
318  // 8. Perform time-integration (looping over the time iterations, ti, with a
319  // time-step dt).
320  bool last_step = false;
321  for (int ti = 1; !last_step; ti++)
322  {
323  double dt_real = min(dt, t_final - t);
324 
325  ode_solver->Step(vx, t, dt_real);
326 
327  last_step = (t >= t_final - 1e-8*dt);
328 
329  if (last_step || (ti % vis_steps) == 0)
330  {
331  double ee = oper.ElasticEnergy(x.GetTrueVector());
332  double ke = oper.KineticEnergy(v.GetTrueVector());
333 
334  cout << "step " << ti << ", t = " << t << ", EE = " << ee << ", KE = "
335  << ke << ", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
336 
337  if (visualization)
338  {
340  visualize(vis_v, mesh, &x, &v);
341  if (vis_w)
342  {
343  oper.GetElasticEnergyDensity(x, w);
344  visualize(vis_w, mesh, &x, &w);
345  }
346  }
347  }
348  }
349 
350  // 9. Save the displaced mesh, the velocity and elastic energy.
351  {
353  GridFunction *nodes = &x;
354  int owns_nodes = 0;
355  mesh->SwapNodes(nodes, owns_nodes);
356  ofstream mesh_ofs("deformed.mesh");
357  mesh_ofs.precision(8);
358  mesh->Print(mesh_ofs);
359  mesh->SwapNodes(nodes, owns_nodes);
360  ofstream velo_ofs("velocity.sol");
361  velo_ofs.precision(8);
362  v.Save(velo_ofs);
363  ofstream ee_ofs("elastic_energy.sol");
364  ee_ofs.precision(8);
365  oper.GetElasticEnergyDensity(x, w);
366  w.Save(ee_ofs);
367  }
368 
369  // 10. Free the used memory.
370  delete ode_solver;
371  delete mesh;
372 
373  return 0;
374 }
375 
376 
377 void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes,
378  GridFunction *field, const char *field_name, bool init_vis)
379 {
380  if (!out)
381  {
382  return;
383  }
384 
385  GridFunction *nodes = deformed_nodes;
386  int owns_nodes = 0;
387 
388  mesh->SwapNodes(nodes, owns_nodes);
389 
390  out << "solution\n" << *mesh << *field;
391 
392  mesh->SwapNodes(nodes, owns_nodes);
393 
394  if (init_vis)
395  {
396  out << "window_size 800 800\n";
397  out << "window_title '" << field_name << "'\n";
398  if (mesh->SpaceDimension() == 2)
399  {
400  out << "view 0 0\n"; // view from top
401  out << "keys jl\n"; // turn off perspective and light
402  }
403  out << "keys cm\n"; // show colorbar and mesh
404  out << "autoscale value\n"; // update value-range; keep mesh-extents fixed
405  out << "pause\n";
406  }
407  out << flush;
408 }
409 
410 
411 ReducedSystemOperator::ReducedSystemOperator(
413  : Operator(M_->Height()), M(M_), S(S_), H(H_), Jacobian(NULL),
414  dt(0.0), v(NULL), x(NULL), w(height), z(height)
415 { }
416 
417 void ReducedSystemOperator::SetParameters(double dt_, const Vector *v_,
418  const Vector *x_)
419 {
420  dt = dt_; v = v_; x = x_;
421 }
422 
423 void ReducedSystemOperator::Mult(const Vector &k, Vector &y) const
424 {
425  // compute: y = H(x + dt*(v + dt*k)) + M*k + S*(v + dt*k)
426  add(*v, dt, k, w);
427  add(*x, dt, w, z);
428  H->Mult(z, y);
429  M->AddMult(k, y);
430  S->AddMult(w, y);
431 }
432 
433 Operator &ReducedSystemOperator::GetGradient(const Vector &k) const
434 {
435  delete Jacobian;
436  Jacobian = Add(1.0, M->SpMat(), dt, S->SpMat());
437  add(*v, dt, k, w);
438  add(*x, dt, w, z);
439  SparseMatrix *grad_H = dynamic_cast<SparseMatrix *>(&H->GetGradient(z));
440  Jacobian->Add(dt*dt, *grad_H);
441  return *Jacobian;
442 }
443 
444 ReducedSystemOperator::~ReducedSystemOperator()
445 {
446  delete Jacobian;
447 }
448 
449 
450 HyperelasticOperator::HyperelasticOperator(FiniteElementSpace &f,
451  Array<int> &ess_bdr, double visc,
452  double mu, double K)
453  : TimeDependentOperator(2*f.GetTrueVSize(), 0.0), fespace(f),
454  M(&fespace), S(&fespace), H(&fespace),
455  viscosity(visc), z(height/2)
456 {
457  const double rel_tol = 1e-8;
458  const int skip_zero_entries = 0;
459 
460  const double ref_density = 1.0; // density in the reference configuration
461  ConstantCoefficient rho0(ref_density);
462  M.AddDomainIntegrator(new VectorMassIntegrator(rho0));
463  M.Assemble(skip_zero_entries);
464  Array<int> ess_tdof_list;
465  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
466  SparseMatrix tmp;
467  M.FormSystemMatrix(ess_tdof_list, tmp);
468 
469  M_solver.iterative_mode = false;
470  M_solver.SetRelTol(rel_tol);
471  M_solver.SetAbsTol(0.0);
472  M_solver.SetMaxIter(30);
473  M_solver.SetPrintLevel(0);
474  M_solver.SetPreconditioner(M_prec);
475  M_solver.SetOperator(M.SpMat());
476 
477  model = new NeoHookeanModel(mu, K);
478  H.AddDomainIntegrator(new HyperelasticNLFIntegrator(model));
479  H.SetEssentialTrueDofs(ess_tdof_list);
480 
481  ConstantCoefficient visc_coeff(viscosity);
482  S.AddDomainIntegrator(new VectorDiffusionIntegrator(visc_coeff));
483  S.Assemble(skip_zero_entries);
484  S.FormSystemMatrix(ess_tdof_list, tmp);
485 
486  reduced_oper = new ReducedSystemOperator(&M, &S, &H);
487 
488 #ifndef MFEM_USE_SUITESPARSE
489  J_prec = new DSmoother(1);
490  MINRESSolver *J_minres = new MINRESSolver;
491  J_minres->SetRelTol(rel_tol);
492  J_minres->SetAbsTol(0.0);
493  J_minres->SetMaxIter(300);
494  J_minres->SetPrintLevel(-1);
495  J_minres->SetPreconditioner(*J_prec);
496  J_solver = J_minres;
497 #else
498  J_solver = new UMFPackSolver;
499  J_prec = NULL;
500 #endif
501 
502  newton_solver.iterative_mode = false;
503  newton_solver.SetSolver(*J_solver);
504  newton_solver.SetOperator(*reduced_oper);
505  newton_solver.SetPrintLevel(1); // print Newton iterations
506  newton_solver.SetRelTol(rel_tol);
507  newton_solver.SetAbsTol(0.0);
508  newton_solver.SetMaxIter(10);
509 }
510 
511 void HyperelasticOperator::Mult(const Vector &vx, Vector &dvx_dt) const
512 {
513  // Create views to the sub-vectors v, x of vx, and dv_dt, dx_dt of dvx_dt
514  int sc = height/2;
515  Vector v(vx.GetData() + 0, sc);
516  Vector x(vx.GetData() + sc, sc);
517  Vector dv_dt(dvx_dt.GetData() + 0, sc);
518  Vector dx_dt(dvx_dt.GetData() + sc, sc);
519 
520  H.Mult(x, z);
521  if (viscosity != 0.0)
522  {
523  S.AddMult(v, z);
524  }
525  z.Neg(); // z = -z
526  M_solver.Mult(z, dv_dt);
527 
528  dx_dt = v;
529 }
530 
531 void HyperelasticOperator::ImplicitSolve(const double dt,
532  const Vector &vx, Vector &dvx_dt)
533 {
534  int sc = height/2;
535  Vector v(vx.GetData() + 0, sc);
536  Vector x(vx.GetData() + sc, sc);
537  Vector dv_dt(dvx_dt.GetData() + 0, sc);
538  Vector dx_dt(dvx_dt.GetData() + sc, sc);
539 
540  // By eliminating kx from the coupled system:
541  // kv = -M^{-1}*[H(x + dt*kx) + S*(v + dt*kv)]
542  // kx = v + dt*kv
543  // we reduce it to a nonlinear equation for kv, represented by the
544  // reduced_oper. This equation is solved with the newton_solver
545  // object (using J_solver and J_prec internally).
546  reduced_oper->SetParameters(dt, &v, &x);
547  Vector zero; // empty vector is interpreted as zero r.h.s. by NewtonSolver
548  newton_solver.Mult(zero, dv_dt);
549  MFEM_VERIFY(newton_solver.GetConverged(), "Newton solver did not converge.");
550  add(v, dt, dv_dt, dx_dt);
551 }
552 
553 double HyperelasticOperator::ElasticEnergy(const Vector &x) const
554 {
555  return H.GetEnergy(x);
556 }
557 
558 double HyperelasticOperator::KineticEnergy(const Vector &v) const
559 {
560  return 0.5*M.InnerProduct(v, v);
561 }
562 
563 void HyperelasticOperator::GetElasticEnergyDensity(
564  const GridFunction &x, GridFunction &w) const
565 {
566  ElasticEnergyCoefficient w_coeff(*model, x);
567  w.ProjectCoefficient(w_coeff);
568 }
569 
570 HyperelasticOperator::~HyperelasticOperator()
571 {
572  delete J_solver;
573  delete J_prec;
574  delete reduced_oper;
575  delete model;
576 }
577 
578 
580  const IntegrationPoint &ip)
581 {
582  model.SetTransformation(T);
583  x.GetVectorGradient(T, J);
584  // return model.EvalW(J); // in reference configuration
585  return model.EvalW(J)/J.Det(); // in deformed configuration
586 }
587 
588 
589 void InitialDeformation(const Vector &x, Vector &y)
590 {
591  // set the initial configuration to be the same as the reference, stress
592  // free, configuration
593  y = x;
594 }
595 
596 void InitialVelocity(const Vector &x, Vector &v)
597 {
598  const int dim = x.Size();
599  const double s = 0.1/64.;
600 
601  v = 0.0;
602  v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
603  v(0) = -s*x(0)*x(0);
604 }
void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
Definition: ex10.cpp:377
double Eval(ElementTransformation &T, const IntegrationPoint &ip, double t)
Evaluate the coefficient in the element described by T at the point ip at time t. ...
Definition: coefficient.hpp:55
void Add(const int i, const int j, const double a)
Definition: sparsemat.cpp:2244
void InitialDeformation(const Vector &x, Vector &y)
Definition: ex10.cpp:589
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1156
Conjugate gradient method.
Definition: solvers.hpp:111
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
Data type for scaled Jacobi-type smoother of sparse matrix.
A class to handle Vectors in a block fashion.
Definition: blockvector.hpp:30
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:139
Subclass constant coefficient.
Definition: coefficient.hpp:67
Base abstract class for time dependent operators.
Definition: operator.hpp:162
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:5875
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
virtual void Step(Vector &x, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:150
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:283
MINRES method.
Definition: solvers.hpp:221
int main(int argc, char *argv[])
Definition: ex1.cpp:58
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:210
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:159
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:238
void InitialVelocity(const Vector &x, Vector &v)
Definition: ex10.cpp:596
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:344
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2627
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:3030
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:133
int dim
Definition: ex3.cpp:48
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:234
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:67
Data type sparse matrix.
Definition: sparsemat.hpp:40
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:272
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7591
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
Definition: gridfunc.cpp:207
void SetMaxIter(int max_it)
Definition: solvers.hpp:63
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:259
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:350
int Dimension() const
Definition: mesh.hpp:713
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:121
int SpaceDimension() const
Definition: mesh.hpp:714
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:145
void SetAbsTol(double atol)
Definition: solvers.hpp:62
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:179
void SetRelTol(double rtol)
Definition: solvers.hpp:61
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:85
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Definition: optparser.hpp:76
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:132
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Definition: gridfunc.cpp:1046
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:223
Class for integration point with weight.
Definition: intrules.hpp:25
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1582
Abstract class for hyperelastic models.
int open(const char hostname[], int port)
Vector data type.
Definition: vector.hpp:48
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5837
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
Base class for solvers.
Definition: operator.hpp:279
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:64
The classical forward Euler method.
Definition: ode.hpp:99
Abstract operator.
Definition: operator.hpp:21
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:84
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:134
void Neg()
(*this) = -(*this)
Definition: vector.cpp:230
bool Good() const
Definition: optparser.hpp:122