MFEM  v4.3.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 //
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 &out, ParMesh *mesh, ParGridFunction *deformed_nodes,
158  ParGridFunction *field, const char *field_name = NULL,
159  bool init_vis = false);
160 
161 
162 int main(int argc, char *argv[])
163 {
164  // 1. Initialize MPI.
165  int num_procs, myid;
166  MPI_Init(&argc, &argv);
167  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
168  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
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  MPI_Finalize();
226  return 1;
227  }
228  if (myid == 0)
229  {
230  args.PrintOptions(cout);
231  }
232 
233  // 3. Read the serial mesh from the given mesh file on all processors. We can
234  // handle triangular, quadrilateral, tetrahedral and hexahedral meshes
235  // with the same code.
236  Mesh *mesh = new Mesh(mesh_file, 1, 1);
237  int dim = mesh->Dimension();
238 
239  // 4. Define the ODE solver used for time integration. Several implicit
240  // singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
241  // explicit Runge-Kutta methods are available.
242  ODESolver *ode_solver;
243  switch (ode_solver_type)
244  {
245  // Implicit L-stable methods
246  case 1: ode_solver = new BackwardEulerSolver; break;
247  case 2: ode_solver = new SDIRK23Solver(2); break;
248  case 3: ode_solver = new SDIRK33Solver; break;
249  // Explicit methods
250  case 11: ode_solver = new ForwardEulerSolver; break;
251  case 12: ode_solver = new RK2Solver(0.5); break; // midpoint method
252  case 13: ode_solver = new RK3SSPSolver; break;
253  case 14: ode_solver = new RK4Solver; break;
254  case 15: ode_solver = new GeneralizedAlphaSolver(0.5); break;
255  // Implicit A-stable methods (not L-stable)
256  case 22: ode_solver = new ImplicitMidpointSolver; break;
257  case 23: ode_solver = new SDIRK23Solver; break;
258  case 24: ode_solver = new SDIRK34Solver; break;
259  default:
260  if (myid == 0)
261  {
262  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
263  }
264  delete mesh;
265  MPI_Finalize();
266  return 3;
267  }
268 
269  // 5. Refine the mesh in serial to increase the resolution. In this example
270  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
271  // a command-line parameter.
272  for (int lev = 0; lev < ser_ref_levels; lev++)
273  {
274  mesh->UniformRefinement();
275  }
276 
277  // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
278  // this mesh further in parallel to increase the resolution. Once the
279  // parallel mesh is defined, the serial mesh can be deleted.
280  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
281  delete mesh;
282  for (int lev = 0; lev < par_ref_levels; lev++)
283  {
284  pmesh->UniformRefinement();
285  }
286 
287  // 7. Define the parallel vector finite element spaces representing the mesh
288  // deformation x_gf, the velocity v_gf, and the initial configuration,
289  // x_ref. Define also the elastic energy density, w_gf, which is in a
290  // discontinuous higher-order space. Since x and v are integrated in time
291  // as a system, we group them together in block vector vx, on the unique
292  // parallel degrees of freedom, with offsets given by array true_offset.
293  H1_FECollection fe_coll(order, dim);
294  ParFiniteElementSpace fespace(pmesh, &fe_coll, dim);
295 
296  HYPRE_BigInt glob_size = fespace.GlobalTrueVSize();
297  if (myid == 0)
298  {
299  cout << "Number of velocity/deformation unknowns: " << glob_size << endl;
300  }
301  int true_size = fespace.TrueVSize();
302  Array<int> true_offset(3);
303  true_offset[0] = 0;
304  true_offset[1] = true_size;
305  true_offset[2] = 2*true_size;
306 
307  BlockVector vx(true_offset);
308  ParGridFunction v_gf, x_gf;
309  v_gf.MakeTRef(&fespace, vx, true_offset[0]);
310  x_gf.MakeTRef(&fespace, vx, true_offset[1]);
311 
312  ParGridFunction x_ref(&fespace);
313  pmesh->GetNodes(x_ref);
314 
315  L2_FECollection w_fec(order + 1, dim);
316  ParFiniteElementSpace w_fespace(pmesh, &w_fec);
317  ParGridFunction w_gf(&w_fespace);
318 
319  // 8. Set the initial conditions for v_gf, x_gf and vx, and define the
320  // boundary conditions on a beam-like mesh (see description above).
322  v_gf.ProjectCoefficient(velo);
323  v_gf.SetTrueVector();
325  x_gf.ProjectCoefficient(deform);
326  x_gf.SetTrueVector();
327 
328  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
329 
330  Array<int> ess_bdr(fespace.GetMesh()->bdr_attributes.Max());
331  ess_bdr = 0;
332  ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed
333 
334  // 9. Initialize the hyperelastic operator, the GLVis visualization and print
335  // the initial energies.
336  HyperelasticOperator oper(fespace, ess_bdr, visc, mu, K);
337 
338  socketstream vis_v, vis_w;
339  if (visualization)
340  {
341  char vishost[] = "localhost";
342  int visport = 19916;
343  vis_v.open(vishost, visport);
344  vis_v.precision(8);
345  visualize(vis_v, pmesh, &x_gf, &v_gf, "Velocity", true);
346  // Make sure all ranks have sent their 'v' solution before initiating
347  // another set of GLVis connections (one from each rank):
348  MPI_Barrier(pmesh->GetComm());
349  vis_w.open(vishost, visport);
350  if (vis_w)
351  {
352  oper.GetElasticEnergyDensity(x_gf, w_gf);
353  vis_w.precision(8);
354  visualize(vis_w, pmesh, &x_gf, &w_gf, "Elastic energy density", true);
355  }
356  }
357 
358  double ee0 = oper.ElasticEnergy(x_gf);
359  double ke0 = oper.KineticEnergy(v_gf);
360  if (myid == 0)
361  {
362  cout << "initial elastic energy (EE) = " << ee0 << endl;
363  cout << "initial kinetic energy (KE) = " << ke0 << endl;
364  cout << "initial total energy (TE) = " << (ee0 + ke0) << endl;
365  }
366 
367  double t = 0.0;
368  oper.SetTime(t);
369  ode_solver->Init(oper);
370 
371  // 10. Perform time-integration
372  // (looping over the time iterations, ti, with a time-step dt).
373  bool last_step = false;
374  for (int ti = 1; !last_step; ti++)
375  {
376  double dt_real = min(dt, t_final - t);
377 
378  ode_solver->Step(vx, t, dt_real);
379 
380  last_step = (t >= t_final - 1e-8*dt);
381 
382  if (last_step || (ti % vis_steps) == 0)
383  {
384  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
385 
386  double ee = oper.ElasticEnergy(x_gf);
387  double ke = oper.KineticEnergy(v_gf);
388 
389  if (myid == 0)
390  {
391  cout << "step " << ti << ", t = " << t << ", EE = " << ee
392  << ", KE = " << ke << ", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
393  }
394 
395  if (visualization)
396  {
397  visualize(vis_v, pmesh, &x_gf, &v_gf);
398  if (vis_w)
399  {
400  oper.GetElasticEnergyDensity(x_gf, w_gf);
401  visualize(vis_w, pmesh, &x_gf, &w_gf);
402  }
403  }
404  }
405  }
406 
407  // 11. Save the displaced mesh, the velocity and elastic energy.
408  {
409  v_gf.SetFromTrueVector(); x_gf.SetFromTrueVector();
410  GridFunction *nodes = &x_gf;
411  int owns_nodes = 0;
412  pmesh->SwapNodes(nodes, owns_nodes);
413 
414  ostringstream mesh_name, velo_name, ee_name;
415  mesh_name << "deformed." << setfill('0') << setw(6) << myid;
416  velo_name << "velocity." << setfill('0') << setw(6) << myid;
417  ee_name << "elastic_energy." << setfill('0') << setw(6) << myid;
418 
419  ofstream mesh_ofs(mesh_name.str().c_str());
420  mesh_ofs.precision(8);
421  pmesh->Print(mesh_ofs);
422  pmesh->SwapNodes(nodes, owns_nodes);
423  ofstream velo_ofs(velo_name.str().c_str());
424  velo_ofs.precision(8);
425  v_gf.Save(velo_ofs);
426  ofstream ee_ofs(ee_name.str().c_str());
427  ee_ofs.precision(8);
428  oper.GetElasticEnergyDensity(x_gf, w_gf);
429  w_gf.Save(ee_ofs);
430  }
431 
432  // 12. Free the used memory.
433  delete ode_solver;
434  delete pmesh;
435 
436  MPI_Finalize();
437 
438  return 0;
439 }
440 
441 void visualize(ostream &out, ParMesh *mesh, ParGridFunction *deformed_nodes,
442  ParGridFunction *field, const char *field_name, bool init_vis)
443 {
444  if (!out)
445  {
446  return;
447  }
448 
449  GridFunction *nodes = deformed_nodes;
450  int owns_nodes = 0;
451 
452  mesh->SwapNodes(nodes, owns_nodes);
453 
454  out << "parallel " << mesh->GetNRanks() << " " << mesh->GetMyRank() << "\n";
455  out << "solution\n" << *mesh << *field;
456 
457  mesh->SwapNodes(nodes, owns_nodes);
458 
459  if (init_vis)
460  {
461  out << "window_size 800 800\n";
462  out << "window_title '" << field_name << "'\n";
463  if (mesh->SpaceDimension() == 2)
464  {
465  out << "view 0 0\n"; // view from top
466  out << "keys jl\n"; // turn off perspective and light
467  }
468  out << "keys cm\n"; // show colorbar and mesh
469  out << "autoscale value\n"; // update value-range; keep mesh-extents fixed
470  out << "pause\n";
471  }
472  out << flush;
473 }
474 
475 
476 ReducedSystemOperator::ReducedSystemOperator(
478  const Array<int> &ess_tdof_list_)
479  : Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
480  Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height),
481  ess_tdof_list(ess_tdof_list_)
482 { }
483 
484 void ReducedSystemOperator::SetParameters(double dt_, const Vector *v_,
485  const Vector *x_)
486 {
487  dt = dt_; v = v_; x = x_;
488 }
489 
490 void ReducedSystemOperator::Mult(const Vector &k, Vector &y) const
491 {
492  // compute: y = H(x + dt*(v + dt*k)) + M*k + S*(v + dt*k)
493  add(*v, dt, k, w);
494  add(*x, dt, w, z);
495  H->Mult(z, y);
496  M->TrueAddMult(k, y);
497  S->TrueAddMult(w, y);
498  y.SetSubVector(ess_tdof_list, 0.0);
499 }
500 
501 Operator &ReducedSystemOperator::GetGradient(const Vector &k) const
502 {
503  delete Jacobian;
504  SparseMatrix *localJ = Add(1.0, M->SpMat(), dt, S->SpMat());
505  add(*v, dt, k, w);
506  add(*x, dt, w, z);
507  localJ->Add(dt*dt, H->GetLocalGradient(z));
508  Jacobian = M->ParallelAssemble(localJ);
509  delete localJ;
510  HypreParMatrix *Je = Jacobian->EliminateRowsCols(ess_tdof_list);
511  delete Je;
512  return *Jacobian;
513 }
514 
515 ReducedSystemOperator::~ReducedSystemOperator()
516 {
517  delete Jacobian;
518 }
519 
520 
521 HyperelasticOperator::HyperelasticOperator(ParFiniteElementSpace &f,
522  Array<int> &ess_bdr, double visc,
523  double mu, double K)
524  : TimeDependentOperator(2*f.TrueVSize(), 0.0), fespace(f),
525  M(&fespace), S(&fespace), H(&fespace),
526  viscosity(visc), M_solver(f.GetComm()), newton_solver(f.GetComm()),
527  z(height/2)
528 {
529  const double rel_tol = 1e-8;
530  const int skip_zero_entries = 0;
531 
532  const double ref_density = 1.0; // density in the reference configuration
533  ConstantCoefficient rho0(ref_density);
534  M.AddDomainIntegrator(new VectorMassIntegrator(rho0));
535  M.Assemble(skip_zero_entries);
536  M.Finalize(skip_zero_entries);
537  Mmat = M.ParallelAssemble();
538  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
539  HypreParMatrix *Me = Mmat->EliminateRowsCols(ess_tdof_list);
540  delete Me;
541 
542  M_solver.iterative_mode = false;
543  M_solver.SetRelTol(rel_tol);
544  M_solver.SetAbsTol(0.0);
545  M_solver.SetMaxIter(30);
546  M_solver.SetPrintLevel(0);
547  M_prec.SetType(HypreSmoother::Jacobi);
548  M_solver.SetPreconditioner(M_prec);
549  M_solver.SetOperator(*Mmat);
550 
551  model = new NeoHookeanModel(mu, K);
552  H.AddDomainIntegrator(new HyperelasticNLFIntegrator(model));
553  H.SetEssentialTrueDofs(ess_tdof_list);
554 
555  ConstantCoefficient visc_coeff(viscosity);
556  S.AddDomainIntegrator(new VectorDiffusionIntegrator(visc_coeff));
557  S.Assemble(skip_zero_entries);
558  S.Finalize(skip_zero_entries);
559 
560  reduced_oper = new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
561 
562  HypreSmoother *J_hypreSmoother = new HypreSmoother;
563  J_hypreSmoother->SetType(HypreSmoother::l1Jacobi);
564  J_hypreSmoother->SetPositiveDiagonal(true);
565  J_prec = J_hypreSmoother;
566 
567  MINRESSolver *J_minres = new MINRESSolver(f.GetComm());
568  J_minres->SetRelTol(rel_tol);
569  J_minres->SetAbsTol(0.0);
570  J_minres->SetMaxIter(300);
571  J_minres->SetPrintLevel(-1);
572  J_minres->SetPreconditioner(*J_prec);
573  J_solver = J_minres;
574 
575  newton_solver.iterative_mode = false;
576  newton_solver.SetSolver(*J_solver);
577  newton_solver.SetOperator(*reduced_oper);
578  newton_solver.SetPrintLevel(1); // print Newton iterations
579  newton_solver.SetRelTol(rel_tol);
580  newton_solver.SetAbsTol(0.0);
581  newton_solver.SetAdaptiveLinRtol(2, 0.5, 0.9);
582  newton_solver.SetMaxIter(10);
583 }
584 
585 void HyperelasticOperator::Mult(const Vector &vx, Vector &dvx_dt) const
586 {
587  // Create views to the sub-vectors v, x of vx, and dv_dt, dx_dt of dvx_dt
588  int sc = height/2;
589  Vector v(vx.GetData() + 0, sc);
590  Vector x(vx.GetData() + sc, sc);
591  Vector dv_dt(dvx_dt.GetData() + 0, sc);
592  Vector dx_dt(dvx_dt.GetData() + sc, sc);
593 
594  H.Mult(x, z);
595  if (viscosity != 0.0)
596  {
597  S.TrueAddMult(v, z);
598  z.SetSubVector(ess_tdof_list, 0.0);
599  }
600  z.Neg(); // z = -z
601  M_solver.Mult(z, dv_dt);
602 
603  dx_dt = v;
604 }
605 
606 void HyperelasticOperator::ImplicitSolve(const double dt,
607  const Vector &vx, Vector &dvx_dt)
608 {
609  int sc = height/2;
610  Vector v(vx.GetData() + 0, sc);
611  Vector x(vx.GetData() + sc, sc);
612  Vector dv_dt(dvx_dt.GetData() + 0, sc);
613  Vector dx_dt(dvx_dt.GetData() + sc, sc);
614 
615  // By eliminating kx from the coupled system:
616  // kv = -M^{-1}*[H(x + dt*kx) + S*(v + dt*kv)]
617  // kx = v + dt*kv
618  // we reduce it to a nonlinear equation for kv, represented by the
619  // reduced_oper. This equation is solved with the newton_solver
620  // object (using J_solver and J_prec internally).
621  reduced_oper->SetParameters(dt, &v, &x);
622  Vector zero; // empty vector is interpreted as zero r.h.s. by NewtonSolver
623  newton_solver.Mult(zero, dv_dt);
624  MFEM_VERIFY(newton_solver.GetConverged(), "Newton solver did not converge.");
625  add(v, dt, dv_dt, dx_dt);
626 }
627 
628 double HyperelasticOperator::ElasticEnergy(const ParGridFunction &x) const
629 {
630  return H.GetEnergy(x);
631 }
632 
633 double HyperelasticOperator::KineticEnergy(const ParGridFunction &v) const
634 {
635  double loc_energy = 0.5*M.InnerProduct(v, v);
636  double energy;
637  MPI_Allreduce(&loc_energy, &energy, 1, MPI_DOUBLE, MPI_SUM,
638  fespace.GetComm());
639  return energy;
640 }
641 
642 void HyperelasticOperator::GetElasticEnergyDensity(
643  const ParGridFunction &x, ParGridFunction &w) const
644 {
645  ElasticEnergyCoefficient w_coeff(*model, x);
646  w.ProjectCoefficient(w_coeff);
647 }
648 
649 HyperelasticOperator::~HyperelasticOperator()
650 {
651  delete J_solver;
652  delete J_prec;
653  delete reduced_oper;
654  delete model;
655  delete Mmat;
656 }
657 
658 
660  const IntegrationPoint &ip)
661 {
662  model.SetTransformation(T);
663  x.GetVectorGradient(T, J);
664  // return model.EvalW(J); // in reference configuration
665  return model.EvalW(J)/J.Det(); // in deformed configuration
666 }
667 
668 
669 void InitialDeformation(const Vector &x, Vector &y)
670 {
671  // set the initial configuration to be the same as the reference, stress
672  // free, configuration
673  y = x;
674 }
675 
676 void InitialVelocity(const Vector &x, Vector &v)
677 {
678  const int dim = x.Size();
679  const double s = 0.1/64.;
680 
681  v = 0.0;
682  v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
683  v(0) = -s*x(0)*x(0);
684 }
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:2052
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:551
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:2577
void InitialDeformation(const Vector &x, Vector &y)
Definition: ex10.cpp:591
Conjugate gradient method.
Definition: solvers.hpp:316
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:143
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition: pfespace.hpp:422
Base abstract class for first order time dependent operators.
Definition: operator.hpp:282
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:7386
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
virtual void Step(Vector &x, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:275
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
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:841
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:493
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:652
MINRES method.
Definition: solvers.hpp:426
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:199
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:291
int GetNRanks() const
Definition: pmesh.hpp:277
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:946
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:1931
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:137
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:439
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
MPI_Comm GetComm() const
Definition: pfespace.hpp:263
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:41
constexpr char vishost[]
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:410
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
constexpr int visport
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
Definition: gridfunc.cpp:219
void SetMaxIter(int max_it)
Definition: solvers.hpp:100
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:464
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4382
Parallel smoothers in hypre.
Definition: hypre.hpp:840
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
A general vector function coefficient.
int SpaceDimension() const
Definition: mesh.hpp:912
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
void SetAbsTol(double atol)
Definition: solvers.hpp:99
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:204
int GetMyRank() const
Definition: pmesh.hpp:278
void SetRelTol(double rtol)
Definition: solvers.hpp:98
MPI_Comm GetComm() const
Definition: pmesh.hpp:276
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
HYPRE_Int HYPRE_BigInt
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Definition: gridfunc.cpp:1661
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:402
Class for integration point with weight.
Definition: intrules.hpp:25
double mu
Definition: ex25.cpp:139
int dim
Definition: ex24.cpp:53
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
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:60
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7343
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
RefCoord s[3]
Base class for solvers.
Definition: operator.hpp:648
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:277
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:3028
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
int main()
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:285
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150