MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex10p.cpp
Go to the documentation of this file.
1 // MFEM Example 10 - Parallel Version
2 // SUNDIALS Modification
3 //
4 // Compile with: make ex10p
5 //
6 // Sample runs:
7 // mpirun -np 4 ex10p -m ../../data/beam-quad.mesh -rp 1 -o 2 -s 12 -dt 0.15 -vs 10
8 // mpirun -np 4 ex10p -m ../../data/beam-tri.mesh -rp 1 -o 2 -s 16 -dt 0.25 -vs 10
9 // mpirun -np 4 ex10p -m ../../data/beam-hex.mesh -rp 0 -o 2 -s 12 -dt 0.15 -vs 10
10 // mpirun -np 4 ex10p -m ../../data/beam-tri.mesh -rp 1 -o 2 -s 2 -dt 3 -nls kinsol
11 // mpirun -np 4 ex10p -m ../../data/beam-quad.mesh -rp 1 -o 2 -s 2 -dt 3 -nls kinsol
12 // mpirun -np 4 ex10p -m ../../data/beam-hex.mesh -rs 1 -o 2 -s 2 -dt 3 -nls kinsol
13 // mpirun -np 4 ex10p -m ../../data/beam-quad.mesh -rp 1 -o 2 -s 14 -dt 0.15 -vs 10
14 // mpirun -np 4 ex10p -m ../../data/beam-tri.mesh -rp 1 -o 2 -s 17 -dt 5e-3 -vs 60
15 // mpirun -np 4 ex10p -m ../../data/beam-hex.mesh -rp 0 -o 2 -s 14 -dt 0.15 -vs 10
16 // mpirun -np 4 ex10p -m ../../data/beam-quad-amr.mesh -rp 1 -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  ParFiniteElementSpace &fespace;
71  Array<int> ess_tdof_list;
72 
73  ParBilinearForm M, S;
75  double viscosity;
76  HyperelasticModel *model;
77 
78  HypreParMatrix *Mmat; // Mass matrix from ParallelAssemble()
79  CGSolver M_solver; // Krylov solver for inverting the mass matrix M
80  HypreSmoother M_prec; // Preconditioner for the mass matrix M
81 
82  /** Nonlinear operator defining the reduced backward Euler equation for the
83  velocity. Used in the implementation of method ImplicitSolve. */
84  ReducedSystemOperator *reduced_oper;
85 
86  /// Newton solver for the reduced backward Euler equation
87  NewtonSolver *newton_solver;
88 
89  /// Solver for the Jacobian solve in the Newton method
90  Solver *J_solver;
91  /// Preconditioner for the Jacobian solve in the Newton method
92  Solver *J_prec;
93 
94  mutable Vector z; // auxiliary vector
95 
96  const SparseMatrix *local_grad_H;
97  HypreParMatrix *Jacobian;
98 
99  double saved_gamma; // saved gamma value from implicit setup
100 
101 public:
102  /// Solver type to use in the ImplicitSolve() method, used by SDIRK methods.
103  enum NonlinearSolverType
104  {
105  NEWTON = 0, ///< Use MFEM's plain NewtonSolver
106  KINSOL = 1 ///< Use SUNDIALS' KINSOL (through MFEM's class KINSolver)
107  };
108 
109  HyperelasticOperator(ParFiniteElementSpace &f, Array<int> &ess_bdr,
110  double visc, double mu, double K,
111  NonlinearSolverType nls_type);
112 
113  /// Compute the right-hand side of the ODE system.
114  virtual void Mult(const Vector &vx, Vector &dvx_dt) const;
115 
116  /** Solve the Backward-Euler equation: k = f(x + dt*k, t), for the unknown k.
117  This is the only requirement for high-order SDIRK implicit integration.*/
118  virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k);
119 
120 
121  /// Custom Jacobian system solver for the SUNDIALS time integrators.
122  /** For the ODE system represented by HyperelasticOperator
123 
124  M dv/dt = -(H(x) + S*v)
125  dx/dt = v,
126 
127  this class facilitates the solution of linear systems of the form
128 
129  (M + γS) yv + γJ yx = M bv, J=(dH/dx)(x)
130  - γ yv + yx = bx
131 
132  for given bv, bx, x, and γ = GetTimeStep(). */
133 
134  /** Linear solve applicable to the SUNDIALS format.
135  Solves (Mass - dt J) y = Mass b, where in our case:
136  Mass = | M 0 | J = | -S -grad_H | y = | v_hat | b = | b_v |
137  | 0 I | | I 0 | | x_hat | | b_x |
138  The result replaces the rhs b.
139  We substitute x_hat = b_x + dt v_hat and solve
140  (M + dt S + dt^2 grad_H) v_hat = M b_v - dt grad_H b_x. */
141 
142  /** Setup the linear system. This method is used by the implicit
143  SUNDIALS solvers. */
144  virtual int SUNImplicitSetup(const Vector &y, const Vector &fy,
145  int jok, int *jcur, double gamma);
146 
147  /** Solve the linear system. This method is used by the implicit
148  SUNDIALS solvers. */
149  virtual int SUNImplicitSolve(const Vector &b, Vector &x, double tol);
150 
151  double ElasticEnergy(const ParGridFunction &x) const;
152  double KineticEnergy(const ParGridFunction &v) const;
153  void GetElasticEnergyDensity(const ParGridFunction &x,
154  ParGridFunction &w) const;
155 
156  virtual ~HyperelasticOperator();
157 };
158 
159 /** Nonlinear operator of the form:
160  k --> (M + dt*S)*k + H(x + dt*v + dt^2*k) + S*v,
161  where M and S are given BilinearForms, H is a given NonlinearForm, v and x
162  are given vectors, and dt is a scalar. */
163 class ReducedSystemOperator : public Operator
164 {
165 private:
166  ParBilinearForm *M, *S;
167  ParNonlinearForm *H;
168  mutable HypreParMatrix *Jacobian;
169  double dt;
170  const Vector *v, *x;
171  mutable Vector w, z;
172  const Array<int> &ess_tdof_list;
173 
174 public:
175  ReducedSystemOperator(ParBilinearForm *M_, ParBilinearForm *S_,
176  ParNonlinearForm *H_, const Array<int> &ess_tdof_list);
177 
178  /// Set current dt, v, x values - needed to compute action and Jacobian.
179  void SetParameters(double dt_, const Vector *v_, const Vector *x_);
180 
181  /// Compute y = H(x + dt (v + dt k)) + M k + S (v + dt k).
182  virtual void Mult(const Vector &k, Vector &y) const;
183 
184  /// Compute J = M + dt S + dt^2 grad_H(x + dt (v + dt k)).
185  virtual Operator &GetGradient(const Vector &k) const;
186 
187  virtual ~ReducedSystemOperator();
188 };
189 
190 
191 /** Function representing the elastic energy density for the given hyperelastic
192  model+deformation. Used in HyperelasticOperator::GetElasticEnergyDensity. */
193 class ElasticEnergyCoefficient : public Coefficient
194 {
195 private:
196  HyperelasticModel &model;
197  const ParGridFunction &x;
198  DenseMatrix J;
199 
200 public:
201  ElasticEnergyCoefficient(HyperelasticModel &m, const ParGridFunction &x_)
202  : model(m), x(x_) { }
203  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
204  virtual ~ElasticEnergyCoefficient() { }
205 };
206 
207 void InitialDeformation(const Vector &x, Vector &y);
208 
209 void InitialVelocity(const Vector &x, Vector &v);
210 
211 void visualize(ostream &out, ParMesh *mesh, ParGridFunction *deformed_nodes,
212  ParGridFunction *field, const char *field_name = NULL,
213  bool init_vis = false);
214 
215 
216 int main(int argc, char *argv[])
217 {
218  // 1. Initialize MPI.
219  int num_procs, myid;
220  MPI_Init(&argc, &argv);
221  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
222  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
223 
224  // 2. Parse command-line options.
225  const char *mesh_file = "../../data/beam-quad.mesh";
226  int ser_ref_levels = 2;
227  int par_ref_levels = 0;
228  int order = 2;
229  int ode_solver_type = 3;
230  double t_final = 300.0;
231  double dt = 3.0;
232  double visc = 1e-2;
233  double mu = 0.25;
234  double K = 5.0;
235  bool visualization = true;
236  const char *nls = "newton";
237  int vis_steps = 1;
238 
239  // Relative and absolute tolerances for CVODE and ARKODE.
240  const double reltol = 1e-1, abstol = 1e-1;
241  // Since this example uses the loose tolerances defined above, it is
242  // necessary to lower the linear solver tolerance for CVODE which is relative
243  // to the above tolerances.
244  const double cvode_eps_lin = 1e-4;
245  // Similarly, the nonlinear tolerance for ARKODE needs to be tightened.
246  const double arkode_eps_nonlin = 1e-6;
247 
248  OptionsParser args(argc, argv);
249  args.AddOption(&mesh_file, "-m", "--mesh",
250  "Mesh file to use.");
251  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
252  "Number of times to refine the mesh uniformly in serial.");
253  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
254  "Number of times to refine the mesh uniformly in parallel.");
255  args.AddOption(&order, "-o", "--order",
256  "Order (degree) of the finite elements.");
257  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
258  "ODE solver:\n\t"
259  "1 - Backward Euler,\n\t"
260  "2 - SDIRK2, L-stable\n\t"
261  "3 - SDIRK3, L-stable\n\t"
262  "4 - Implicit Midpoint,\n\t"
263  "5 - SDIRK2, A-stable,\n\t"
264  "6 - SDIRK3, A-stable,\n\t"
265  "7 - Forward Euler,\n\t"
266  "8 - RK2,\n\t"
267  "9 - RK3 SSP,\n\t"
268  "10 - RK4,\n\t"
269  "11 - CVODE implicit BDF, approximate Jacobian,\n\t"
270  "12 - CVODE implicit BDF, specified Jacobian,\n\t"
271  "13 - CVODE implicit ADAMS, approximate Jacobian,\n\t"
272  "14 - CVODE implicit ADAMS, specified Jacobian,\n\t"
273  "15 - ARKODE implicit, approximate Jacobian,\n\t"
274  "16 - ARKODE implicit, specified Jacobian,\n\t"
275  "17 - ARKODE explicit, 4th order.");
276  args.AddOption(&nls, "-nls", "--nonlinear-solver",
277  "Nonlinear systems solver: "
278  "\"newton\" (plain Newton) or \"kinsol\" (KINSOL).");
279  args.AddOption(&t_final, "-tf", "--t-final",
280  "Final time; start time is 0.");
281  args.AddOption(&dt, "-dt", "--time-step",
282  "Time step.");
283  args.AddOption(&visc, "-v", "--viscosity",
284  "Viscosity coefficient.");
285  args.AddOption(&mu, "-mu", "--shear-modulus",
286  "Shear modulus in the Neo-Hookean hyperelastic model.");
287  args.AddOption(&K, "-K", "--bulk-modulus",
288  "Bulk modulus in the Neo-Hookean hyperelastic model.");
289  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
290  "--no-visualization",
291  "Enable or disable GLVis visualization.");
292  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
293  "Visualize every n-th timestep.");
294  args.Parse();
295  if (!args.Good())
296  {
297  if (myid == 0)
298  {
299  args.PrintUsage(cout);
300  }
301  MPI_Finalize();
302  return 1;
303  }
304  if (myid == 0)
305  {
306  args.PrintOptions(cout);
307  }
308 
309  // check for valid ODE solver option
310  if (ode_solver_type < 1 || ode_solver_type > 17)
311  {
312  if (myid == 0)
313  {
314  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
315  }
316  MPI_Finalize();
317  return 1;
318  }
319 
320  // 3. Read the serial mesh from the given mesh file on all processors. We can
321  // handle triangular, quadrilateral, tetrahedral and hexahedral meshes
322  // with the same code.
323  Mesh *mesh = new Mesh(mesh_file, 1, 1);
324  int dim = mesh->Dimension();
325 
326  // 4. Nonlinear solver
327  map<string,HyperelasticOperator::NonlinearSolverType> nls_map;
328  nls_map["newton"] = HyperelasticOperator::NEWTON;
329  nls_map["kinsol"] = HyperelasticOperator::KINSOL;
330  if (nls_map.find(nls) == nls_map.end())
331  {
332  if (myid == 0)
333  {
334  cout << "Unknown type of nonlinear solver: " << nls << endl;
335  }
336  delete mesh;
337  MPI_Finalize();
338  return 4;
339  }
340 
341  // 5. Refine the mesh in serial to increase the resolution. In this example
342  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
343  // a command-line parameter.
344  for (int lev = 0; lev < ser_ref_levels; lev++)
345  {
346  mesh->UniformRefinement();
347  }
348 
349  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
350  // this mesh further in parallel to increase the resolution. Once the
351  // parallel mesh is defined, the serial mesh can be deleted.
352  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
353  delete mesh;
354  for (int lev = 0; lev < par_ref_levels; lev++)
355  {
356  pmesh->UniformRefinement();
357  }
358 
359  // 7. Define the parallel vector finite element spaces representing the mesh
360  // deformation x_gf, the velocity v_gf, and the initial configuration,
361  // x_ref. Define also the elastic energy density, w_gf, which is in a
362  // discontinuous higher-order space. Since x and v are integrated in time
363  // as a system, we group them together in block vector vx, on the unique
364  // parallel degrees of freedom, with offsets given by array true_offset.
365  H1_FECollection fe_coll(order, dim);
366  ParFiniteElementSpace fespace(pmesh, &fe_coll, dim);
367 
368  HYPRE_Int glob_size = fespace.GlobalTrueVSize();
369  if (myid == 0)
370  {
371  cout << "Number of velocity/deformation unknowns: " << glob_size << endl;
372  }
373  int true_size = fespace.TrueVSize();
374  Array<int> true_offset(3);
375  true_offset[0] = 0;
376  true_offset[1] = true_size;
377  true_offset[2] = 2*true_size;
378 
379  BlockVector vx(true_offset);
380  ParGridFunction v_gf, x_gf;
381  v_gf.MakeTRef(&fespace, vx, true_offset[0]);
382  x_gf.MakeTRef(&fespace, vx, true_offset[1]);
383 
384  ParGridFunction x_ref(&fespace);
385  pmesh->GetNodes(x_ref);
386 
387  L2_FECollection w_fec(order + 1, dim);
388  ParFiniteElementSpace w_fespace(pmesh, &w_fec);
389  ParGridFunction w_gf(&w_fespace);
390 
391  // 8. Set the initial conditions for v_gf, x_gf and vx, and define the
392  // boundary conditions on a beam-like mesh (see description above).
394  v_gf.ProjectCoefficient(velo);
395  v_gf.SetTrueVector();
397  x_gf.ProjectCoefficient(deform);
398  x_gf.SetTrueVector();
399 
400  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
401 
402  Array<int> ess_bdr(fespace.GetMesh()->bdr_attributes.Max());
403  ess_bdr = 0;
404  ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed
405 
406  // 9. Initialize the hyperelastic operator, the GLVis visualization and print
407  // the initial energies.
408  HyperelasticOperator oper(fespace, ess_bdr, visc, mu, K, nls_map[nls]);
409 
410  socketstream vis_v, vis_w;
411  if (visualization)
412  {
413  char vishost[] = "localhost";
414  int visport = 19916;
415  vis_v.open(vishost, visport);
416  vis_v.precision(8);
417  visualize(vis_v, pmesh, &x_gf, &v_gf, "Velocity", true);
418  // Make sure all ranks have sent their 'v' solution before initiating
419  // another set of GLVis connections (one from each rank):
420  MPI_Barrier(pmesh->GetComm());
421  vis_w.open(vishost, visport);
422  if (vis_w)
423  {
424  oper.GetElasticEnergyDensity(x_gf, w_gf);
425  vis_w.precision(8);
426  visualize(vis_w, pmesh, &x_gf, &w_gf, "Elastic energy density", true);
427  }
428  }
429 
430  double ee0 = oper.ElasticEnergy(x_gf);
431  double ke0 = oper.KineticEnergy(v_gf);
432  if (myid == 0)
433  {
434  cout << "initial elastic energy (EE) = " << ee0 << endl;
435  cout << "initial kinetic energy (KE) = " << ke0 << endl;
436  cout << "initial total energy (TE) = " << (ee0 + ke0) << endl;
437  }
438 
439  // 10. Define the ODE solver used for time integration. Several implicit
440  // singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
441  // explicit Runge-Kutta methods are available.
442  double t = 0.0;
443  oper.SetTime(t);
444 
445  ODESolver *ode_solver = NULL;
446  CVODESolver *cvode = NULL;
447  ARKStepSolver *arkode = NULL;
448  switch (ode_solver_type)
449  {
450  // Implicit L-stable methods
451  case 1: ode_solver = new BackwardEulerSolver; break;
452  case 2: ode_solver = new SDIRK23Solver(2); break;
453  case 3: ode_solver = new SDIRK33Solver; break;
454  // Implicit A-stable methods (not L-stable)
455  case 4: ode_solver = new ImplicitMidpointSolver; break;
456  case 5: ode_solver = new SDIRK23Solver; break;
457  case 6: ode_solver = new SDIRK34Solver; break;
458  // Explicit methods
459  case 7: ode_solver = new ForwardEulerSolver; break;
460  case 8: ode_solver = new RK2Solver(0.5); break; // midpoint method
461  case 9: ode_solver = new RK3SSPSolver; break;
462  case 10: ode_solver = new RK4Solver; break;
463  // CVODE BDF
464  case 11:
465  case 12:
466  cvode = new CVODESolver(MPI_COMM_WORLD, CV_BDF);
467  cvode->Init(oper);
468  cvode->SetSStolerances(reltol, abstol);
469  CVodeSetEpsLin(cvode->GetMem(), cvode_eps_lin);
470  cvode->SetMaxStep(dt);
471  if (ode_solver_type == 11)
472  {
473  cvode->UseSundialsLinearSolver();
474  }
475  ode_solver = cvode; break;
476  // CVODE Adams
477  case 13:
478  case 14:
479  cvode = new CVODESolver(MPI_COMM_WORLD, CV_ADAMS);
480  cvode->Init(oper);
481  cvode->SetSStolerances(reltol, abstol);
482  CVodeSetEpsLin(cvode->GetMem(), cvode_eps_lin);
483  cvode->SetMaxStep(dt);
484  if (ode_solver_type == 13)
485  {
486  cvode->UseSundialsLinearSolver();
487  }
488  ode_solver = cvode; break;
489  // ARKStep Implicit methods
490  case 15:
491  case 16:
492  arkode = new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::IMPLICIT);
493  arkode->Init(oper);
494  arkode->SetSStolerances(reltol, abstol);
495  ARKStepSetNonlinConvCoef(arkode->GetMem(), arkode_eps_nonlin);
496  arkode->SetMaxStep(dt);
497  if (ode_solver_type == 15)
498  {
499  arkode->UseSundialsLinearSolver();
500  }
501  ode_solver = arkode; break;
502  // ARKStep Explicit methods
503  case 17:
504  arkode = new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::EXPLICIT);
505  arkode->Init(oper);
506  arkode->SetSStolerances(reltol, abstol);
507  arkode->SetMaxStep(dt);
508  ode_solver = arkode; break;
509  }
510 
511  // Initialize MFEM integrators, SUNDIALS integrators are initialized above
512  if (ode_solver_type < 11) { ode_solver->Init(oper); }
513 
514  // 11. Perform time-integration
515  // (looping over the time iterations, ti, with a time-step dt).
516  bool last_step = false;
517  for (int ti = 1; !last_step; ti++)
518  {
519  double dt_real = min(dt, t_final - t);
520 
521  ode_solver->Step(vx, t, dt_real);
522 
523  last_step = (t >= t_final - 1e-8*dt);
524 
525  if (last_step || (ti % vis_steps) == 0)
526  {
527  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
528 
529  double ee = oper.ElasticEnergy(x_gf);
530  double ke = oper.KineticEnergy(v_gf);
531 
532  if (myid == 0)
533  {
534  cout << "step " << ti << ", t = " << t << ", EE = " << ee
535  << ", KE = " << ke << ", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
536 
537  if (cvode) { cvode->PrintInfo(); }
538  else if (arkode) { arkode->PrintInfo(); }
539  }
540 
541  if (visualization)
542  {
543  visualize(vis_v, pmesh, &x_gf, &v_gf);
544  if (vis_w)
545  {
546  oper.GetElasticEnergyDensity(x_gf, w_gf);
547  visualize(vis_w, pmesh, &x_gf, &w_gf);
548  }
549  }
550  }
551  }
552 
553  // 12. Save the displaced mesh, the velocity and elastic energy.
554  {
555  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
556  GridFunction *nodes = &x_gf;
557  int owns_nodes = 0;
558  pmesh->SwapNodes(nodes, owns_nodes);
559 
560  ostringstream mesh_name, velo_name, ee_name;
561  mesh_name << "deformed." << setfill('0') << setw(6) << myid;
562  velo_name << "velocity." << setfill('0') << setw(6) << myid;
563  ee_name << "elastic_energy." << setfill('0') << setw(6) << myid;
564 
565  ofstream mesh_ofs(mesh_name.str().c_str());
566  mesh_ofs.precision(8);
567  pmesh->Print(mesh_ofs);
568  pmesh->SwapNodes(nodes, owns_nodes);
569  ofstream velo_ofs(velo_name.str().c_str());
570  velo_ofs.precision(8);
571  v_gf.Save(velo_ofs);
572  ofstream ee_ofs(ee_name.str().c_str());
573  ee_ofs.precision(8);
574  oper.GetElasticEnergyDensity(x_gf, w_gf);
575  w_gf.Save(ee_ofs);
576  }
577 
578  // 13. Free the used memory.
579  delete ode_solver;
580  delete pmesh;
581 
582  MPI_Finalize();
583 
584  return 0;
585 }
586 
587 void visualize(ostream &out, ParMesh *mesh, ParGridFunction *deformed_nodes,
588  ParGridFunction *field, const char *field_name, bool init_vis)
589 {
590  if (!out)
591  {
592  return;
593  }
594 
595  GridFunction *nodes = deformed_nodes;
596  int owns_nodes = 0;
597 
598  mesh->SwapNodes(nodes, owns_nodes);
599 
600  out << "parallel " << mesh->GetNRanks() << " " << mesh->GetMyRank() << "\n";
601  out << "solution\n" << *mesh << *field;
602 
603  mesh->SwapNodes(nodes, owns_nodes);
604 
605  if (init_vis)
606  {
607  out << "window_size 800 800\n";
608  out << "window_title '" << field_name << "'\n";
609  if (mesh->SpaceDimension() == 2)
610  {
611  out << "view 0 0\n"; // view from top
612  out << "keys jl\n"; // turn off perspective and light
613  }
614  out << "keys cm\n"; // show colorbar and mesh
615  out << "autoscale value\n"; // update value-range; keep mesh-extents fixed
616  out << "pause\n";
617  }
618  out << flush;
619 }
620 
621 
622 ReducedSystemOperator::ReducedSystemOperator(
624  const Array<int> &ess_tdof_list_)
625  : Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
626  Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height),
627  ess_tdof_list(ess_tdof_list_)
628 { }
629 
630 void ReducedSystemOperator::SetParameters(double dt_, const Vector *v_,
631  const Vector *x_)
632 {
633  dt = dt_; v = v_; x = x_;
634 }
635 
636 void ReducedSystemOperator::Mult(const Vector &k, Vector &y) const
637 {
638  // compute: y = H(x + dt*(v + dt*k)) + M*k + S*(v + dt*k)
639  add(*v, dt, k, w);
640  add(*x, dt, w, z);
641  H->Mult(z, y);
642  M->TrueAddMult(k, y);
643  S->TrueAddMult(w, y);
644  y.SetSubVector(ess_tdof_list, 0.0);
645 }
646 
647 Operator &ReducedSystemOperator::GetGradient(const Vector &k) const
648 {
649  delete Jacobian;
650  SparseMatrix *localJ = Add(1.0, M->SpMat(), dt, S->SpMat());
651  add(*v, dt, k, w);
652  add(*x, dt, w, z);
653  localJ->Add(dt*dt, H->GetLocalGradient(z));
654  Jacobian = M->ParallelAssemble(localJ);
655  delete localJ;
656  HypreParMatrix *Je = Jacobian->EliminateRowsCols(ess_tdof_list);
657  delete Je;
658  return *Jacobian;
659 }
660 
661 ReducedSystemOperator::~ReducedSystemOperator()
662 {
663  delete Jacobian;
664 }
665 
666 
667 HyperelasticOperator::HyperelasticOperator(ParFiniteElementSpace &f,
668  Array<int> &ess_bdr, double visc,
669  double mu, double K,
670  NonlinearSolverType nls_type)
671  : TimeDependentOperator(2*f.TrueVSize(), 0.0), fespace(f),
672  M(&fespace), S(&fespace), H(&fespace),
673  viscosity(visc), M_solver(f.GetComm()), z(height/2),
674  local_grad_H(NULL), Jacobian(NULL)
675 {
676  const double rel_tol = 1e-8;
677  const int skip_zero_entries = 0;
678 
679  const double ref_density = 1.0; // density in the reference configuration
680  ConstantCoefficient rho0(ref_density);
681  M.AddDomainIntegrator(new VectorMassIntegrator(rho0));
682  M.Assemble(skip_zero_entries);
683  M.Finalize(skip_zero_entries);
684  Mmat = M.ParallelAssemble();
685  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
686  HypreParMatrix *Me = Mmat->EliminateRowsCols(ess_tdof_list);
687  delete Me;
688 
689  M_solver.iterative_mode = false;
690  M_solver.SetRelTol(rel_tol);
691  M_solver.SetAbsTol(0.0);
692  M_solver.SetMaxIter(30);
693  M_solver.SetPrintLevel(0);
694  M_prec.SetType(HypreSmoother::Jacobi);
695  M_solver.SetPreconditioner(M_prec);
696  M_solver.SetOperator(*Mmat);
697 
698  model = new NeoHookeanModel(mu, K);
699  H.AddDomainIntegrator(new HyperelasticNLFIntegrator(model));
700  H.SetEssentialTrueDofs(ess_tdof_list);
701 
702  ConstantCoefficient visc_coeff(viscosity);
703  S.AddDomainIntegrator(new VectorDiffusionIntegrator(visc_coeff));
704  S.Assemble(skip_zero_entries);
705  S.Finalize(skip_zero_entries);
706 
707  reduced_oper = new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
708 
709  HypreSmoother *J_hypreSmoother = new HypreSmoother;
710  J_hypreSmoother->SetType(HypreSmoother::l1Jacobi);
711  J_hypreSmoother->SetPositiveDiagonal(true);
712  J_prec = J_hypreSmoother;
713 
714  MINRESSolver *J_minres = new MINRESSolver(f.GetComm());
715  J_minres->SetRelTol(rel_tol);
716  J_minres->SetAbsTol(0.0);
717  J_minres->SetMaxIter(300);
718  J_minres->SetPrintLevel(-1);
719  J_minres->SetPreconditioner(*J_prec);
720  J_solver = J_minres;
721 
722  if (nls_type == KINSOL)
723  {
724  KINSolver *kinsolver = new KINSolver(f.GetComm(), KIN_LINESEARCH, true);
725  kinsolver->SetJFNK(true);
726  kinsolver->SetLSMaxIter(100);
727  newton_solver = kinsolver;
728  newton_solver->SetOperator(*reduced_oper);
729  newton_solver->SetMaxIter(200);
730  newton_solver->SetRelTol(rel_tol);
731  newton_solver->SetPrintLevel(1);
732  kinsolver->SetMaxSetupCalls(4);
733  }
734  else
735  {
736  newton_solver = new NewtonSolver(f.GetComm());
737  newton_solver->SetOperator(*reduced_oper);
738  newton_solver->SetMaxIter(10);
739  newton_solver->SetRelTol(rel_tol);
740  newton_solver->SetPrintLevel(-1);
741  }
742  newton_solver->SetSolver(*J_solver);
743  newton_solver->iterative_mode = false;
744 }
745 
746 void HyperelasticOperator::Mult(const Vector &vx, Vector &dvx_dt) const
747 {
748  // Create views to the sub-vectors v, x of vx, and dv_dt, dx_dt of dvx_dt
749  int sc = height/2;
750  Vector v(vx.GetData() + 0, sc);
751  Vector x(vx.GetData() + sc, sc);
752  Vector dv_dt(dvx_dt.GetData() + 0, sc);
753  Vector dx_dt(dvx_dt.GetData() + sc, sc);
754 
755  H.Mult(x, z);
756  if (viscosity != 0.0)
757  {
758  S.TrueAddMult(v, z);
759  z.SetSubVector(ess_tdof_list, 0.0);
760  }
761  z.Neg(); // z = -z
762  M_solver.Mult(z, dv_dt);
763 
764  dx_dt = v;
765 }
766 
767 void HyperelasticOperator::ImplicitSolve(const double dt,
768  const Vector &vx, Vector &dvx_dt)
769 {
770  int sc = height/2;
771  Vector v(vx.GetData() + 0, sc);
772  Vector x(vx.GetData() + sc, sc);
773  Vector dv_dt(dvx_dt.GetData() + 0, sc);
774  Vector dx_dt(dvx_dt.GetData() + sc, sc);
775 
776  // By eliminating kx from the coupled system:
777  // kv = -M^{-1}*[H(x + dt*kx) + S*(v + dt*kv)]
778  // kx = v + dt*kv
779  // we reduce it to a nonlinear equation for kv, represented by the
780  // reduced_oper. This equation is solved with the newton_solver
781  // object (using J_solver and J_prec internally).
782  reduced_oper->SetParameters(dt, &v, &x);
783  Vector zero; // empty vector is interpreted as zero r.h.s. by NewtonSolver
784  newton_solver->Mult(zero, dv_dt);
785  MFEM_VERIFY(newton_solver->GetConverged(),
786  "Nonlinear solver did not converge.");
787 #ifdef MFEM_DEBUG
788  if (fespace.GetMyRank() == 0)
789  {
790  cout << " num nonlin sol iters = " << newton_solver->GetNumIterations()
791  << ", final norm = " << newton_solver->GetFinalNorm() << '\n';
792  }
793 #endif
794  add(v, dt, dv_dt, dx_dt);
795 }
796 
797 int HyperelasticOperator::SUNImplicitSetup(const Vector &y,
798  const Vector &fy, int jok, int *jcur,
799  double gamma)
800 {
801  int sc = y.Size() / 2;
802  const Vector x(y.GetData() + sc, sc);
803 
804  // J = M + dt*(S + dt*grad(H))
805  if (Jacobian) { delete Jacobian; }
806  SparseMatrix *localJ = Add(1.0, M.SpMat(), gamma, S.SpMat());
807  local_grad_H = &H.GetLocalGradient(x);
808  localJ->Add(gamma*gamma, *local_grad_H);
809  Jacobian = M.ParallelAssemble(localJ);
810  delete localJ;
811  HypreParMatrix *Je = Jacobian->EliminateRowsCols(ess_tdof_list);
812  delete Je;
813 
814  // Set Jacobian solve operator
815  J_solver->SetOperator(*Jacobian);
816 
817  // Indicate that the Jacobian was updated
818  *jcur = 1;
819 
820  // Save gamma for use in solve
821  saved_gamma = gamma;
822 
823  // Return success
824  return 0;
825 }
826 
827 int HyperelasticOperator::SUNImplicitSolve(const Vector &b, Vector &x,
828  double tol)
829 {
830  int sc = b.Size() / 2;
831  ParFiniteElementSpace *fes = H.ParFESpace();
832  Vector b_v(b.GetData() + 0, sc);
833  Vector b_x(b.GetData() + sc, sc);
834  Vector x_v(x.GetData() + 0, sc);
835  Vector x_x(x.GetData() + sc, sc);
836  Vector rhs(sc);
837 
838  // We can assume that b_v and b_x have zeros at essential tdofs.
839 
840  // rhs = M b_v - dt*grad(H) b_x
841  ParGridFunction lb_x(fes), lrhs(fes);
842  lb_x.Distribute(b_x);
843  local_grad_H->Mult(lb_x, lrhs);
844  lrhs.ParallelAssemble(rhs);
845  rhs *= -saved_gamma;
846  M.TrueAddMult(b_v, rhs);
847  rhs.SetSubVector(ess_tdof_list, 0.0);
848 
849  J_solver->iterative_mode = false;
850  J_solver->Mult(rhs, x_v);
851 
852  add(b_x, saved_gamma, x_v, x_x);
853 
854  return 0;
855 }
856 
857 double HyperelasticOperator::ElasticEnergy(const ParGridFunction &x) const
858 {
859  return H.GetEnergy(x);
860 }
861 
862 double HyperelasticOperator::KineticEnergy(const ParGridFunction &v) const
863 {
864  double loc_energy = 0.5*M.InnerProduct(v, v);
865  double energy;
866  MPI_Allreduce(&loc_energy, &energy, 1, MPI_DOUBLE, MPI_SUM,
867  fespace.GetComm());
868  return energy;
869 }
870 
871 void HyperelasticOperator::GetElasticEnergyDensity(
872  const ParGridFunction &x, ParGridFunction &w) const
873 {
874  ElasticEnergyCoefficient w_coeff(*model, x);
875  w.ProjectCoefficient(w_coeff);
876 }
877 
878 HyperelasticOperator::~HyperelasticOperator()
879 {
880  delete Jacobian;
881  delete newton_solver;
882  delete J_solver;
883  delete J_prec;
884  delete reduced_oper;
885  delete model;
886  delete Mmat;
887 }
888 
889 
891  const IntegrationPoint &ip)
892 {
893  model.SetTransformation(T);
894  x.GetVectorGradient(T, J);
895  // return model.EvalW(J); // in reference configuration
896  return model.EvalW(J)/J.Det(); // in deformed configuration
897 }
898 
899 
900 void InitialDeformation(const Vector &x, Vector &y)
901 {
902  // set the initial configuration to be the same as the reference, stress
903  // free, configuration
904  y = x;
905 }
906 
907 void InitialVelocity(const Vector &x, Vector &v)
908 {
909  const int dim = x.Size();
910  const double s = 0.1/64.;
911 
912  v = 0.0;
913  v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
914  v(0) = -s*x(0)*x(0);
915 }
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
Definition: sundials.cpp:532
void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
Definition: ex10.cpp:379
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:1368
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:521
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:66
void Add(const int i, const int j, const double a)
Definition: sparsemat.cpp:2493
void InitialDeformation(const Vector &x, Vector &y)
Definition: ex10.cpp:591
Conjugate gradient method.
Definition: solvers.hpp:258
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:698
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:142
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition: pfespace.hpp:382
Base abstract class for first order time dependent operators.
Definition: operator.hpp:267
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:6641
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
void * GetMem() const
Access the SUNDIALS memory structure.
Definition: sundials.hpp:240
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
Parallel non-linear operator on the true dofs.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:819
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:1540
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:474
MINRES method.
Definition: solvers.hpp:368
int main(int argc, char *argv[])
Definition: ex1.cpp:66
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:389
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:169
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:261
void SetMaxSetupCalls(int max_calls)
Set maximum number of nonlinear iterations without a Jacobian update.
Definition: sundials.cpp:1961
int GetNRanks() const
Definition: pmesh.hpp:237
void InitialVelocity(const Vector &x, Vector &v)
Definition: ex10.cpp:598
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
Definition: hypre.hpp:685
Interface to ARKode&#39;s ARKStep module – additive Runge-Kutta methods.
Definition: sundials.hpp:539
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:1928
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:136
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:381
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
MPI_Comm GetComm() const
Definition: pfespace.hpp:239
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
Data type sparse matrix.
Definition: sparsemat.hpp:46
constexpr char vishost[]
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:314
Interface to the CVODE library – linear multi-step methods.
Definition: sundials.hpp:252
double b
Definition: lissajous.cpp:42
Interface to the KINSOL library – nonlinear solver methods.
Definition: sundials.hpp:713
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
constexpr int visport
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
Definition: gridfunc.cpp:206
void SetMaxIter(int max_it)
Definition: solvers.hpp:98
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:716
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:406
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4169
void SetJFNK(bool use_jfnk)
Set the Jacobian Free Newton Krylov flag. The default is false.
Definition: sundials.hpp:826
Parallel smoothers in hypre.
Definition: hypre.hpp:592
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
A general vector function coefficient.
int SpaceDimension() const
Definition: mesh.hpp:789
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
void SetAbsTol(double atol)
Definition: solvers.hpp:97
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:201
int GetMyRank() const
Definition: pmesh.hpp:238
void SetRelTol(double rtol)
Definition: solvers.hpp:96
MPI_Comm GetComm() const
Definition: pmesh.hpp:236
virtual void SetOperator(const Operator &op)
Set the nonlinear Operator of the system and initialize KINSOL.
Definition: sundials.cpp:1776
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
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 SetLSMaxIter(int m)
Set the maximum number of linear solver iterations.
Definition: sundials.hpp:830
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Definition: gridfunc.cpp:1642
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:402
Class for integration point with weight.
Definition: intrules.hpp:25
void PrintInfo() const
Print various ARKStep statistics.
Definition: sundials.cpp:1576
double mu
Definition: ex25.cpp:135
int dim
Definition: ex24.cpp:53
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
Class for parallel bilinear form.
Abstract class for hyperelastic models.
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
Vector data type.
Definition: vector.hpp:51
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:6603
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
void PrintInfo() const
Print various CVODE statistics.
Definition: sundials.cpp:742
Base class for solvers.
Definition: operator.hpp:634
Class for parallel grid function.
Definition: pgridfunc.hpp:32
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:66
The classical forward Euler method.
Definition: ode.hpp:116
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:181
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
Definition: sundials.cpp:1297
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:1534
void UseSundialsLinearSolver()
Attach a SUNDIALS GMRES linear solver to ARKode.
Definition: sundials.cpp:1462
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:2187
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
Definition: sundials.cpp:678
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1552
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:221
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145