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