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