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