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