MFEM  v4.5.1
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 // SUNDIALS Modification
3 //
4 // Compile with: make ex10
5 //
6 // Sample runs:
7 // ex10 -m ../../data/beam-quad.mesh -r 2 -o 2 -s 12 -dt 0.15 -vs 10
8 // ex10 -m ../../data/beam-tri.mesh -r 2 -o 2 -s 16 -dt 0.3 -vs 5
9 // ex10 -m ../../data/beam-hex.mesh -r 1 -o 2 -s 12 -dt 0.2 -vs 5
10 // ex10 -m ../../data/beam-tri.mesh -r 2 -o 2 -s 2 -dt 3 -nls kinsol
11 // ex10 -m ../../data/beam-quad.mesh -r 2 -o 2 -s 2 -dt 3 -nls kinsol
12 // ex10 -m ../../data/beam-hex.mesh -r 1 -o 2 -s 2 -dt 3 -nls kinsol
13 // ex10 -m ../../data/beam-quad.mesh -r 2 -o 2 -s 14 -dt 0.15 -vs 10
14 // ex10 -m ../../data/beam-tri.mesh -r 2 -o 2 -s 17 -dt 0.01 -vs 30
15 // ex10 -m ../../data/beam-hex.mesh -r 1 -o 2 -s 14 -dt 0.15 -vs 10
16 // ex10 -m ../../data/beam-quad-amr.mesh -r 2 -o 2 -s 12 -dt 0.15 -vs 10
17 //
18 // Description: This examples solves a time dependent nonlinear elasticity
19 // problem of the form dv/dt = H(x) + S v, dx/dt = v, where H is a
20 // hyperelastic model and S is a viscosity operator of Laplacian
21 // type. The geometry of the domain is assumed to be as follows:
22 //
23 // +---------------------+
24 // boundary --->| |
25 // attribute 1 | |
26 // (fixed) +---------------------+
27 //
28 // The example demonstrates the use of nonlinear operators (the
29 // class HyperelasticOperator defining H(x)), as well as their
30 // implicit time integration using a Newton method for solving an
31 // associated reduced backward-Euler type nonlinear equation
32 // (class ReducedSystemOperator). Each Newton step requires the
33 // inversion of a Jacobian matrix, which is done through a
34 // (preconditioned) inner solver. Note that implementing the
35 // method HyperelasticOperator::ImplicitSolve is the only
36 // requirement for high-order implicit (SDIRK) time integration.
37 //
38 // We recommend viewing examples 2 and 9 before viewing this
39 // example.
40 
41 #include "mfem.hpp"
42 #include <memory>
43 #include <iostream>
44 #include <fstream>
45 #include <string>
46 #include <map>
47 
48 #ifndef MFEM_USE_SUNDIALS
49 #error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES
50 #endif
51 
52 using namespace std;
53 using namespace mfem;
54 
55 class ReducedSystemOperator;
56 
57 /** After spatial discretization, the hyperelastic model can be written as a
58  * system of ODEs:
59  * dv/dt = -M^{-1}*(H(x) + S*v)
60  * dx/dt = v,
61  * where x is the vector representing the deformation, v is the velocity field,
62  * M is the mass matrix, S is the viscosity matrix, and H(x) is the nonlinear
63  * hyperelastic operator.
64  *
65  * Class HyperelasticOperator represents the right-hand side of the above
66  * system of ODEs. */
67 class HyperelasticOperator : public TimeDependentOperator
68 {
69 protected:
70  FiniteElementSpace &fespace;
71 
72  BilinearForm M, S;
73  NonlinearForm H;
74  double viscosity;
75  HyperelasticModel *model;
76 
77  CGSolver M_solver; // Krylov solver for inverting the mass matrix M
78  DSmoother M_prec; // Preconditioner for the mass matrix M
79 
80  /** Nonlinear operator defining the reduced backward Euler equation for the
81  velocity. Used in the implementation of method ImplicitSolve. */
82  ReducedSystemOperator *reduced_oper;
83 
84  /// Newton solver for the reduced backward Euler equation
85  NewtonSolver *newton_solver;
86 
87  /// Solver for the Jacobian solve in the Newton method
88  Solver *J_solver;
89  /// Preconditioner for the Jacobian solve in the Newton method
90  Solver *J_prec;
91 
92  mutable Vector z; // auxiliary vector
93 
94  SparseMatrix *grad_H;
95  SparseMatrix *Jacobian;
96 
97  double saved_gamma; // saved gamma value from implicit setup
98 
99 public:
100  /// Solver type to use in the ImplicitSolve() method, used by SDIRK methods.
101  enum NonlinearSolverType
102  {
103  NEWTON = 0, ///< Use MFEM's plain NewtonSolver
104  KINSOL = 1 ///< Use SUNDIALS' KINSOL (through MFEM's class KINSolver)
105  };
106 
107  HyperelasticOperator(FiniteElementSpace &f, Array<int> &ess_bdr,
108  double visc, double mu, double K,
109  NonlinearSolverType nls_type);
110 
111  /// Compute the right-hand side of the ODE system.
112  virtual void Mult(const Vector &vx, Vector &dvx_dt) const;
113 
114  /** Solve the Backward-Euler equation: k = f(x + dt*k, t), for the unknown k.
115  This is the only requirement for high-order SDIRK implicit integration.*/
116  virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k);
117 
118 
119  /// Custom Jacobian system solver for the SUNDIALS time integrators.
120  /** For the ODE system represented by HyperelasticOperator
121 
122  M dv/dt = -(H(x) + S*v)
123  dx/dt = v,
124 
125  this class facilitates the solution of linear systems of the form
126 
127  (M + γS) yv + γJ yx = M bv, J=(dH/dx)(x)
128  - γ yv + yx = bx
129 
130  for given bv, bx, x, and γ = GetTimeStep(). */
131 
132  /** Linear solve applicable to the SUNDIALS format.
133  Solves (Mass - dt J) y = Mass b, where in our case:
134  Mass = | M 0 | J = | -S -grad_H | y = | v_hat | b = | b_v |
135  | 0 I | | I 0 | | x_hat | | b_x |
136  The result replaces the rhs b.
137  We substitute x_hat = b_x + dt v_hat and solve
138  (M + dt S + dt^2 grad_H) v_hat = M b_v - dt grad_H b_x. */
139 
140  /** Setup the linear system. This method is used by the implicit
141  SUNDIALS solvers. */
142  virtual int SUNImplicitSetup(const Vector &y, const Vector &fy,
143  int jok, int *jcur, double gamma);
144 
145  /** Solve the linear system. This method is used by the implicit
146  SUNDIALS solvers. */
147  virtual int SUNImplicitSolve(const Vector &b, Vector &x, double tol);
148 
149  double ElasticEnergy(const Vector &x) const;
150  double KineticEnergy(const Vector &v) const;
151  void GetElasticEnergyDensity(const GridFunction &x, GridFunction &w) const;
152 
153  virtual ~HyperelasticOperator();
154 };
155 
156 /** Nonlinear operator of the form:
157  k --> (M + dt*S)*k + H(x + dt*v + dt^2*k) + S*v,
158  where M and S are given BilinearForms, H is a given NonlinearForm, v and x
159  are given vectors, and dt is a scalar. */
160 class ReducedSystemOperator : public Operator
161 {
162 private:
163  BilinearForm *M, *S;
164  NonlinearForm *H;
165  mutable SparseMatrix *Jacobian;
166  double dt;
167  const Vector *v, *x;
168  mutable Vector w, z;
169 
170 public:
171  ReducedSystemOperator(BilinearForm *M_, BilinearForm *S_, NonlinearForm *H_);
172 
173  /// Set current dt, v, x values - needed to compute action and Jacobian.
174  void SetParameters(double dt_, const Vector *v_, const Vector *x_);
175 
176  /// Compute y = H(x + dt (v + dt k)) + M k + S (v + dt k).
177  virtual void Mult(const Vector &k, Vector &y) const;
178 
179  /// Compute J = M + dt S + dt^2 grad_H(x + dt (v + dt k)).
180  virtual Operator &GetGradient(const Vector &k) const;
181 
182  virtual ~ReducedSystemOperator();
183 };
184 
185 
186 /** Function representing the elastic energy density for the given hyperelastic
187  model+deformation. Used in HyperelasticOperator::GetElasticEnergyDensity. */
188 class ElasticEnergyCoefficient : public Coefficient
189 {
190 private:
191  HyperelasticModel &model;
192  const GridFunction &x;
193  DenseMatrix J;
194 
195 public:
196  ElasticEnergyCoefficient(HyperelasticModel &m, const GridFunction &x_)
197  : model(m), x(x_) { }
198  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
199  virtual ~ElasticEnergyCoefficient() { }
200 };
201 
202 void InitialDeformation(const Vector &x, Vector &y);
203 
204 void InitialVelocity(const Vector &x, Vector &v);
205 
206 void visualize(ostream &os, Mesh *mesh, GridFunction *deformed_nodes,
207  GridFunction *field, const char *field_name = NULL,
208  bool init_vis = false);
209 
210 
211 int main(int argc, char *argv[])
212 {
213  // 0. Initialize SUNDIALS.
214  Sundials::Init();
215 
216  // 1. Parse command-line options.
217  const char *mesh_file = "../../data/beam-quad.mesh";
218  int ref_levels = 2;
219  int order = 2;
220  int ode_solver_type = 3;
221  double t_final = 300.0;
222  double dt = 3.0;
223  double visc = 1e-2;
224  double mu = 0.25;
225  double K = 5.0;
226  bool visualization = true;
227  const char *nls = "newton";
228  int vis_steps = 1;
229 
230  // Relative and absolute tolerances for CVODE and ARKODE.
231  const double reltol = 1e-1, abstol = 1e-1;
232  // Since this example uses the loose tolerances defined above, it is
233  // necessary to lower the linear solver tolerance for CVODE which is relative
234  // to the above tolerances.
235  const double cvode_eps_lin = 1e-4;
236  // Similarly, the nonlinear tolerance for ARKODE needs to be tightened.
237  const double arkode_eps_nonlin = 1e-6;
238 
239  OptionsParser args(argc, argv);
240  args.AddOption(&mesh_file, "-m", "--mesh",
241  "Mesh file to use.");
242  args.AddOption(&ref_levels, "-r", "--refine",
243  "Number of times to refine the mesh uniformly.");
244  args.AddOption(&order, "-o", "--order",
245  "Order (degree) of the finite elements.");
246  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
247  "ODE solver:\n\t"
248  "1 - Backward Euler,\n\t"
249  "2 - SDIRK2, L-stable\n\t"
250  "3 - SDIRK3, L-stable\n\t"
251  "4 - Implicit Midpoint,\n\t"
252  "5 - SDIRK2, A-stable,\n\t"
253  "6 - SDIRK3, A-stable,\n\t"
254  "7 - Forward Euler,\n\t"
255  "8 - RK2,\n\t"
256  "9 - RK3 SSP,\n\t"
257  "10 - RK4,\n\t"
258  "11 - CVODE implicit BDF, approximate Jacobian,\n\t"
259  "12 - CVODE implicit BDF, specified Jacobian,\n\t"
260  "13 - CVODE implicit ADAMS, approximate Jacobian,\n\t"
261  "14 - CVODE implicit ADAMS, specified Jacobian,\n\t"
262  "15 - ARKODE implicit, approximate Jacobian,\n\t"
263  "16 - ARKODE implicit, specified Jacobian,\n\t"
264  "17 - ARKODE explicit, 4th order.");
265  args.AddOption(&nls, "-nls", "--nonlinear-solver",
266  "Nonlinear systems solver: "
267  "\"newton\" (plain Newton) or \"kinsol\" (KINSOL).");
268  args.AddOption(&t_final, "-tf", "--t-final",
269  "Final time; start time is 0.");
270  args.AddOption(&dt, "-dt", "--time-step",
271  "Time step.");
272  args.AddOption(&visc, "-v", "--viscosity",
273  "Viscosity coefficient.");
274  args.AddOption(&mu, "-mu", "--shear-modulus",
275  "Shear modulus in the Neo-Hookean hyperelastic model.");
276  args.AddOption(&K, "-K", "--bulk-modulus",
277  "Bulk modulus in the Neo-Hookean hyperelastic model.");
278  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
279  "--no-visualization",
280  "Enable or disable GLVis visualization.");
281  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
282  "Visualize every n-th timestep.");
283  args.Parse();
284  if (!args.Good())
285  {
286  args.PrintUsage(cout);
287  return 1;
288  }
289  args.PrintOptions(cout);
290 
291  // check for valid ODE solver option
292  if (ode_solver_type < 1 || ode_solver_type > 17)
293  {
294  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
295  return 1;
296  }
297 
298  // 2. Read the mesh from the given mesh file. We can handle triangular,
299  // quadrilateral, tetrahedral and hexahedral meshes with the same code.
300  Mesh *mesh = new Mesh(mesh_file, 1, 1);
301  int dim = mesh->Dimension();
302 
303  // 3. Setup the nonlinear solver
304  map<string,HyperelasticOperator::NonlinearSolverType> nls_map;
305  nls_map["newton"] = HyperelasticOperator::NEWTON;
306  nls_map["kinsol"] = HyperelasticOperator::KINSOL;
307  if (nls_map.find(nls) == nls_map.end())
308  {
309  cout << "Unknown type of nonlinear solver: " << nls << endl;
310  return 4;
311  }
312 
313  // 4. Refine the mesh to increase the resolution. In this example we do
314  // 'ref_levels' of uniform refinement, where 'ref_levels' is a
315  // command-line parameter.
316  for (int lev = 0; lev < ref_levels; lev++)
317  {
318  mesh->UniformRefinement();
319  }
320 
321  // 5. Define the vector finite element spaces representing the mesh
322  // deformation x, the velocity v, and the initial configuration, x_ref.
323  // Define also the elastic energy density, w, which is in a discontinuous
324  // higher-order space. Since x and v are integrated in time as a system,
325  // we group them together in block vector vx, with offsets given by the
326  // fe_offset array.
327  H1_FECollection fe_coll(order, dim);
328  FiniteElementSpace fespace(mesh, &fe_coll, dim);
329 
330  int fe_size = fespace.GetTrueVSize();
331  cout << "Number of velocity/deformation unknowns: " << fe_size << endl;
332  Array<int> fe_offset(3);
333  fe_offset[0] = 0;
334  fe_offset[1] = fe_size;
335  fe_offset[2] = 2*fe_size;
336 
337  BlockVector vx(fe_offset);
338  GridFunction v, x;
339  v.MakeTRef(&fespace, vx.GetBlock(0), 0);
340  x.MakeTRef(&fespace, vx.GetBlock(1), 0);
341 
342  GridFunction x_ref(&fespace);
343  mesh->GetNodes(x_ref);
344 
345  L2_FECollection w_fec(order + 1, dim);
346  FiniteElementSpace w_fespace(mesh, &w_fec);
347  GridFunction w(&w_fespace);
348 
349  // 6. Set the initial conditions for v and x, and the boundary conditions on
350  // a beam-like mesh (see description above).
352  v.ProjectCoefficient(velo);
353  v.SetTrueVector();
355  x.ProjectCoefficient(deform);
356  x.SetTrueVector();
357 
358  Array<int> ess_bdr(fespace.GetMesh()->bdr_attributes.Max());
359  ess_bdr = 0;
360  ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed
361 
362  // 7. Initialize the hyperelastic operator, the GLVis visualization and print
363  // the initial energies.
364  HyperelasticOperator oper(fespace, ess_bdr, visc, mu, K, nls_map[nls]);
365 
366  socketstream vis_v, vis_w;
367  if (visualization)
368  {
369  char vishost[] = "localhost";
370  int visport = 19916;
371  vis_v.open(vishost, visport);
372  vis_v.precision(8);
374  visualize(vis_v, mesh, &x, &v, "Velocity", true);
375  vis_w.open(vishost, visport);
376  if (vis_w)
377  {
378  oper.GetElasticEnergyDensity(x, w);
379  vis_w.precision(8);
380  visualize(vis_w, mesh, &x, &w, "Elastic energy density", true);
381  }
382  }
383 
384  double ee0 = oper.ElasticEnergy(x.GetTrueVector());
385  double ke0 = oper.KineticEnergy(v.GetTrueVector());
386  cout << "initial elastic energy (EE) = " << ee0 << endl;
387  cout << "initial kinetic energy (KE) = " << ke0 << endl;
388  cout << "initial total energy (TE) = " << (ee0 + ke0) << endl;
389 
390  // 8. Define the ODE solver used for time integration. Several implicit
391  // singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
392  // explicit Runge-Kutta methods are available.
393  double t = 0.0;
394  oper.SetTime(t);
395 
396  ODESolver *ode_solver = NULL;
397  CVODESolver *cvode = NULL;
398  ARKStepSolver *arkode = NULL;
399  switch (ode_solver_type)
400  {
401  // Implicit L-stable methods
402  case 1: ode_solver = new BackwardEulerSolver; break;
403  case 2: ode_solver = new SDIRK23Solver(2); break;
404  case 3: ode_solver = new SDIRK33Solver; break;
405  // Implicit A-stable methods (not L-stable)
406  case 4: ode_solver = new ImplicitMidpointSolver; break;
407  case 5: ode_solver = new SDIRK23Solver; break;
408  case 6: ode_solver = new SDIRK34Solver; break;
409  // Explicit methods
410  case 7: ode_solver = new ForwardEulerSolver; break;
411  case 8: ode_solver = new RK2Solver(0.5); break; // midpoint method
412  case 9: ode_solver = new RK3SSPSolver; break;
413  case 10: ode_solver = new RK4Solver; break;
414  // CVODE BDF
415  case 11:
416  case 12:
417  cvode = new CVODESolver(CV_BDF);
418  cvode->Init(oper);
419  cvode->SetSStolerances(reltol, abstol);
420  CVodeSetEpsLin(cvode->GetMem(), cvode_eps_lin);
421  cvode->SetMaxStep(dt);
422  if (ode_solver_type == 11)
423  {
424  cvode->UseSundialsLinearSolver();
425  }
426  ode_solver = cvode; break;
427  // CVODE Adams
428  case 13:
429  case 14:
430  cvode = new CVODESolver(CV_ADAMS);
431  cvode->Init(oper);
432  cvode->SetSStolerances(reltol, abstol);
433  CVodeSetEpsLin(cvode->GetMem(), cvode_eps_lin);
434  cvode->SetMaxStep(dt);
435  if (ode_solver_type == 13)
436  {
437  cvode->UseSundialsLinearSolver();
438  }
439  ode_solver = cvode; break;
440  // ARKStep Implicit methods
441  case 15:
442  case 16:
443  arkode = new ARKStepSolver(ARKStepSolver::IMPLICIT);
444  arkode->Init(oper);
445  arkode->SetSStolerances(reltol, abstol);
446  ARKStepSetNonlinConvCoef(arkode->GetMem(), arkode_eps_nonlin);
447  arkode->SetMaxStep(dt);
448  if (ode_solver_type == 15)
449  {
450  arkode->UseSundialsLinearSolver();
451  }
452  ode_solver = arkode; break;
453  // ARKStep Explicit methods
454  case 17:
455  arkode = new ARKStepSolver(ARKStepSolver::EXPLICIT);
456  arkode->Init(oper);
457  arkode->SetSStolerances(reltol, abstol);
458  arkode->SetMaxStep(dt);
459  ode_solver = arkode; break;
460  }
461 
462  // Initialize MFEM integrators, SUNDIALS integrators are initialized above
463  if (ode_solver_type < 11) { ode_solver->Init(oper); }
464 
465  // 9. Perform time-integration (looping over the time iterations, ti, with a
466  // time-step dt).
467  bool last_step = false;
468  for (int ti = 1; !last_step; ti++)
469  {
470  double dt_real = min(dt, t_final - t);
471 
472  ode_solver->Step(vx, t, dt_real);
473 
474  last_step = (t >= t_final - 1e-8*dt);
475 
476  if (last_step || (ti % vis_steps) == 0)
477  {
478  double ee = oper.ElasticEnergy(x.GetTrueVector());
479  double ke = oper.KineticEnergy(v.GetTrueVector());
480 
481  cout << "step " << ti << ", t = " << t << ", EE = " << ee << ", KE = "
482  << ke << ", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
483 
484  if (cvode) { cvode->PrintInfo(); }
485  else if (arkode) { arkode->PrintInfo(); }
486 
487  if (visualization)
488  {
490  visualize(vis_v, mesh, &x, &v);
491  if (vis_w)
492  {
493  oper.GetElasticEnergyDensity(x, w);
494  visualize(vis_w, mesh, &x, &w);
495  }
496  }
497  }
498  }
499 
500  // 10. Save the displaced mesh, the velocity and elastic energy.
501  {
503  GridFunction *nodes = &x;
504  int owns_nodes = 0;
505  mesh->SwapNodes(nodes, owns_nodes);
506  ofstream mesh_ofs("deformed.mesh");
507  mesh_ofs.precision(8);
508  mesh->Print(mesh_ofs);
509  mesh->SwapNodes(nodes, owns_nodes);
510  ofstream velo_ofs("velocity.sol");
511  velo_ofs.precision(8);
512  v.Save(velo_ofs);
513  ofstream ee_ofs("elastic_energy.sol");
514  ee_ofs.precision(8);
515  oper.GetElasticEnergyDensity(x, w);
516  w.Save(ee_ofs);
517  }
518 
519  // 11. Free the used memory.
520  delete ode_solver;
521  delete mesh;
522 
523  return 0;
524 }
525 
526 
527 void visualize(ostream &os, Mesh *mesh, GridFunction *deformed_nodes,
528  GridFunction *field, const char *field_name, bool init_vis)
529 {
530  if (!os)
531  {
532  return;
533  }
534 
535  GridFunction *nodes = deformed_nodes;
536  int owns_nodes = 0;
537 
538  mesh->SwapNodes(nodes, owns_nodes);
539 
540  os << "solution\n" << *mesh << *field;
541 
542  mesh->SwapNodes(nodes, owns_nodes);
543 
544  if (init_vis)
545  {
546  os << "window_size 800 800\n";
547  os << "window_title '" << field_name << "'\n";
548  if (mesh->SpaceDimension() == 2)
549  {
550  os << "view 0 0\n"; // view from top
551  os << "keys jl\n"; // turn off perspective and light
552  }
553  os << "keys cm\n"; // show colorbar and mesh
554  os << "autoscale value\n"; // update value-range; keep mesh-extents fixed
555  os << "pause\n";
556  }
557  os << flush;
558 }
559 
560 
561 ReducedSystemOperator::ReducedSystemOperator(
563  : Operator(M_->Height()), M(M_), S(S_), H(H_), Jacobian(NULL),
564  dt(0.0), v(NULL), x(NULL), w(height), z(height)
565 { }
566 
567 void ReducedSystemOperator::SetParameters(double dt_, const Vector *v_,
568  const Vector *x_)
569 {
570  dt = dt_; v = v_; x = x_;
571 }
572 
573 void ReducedSystemOperator::Mult(const Vector &k, Vector &y) const
574 {
575  // compute: y = H(x + dt*(v + dt*k)) + M*k + S*(v + dt*k)
576  add(*v, dt, k, w);
577  add(*x, dt, w, z);
578  H->Mult(z, y);
579  M->AddMult(k, y);
580  S->AddMult(w, y);
581 }
582 
583 Operator &ReducedSystemOperator::GetGradient(const Vector &k) const
584 {
585  delete Jacobian;
586  Jacobian = Add(1.0, M->SpMat(), dt, S->SpMat());
587  add(*v, dt, k, w);
588  add(*x, dt, w, z);
589  SparseMatrix *grad_H = dynamic_cast<SparseMatrix *>(&H->GetGradient(z));
590  Jacobian->Add(dt*dt, *grad_H);
591  return *Jacobian;
592 }
593 
594 ReducedSystemOperator::~ReducedSystemOperator()
595 {
596  delete Jacobian;
597 }
598 
599 
600 HyperelasticOperator::HyperelasticOperator(FiniteElementSpace &f,
601  Array<int> &ess_bdr, double visc,
602  double mu, double K,
603  NonlinearSolverType nls_type)
604  : TimeDependentOperator(2*f.GetTrueVSize(), 0.0), fespace(f),
605  M(&fespace), S(&fespace), H(&fespace),
606  viscosity(visc), z(height/2),
607  grad_H(NULL), Jacobian(NULL)
608 {
609  const double rel_tol = 1e-8;
610  const int skip_zero_entries = 0;
611 
612  const double ref_density = 1.0; // density in the reference configuration
613  ConstantCoefficient rho0(ref_density);
614  M.AddDomainIntegrator(new VectorMassIntegrator(rho0));
615  M.Assemble(skip_zero_entries);
616  Array<int> ess_tdof_list;
617  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
618  SparseMatrix tmp;
619  M.FormSystemMatrix(ess_tdof_list, tmp);
620 
621  M_solver.iterative_mode = false;
622  M_solver.SetRelTol(rel_tol);
623  M_solver.SetAbsTol(0.0);
624  M_solver.SetMaxIter(30);
625  M_solver.SetPrintLevel(0);
626  M_solver.SetPreconditioner(M_prec);
627  M_solver.SetOperator(M.SpMat());
628 
629  model = new NeoHookeanModel(mu, K);
630  H.AddDomainIntegrator(new HyperelasticNLFIntegrator(model));
631  H.SetEssentialTrueDofs(ess_tdof_list);
632 
633  ConstantCoefficient visc_coeff(viscosity);
634  S.AddDomainIntegrator(new VectorDiffusionIntegrator(visc_coeff));
635  S.Assemble(skip_zero_entries);
636  S.FormSystemMatrix(ess_tdof_list, tmp);
637 
638  reduced_oper = new ReducedSystemOperator(&M, &S, &H);
639 
640 #ifndef MFEM_USE_SUITESPARSE
641  J_prec = new DSmoother(1);
642  MINRESSolver *J_minres = new MINRESSolver;
643  J_minres->SetRelTol(rel_tol);
644  J_minres->SetAbsTol(0.0);
645  J_minres->SetMaxIter(300);
646  J_minres->SetPrintLevel(-1);
647  J_minres->SetPreconditioner(*J_prec);
648  J_solver = J_minres;
649 #else
650  J_solver = new UMFPackSolver;
651  J_prec = NULL;
652 #endif
653 
654  if (nls_type == KINSOL)
655  {
656  KINSolver *kinsolver = new KINSolver(KIN_NONE, true);
657  newton_solver = kinsolver;
658  newton_solver->SetOperator(*reduced_oper);
659  newton_solver->SetMaxIter(200);
660  newton_solver->SetRelTol(rel_tol);
661  newton_solver->SetPrintLevel(0);
662  kinsolver->SetMaxSetupCalls(4);
663  }
664  else
665  {
666  newton_solver = new NewtonSolver();
667  newton_solver->SetOperator(*reduced_oper);
668  newton_solver->SetMaxIter(10);
669  newton_solver->SetRelTol(rel_tol);
670  newton_solver->SetPrintLevel(-1);
671  }
672  newton_solver->SetSolver(*J_solver);
673  newton_solver->iterative_mode = false;
674 }
675 
676 void HyperelasticOperator::Mult(const Vector &vx, Vector &dvx_dt) const
677 {
678  // Create views to the sub-vectors v, x of vx, and dv_dt, dx_dt of dvx_dt
679  int sc = height/2;
680  Vector v(vx.GetData() + 0, sc);
681  Vector x(vx.GetData() + sc, sc);
682  Vector dv_dt(dvx_dt.GetData() + 0, sc);
683  Vector dx_dt(dvx_dt.GetData() + sc, sc);
684 
685  H.Mult(x, z);
686  if (viscosity != 0.0)
687  {
688  S.AddMult(v, z);
689  }
690  z.Neg(); // z = -z
691  M_solver.Mult(z, dv_dt);
692 
693  dx_dt = v;
694 }
695 
696 void HyperelasticOperator::ImplicitSolve(const double dt,
697  const Vector &vx, Vector &dvx_dt)
698 {
699  int sc = height/2;
700  Vector v(vx.GetData() + 0, sc);
701  Vector x(vx.GetData() + sc, sc);
702  Vector dv_dt(dvx_dt.GetData() + 0, sc);
703  Vector dx_dt(dvx_dt.GetData() + sc, sc);
704 
705  // By eliminating kx from the coupled system:
706  // kv = -M^{-1}*[H(x + dt*kx) + S*(v + dt*kv)]
707  // kx = v + dt*kv
708  // we reduce it to a nonlinear equation for kv, represented by the
709  // reduced_oper. This equation is solved with the newton_solver
710  // object (using J_solver and J_prec internally).
711  reduced_oper->SetParameters(dt, &v, &x);
712  Vector zero; // empty vector is interpreted as zero r.h.s. by NewtonSolver
713  newton_solver->Mult(zero, dv_dt);
714  MFEM_VERIFY(newton_solver->GetConverged(),
715  "Nonlinear solver did not converge.");
716 #ifdef MFEM_DEBUG
717  cout << " num nonlin sol iters = " << newton_solver->GetNumIterations()
718  << ", final norm = " << newton_solver->GetFinalNorm() << '\n';
719 #endif
720  add(v, dt, dv_dt, dx_dt);
721 }
722 
723 int HyperelasticOperator::SUNImplicitSetup(const Vector &y,
724  const Vector &fy, int jok, int *jcur,
725  double gamma)
726 {
727  int sc = y.Size() / 2;
728  const Vector x(y.GetData() + sc, sc);
729 
730  // J = M + dt*(S + dt*grad(H))
731  if (Jacobian) { delete Jacobian; }
732  Jacobian = Add(1.0, M.SpMat(), gamma, S.SpMat());
733  grad_H = dynamic_cast<SparseMatrix *>(&H.GetGradient(x));
734  Jacobian->Add(gamma * gamma, *grad_H);
735 
736  // Set Jacobian solve operator
737  J_solver->SetOperator(*Jacobian);
738 
739  // Indicate that the Jacobian was updated
740  *jcur = 1;
741 
742  // Save gamma for use in solve
743  saved_gamma = gamma;
744 
745  // Return success
746  return 0;
747 }
748 
749 int HyperelasticOperator::SUNImplicitSolve(const Vector &b, Vector &x,
750  double tol)
751 {
752  int sc = b.Size() / 2;
753  Vector b_v(b.GetData() + 0, sc);
754  Vector b_x(b.GetData() + sc, sc);
755  Vector x_v(x.GetData() + 0, sc);
756  Vector x_x(x.GetData() + sc, sc);
757  Vector rhs(sc);
758 
759  // rhs = M b_v - dt*grad(H) b_x
760  grad_H->Mult(b_x, rhs);
761  rhs *= -saved_gamma;
762  M.AddMult(b_v, rhs);
763 
764  J_solver->iterative_mode = false;
765  J_solver->Mult(rhs, x_v);
766 
767  add(b_x, saved_gamma, x_v, x_x);
768 
769  return 0;
770 }
771 
772 double HyperelasticOperator::ElasticEnergy(const Vector &x) const
773 {
774  return H.GetEnergy(x);
775 }
776 
777 double HyperelasticOperator::KineticEnergy(const Vector &v) const
778 {
779  return 0.5*M.InnerProduct(v, v);
780 }
781 
782 void HyperelasticOperator::GetElasticEnergyDensity(
783  const GridFunction &x, GridFunction &w) const
784 {
785  ElasticEnergyCoefficient w_coeff(*model, x);
786  w.ProjectCoefficient(w_coeff);
787 }
788 
789 HyperelasticOperator::~HyperelasticOperator()
790 {
791  delete Jacobian;
792  delete newton_solver;
793  delete J_solver;
794  delete J_prec;
795  delete reduced_oper;
796  delete model;
797 }
798 
799 
801  const IntegrationPoint &ip)
802 {
803  model.SetTransformation(T);
804  x.GetVectorGradient(T, J);
805  // return model.EvalW(J); // in reference configuration
806  return model.EvalW(J)/J.Det(); // in deformed configuration
807 }
808 
809 
810 void InitialDeformation(const Vector &x, Vector &y)
811 {
812  // set the initial configuration to be the same as the reference, stress
813  // free, configuration
814  y = x;
815 }
816 
817 void InitialVelocity(const Vector &x, Vector &v)
818 {
819  const int dim = x.Size();
820  const double s = 0.1/64.;
821 
822  v = 0.0;
823  v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
824  v(0) = -s*x(0)*x(0);
825 }
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
Definition: sundials.cpp:696
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:68
void InitialDeformation(const Vector &x, Vector &y)
Definition: ex10.cpp:594
Conjugate gradient method.
Definition: solvers.hpp:465
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
Data type for scaled Jacobi-type smoother of sparse matrix.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:862
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:152
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
Base abstract class for first order time dependent operators.
Definition: operator.hpp:289
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:7957
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
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]...
void * GetMem() const
Access the SUNDIALS memory structure.
Definition: sundials.hpp:362
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:1706
MINRES method.
Definition: solvers.hpp:575
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:391
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:209
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:313
void SetMaxSetupCalls(int max_calls)
Set maximum number of nonlinear iterations without a Jacobian update.
Definition: sundials.cpp:2129
void InitialVelocity(const Vector &x, Vector &v)
Definition: ex10.cpp:601
Interface to ARKode&#39;s ARKStep module – additive Runge-Kutta methods.
Definition: sundials.hpp:661
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:1025
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3676
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2282
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:146
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:588
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
Data type sparse matrix.
Definition: sparsemat.hpp:50
constexpr char vishost[]
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:441
Interface to the CVODE library – linear multi-step methods.
Definition: sundials.hpp:374
double b
Definition: lissajous.cpp:42
Interface to the KINSOL library – nonlinear solver methods.
Definition: sundials.hpp:835
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9816
constexpr int visport
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
Definition: gridfunc.cpp:222
void SetMaxIter(int max_it)
Definition: solvers.hpp:200
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:880
void visualize(ostream &os, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
Definition: ex10.cpp:381
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:613
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:590
int Dimension() const
Definition: mesh.hpp:1006
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
void Add(const int i, const int j, const double val)
Definition: sparsemat.cpp:2783
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:132
A general vector function coefficient.
int SpaceDimension() const
Definition: mesh.hpp:1007
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
void SetAbsTol(double atol)
Definition: solvers.hpp:199
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
void SetRelTol(double rtol)
Definition: solvers.hpp:198
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
virtual void SetOperator(const Operator &op)
Set the nonlinear Operator of the system and initialize KINSOL.
Definition: sundials.cpp:1943
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
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)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:149
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Definition: gridfunc.cpp:1742
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:404
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
Class for integration point with weight.
Definition: intrules.hpp:25
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:1671
void PrintInfo() const
Print various ARKStep statistics.
Definition: sundials.cpp:1743
double mu
Definition: ex25.cpp:139
int dim
Definition: ex24.cpp:53
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2399
Abstract class for hyperelastic models.
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
RefCoord t[3]
Vector data type.
Definition: vector.hpp:60
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7908
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:220
void PrintInfo() const
Print various CVODE statistics.
Definition: sundials.cpp:906
RefCoord s[3]
Base class for solvers.
Definition: operator.hpp:655
The classical forward Euler method.
Definition: ode.hpp:116
Abstract operator.
Definition: operator.hpp:24
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
Definition: sundials.cpp:1461
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:1700
void UseSundialsLinearSolver()
Attach a SUNDIALS GMRES linear solver to ARKode.
Definition: sundials.cpp:1628
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
int main()
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
Definition: sundials.cpp:842
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:90
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:288
double f(const Vector &p)
void Neg()
(*this) = -(*this)
Definition: vector.cpp:305
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150