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