MFEM  v4.6.0
Finite element discretization library
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 &os, 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 and HYPRE.
189  Mpi::Init(argc, argv);
190  int myid = Mpi::WorldRank();
191  Hypre::Init();
192 
193  // 2. Parse command-line options.
194  const char *mesh_file = "../../data/beam-quad.mesh";
195  int ser_ref_levels = 2;
196  int par_ref_levels = 0;
197  int order = 2;
198  int ode_solver_type = 3;
199  double t_final = 300.0;
200  double dt = 3.0;
201  double visc = 1e-2;
202  double mu = 0.25;
203  double K = 5.0;
204  bool visualization = true;
205  int vis_steps = 1;
206  bool use_petsc = true;
207  const char *petscrc_file = "";
208  bool petsc_use_jfnk = false;
209 
210  OptionsParser args(argc, argv);
211  args.AddOption(&mesh_file, "-m", "--mesh",
212  "Mesh file to use.");
213  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
214  "Number of times to refine the mesh uniformly in serial.");
215  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
216  "Number of times to refine the mesh uniformly in parallel.");
217  args.AddOption(&order, "-o", "--order",
218  "Order (degree) of the finite elements.");
219  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
220  "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
221  " 11 - Forward Euler, 12 - RK2,\n\t"
222  " 13 - RK3 SSP, 14 - RK4.");
223  args.AddOption(&t_final, "-tf", "--t-final",
224  "Final time; start time is 0.");
225  args.AddOption(&dt, "-dt", "--time-step",
226  "Time step.");
227  args.AddOption(&visc, "-v", "--viscosity",
228  "Viscosity coefficient.");
229  args.AddOption(&mu, "-mu", "--shear-modulus",
230  "Shear modulus in the Neo-Hookean hyperelastic model.");
231  args.AddOption(&K, "-K", "--bulk-modulus",
232  "Bulk modulus in the Neo-Hookean hyperelastic model.");
233  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
234  "--no-visualization",
235  "Enable or disable GLVis visualization.");
236  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
237  "Visualize every n-th timestep.");
238  args.AddOption(&use_petsc, "-usepetsc", "--usepetsc", "-no-petsc",
239  "--no-petsc",
240  "Use or not PETSc to solve the nonlinear system.");
241  args.AddOption(&petscrc_file, "-petscopts", "--petscopts",
242  "PetscOptions file to use.");
243  args.AddOption(&petsc_use_jfnk, "-jfnk", "--jfnk", "-no-jfnk",
244  "--no-jfnk",
245  "Use JFNK with user-defined preconditioner factory.");
246  args.Parse();
247  if (!args.Good())
248  {
249  if (myid == 0)
250  {
251  args.PrintUsage(cout);
252  }
253  return 1;
254  }
255  if (myid == 0)
256  {
257  args.PrintOptions(cout);
258  }
259 
260  // 2b. We initialize PETSc
261  if (use_petsc)
262  {
263  MFEMInitializePetsc(NULL,NULL,petscrc_file,NULL);
264  }
265 
266  // 3. Read the serial mesh from the given mesh file on all processors. We can
267  // handle triangular, quadrilateral, tetrahedral and hexahedral meshes
268  // with the same code.
269  Mesh *mesh = new Mesh(mesh_file, 1, 1);
270  int dim = mesh->Dimension();
271 
272  // 4. Define the ODE solver used for time integration. Several implicit
273  // singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
274  // explicit Runge-Kutta methods are available.
275  ODESolver *ode_solver;
276  switch (ode_solver_type)
277  {
278  // Implicit L-stable methods
279  case 1: ode_solver = new BackwardEulerSolver; break;
280  case 2: ode_solver = new SDIRK23Solver(2); break;
281  case 3: ode_solver = new SDIRK33Solver; break;
282  // Explicit methods
283  case 11: ode_solver = new ForwardEulerSolver; break;
284  case 12: ode_solver = new RK2Solver(0.5); break; // midpoint method
285  case 13: ode_solver = new RK3SSPSolver; break;
286  case 14: ode_solver = new RK4Solver; break;
287  // Implicit A-stable methods (not L-stable)
288  case 22: ode_solver = new ImplicitMidpointSolver; break;
289  case 23: ode_solver = new SDIRK23Solver; break;
290  case 24: ode_solver = new SDIRK34Solver; break;
291  default:
292  if (myid == 0)
293  {
294  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
295  }
296  return 3;
297  }
298 
299  // 5. Refine the mesh in serial to increase the resolution. In this example
300  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
301  // a command-line parameter.
302  for (int lev = 0; lev < ser_ref_levels; lev++)
303  {
304  mesh->UniformRefinement();
305  }
306 
307  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
308  // this mesh further in parallel to increase the resolution. Once the
309  // parallel mesh is defined, the serial mesh can be deleted.
310  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
311  delete mesh;
312  for (int lev = 0; lev < par_ref_levels; lev++)
313  {
314  pmesh->UniformRefinement();
315  }
316 
317  // 7. Define the parallel vector finite element spaces representing the mesh
318  // deformation x_gf, the velocity v_gf, and the initial configuration,
319  // x_ref. Define also the elastic energy density, w_gf, which is in a
320  // discontinuous higher-order space. Since x and v are integrated in time
321  // as a system, we group them together in block vector vx, on the unique
322  // parallel degrees of freedom, with offsets given by array true_offset.
323  H1_FECollection fe_coll(order, dim);
324  ParFiniteElementSpace fespace(pmesh, &fe_coll, dim);
325 
326  HYPRE_BigInt glob_size = fespace.GlobalTrueVSize();
327  if (myid == 0)
328  {
329  cout << "Number of velocity/deformation unknowns: " << glob_size << endl;
330  }
331  int true_size = fespace.TrueVSize();
332  Array<int> true_offset(3);
333  true_offset[0] = 0;
334  true_offset[1] = true_size;
335  true_offset[2] = 2*true_size;
336 
337  BlockVector vx(true_offset);
338  ParGridFunction v_gf, x_gf;
339  v_gf.MakeTRef(&fespace, vx, true_offset[0]);
340  x_gf.MakeTRef(&fespace, vx, true_offset[1]);
341 
342  ParGridFunction x_ref(&fespace);
343  pmesh->GetNodes(x_ref);
344 
345  L2_FECollection w_fec(order + 1, dim);
346  ParFiniteElementSpace w_fespace(pmesh, &w_fec);
347  ParGridFunction w_gf(&w_fespace);
348 
349  // 8. Set the initial conditions for v_gf, x_gf and vx, and define the
350  // boundary conditions on a beam-like mesh (see description above).
352  v_gf.ProjectCoefficient(velo);
353  v_gf.SetTrueVector();
355  x_gf.ProjectCoefficient(deform);
356  x_gf.SetTrueVector();
357 
358  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
359 
360  Array<int> ess_bdr(fespace.GetMesh()->bdr_attributes.Max());
361  ess_bdr = 0;
362  ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed
363 
364  // 9. Initialize the hyperelastic operator, the GLVis visualization and print
365  // the initial energies.
366  HyperelasticOperator *oper = new HyperelasticOperator(fespace, ess_bdr, visc,
367  mu, K, use_petsc,
368  petsc_use_jfnk);
369 
370  socketstream vis_v, vis_w;
371  if (visualization)
372  {
373  char vishost[] = "localhost";
374  int visport = 19916;
375  vis_v.open(vishost, visport);
376  vis_v.precision(8);
377  visualize(vis_v, pmesh, &x_gf, &v_gf, "Velocity", true);
378  // Make sure all ranks have sent their 'v' solution before initiating
379  // another set of GLVis connections (one from each rank):
380  MPI_Barrier(pmesh->GetComm());
381  vis_w.open(vishost, visport);
382  if (vis_w)
383  {
384  oper->GetElasticEnergyDensity(x_gf, w_gf);
385  vis_w.precision(8);
386  visualize(vis_w, pmesh, &x_gf, &w_gf, "Elastic energy density", true);
387  }
388  }
389 
390  double ee0 = oper->ElasticEnergy(x_gf);
391  double ke0 = oper->KineticEnergy(v_gf);
392  if (myid == 0)
393  {
394  cout << "initial elastic energy (EE) = " << ee0 << endl;
395  cout << "initial kinetic energy (KE) = " << ke0 << endl;
396  cout << "initial total energy (TE) = " << (ee0 + ke0) << endl;
397  }
398 
399  double t = 0.0;
400  oper->SetTime(t);
401  ode_solver->Init(*oper);
402 
403  // 10. Perform time-integration
404  // (looping over the time iterations, ti, with a time-step dt).
405  bool last_step = false;
406  for (int ti = 1; !last_step; ti++)
407  {
408  double dt_real = min(dt, t_final - t);
409 
410  ode_solver->Step(vx, t, dt_real);
411 
412  last_step = (t >= t_final - 1e-8*dt);
413 
414  if (last_step || (ti % vis_steps) == 0)
415  {
416  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
417 
418  double ee = oper->ElasticEnergy(x_gf);
419  double ke = oper->KineticEnergy(v_gf);
420 
421  if (myid == 0)
422  {
423  cout << "step " << ti << ", t = " << t << ", EE = " << ee
424  << ", KE = " << ke << ", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
425  }
426 
427  if (visualization)
428  {
429  visualize(vis_v, pmesh, &x_gf, &v_gf);
430  if (vis_w)
431  {
432  oper->GetElasticEnergyDensity(x_gf, w_gf);
433  visualize(vis_w, pmesh, &x_gf, &w_gf);
434  }
435  }
436  }
437  }
438 
439  // 11. Save the displaced mesh, the velocity and elastic energy.
440  {
441  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
442  GridFunction *nodes = &x_gf;
443  int owns_nodes = 0;
444  pmesh->SwapNodes(nodes, owns_nodes);
445 
446  ostringstream mesh_name, velo_name, ee_name;
447  mesh_name << "deformed." << setfill('0') << setw(6) << myid;
448  velo_name << "velocity." << setfill('0') << setw(6) << myid;
449  ee_name << "elastic_energy." << setfill('0') << setw(6) << myid;
450 
451  ofstream mesh_ofs(mesh_name.str().c_str());
452  mesh_ofs.precision(8);
453  pmesh->Print(mesh_ofs);
454  pmesh->SwapNodes(nodes, owns_nodes);
455  ofstream velo_ofs(velo_name.str().c_str());
456  velo_ofs.precision(8);
457  v_gf.Save(velo_ofs);
458  ofstream ee_ofs(ee_name.str().c_str());
459  ee_ofs.precision(8);
460  oper->GetElasticEnergyDensity(x_gf, w_gf);
461  w_gf.Save(ee_ofs);
462  }
463 
464  // 12. Free the used memory.
465  delete ode_solver;
466  delete pmesh;
467  delete oper;
468 
469  // We finalize PETSc
470  if (use_petsc) { MFEMFinalizePetsc(); }
471 
472  return 0;
473 }
474 
475 void visualize(ostream &os, ParMesh *mesh, ParGridFunction *deformed_nodes,
476  ParGridFunction *field, const char *field_name, bool init_vis)
477 {
478  if (!os)
479  {
480  return;
481  }
482 
483  GridFunction *nodes = deformed_nodes;
484  int owns_nodes = 0;
485 
486  mesh->SwapNodes(nodes, owns_nodes);
487 
488  os << "parallel " << mesh->GetNRanks() << " " << mesh->GetMyRank() << "\n";
489  os << "solution\n" << *mesh << *field;
490 
491  mesh->SwapNodes(nodes, owns_nodes);
492 
493  if (init_vis)
494  {
495  os << "window_size 800 800\n";
496  os << "window_title '" << field_name << "'\n";
497  if (mesh->SpaceDimension() == 2)
498  {
499  os << "view 0 0\n"; // view from top
500  os << "keys jl\n"; // turn off perspective and light
501  }
502  os << "keys cm\n"; // show colorbar and mesh
503  os << "autoscale value\n"; // update value-range; keep mesh-extents fixed
504  os << "pause\n";
505  }
506  os << flush;
507 }
508 
509 
510 ReducedSystemOperator::ReducedSystemOperator(
512  const Array<int> &ess_tdof_list_)
513  : Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
514  Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height),
515  ess_tdof_list(ess_tdof_list_)
516 { }
517 
518 void ReducedSystemOperator::SetParameters(double dt_, const Vector *v_,
519  const Vector *x_)
520 {
521  dt = dt_; v = v_; x = x_;
522 }
523 
524 void ReducedSystemOperator::Mult(const Vector &k, Vector &y) const
525 {
526  // compute: y = H(x + dt*(v + dt*k)) + M*k + S*(v + dt*k)
527  add(*v, dt, k, w);
528  add(*x, dt, w, z);
529  H->Mult(z, y);
530  M->TrueAddMult(k, y);
531  S->TrueAddMult(w, y);
532  y.SetSubVector(ess_tdof_list, 0.0);
533 }
534 
535 Operator &ReducedSystemOperator::GetGradient(const Vector &k) const
536 {
537  delete Jacobian;
538  SparseMatrix *localJ = Add(1.0, M->SpMat(), dt, S->SpMat());
539  add(*v, dt, k, w);
540  add(*x, dt, w, z);
541  localJ->Add(dt*dt, H->GetLocalGradient(z));
542  // if we are using PETSc, the HypreParCSR Jacobian will be converted to
543  // PETSc's AIJ on the fly
544  Jacobian = M->ParallelAssemble(localJ);
545  delete localJ;
546  HypreParMatrix *Je = Jacobian->EliminateRowsCols(ess_tdof_list);
547  delete Je;
548  return *Jacobian;
549 }
550 
551 ReducedSystemOperator::~ReducedSystemOperator()
552 {
553  delete Jacobian;
554 }
555 
556 
557 HyperelasticOperator::HyperelasticOperator(ParFiniteElementSpace &f,
558  Array<int> &ess_bdr, double visc,
559  double mu, double K, bool use_petsc,
560  bool use_petsc_factory)
561  : TimeDependentOperator(2*f.TrueVSize(), 0.0), fespace(f),
562  M(&fespace), S(&fespace), H(&fespace),
563  viscosity(visc), M_solver(f.GetComm()),
564  newton_solver(f.GetComm()), pnewton_solver(NULL), z(height/2)
565 {
566  const double rel_tol = 1e-8;
567  const int skip_zero_entries = 0;
568 
569  const double ref_density = 1.0; // density in the reference configuration
570  ConstantCoefficient rho0(ref_density);
571  M.AddDomainIntegrator(new VectorMassIntegrator(rho0));
572  M.Assemble(skip_zero_entries);
573  M.Finalize(skip_zero_entries);
574  Mmat = M.ParallelAssemble();
575  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
576  HypreParMatrix *Me = Mmat->EliminateRowsCols(ess_tdof_list);
577  delete Me;
578 
579  M_solver.iterative_mode = false;
580  M_solver.SetRelTol(rel_tol);
581  M_solver.SetAbsTol(0.0);
582  M_solver.SetMaxIter(30);
583  M_solver.SetPrintLevel(0);
584  M_prec.SetType(HypreSmoother::Jacobi);
585  M_solver.SetPreconditioner(M_prec);
586  M_solver.SetOperator(*Mmat);
587 
588  model = new NeoHookeanModel(mu, K);
589  H.AddDomainIntegrator(new HyperelasticNLFIntegrator(model));
590  H.SetEssentialTrueDofs(ess_tdof_list);
591 
592  ConstantCoefficient visc_coeff(viscosity);
593  S.AddDomainIntegrator(new VectorDiffusionIntegrator(visc_coeff));
594  S.Assemble(skip_zero_entries);
595  S.Finalize(skip_zero_entries);
596 
597  reduced_oper = new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
598  if (!use_petsc)
599  {
600  HypreSmoother *J_hypreSmoother = new HypreSmoother;
601  J_hypreSmoother->SetType(HypreSmoother::l1Jacobi);
602  J_hypreSmoother->SetPositiveDiagonal(true);
603  J_prec = J_hypreSmoother;
604 
605  MINRESSolver *J_minres = new MINRESSolver(f.GetComm());
606  J_minres->SetRelTol(rel_tol);
607  J_minres->SetAbsTol(0.0);
608  J_minres->SetMaxIter(300);
609  J_minres->SetPrintLevel(-1);
610  J_minres->SetPreconditioner(*J_prec);
611  J_solver = J_minres;
612 
613  J_factory = NULL;
614 
615  newton_solver.iterative_mode = false;
616  newton_solver.SetSolver(*J_solver);
617  newton_solver.SetOperator(*reduced_oper);
618  newton_solver.SetPrintLevel(1); // print Newton iterations
619  newton_solver.SetRelTol(rel_tol);
620  newton_solver.SetAbsTol(0.0);
621  newton_solver.SetMaxIter(10);
622  }
623  else
624  {
625  // if using PETSc, we create the same solver (Newton + MINRES + Jacobi)
626  // by command line options (see rc_ex10p)
627  J_solver = NULL;
628  J_prec = NULL;
629  J_factory = NULL;
630  pnewton_solver = new PetscNonlinearSolver(f.GetComm(),
631  *reduced_oper);
632 
633  // we can setup a factory to construct a "physics-based" preconditioner
634  if (use_petsc_factory)
635  {
636  J_factory = new PreconditionerFactory(*reduced_oper, "JFNK preconditioner");
637  pnewton_solver->SetPreconditionerFactory(J_factory);
638  }
639  pnewton_solver->SetPrintLevel(1); // print Newton iterations
640  pnewton_solver->SetRelTol(rel_tol);
641  pnewton_solver->SetAbsTol(0.0);
642  pnewton_solver->SetMaxIter(10);
643  }
644 }
645 
646 void HyperelasticOperator::Mult(const Vector &vx, Vector &dvx_dt) const
647 {
648  // Create views to the sub-vectors v, x of vx, and dv_dt, dx_dt of dvx_dt
649  int sc = height/2;
650  Vector v(vx.GetData() + 0, sc);
651  Vector x(vx.GetData() + sc, sc);
652  Vector dv_dt(dvx_dt.GetData() + 0, sc);
653  Vector dx_dt(dvx_dt.GetData() + sc, sc);
654 
655  H.Mult(x, z);
656  if (viscosity != 0.0)
657  {
658  S.TrueAddMult(v, z);
659  z.SetSubVector(ess_tdof_list, 0.0);
660  }
661  z.Neg(); // z = -z
662  M_solver.Mult(z, dv_dt);
663 
664  dx_dt = v;
665 }
666 
667 void HyperelasticOperator::ImplicitSolve(const double dt,
668  const Vector &vx, Vector &dvx_dt)
669 {
670  int sc = height/2;
671  Vector v(vx.GetData() + 0, sc);
672  Vector x(vx.GetData() + sc, sc);
673  Vector dv_dt(dvx_dt.GetData() + 0, sc);
674  Vector dx_dt(dvx_dt.GetData() + sc, sc);
675 
676  // By eliminating kx from the coupled system:
677  // kv = -M^{-1}*[H(x + dt*kx) + S*(v + dt*kv)]
678  // kx = v + dt*kv
679  // we reduce it to a nonlinear equation for kv, represented by the
680  // reduced_oper. This equation is solved with the newton_solver
681  // object (using J_solver and J_prec internally).
682  reduced_oper->SetParameters(dt, &v, &x);
683  Vector zero; // empty vector is interpreted as zero r.h.s. by NewtonSolver
684  if (!pnewton_solver)
685  {
686  newton_solver.Mult(zero, dv_dt);
687  MFEM_VERIFY(newton_solver.GetConverged(),
688  "Newton solver did not converge.");
689  }
690  else
691  {
692  pnewton_solver->Mult(zero, dv_dt);
693  MFEM_VERIFY(pnewton_solver->GetConverged(),
694  "Newton solver did not converge.");
695  }
696  add(v, dt, dv_dt, dx_dt);
697 }
698 
699 double HyperelasticOperator::ElasticEnergy(const ParGridFunction &x) const
700 {
701  return H.GetEnergy(x);
702 }
703 
704 double HyperelasticOperator::KineticEnergy(const ParGridFunction &v) const
705 {
706  double loc_energy = 0.5*M.InnerProduct(v, v);
707  double energy;
708  MPI_Allreduce(&loc_energy, &energy, 1, MPI_DOUBLE, MPI_SUM,
709  fespace.GetComm());
710  return energy;
711 }
712 
713 void HyperelasticOperator::GetElasticEnergyDensity(
714  const ParGridFunction &x, ParGridFunction &w) const
715 {
716  ElasticEnergyCoefficient w_coeff(*model, x);
717  w.ProjectCoefficient(w_coeff);
718 }
719 
720 HyperelasticOperator::~HyperelasticOperator()
721 {
722  delete J_solver;
723  delete J_prec;
724  delete J_factory;
725  delete reduced_oper;
726  delete model;
727  delete Mmat;
728  delete pnewton_solver;
729 }
730 
731 // This method gets called every time we need a preconditioner "oh"
732 // contains the PetscParMatrix that wraps the operator constructed in
733 // the GetGradient() method (see also PetscSolver::SetJacobianType()).
734 // In this example, we just return a customizable PetscPreconditioner
735 // using that matrix. However, the OperatorHandle argument can be
736 // ignored, and any "physics-based" solver can be constructed since we
737 // have access to the HyperElasticOperator class.
738 Solver* PreconditionerFactory::NewPreconditioner(const mfem::OperatorHandle& oh)
739 {
740  PetscParMatrix *pP;
741  oh.Get(pP);
742  return new PetscPreconditioner(*pP,"jfnk_");
743 }
744 
746  const IntegrationPoint &ip)
747 {
748  model.SetTransformation(T);
749  x.GetVectorGradient(T, J);
750  // return model.EvalW(J); // in reference configuration
751  return model.EvalW(J)/J.Det(); // in deformed configuration
752 }
753 
754 
755 void InitialDeformation(const Vector &x, Vector &y)
756 {
757  // set the initial configuration to be the same as the reference,
758  // stress free, configuration
759  y = x;
760 }
761 
762 void InitialVelocity(const Vector &x, Vector &v)
763 {
764  const int dim = x.Size();
765  const double s = 0.1/64.;
766 
767  v = 0.0;
768  v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
769  v(0) = -s*x(0)*x(0);
770 }
void InitialDeformation(const Vector &x, Vector &y)
Definition: ex10p.cpp:673
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition: hypre.cpp:2276
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:605
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
int visport
Conjugate gradient method.
Definition: solvers.hpp:493
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
Abstract class for PETSc&#39;s preconditioners.
Definition: petsc.hpp:801
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:150
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition: pfespace.hpp:424
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
Wrapper for PETSc&#39;s matrix class.
Definition: petsc.hpp:315
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Compute the vector gradient with respect to the physical element variable.
Definition: gridfunc.cpp:1736
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
Base abstract class for first order time dependent operators.
Definition: operator.hpp:316
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Swap the internal node GridFunction pointer and ownership flag members with the given ones...
Definition: mesh.cpp:8351
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
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]...
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Parallel non-linear operator on the true dofs.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:23
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
Definition: handle.hpp:114
int main(int argc, char *argv[])
Definition: ex10p.cpp:163
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:561
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:686
STL namespace.
MINRES method.
Definition: solvers.hpp:603
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:411
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:317
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
Definition: hypre.hpp:1077
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2320
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:144
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:616
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
char vishost[]
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
Definition: gridfunc.cpp:221
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:641
Parallel smoothers in hypre.
Definition: hypre.hpp:971
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:281
MPI_Comm GetComm() const
Definition: pmesh.hpp:351
void Add(const int i, const int j, const double val)
Definition: sparsemat.cpp:2787
Abstract class for PETSc&#39;s nonlinear solvers.
Definition: petsc.hpp:903
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:206
A general vector function coefficient.
int GetMyRank() const
Definition: pmesh.hpp:353
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:555
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:163
void SetAbsTol(double atol)
Definition: solvers.hpp:200
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
void SetRelTol(double rtol)
Definition: solvers.hpp:199
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
void visualize(ostream &os, ParMesh *mesh, ParGridFunction *deformed_nodes, ParGridFunction *field, const char *field_name=NULL, bool init_vis=false)
Definition: ex10p.cpp:442
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:150
HYPRE_Int HYPRE_BigInt
void InitialVelocity(const Vector &x, Vector &v)
Definition: ex10p.cpp:680
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition: mesh.hpp:1023
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:909
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:424
void MFEMFinalizePetsc()
Definition: petsc.cpp:241
Class for integration point with weight.
Definition: intrules.hpp:31
int GetNRanks() const
Definition: pmesh.hpp:352
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Definition: petsc.cpp:201
double mu
Definition: ex25.cpp:140
int dim
Definition: ex24.cpp:53
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;.
RefCoord t[3]
Vector data type.
Definition: vector.hpp:58
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:259
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:8302
RefCoord s[3]
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4825
Base class for solvers.
Definition: operator.hpp:682
Class for parallel grid function.
Definition: pgridfunc.hpp:32
The classical forward Euler method.
Definition: ode.hpp:117
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
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:3429
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition: fe_coll.hpp:327
double f(const Vector &p)