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