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 // PETSc Modification
3 //
4 // Compile with: make ex10p
5 //
6 // Sample runs:
7 // mpirun -np 4 ex10p -m ../../data/beam-quad.mesh --petscopts rc_ex10p -s 3 -rs 2 -dt 3
8 // mpirun -np 4 ex10p -m ../../data/beam-quad-amr.mesh --petscopts rc_ex10p -s 3 -rs 2 -dt 3
9 //
10 // Description: This examples solves a time dependent nonlinear elasticity
11 // problem of the form dv/dt = H(x) + S v, dx/dt = v, where H is a
12 // hyperelastic model and S is a viscosity operator of Laplacian
13 // type. The geometry of the domain is assumed to be as follows:
14 //
15 // +---------------------+
16 // boundary --->| |
17 // attribute 1 | |
18 // (fixed) +---------------------+
19 //
20 // The example demonstrates the use of nonlinear operators (the
21 // class HyperelasticOperator defining H(x)), as well as their
22 // implicit time integration using a Newton method for solving an
23 // associated reduced backward-Euler type nonlinear equation
24 // (class ReducedSystemOperator). Each Newton step requires the
25 // inversion of a Jacobian matrix, which is done through a
26 // (preconditioned) inner solver. Note that implementing the
27 // method HyperelasticOperator::ImplicitSolve is the only
28 // requirement for high-order implicit (SDIRK) time integration.
29 // If using PETSc to solve the nonlinear problem, use the option
30 // files provided (see rc_ex10p, rc_ex10p_mf, rc_ex10p_mfop) that
31 // customize the Newton-Krylov method.
32 // When option --jfnk is used, PETSc will use a Jacobian-free
33 // Newton-Krylov method, using a user-defined preconditioner
34 // constructed with the PetscPreconditionerFactory class.
35 //
36 // We recommend viewing examples 2 and 9 before viewing this
37 // example.
38 
39 #include "mfem.hpp"
40 #include <memory>
41 #include <iostream>
42 #include <fstream>
43 
44 #ifndef MFEM_USE_PETSC
45 #error This example requires that MFEM is built with MFEM_USE_PETSC=YES
46 #endif
47 
48 using namespace std;
49 using namespace mfem;
50 
51 class ReducedSystemOperator;
52 
53 /** After spatial discretization, the hyperelastic model can be written as a
54  * system of ODEs:
55  * dv/dt = -M^{-1}*(H(x) + S*v)
56  * dx/dt = v,
57  * where x is the vector representing the deformation, v is the velocity field,
58  * M is the mass matrix, S is the viscosity matrix, and H(x) is the nonlinear
59  * hyperelastic operator.
60  *
61  * Class HyperelasticOperator represents the right-hand side of the above
62  * system of ODEs. */
63 class HyperelasticOperator : public TimeDependentOperator
64 {
65 protected:
66  ParFiniteElementSpace &fespace;
67  Array<int> ess_tdof_list;
68 
69  ParBilinearForm M, S;
71  double viscosity;
72  HyperelasticModel *model;
73 
74  HypreParMatrix *Mmat; // Mass matrix from ParallelAssemble()
75  CGSolver M_solver; // Krylov solver for inverting the mass matrix M
76  HypreSmoother M_prec; // Preconditioner for the mass matrix M
77 
78  /** Nonlinear operator defining the reduced backward Euler equation for the
79  velocity. Used in the implementation of method ImplicitSolve. */
80  ReducedSystemOperator *reduced_oper;
81 
82  /// Newton solver for the reduced backward Euler equation
83  NewtonSolver newton_solver;
84 
85  /// Newton solver for the reduced backward Euler equation (PETSc SNES)
86  PetscNonlinearSolver* pnewton_solver;
87 
88  /// Solver for the Jacobian solve in the Newton method
89  Solver *J_solver;
90  /// Preconditioner for the Jacobian solve in the Newton method
91  Solver *J_prec;
92  /// Preconditioner factory for JFNK
93  PetscPreconditionerFactory *J_factory;
94 
95  mutable Vector z; // auxiliary vector
96 
97 public:
98  HyperelasticOperator(ParFiniteElementSpace &f, Array<int> &ess_bdr,
99  double visc, double mu, double K,
100  bool use_petsc, bool petsc_use_jfnk);
101 
102  /// Compute the right-hand side of the ODE system.
103  virtual void Mult(const Vector &vx, Vector &dvx_dt) const;
104  /** Solve the Backward-Euler equation: k = f(x + dt*k, t), for the unknown k.
105  This is the only requirement for high-order SDIRK implicit integration.*/
106  virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k);
107 
108  double ElasticEnergy(const ParGridFunction &x) const;
109  double KineticEnergy(const ParGridFunction &v) const;
110  void GetElasticEnergyDensity(const ParGridFunction &x,
111  ParGridFunction &w) const;
112 
113  virtual ~HyperelasticOperator();
114 };
115 
116 /** Nonlinear operator of the form:
117  k --> (M + dt*S)*k + H(x + dt*v + dt^2*k) + S*v,
118  where M and S are given BilinearForms, H is a given NonlinearForm, v and x
119  are given vectors, and dt is a scalar. */
120 class ReducedSystemOperator : public Operator
121 {
122 private:
123  ParBilinearForm *M, *S;
124  ParNonlinearForm *H;
125  mutable HypreParMatrix *Jacobian;
126  double dt;
127  const Vector *v, *x;
128  mutable Vector w, z;
129  const Array<int> &ess_tdof_list;
130 
131 public:
132  ReducedSystemOperator(ParBilinearForm *M_, ParBilinearForm *S_,
133  ParNonlinearForm *H_, const Array<int> &ess_tdof_list);
134 
135  /// Set current dt, v, x values - needed to compute action and Jacobian.
136  void SetParameters(double dt_, const Vector *v_, const Vector *x_);
137 
138  /// Compute y = H(x + dt (v + dt k)) + M k + S (v + dt k).
139  virtual void Mult(const Vector &k, Vector &y) const;
140 
141  /// Compute J = M + dt S + dt^2 grad_H(x + dt (v + dt k)).
142  virtual Operator &GetGradient(const Vector &k) const;
143 
144  virtual ~ReducedSystemOperator();
145 
146 };
147 
148 /** Auxiliary class to provide preconditioners for matrix-free methods */
149 class PreconditionerFactory : public PetscPreconditionerFactory
150 {
151 private:
152  // const ReducedSystemOperator& op; // unused for now (generates warning)
153 
154 public:
155  PreconditionerFactory(const ReducedSystemOperator& op_, const string& name_)
156  : PetscPreconditionerFactory(name_) /* , op(op_) */ {}
157  virtual mfem::Solver* NewPreconditioner(const mfem::OperatorHandle&);
158  virtual ~PreconditionerFactory() {}
159 };
160 
161 /** Function representing the elastic energy density for the given hyperelastic
162  model+deformation. Used in HyperelasticOperator::GetElasticEnergyDensity. */
163 class ElasticEnergyCoefficient : public Coefficient
164 {
165 private:
166  HyperelasticModel &model;
167  const ParGridFunction &x;
168  DenseMatrix J;
169 
170 public:
171  ElasticEnergyCoefficient(HyperelasticModel &m, const ParGridFunction &x_)
172  : model(m), x(x_) { }
173  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
174  virtual ~ElasticEnergyCoefficient() { }
175 };
176 
177 void InitialDeformation(const Vector &x, Vector &y);
178 
179 void InitialVelocity(const Vector &x, Vector &v);
180 
181 void visualize(ostream &out, ParMesh *mesh, ParGridFunction *deformed_nodes,
182  ParGridFunction *field, const char *field_name = NULL,
183  bool init_vis = false);
184 
185 
186 int main(int argc, char *argv[])
187 {
188  // 1. Initialize MPI.
189  int num_procs, myid;
190  MPI_Init(&argc, &argv);
191  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
192  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
193 
194  // 2. Parse command-line options.
195  const char *mesh_file = "../../data/beam-quad.mesh";
196  int ser_ref_levels = 2;
197  int par_ref_levels = 0;
198  int order = 2;
199  int ode_solver_type = 3;
200  double t_final = 300.0;
201  double dt = 3.0;
202  double visc = 1e-2;
203  double mu = 0.25;
204  double K = 5.0;
205  bool visualization = true;
206  int vis_steps = 1;
207  bool use_petsc = true;
208  const char *petscrc_file = "";
209  bool petsc_use_jfnk = false;
210 
211  OptionsParser args(argc, argv);
212  args.AddOption(&mesh_file, "-m", "--mesh",
213  "Mesh file to use.");
214  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
215  "Number of times to refine the mesh uniformly in serial.");
216  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
217  "Number of times to refine the mesh uniformly in parallel.");
218  args.AddOption(&order, "-o", "--order",
219  "Order (degree) of the finite elements.");
220  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
221  "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
222  " 11 - Forward Euler, 12 - RK2,\n\t"
223  " 13 - RK3 SSP, 14 - RK4.");
224  args.AddOption(&t_final, "-tf", "--t-final",
225  "Final time; start time is 0.");
226  args.AddOption(&dt, "-dt", "--time-step",
227  "Time step.");
228  args.AddOption(&visc, "-v", "--viscosity",
229  "Viscosity coefficient.");
230  args.AddOption(&mu, "-mu", "--shear-modulus",
231  "Shear modulus in the Neo-Hookean hyperelastic model.");
232  args.AddOption(&K, "-K", "--bulk-modulus",
233  "Bulk modulus in the Neo-Hookean hyperelastic model.");
234  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
235  "--no-visualization",
236  "Enable or disable GLVis visualization.");
237  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
238  "Visualize every n-th timestep.");
239  args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
240  "--no-petsc",
241  "Use or not PETSc to solve the nonlinear system.");
242  args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
243  "PetscOptions file to use.");
244  args.AddOption(&petsc_use_jfnk, "-jfnk", "--jfnk", "-no-jfnk",
245  "--no-jfnk",
246  "Use JFNK with user-defined preconditioner factory.");
247  args.Parse();
248  if (!args.Good())
249  {
250  if (myid == 0)
251  {
252  args.PrintUsage(cout);
253  }
254  MPI_Finalize();
255  return 1;
256  }
257  if (myid == 0)
258  {
259  args.PrintOptions(cout);
260  }
261 
262  // 2b. We initialize PETSc
263  if (use_petsc)
264  {
265  MFEMInitializePetsc(NULL,NULL,petscrc_file,NULL);
266  }
267 
268  // 3. Read the serial mesh from the given mesh file on all processors. We can
269  // handle triangular, quadrilateral, tetrahedral and hexahedral meshes
270  // with the same code.
271  Mesh *mesh = new Mesh(mesh_file, 1, 1);
272  int dim = mesh->Dimension();
273 
274  // 4. Define the ODE solver used for time integration. Several implicit
275  // singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
276  // explicit Runge-Kutta methods are available.
277  ODESolver *ode_solver;
278  switch (ode_solver_type)
279  {
280  // Implicit L-stable methods
281  case 1: ode_solver = new BackwardEulerSolver; break;
282  case 2: ode_solver = new SDIRK23Solver(2); break;
283  case 3: ode_solver = new SDIRK33Solver; break;
284  // Explicit methods
285  case 11: ode_solver = new ForwardEulerSolver; break;
286  case 12: ode_solver = new RK2Solver(0.5); break; // midpoint method
287  case 13: ode_solver = new RK3SSPSolver; break;
288  case 14: ode_solver = new RK4Solver; break;
289  // Implicit A-stable methods (not L-stable)
290  case 22: ode_solver = new ImplicitMidpointSolver; break;
291  case 23: ode_solver = new SDIRK23Solver; break;
292  case 24: ode_solver = new SDIRK34Solver; break;
293  default:
294  if (myid == 0)
295  {
296  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
297  }
298  MPI_Finalize();
299  return 3;
300  }
301 
302  // 5. Refine the mesh in serial to increase the resolution. In this example
303  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
304  // a command-line parameter.
305  for (int lev = 0; lev < ser_ref_levels; lev++)
306  {
307  mesh->UniformRefinement();
308  }
309 
310  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
311  // this mesh further in parallel to increase the resolution. Once the
312  // parallel mesh is defined, the serial mesh can be deleted.
313  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
314  delete mesh;
315  for (int lev = 0; lev < par_ref_levels; lev++)
316  {
317  pmesh->UniformRefinement();
318  }
319 
320  // 7. Define the parallel vector finite element spaces representing the mesh
321  // deformation x_gf, the velocity v_gf, and the initial configuration,
322  // x_ref. Define also the elastic energy density, w_gf, which is in a
323  // discontinuous higher-order space. Since x and v are integrated in time
324  // as a system, we group them together in block vector vx, on the unique
325  // parallel degrees of freedom, with offsets given by array true_offset.
326  H1_FECollection fe_coll(order, dim);
327  ParFiniteElementSpace fespace(pmesh, &fe_coll, dim);
328 
329  HYPRE_Int glob_size = fespace.GlobalTrueVSize();
330  if (myid == 0)
331  {
332  cout << "Number of velocity/deformation unknowns: " << glob_size << endl;
333  }
334  int true_size = fespace.TrueVSize();
335  Array<int> true_offset(3);
336  true_offset[0] = 0;
337  true_offset[1] = true_size;
338  true_offset[2] = 2*true_size;
339 
340  BlockVector vx(true_offset);
341  ParGridFunction v_gf, x_gf;
342  v_gf.MakeTRef(&fespace, vx, true_offset[0]);
343  x_gf.MakeTRef(&fespace, vx, true_offset[1]);
344 
345  ParGridFunction x_ref(&fespace);
346  pmesh->GetNodes(x_ref);
347 
348  L2_FECollection w_fec(order + 1, dim);
349  ParFiniteElementSpace w_fespace(pmesh, &w_fec);
350  ParGridFunction w_gf(&w_fespace);
351 
352  // 8. Set the initial conditions for v_gf, x_gf and vx, and define the
353  // boundary conditions on a beam-like mesh (see description above).
355  v_gf.ProjectCoefficient(velo);
356  v_gf.SetTrueVector();
358  x_gf.ProjectCoefficient(deform);
359  x_gf.SetTrueVector();
360 
361  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
362 
363  Array<int> ess_bdr(fespace.GetMesh()->bdr_attributes.Max());
364  ess_bdr = 0;
365  ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed
366 
367  // 9. Initialize the hyperelastic operator, the GLVis visualization and print
368  // the initial energies.
369  HyperelasticOperator *oper = new HyperelasticOperator(fespace, ess_bdr, visc,
370  mu, K, use_petsc,
371  petsc_use_jfnk);
372 
373  socketstream vis_v, vis_w;
374  if (visualization)
375  {
376  char vishost[] = "localhost";
377  int visport = 19916;
378  vis_v.open(vishost, visport);
379  vis_v.precision(8);
380  visualize(vis_v, pmesh, &x_gf, &v_gf, "Velocity", true);
381  // Make sure all ranks have sent their 'v' solution before initiating
382  // another set of GLVis connections (one from each rank):
383  MPI_Barrier(pmesh->GetComm());
384  vis_w.open(vishost, visport);
385  if (vis_w)
386  {
387  oper->GetElasticEnergyDensity(x_gf, w_gf);
388  vis_w.precision(8);
389  visualize(vis_w, pmesh, &x_gf, &w_gf, "Elastic energy density", true);
390  }
391  }
392 
393  double ee0 = oper->ElasticEnergy(x_gf);
394  double ke0 = oper->KineticEnergy(v_gf);
395  if (myid == 0)
396  {
397  cout << "initial elastic energy (EE) = " << ee0 << endl;
398  cout << "initial kinetic energy (KE) = " << ke0 << endl;
399  cout << "initial total energy (TE) = " << (ee0 + ke0) << endl;
400  }
401 
402  double t = 0.0;
403  oper->SetTime(t);
404  ode_solver->Init(*oper);
405 
406  // 10. Perform time-integration
407  // (looping over the time iterations, ti, with a time-step dt).
408  bool last_step = false;
409  for (int ti = 1; !last_step; ti++)
410  {
411  double dt_real = min(dt, t_final - t);
412 
413  ode_solver->Step(vx, t, dt_real);
414 
415  last_step = (t >= t_final - 1e-8*dt);
416 
417  if (last_step || (ti % vis_steps) == 0)
418  {
419  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
420 
421  double ee = oper->ElasticEnergy(x_gf);
422  double ke = oper->KineticEnergy(v_gf);
423 
424  if (myid == 0)
425  {
426  cout << "step " << ti << ", t = " << t << ", EE = " << ee
427  << ", KE = " << ke << ", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
428  }
429 
430  if (visualization)
431  {
432  visualize(vis_v, pmesh, &x_gf, &v_gf);
433  if (vis_w)
434  {
435  oper->GetElasticEnergyDensity(x_gf, w_gf);
436  visualize(vis_w, pmesh, &x_gf, &w_gf);
437  }
438  }
439  }
440  }
441 
442  // 11. Save the displaced mesh, the velocity and elastic energy.
443  {
444  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
445  GridFunction *nodes = &x_gf;
446  int owns_nodes = 0;
447  pmesh->SwapNodes(nodes, owns_nodes);
448 
449  ostringstream mesh_name, velo_name, ee_name;
450  mesh_name << "deformed." << setfill('0') << setw(6) << myid;
451  velo_name << "velocity." << setfill('0') << setw(6) << myid;
452  ee_name << "elastic_energy." << setfill('0') << setw(6) << myid;
453 
454  ofstream mesh_ofs(mesh_name.str().c_str());
455  mesh_ofs.precision(8);
456  pmesh->Print(mesh_ofs);
457  pmesh->SwapNodes(nodes, owns_nodes);
458  ofstream velo_ofs(velo_name.str().c_str());
459  velo_ofs.precision(8);
460  v_gf.Save(velo_ofs);
461  ofstream ee_ofs(ee_name.str().c_str());
462  ee_ofs.precision(8);
463  oper->GetElasticEnergyDensity(x_gf, w_gf);
464  w_gf.Save(ee_ofs);
465  }
466 
467  // 12. Free the used memory.
468  delete ode_solver;
469  delete pmesh;
470  delete oper;
471 
472  // We finalize PETSc
473  if (use_petsc) { MFEMFinalizePetsc(); }
474 
475  MPI_Finalize();
476 
477  return 0;
478 }
479 
480 void visualize(ostream &out, ParMesh *mesh, ParGridFunction *deformed_nodes,
481  ParGridFunction *field, const char *field_name, bool init_vis)
482 {
483  if (!out)
484  {
485  return;
486  }
487 
488  GridFunction *nodes = deformed_nodes;
489  int owns_nodes = 0;
490 
491  mesh->SwapNodes(nodes, owns_nodes);
492 
493  out << "parallel " << mesh->GetNRanks() << " " << mesh->GetMyRank() << "\n";
494  out << "solution\n" << *mesh << *field;
495 
496  mesh->SwapNodes(nodes, owns_nodes);
497 
498  if (init_vis)
499  {
500  out << "window_size 800 800\n";
501  out << "window_title '" << field_name << "'\n";
502  if (mesh->SpaceDimension() == 2)
503  {
504  out << "view 0 0\n"; // view from top
505  out << "keys jl\n"; // turn off perspective and light
506  }
507  out << "keys cm\n"; // show colorbar and mesh
508  out << "autoscale value\n"; // update value-range; keep mesh-extents fixed
509  out << "pause\n";
510  }
511  out << flush;
512 }
513 
514 
515 ReducedSystemOperator::ReducedSystemOperator(
517  const Array<int> &ess_tdof_list_)
518  : Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
519  Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height),
520  ess_tdof_list(ess_tdof_list_)
521 { }
522 
523 void ReducedSystemOperator::SetParameters(double dt_, const Vector *v_,
524  const Vector *x_)
525 {
526  dt = dt_; v = v_; x = x_;
527 }
528 
529 void ReducedSystemOperator::Mult(const Vector &k, Vector &y) const
530 {
531  // compute: y = H(x + dt*(v + dt*k)) + M*k + S*(v + dt*k)
532  add(*v, dt, k, w);
533  add(*x, dt, w, z);
534  H->Mult(z, y);
535  M->TrueAddMult(k, y);
536  S->TrueAddMult(w, y);
537  y.SetSubVector(ess_tdof_list, 0.0);
538 }
539 
540 Operator &ReducedSystemOperator::GetGradient(const Vector &k) const
541 {
542  delete Jacobian;
543  SparseMatrix *localJ = Add(1.0, M->SpMat(), dt, S->SpMat());
544  add(*v, dt, k, w);
545  add(*x, dt, w, z);
546  localJ->Add(dt*dt, H->GetLocalGradient(z));
547  // if we are using PETSc, the HypreParCSR Jacobian will be converted to
548  // PETSc's AIJ on the fly
549  Jacobian = M->ParallelAssemble(localJ);
550  delete localJ;
551  HypreParMatrix *Je = Jacobian->EliminateRowsCols(ess_tdof_list);
552  delete Je;
553  return *Jacobian;
554 }
555 
556 ReducedSystemOperator::~ReducedSystemOperator()
557 {
558  delete Jacobian;
559 }
560 
561 
562 HyperelasticOperator::HyperelasticOperator(ParFiniteElementSpace &f,
563  Array<int> &ess_bdr, double visc,
564  double mu, double K, bool use_petsc,
565  bool use_petsc_factory)
566  : TimeDependentOperator(2*f.TrueVSize(), 0.0), fespace(f),
567  M(&fespace), S(&fespace), H(&fespace),
568  viscosity(visc), M_solver(f.GetComm()),
569  newton_solver(f.GetComm()), pnewton_solver(NULL), z(height/2)
570 {
571  const double rel_tol = 1e-8;
572  const int skip_zero_entries = 0;
573 
574  const double ref_density = 1.0; // density in the reference configuration
575  ConstantCoefficient rho0(ref_density);
576  M.AddDomainIntegrator(new VectorMassIntegrator(rho0));
577  M.Assemble(skip_zero_entries);
578  M.Finalize(skip_zero_entries);
579  Mmat = M.ParallelAssemble();
580  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
581  HypreParMatrix *Me = Mmat->EliminateRowsCols(ess_tdof_list);
582  delete Me;
583 
584  M_solver.iterative_mode = false;
585  M_solver.SetRelTol(rel_tol);
586  M_solver.SetAbsTol(0.0);
587  M_solver.SetMaxIter(30);
588  M_solver.SetPrintLevel(0);
589  M_prec.SetType(HypreSmoother::Jacobi);
590  M_solver.SetPreconditioner(M_prec);
591  M_solver.SetOperator(*Mmat);
592 
593  model = new NeoHookeanModel(mu, K);
594  H.AddDomainIntegrator(new HyperelasticNLFIntegrator(model));
595  H.SetEssentialTrueDofs(ess_tdof_list);
596 
597  ConstantCoefficient visc_coeff(viscosity);
598  S.AddDomainIntegrator(new VectorDiffusionIntegrator(visc_coeff));
599  S.Assemble(skip_zero_entries);
600  S.Finalize(skip_zero_entries);
601 
602  reduced_oper = new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
603  if (!use_petsc)
604  {
605  HypreSmoother *J_hypreSmoother = new HypreSmoother;
606  J_hypreSmoother->SetType(HypreSmoother::l1Jacobi);
607  J_hypreSmoother->SetPositiveDiagonal(true);
608  J_prec = J_hypreSmoother;
609 
610  MINRESSolver *J_minres = new MINRESSolver(f.GetComm());
611  J_minres->SetRelTol(rel_tol);
612  J_minres->SetAbsTol(0.0);
613  J_minres->SetMaxIter(300);
614  J_minres->SetPrintLevel(-1);
615  J_minres->SetPreconditioner(*J_prec);
616  J_solver = J_minres;
617 
618  J_factory = NULL;
619 
620  newton_solver.iterative_mode = false;
621  newton_solver.SetSolver(*J_solver);
622  newton_solver.SetOperator(*reduced_oper);
623  newton_solver.SetPrintLevel(1); // print Newton iterations
624  newton_solver.SetRelTol(rel_tol);
625  newton_solver.SetAbsTol(0.0);
626  newton_solver.SetMaxIter(10);
627  }
628  else
629  {
630  // if using PETSc, we create the same solver (Newton + MINRES + Jacobi)
631  // by command line options (see rc_ex10p)
632  J_solver = NULL;
633  J_prec = NULL;
634  J_factory = NULL;
635  pnewton_solver = new PetscNonlinearSolver(f.GetComm(),
636  *reduced_oper);
637 
638  // we can setup a factory to construct a "physics-based" preconditioner
639  if (use_petsc_factory)
640  {
641  J_factory = new PreconditionerFactory(*reduced_oper, "JFNK preconditioner");
642  pnewton_solver->SetPreconditionerFactory(J_factory);
643  }
644  pnewton_solver->SetPrintLevel(1); // print Newton iterations
645  pnewton_solver->SetRelTol(rel_tol);
646  pnewton_solver->SetAbsTol(0.0);
647  pnewton_solver->SetMaxIter(10);
648  }
649 }
650 
651 void HyperelasticOperator::Mult(const Vector &vx, Vector &dvx_dt) const
652 {
653  // Create views to the sub-vectors v, x of vx, and dv_dt, dx_dt of dvx_dt
654  int sc = height/2;
655  Vector v(vx.GetData() + 0, sc);
656  Vector x(vx.GetData() + sc, sc);
657  Vector dv_dt(dvx_dt.GetData() + 0, sc);
658  Vector dx_dt(dvx_dt.GetData() + sc, sc);
659 
660  H.Mult(x, z);
661  if (viscosity != 0.0)
662  {
663  S.TrueAddMult(v, z);
664  z.SetSubVector(ess_tdof_list, 0.0);
665  }
666  z.Neg(); // z = -z
667  M_solver.Mult(z, dv_dt);
668 
669  dx_dt = v;
670 }
671 
672 void HyperelasticOperator::ImplicitSolve(const double dt,
673  const Vector &vx, Vector &dvx_dt)
674 {
675  int sc = height/2;
676  Vector v(vx.GetData() + 0, sc);
677  Vector x(vx.GetData() + sc, sc);
678  Vector dv_dt(dvx_dt.GetData() + 0, sc);
679  Vector dx_dt(dvx_dt.GetData() + sc, sc);
680 
681  // By eliminating kx from the coupled system:
682  // kv = -M^{-1}*[H(x + dt*kx) + S*(v + dt*kv)]
683  // kx = v + dt*kv
684  // we reduce it to a nonlinear equation for kv, represented by the
685  // reduced_oper. This equation is solved with the newton_solver
686  // object (using J_solver and J_prec internally).
687  reduced_oper->SetParameters(dt, &v, &x);
688  Vector zero; // empty vector is interpreted as zero r.h.s. by NewtonSolver
689  if (!pnewton_solver)
690  {
691  newton_solver.Mult(zero, dv_dt);
692  MFEM_VERIFY(newton_solver.GetConverged(),
693  "Newton solver did not converge.");
694  }
695  else
696  {
697  pnewton_solver->Mult(zero, dv_dt);
698  MFEM_VERIFY(pnewton_solver->GetConverged(),
699  "Newton solver did not converge.");
700  }
701  add(v, dt, dv_dt, dx_dt);
702 }
703 
704 double HyperelasticOperator::ElasticEnergy(const ParGridFunction &x) const
705 {
706  return H.GetEnergy(x);
707 }
708 
709 double HyperelasticOperator::KineticEnergy(const ParGridFunction &v) const
710 {
711  double loc_energy = 0.5*M.InnerProduct(v, v);
712  double energy;
713  MPI_Allreduce(&loc_energy, &energy, 1, MPI_DOUBLE, MPI_SUM,
714  fespace.GetComm());
715  return energy;
716 }
717 
718 void HyperelasticOperator::GetElasticEnergyDensity(
719  const ParGridFunction &x, ParGridFunction &w) const
720 {
721  ElasticEnergyCoefficient w_coeff(*model, x);
722  w.ProjectCoefficient(w_coeff);
723 }
724 
725 HyperelasticOperator::~HyperelasticOperator()
726 {
727  delete J_solver;
728  delete J_prec;
729  delete J_factory;
730  delete reduced_oper;
731  delete model;
732  delete Mmat;
733  delete pnewton_solver;
734 }
735 
736 // This method gets called every time we need a preconditioner "oh"
737 // contains the PetscParMatrix that wraps the operator constructed in
738 // the GetGradient() method (see also PetscSolver::SetJacobianType()).
739 // In this example, we just return a customizable PetscPreconditioner
740 // using that matrix. However, the OperatorHandle argument can be
741 // ignored, and any "physics-based" solver can be constructed since we
742 // have access to the HyperElasticOperator class.
743 Solver* PreconditionerFactory::NewPreconditioner(const mfem::OperatorHandle& oh)
744 {
745  PetscParMatrix *pP;
746  oh.Get(pP);
747  return new PetscPreconditioner(*pP,"jfnk_");
748 }
749 
751  const IntegrationPoint &ip)
752 {
753  model.SetTransformation(T);
754  x.GetVectorGradient(T, J);
755  // return model.EvalW(J); // in reference configuration
756  return model.EvalW(J)/J.Det(); // in deformed configuration
757 }
758 
759 
760 void InitialDeformation(const Vector &x, Vector &y)
761 {
762  // set the initial configuration to be the same as the reference,
763  // stress free, configuration
764  y = x;
765 }
766 
767 void InitialVelocity(const Vector &x, Vector &v)
768 {
769  const int dim = x.Size();
770  const double s = 0.1/64.;
771 
772  v = 0.0;
773  v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
774  v(0) = -s*x(0)*x(0);
775 }
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
Abstract class for PETSc&#39;s preconditioners.
Definition: petsc.hpp:663
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
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:197
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
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
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: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
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:474
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:638
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
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
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
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
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
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
Abstract class for PETSc&#39;s nonlinear solvers.
Definition: petsc.hpp:754
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
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 GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Definition: gridfunc.cpp:1642
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:402
void MFEMFinalizePetsc()
Definition: petsc.cpp:164
Class for integration point with weight.
Definition: intrules.hpp:25
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Definition: petsc.cpp:147
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;.
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
Definition: handle.hpp:106
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
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
Class for parallel meshes.
Definition: pmesh.hpp:32
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition: hypre.cpp:2187
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