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