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