MFEM  v3.3
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex10.cpp
Go to the documentation of this file.
1 // MFEM Example 10
2 //
3 // Compile with: make ex10
4 //
5 // Sample runs:
6 // ex10 -m ../data/beam-quad.mesh -s 3 -r 2 -o 2 -dt 3
7 // ex10 -m ../data/beam-tri.mesh -s 3 -r 2 -o 2 -dt 3
8 // ex10 -m ../data/beam-hex.mesh -s 2 -r 1 -o 2 -dt 3
9 // ex10 -m ../data/beam-tet.mesh -s 2 -r 1 -o 2 -dt 3
10 // ex10 -m ../data/beam-quad.mesh -s 14 -r 2 -o 2 -dt 0.03 -vs 20
11 // ex10 -m ../data/beam-hex.mesh -s 14 -r 1 -o 2 -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 ReducedSystemOperator). 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 #include "mfem.hpp"
37 #include <memory>
38 #include <iostream>
39 #include <fstream>
40 
41 using namespace std;
42 using namespace mfem;
43 
44 class ReducedSystemOperator;
45 
46 /** After spatial discretization, the hyperelastic model can be written as a
47  * system of ODEs:
48  * dv/dt = -M^{-1}*(H(x) + S*v)
49  * dx/dt = v,
50  * where x is the vector representing the deformation, v is the velocity field,
51  * M is the mass matrix, S is the viscosity matrix, and H(x) is the nonlinear
52  * hyperelastic operator.
53  *
54  * Class HyperelasticOperator represents the right-hand side of the above
55  * system of ODEs. */
56 class HyperelasticOperator : public TimeDependentOperator
57 {
58 protected:
59  FiniteElementSpace &fespace;
60 
61  BilinearForm M, S;
62  NonlinearForm H;
63  double viscosity;
64  HyperelasticModel *model;
65 
66  CGSolver M_solver; // Krylov solver for inverting the mass matrix M
67  DSmoother M_prec; // Preconditioner for the mass matrix M
68 
69  /** Nonlinear operator defining the reduced backward Euler equation for the
70  velocity. Used in the implementation of method ImplicitSolve. */
71  ReducedSystemOperator *reduced_oper;
72 
73  /// Newton solver for the reduced backward Euler equation
74  NewtonSolver newton_solver;
75 
76  /// Solver for the Jacobian solve in the Newton method
77  Solver *J_solver;
78  /// Preconditioner for the Jacobian solve in the Newton method
79  Solver *J_prec;
80 
81  mutable Vector z; // auxiliary vector
82 
83 public:
84  HyperelasticOperator(FiniteElementSpace &f, Array<int> &ess_bdr,
85  double visc, double mu, double K);
86 
87  /// Compute the right-hand side of the ODE system.
88  virtual void Mult(const Vector &vx, Vector &dvx_dt) const;
89  /** Solve the Backward-Euler equation: k = f(x + dt*k, t), for the unknown k.
90  This is the only requirement for high-order SDIRK implicit integration.*/
91  virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k);
92 
93  double ElasticEnergy(Vector &x) const;
94  double KineticEnergy(Vector &v) const;
95  void GetElasticEnergyDensity(GridFunction &x, GridFunction &w) const;
96 
97  virtual ~HyperelasticOperator();
98 };
99 
100 /** Nonlinear operator of the form:
101  k --> (M + dt*S)*k + H(x + dt*v + dt^2*k) + S*v,
102  where M and S are given BilinearForms, H is a given NonlinearForm, v and x
103  are given vectors, and dt is a scalar. */
104 class ReducedSystemOperator : public Operator
105 {
106 private:
107  BilinearForm *M, *S;
108  NonlinearForm *H;
109  mutable SparseMatrix *Jacobian;
110  double dt;
111  const Vector *v, *x;
112  mutable Vector w, z;
113 
114 public:
115  ReducedSystemOperator(BilinearForm *M_, BilinearForm *S_, NonlinearForm *H_);
116 
117  /// Set current dt, v, x values - needed to compute action and Jacobian.
118  void SetParameters(double dt_, const Vector *v_, const Vector *x_);
119 
120  /// Compute y = H(x + dt (v + dt k)) + M k + S (v + dt k).
121  virtual void Mult(const Vector &k, Vector &y) const;
122 
123  /// Compute J = M + dt S + dt^2 grad_H(x + dt (v + dt k)).
124  virtual Operator &GetGradient(const Vector &k) const;
125 
126  virtual ~ReducedSystemOperator();
127 };
128 
129 
130 /** Function representing the elastic energy density for the given hyperelastic
131  model+deformation. Used in HyperelasticOperator::GetElasticEnergyDensity. */
132 class ElasticEnergyCoefficient : public Coefficient
133 {
134 private:
135  HyperelasticModel &model;
136  GridFunction &x;
137  DenseMatrix J;
138 
139 public:
140  ElasticEnergyCoefficient(HyperelasticModel &m, GridFunction &x_)
141  : model(m), x(x_) { }
142  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
143  virtual ~ElasticEnergyCoefficient() { }
144 };
145 
146 void InitialDeformation(const Vector &x, Vector &y);
147 
148 void InitialVelocity(const Vector &x, Vector &v);
149 
150 void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes,
151  GridFunction *field, const char *field_name = NULL,
152  bool init_vis = false);
153 
154 
155 int main(int argc, char *argv[])
156 {
157  // 1. Parse command-line options.
158  const char *mesh_file = "../data/beam-quad.mesh";
159  int ref_levels = 2;
160  int order = 2;
161  int ode_solver_type = 3;
162  double t_final = 300.0;
163  double dt = 3.0;
164  double visc = 1e-2;
165  double mu = 0.25;
166  double K = 5.0;
167  bool visualization = true;
168  int vis_steps = 1;
169 
170  OptionsParser args(argc, argv);
171  args.AddOption(&mesh_file, "-m", "--mesh",
172  "Mesh file to use.");
173  args.AddOption(&ref_levels, "-r", "--refine",
174  "Number of times to refine the mesh uniformly.");
175  args.AddOption(&order, "-o", "--order",
176  "Order (degree) of the finite elements.");
177  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
178  "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
179  " 11 - Forward Euler, 12 - RK2,\n\t"
180  " 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  args.PrintUsage(cout);
200  return 1;
201  }
202  args.PrintOptions(cout);
203 
204  // 2. Read the mesh from the given mesh file. We can handle triangular,
205  // quadrilateral, tetrahedral and hexahedral meshes with the same code.
206  Mesh *mesh = new Mesh(mesh_file, 1, 1);
207  int dim = mesh->Dimension();
208 
209  // 3. Define the ODE solver used for time integration. Several implicit
210  // singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
211  // explicit Runge-Kutta methods are available.
212  ODESolver *ode_solver;
213  switch (ode_solver_type)
214  {
215  // Implicit L-stable methods
216  case 1: ode_solver = new BackwardEulerSolver; break;
217  case 2: ode_solver = new SDIRK23Solver(2); break;
218  case 3: ode_solver = new SDIRK33Solver; break;
219  // Explicit methods
220  case 11: ode_solver = new ForwardEulerSolver; break;
221  case 12: ode_solver = new RK2Solver(0.5); break; // midpoint method
222  case 13: ode_solver = new RK3SSPSolver; break;
223  case 14: ode_solver = new RK4Solver; break;
224  // Implicit A-stable methods (not L-stable)
225  case 22: ode_solver = new ImplicitMidpointSolver; break;
226  case 23: ode_solver = new SDIRK23Solver; break;
227  case 24: ode_solver = new SDIRK34Solver; break;
228  default:
229  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
230  return 3;
231  }
232 
233  // 4. Refine the mesh to increase the resolution. In this example we do
234  // 'ref_levels' of uniform refinement, where 'ref_levels' is a
235  // command-line parameter.
236  for (int lev = 0; lev < ref_levels; lev++)
237  {
238  mesh->UniformRefinement();
239  }
240 
241  // 5. Define the vector finite element spaces representing the mesh
242  // deformation x, the velocity v, and the initial configuration, x_ref.
243  // Define also the elastic energy density, w, which is in a discontinuous
244  // higher-order space. Since x and v are integrated in time as a system,
245  // we group them together in block vector vx, with offsets given by the
246  // fe_offset array.
247  H1_FECollection fe_coll(order, dim);
248  FiniteElementSpace fespace(mesh, &fe_coll, dim);
249 
250  int fe_size = fespace.GetVSize();
251  cout << "Number of velocity/deformation unknowns: " << fe_size << endl;
252  Array<int> fe_offset(3);
253  fe_offset[0] = 0;
254  fe_offset[1] = fe_size;
255  fe_offset[2] = 2*fe_size;
256 
257  BlockVector vx(fe_offset);
258  GridFunction v, x;
259  v.MakeRef(&fespace, vx.GetBlock(0), 0);
260  x.MakeRef(&fespace, vx.GetBlock(1), 0);
261 
262  GridFunction x_ref(&fespace);
263  mesh->GetNodes(x_ref);
264 
265  L2_FECollection w_fec(order + 1, dim);
266  FiniteElementSpace w_fespace(mesh, &w_fec);
267  GridFunction w(&w_fespace);
268 
269  // 6. Set the initial conditions for v and x, and the boundary conditions on
270  // a beam-like mesh (see description above).
272  v.ProjectCoefficient(velo);
274  x.ProjectCoefficient(deform);
275 
276  Array<int> ess_bdr(fespace.GetMesh()->bdr_attributes.Max());
277  ess_bdr = 0;
278  ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed
279 
280  // 7. Initialize the hyperelastic operator, the GLVis visualization and print
281  // the initial energies.
282  HyperelasticOperator oper(fespace, ess_bdr, visc, mu, K);
283 
284  socketstream vis_v, vis_w;
285  if (visualization)
286  {
287  char vishost[] = "localhost";
288  int visport = 19916;
289  vis_v.open(vishost, visport);
290  vis_v.precision(8);
291  visualize(vis_v, mesh, &x, &v, "Velocity", true);
292  vis_w.open(vishost, visport);
293  if (vis_w)
294  {
295  oper.GetElasticEnergyDensity(x, w);
296  vis_w.precision(8);
297  visualize(vis_w, mesh, &x, &w, "Elastic energy density", true);
298  }
299  }
300 
301  double ee0 = oper.ElasticEnergy(x);
302  double ke0 = oper.KineticEnergy(v);
303  cout << "initial elastic energy (EE) = " << ee0 << endl;
304  cout << "initial kinetic energy (KE) = " << ke0 << endl;
305  cout << "initial total energy (TE) = " << (ee0 + ke0) << endl;
306 
307  double t = 0.0;
308  oper.SetTime(t);
309  ode_solver->Init(oper);
310 
311  // 8. Perform time-integration (looping over the time iterations, ti, with a
312  // time-step dt).
313  bool last_step = false;
314  for (int ti = 1; !last_step; ti++)
315  {
316  double dt_real = min(dt, t_final - t);
317 
318  ode_solver->Step(vx, t, dt_real);
319 
320  last_step = (t >= t_final - 1e-8*dt);
321 
322  if (last_step || (ti % vis_steps) == 0)
323  {
324  double ee = oper.ElasticEnergy(x);
325  double ke = oper.KineticEnergy(v);
326 
327  cout << "step " << ti << ", t = " << t << ", EE = " << ee << ", KE = "
328  << ke << ", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
329 
330  if (visualization)
331  {
332  visualize(vis_v, mesh, &x, &v);
333  if (vis_w)
334  {
335  oper.GetElasticEnergyDensity(x, w);
336  visualize(vis_w, mesh, &x, &w);
337  }
338  }
339  }
340  }
341 
342  // 9. Save the displaced mesh, the velocity and elastic energy.
343  {
344  GridFunction *nodes = &x;
345  int owns_nodes = 0;
346  mesh->SwapNodes(nodes, owns_nodes);
347  ofstream mesh_ofs("deformed.mesh");
348  mesh_ofs.precision(8);
349  mesh->Print(mesh_ofs);
350  mesh->SwapNodes(nodes, owns_nodes);
351  ofstream velo_ofs("velocity.sol");
352  velo_ofs.precision(8);
353  v.Save(velo_ofs);
354  ofstream ee_ofs("elastic_energy.sol");
355  ee_ofs.precision(8);
356  oper.GetElasticEnergyDensity(x, w);
357  w.Save(ee_ofs);
358  }
359 
360  // 10. Free the used memory.
361  delete ode_solver;
362  delete mesh;
363 
364  return 0;
365 }
366 
367 
368 void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes,
369  GridFunction *field, const char *field_name, bool init_vis)
370 {
371  if (!out)
372  {
373  return;
374  }
375 
376  GridFunction *nodes = deformed_nodes;
377  int owns_nodes = 0;
378 
379  mesh->SwapNodes(nodes, owns_nodes);
380 
381  out << "solution\n" << *mesh << *field;
382 
383  mesh->SwapNodes(nodes, owns_nodes);
384 
385  if (init_vis)
386  {
387  out << "window_size 800 800\n";
388  out << "window_title '" << field_name << "'\n";
389  if (mesh->SpaceDimension() == 2)
390  {
391  out << "view 0 0\n"; // view from top
392  out << "keys jl\n"; // turn off perspective and light
393  }
394  out << "keys cm\n"; // show colorbar and mesh
395  out << "autoscale value\n"; // update value-range; keep mesh-extents fixed
396  out << "pause\n";
397  }
398  out << flush;
399 }
400 
401 
402 ReducedSystemOperator::ReducedSystemOperator(
404  : Operator(M_->Height()), M(M_), S(S_), H(H_), Jacobian(NULL),
405  dt(0.0), v(NULL), x(NULL), w(height), z(height)
406 { }
407 
408 void ReducedSystemOperator::SetParameters(double dt_, const Vector *v_,
409  const Vector *x_)
410 {
411  dt = dt_; v = v_; x = x_;
412 }
413 
414 void ReducedSystemOperator::Mult(const Vector &k, Vector &y) const
415 {
416  // compute: y = H(x + dt*(v + dt*k)) + M*k + S*(v + dt*k)
417  add(*v, dt, k, w);
418  add(*x, dt, w, z);
419  H->Mult(z, y);
420  M->AddMult(k, y);
421  S->AddMult(w, y);
422 }
423 
424 Operator &ReducedSystemOperator::GetGradient(const Vector &k) const
425 {
426  delete Jacobian;
427  Jacobian = Add(1.0, M->SpMat(), dt, S->SpMat());
428  add(*v, dt, k, w);
429  add(*x, dt, w, z);
430  SparseMatrix *grad_H = dynamic_cast<SparseMatrix *>(&H->GetGradient(z));
431  Jacobian->Add(dt*dt, *grad_H);
432  return *Jacobian;
433 }
434 
435 ReducedSystemOperator::~ReducedSystemOperator()
436 {
437  delete Jacobian;
438 }
439 
440 
441 HyperelasticOperator::HyperelasticOperator(FiniteElementSpace &f,
442  Array<int> &ess_bdr, double visc,
443  double mu, double K)
444  : TimeDependentOperator(2*f.GetVSize(), 0.0), fespace(f),
445  M(&fespace), S(&fespace), H(&fespace),
446  viscosity(visc), z(height/2)
447 {
448  const double rel_tol = 1e-8;
449  const int skip_zero_entries = 0;
450 
451  const double ref_density = 1.0; // density in the reference configuration
452  ConstantCoefficient rho0(ref_density);
453  M.AddDomainIntegrator(new VectorMassIntegrator(rho0));
454  M.Assemble(skip_zero_entries);
455  M.EliminateEssentialBC(ess_bdr);
456  M.Finalize(skip_zero_entries);
457 
458  M_solver.iterative_mode = false;
459  M_solver.SetRelTol(rel_tol);
460  M_solver.SetAbsTol(0.0);
461  M_solver.SetMaxIter(30);
462  M_solver.SetPrintLevel(0);
463  M_solver.SetPreconditioner(M_prec);
464  M_solver.SetOperator(M.SpMat());
465 
466  model = new NeoHookeanModel(mu, K);
467  H.AddDomainIntegrator(new HyperelasticNLFIntegrator(model));
468  H.SetEssentialBC(ess_bdr);
469 
470  ConstantCoefficient visc_coeff(viscosity);
471  S.AddDomainIntegrator(new VectorDiffusionIntegrator(visc_coeff));
472  S.Assemble(skip_zero_entries);
473  S.EliminateEssentialBC(ess_bdr);
474  S.Finalize(skip_zero_entries);
475 
476  reduced_oper = new ReducedSystemOperator(&M, &S, &H);
477 
478 #ifndef MFEM_USE_SUITESPARSE
479  J_prec = new DSmoother(1);
480  MINRESSolver *J_minres = new MINRESSolver;
481  J_minres->SetRelTol(rel_tol);
482  J_minres->SetAbsTol(0.0);
483  J_minres->SetMaxIter(300);
484  J_minres->SetPrintLevel(-1);
485  J_minres->SetPreconditioner(*J_prec);
486  J_solver = J_minres;
487 #else
488  J_solver = new UMFPackSolver;
489  J_prec = NULL;
490 #endif
491 
492  newton_solver.iterative_mode = false;
493  newton_solver.SetSolver(*J_solver);
494  newton_solver.SetOperator(*reduced_oper);
495  newton_solver.SetPrintLevel(1); // print Newton iterations
496  newton_solver.SetRelTol(rel_tol);
497  newton_solver.SetAbsTol(0.0);
498  newton_solver.SetMaxIter(10);
499 }
500 
501 void HyperelasticOperator::Mult(const Vector &vx, Vector &dvx_dt) const
502 {
503  // Create views to the sub-vectors v, x of vx, and dv_dt, dx_dt of dvx_dt
504  int sc = height/2;
505  Vector v(vx.GetData() + 0, sc);
506  Vector x(vx.GetData() + sc, sc);
507  Vector dv_dt(dvx_dt.GetData() + 0, sc);
508  Vector dx_dt(dvx_dt.GetData() + sc, sc);
509 
510  H.Mult(x, z);
511  if (viscosity != 0.0)
512  {
513  S.AddMult(v, z);
514  }
515  z.Neg(); // z = -z
516  M_solver.Mult(z, dv_dt);
517 
518  dx_dt = v;
519 }
520 
521 void HyperelasticOperator::ImplicitSolve(const double dt,
522  const Vector &vx, Vector &dvx_dt)
523 {
524  int sc = height/2;
525  Vector v(vx.GetData() + 0, sc);
526  Vector x(vx.GetData() + sc, sc);
527  Vector dv_dt(dvx_dt.GetData() + 0, sc);
528  Vector dx_dt(dvx_dt.GetData() + sc, sc);
529 
530  // By eliminating kx from the coupled system:
531  // kv = -M^{-1}*[H(x + dt*kx) + S*(v + dt*kv)]
532  // kx = v + dt*kv
533  // we reduce it to a nonlinear equation for kv, represented by the
534  // reduced_oper. This equation is solved with the newton_solver
535  // object (using J_solver and J_prec internally).
536  reduced_oper->SetParameters(dt, &v, &x);
537  Vector zero; // empty vector is interpreted as zero r.h.s. by NewtonSolver
538  newton_solver.Mult(zero, dv_dt);
539  MFEM_VERIFY(newton_solver.GetConverged(), "Newton solver did not converge.");
540  add(v, dt, dv_dt, dx_dt);
541 }
542 
543 double HyperelasticOperator::ElasticEnergy(Vector &x) const
544 {
545  return H.GetEnergy(x);
546 }
547 
548 double HyperelasticOperator::KineticEnergy(Vector &v) const
549 {
550  return 0.5*M.InnerProduct(v, v);
551 }
552 
553 void HyperelasticOperator::GetElasticEnergyDensity(
554  GridFunction &x, GridFunction &w) const
555 {
556  ElasticEnergyCoefficient w_coeff(*model, x);
557  w.ProjectCoefficient(w_coeff);
558 }
559 
560 HyperelasticOperator::~HyperelasticOperator()
561 {
562  delete J_solver;
563  delete J_prec;
564  delete reduced_oper;
565  delete model;
566 }
567 
568 
570  const IntegrationPoint &ip)
571 {
572  model.SetTransformation(T);
573  x.GetVectorGradient(T, J);
574  // return model.EvalW(J); // in reference configuration
575  return model.EvalW(J)/J.Det(); // in deformed configuration
576 }
577 
578 
579 void InitialDeformation(const Vector &x, Vector &y)
580 {
581  // set the initial configuration to be the same as the reference, stress
582  // free, configuration
583  y = x;
584 }
585 
586 void InitialVelocity(const Vector &x, Vector &v)
587 {
588  const int dim = x.Size();
589  const double s = 0.1/64.;
590 
591  v = 0.0;
592  v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
593  v(0) = -s*x(0)*x(0);
594 }
void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
Definition: ex10.cpp:368
int GetVSize() const
Definition: fespace.hpp:163
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
Data type for scaled Jacobi-type smoother of sparse matrix.
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:140
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:5406
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:22
int Size() const
Returns the size of the vector.
Definition: vector.hpp:106
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:261
MINRES method.
Definition: solvers.hpp:220
Hyperelastic integrator for any given HyperelasticModel.
int main(int argc, char *argv[])
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:212
double * GetData() const
Definition: vector.hpp:114
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:264
void InitialVelocity(const Vector &x, Vector &v)
Definition: ex10.cpp:586
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:337
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2185
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2807
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:71
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:6684
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
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:258
void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
Definition: gridfunc.cpp:179
int Dimension() const
Definition: mesh.hpp:611
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
int SpaceDimension() const
Definition: mesh.hpp:612
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:142
void SetRelTol(double rtol)
Definition: solvers.hpp:61
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
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
void ProjectCoefficient(Coefficient &coeff)
Definition: gridfunc.cpp:1232
Abstract class for hyperelastic models.
Definition: nonlininteg.hpp:49
int open(const char hostname[], int port)
Vector data type.
Definition: vector.hpp:36
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5368
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:146
Base class for solvers.
Definition: operator.hpp:257
The classical forward Euler method.
Definition: ode.hpp:101
Abstract operator.
Definition: operator.hpp:21
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad)
Definition: gridfunc.cpp:1008
Vector & GetBlock(int i)
Get the i-th vector in the block.
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:195
void Neg()
(*this) = -(*this)
Definition: vector.cpp:256
virtual void Print(std::ostream &out=std::cout) const
Definition: mesh.hpp:988
bool Good() const
Definition: optparser.hpp:120