MFEM  v4.3.0
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-wedge.mesh -s 2 -r 1 -o 2 -dt 3
11 // ex10 -m ../data/beam-quad.mesh -s 14 -r 2 -o 2 -dt 0.03 -vs 20
12 // ex10 -m ../data/beam-hex.mesh -s 14 -r 1 -o 2 -dt 0.05 -vs 20
13 // ex10 -m ../data/beam-quad-amr.mesh -s 3 -r 2 -o 2 -dt 3
14 //
15 // Description: This examples solves a time dependent nonlinear elasticity
16 // problem of the form dv/dt = H(x) + S v, dx/dt = v, where H is a
17 // hyperelastic model and S is a viscosity operator of Laplacian
18 // type. The geometry of the domain is assumed to be as follows:
19 //
20 // +---------------------+
21 // boundary --->| |
22 // attribute 1 | |
23 // (fixed) +---------------------+
24 //
25 // The example demonstrates the use of nonlinear operators (the
26 // class HyperelasticOperator defining H(x)), as well as their
27 // implicit time integration using a Newton method for solving an
28 // associated reduced backward-Euler type nonlinear equation
29 // (class ReducedSystemOperator). Each Newton step requires the
30 // inversion of a Jacobian matrix, which is done through a
31 // (preconditioned) inner solver. Note that implementing the
32 // method HyperelasticOperator::ImplicitSolve is the only
33 // requirement for high-order implicit (SDIRK) time integration.
34 //
35 // We recommend viewing examples 2 and 9 before viewing this
36 // example.
37 
38 #include "mfem.hpp"
39 #include <memory>
40 #include <iostream>
41 #include <fstream>
42 
43 using namespace std;
44 using namespace mfem;
45 
46 class ReducedSystemOperator;
47 
48 /** After spatial discretization, the hyperelastic model can be written as a
49  * system of ODEs:
50  * dv/dt = -M^{-1}*(H(x) + S*v)
51  * dx/dt = v,
52  * where x is the vector representing the deformation, v is the velocity field,
53  * M is the mass matrix, S is the viscosity matrix, and H(x) is the nonlinear
54  * hyperelastic operator.
55  *
56  * Class HyperelasticOperator represents the right-hand side of the above
57  * system of ODEs. */
58 class HyperelasticOperator : public TimeDependentOperator
59 {
60 protected:
61  FiniteElementSpace &fespace;
62 
63  BilinearForm M, S;
64  NonlinearForm H;
65  double viscosity;
66  HyperelasticModel *model;
67 
68  CGSolver M_solver; // Krylov solver for inverting the mass matrix M
69  DSmoother M_prec; // Preconditioner for the mass matrix M
70 
71  /** Nonlinear operator defining the reduced backward Euler equation for the
72  velocity. Used in the implementation of method ImplicitSolve. */
73  ReducedSystemOperator *reduced_oper;
74 
75  /// Newton solver for the reduced backward Euler equation
76  NewtonSolver newton_solver;
77 
78  /// Solver for the Jacobian solve in the Newton method
79  Solver *J_solver;
80  /// Preconditioner for the Jacobian solve in the Newton method
81  Solver *J_prec;
82 
83  mutable Vector z; // auxiliary vector
84 
85 public:
86  HyperelasticOperator(FiniteElementSpace &f, Array<int> &ess_bdr,
87  double visc, double mu, double K);
88 
89  /// Compute the right-hand side of the ODE system.
90  virtual void Mult(const Vector &vx, Vector &dvx_dt) const;
91  /** Solve the Backward-Euler equation: k = f(x + dt*k, t), for the unknown k.
92  This is the only requirement for high-order SDIRK implicit integration.*/
93  virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k);
94 
95  double ElasticEnergy(const Vector &x) const;
96  double KineticEnergy(const Vector &v) const;
97  void GetElasticEnergyDensity(const GridFunction &x, GridFunction &w) const;
98 
99  virtual ~HyperelasticOperator();
100 };
101 
102 /** Nonlinear operator of the form:
103  k --> (M + dt*S)*k + H(x + dt*v + dt^2*k) + S*v,
104  where M and S are given BilinearForms, H is a given NonlinearForm, v and x
105  are given vectors, and dt is a scalar. */
106 class ReducedSystemOperator : public Operator
107 {
108 private:
109  BilinearForm *M, *S;
110  NonlinearForm *H;
111  mutable SparseMatrix *Jacobian;
112  double dt;
113  const Vector *v, *x;
114  mutable Vector w, z;
115 
116 public:
117  ReducedSystemOperator(BilinearForm *M_, BilinearForm *S_, NonlinearForm *H_);
118 
119  /// Set current dt, v, x values - needed to compute action and Jacobian.
120  void SetParameters(double dt_, const Vector *v_, const Vector *x_);
121 
122  /// Compute y = H(x + dt (v + dt k)) + M k + S (v + dt k).
123  virtual void Mult(const Vector &k, Vector &y) const;
124 
125  /// Compute J = M + dt S + dt^2 grad_H(x + dt (v + dt k)).
126  virtual Operator &GetGradient(const Vector &k) const;
127 
128  virtual ~ReducedSystemOperator();
129 };
130 
131 
132 /** Function representing the elastic energy density for the given hyperelastic
133  model+deformation. Used in HyperelasticOperator::GetElasticEnergyDensity. */
134 class ElasticEnergyCoefficient : public Coefficient
135 {
136 private:
137  HyperelasticModel &model;
138  const GridFunction &x;
139  DenseMatrix J;
140 
141 public:
142  ElasticEnergyCoefficient(HyperelasticModel &m, const GridFunction &x_)
143  : model(m), x(x_) { }
144  virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
145  virtual ~ElasticEnergyCoefficient() { }
146 };
147 
148 void InitialDeformation(const Vector &x, Vector &y);
149 
150 void InitialVelocity(const Vector &x, Vector &v);
151 
152 void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes,
153  GridFunction *field, const char *field_name = NULL,
154  bool init_vis = false);
155 
156 
157 int main(int argc, char *argv[])
158 {
159  // 1. Parse command-line options.
160  const char *mesh_file = "../data/beam-quad.mesh";
161  int ref_levels = 2;
162  int order = 2;
163  int ode_solver_type = 3;
164  double t_final = 300.0;
165  double dt = 3.0;
166  double visc = 1e-2;
167  double mu = 0.25;
168  double K = 5.0;
169  bool visualization = true;
170  int vis_steps = 1;
171 
172  OptionsParser args(argc, argv);
173  args.AddOption(&mesh_file, "-m", "--mesh",
174  "Mesh file to use.");
175  args.AddOption(&ref_levels, "-r", "--refine",
176  "Number of times to refine the mesh uniformly.");
177  args.AddOption(&order, "-o", "--order",
178  "Order (degree) of the finite elements.");
179  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
180  "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
181  " 11 - Forward Euler, 12 - RK2,\n\t"
182  " 13 - RK3 SSP, 14 - RK4."
183  " 22 - Implicit Midpoint Method,\n\t"
184  " 23 - SDIRK23 (A-stable), 24 - SDIRK34");
185  args.AddOption(&t_final, "-tf", "--t-final",
186  "Final time; start time is 0.");
187  args.AddOption(&dt, "-dt", "--time-step",
188  "Time step.");
189  args.AddOption(&visc, "-v", "--viscosity",
190  "Viscosity coefficient.");
191  args.AddOption(&mu, "-mu", "--shear-modulus",
192  "Shear modulus in the Neo-Hookean hyperelastic model.");
193  args.AddOption(&K, "-K", "--bulk-modulus",
194  "Bulk modulus in the Neo-Hookean hyperelastic model.");
195  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
196  "--no-visualization",
197  "Enable or disable GLVis visualization.");
198  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
199  "Visualize every n-th timestep.");
200  args.Parse();
201  if (!args.Good())
202  {
203  args.PrintUsage(cout);
204  return 1;
205  }
206  args.PrintOptions(cout);
207 
208  // 2. Read the mesh from the given mesh file. We can handle triangular,
209  // quadrilateral, tetrahedral and hexahedral meshes with the same code.
210  Mesh *mesh = new Mesh(mesh_file, 1, 1);
211  int dim = mesh->Dimension();
212 
213  // 3. Define the ODE solver used for time integration. Several implicit
214  // singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
215  // explicit Runge-Kutta methods are available.
216  ODESolver *ode_solver;
217  switch (ode_solver_type)
218  {
219  // Implicit L-stable methods
220  case 1: ode_solver = new BackwardEulerSolver; break;
221  case 2: ode_solver = new SDIRK23Solver(2); break;
222  case 3: ode_solver = new SDIRK33Solver; break;
223  // Explicit methods
224  case 11: ode_solver = new ForwardEulerSolver; break;
225  case 12: ode_solver = new RK2Solver(0.5); break; // midpoint method
226  case 13: ode_solver = new RK3SSPSolver; break;
227  case 14: ode_solver = new RK4Solver; break;
228  case 15: ode_solver = new GeneralizedAlphaSolver(0.5); break;
229  // Implicit A-stable methods (not L-stable)
230  case 22: ode_solver = new ImplicitMidpointSolver; break;
231  case 23: ode_solver = new SDIRK23Solver; break;
232  case 24: ode_solver = new SDIRK34Solver; break;
233  default:
234  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
235  delete mesh;
236  return 3;
237  }
238 
239  // 4. Refine the mesh to increase the resolution. In this example we do
240  // 'ref_levels' of uniform refinement, where 'ref_levels' is a
241  // command-line parameter.
242  for (int lev = 0; lev < ref_levels; lev++)
243  {
244  mesh->UniformRefinement();
245  }
246 
247  // 5. Define the vector finite element spaces representing the mesh
248  // deformation x, the velocity v, and the initial configuration, x_ref.
249  // Define also the elastic energy density, w, which is in a discontinuous
250  // higher-order space. Since x and v are integrated in time as a system,
251  // we group them together in block vector vx, with offsets given by the
252  // fe_offset array.
253  H1_FECollection fe_coll(order, dim);
254  FiniteElementSpace fespace(mesh, &fe_coll, dim);
255 
256  int fe_size = fespace.GetTrueVSize();
257  cout << "Number of velocity/deformation unknowns: " << fe_size << endl;
258  Array<int> fe_offset(3);
259  fe_offset[0] = 0;
260  fe_offset[1] = fe_size;
261  fe_offset[2] = 2*fe_size;
262 
263  BlockVector vx(fe_offset);
264  GridFunction v, x;
265  v.MakeTRef(&fespace, vx.GetBlock(0), 0);
266  x.MakeTRef(&fespace, vx.GetBlock(1), 0);
267 
268  GridFunction x_ref(&fespace);
269  mesh->GetNodes(x_ref);
270 
271  L2_FECollection w_fec(order + 1, dim);
272  FiniteElementSpace w_fespace(mesh, &w_fec);
273  GridFunction w(&w_fespace);
274 
275  // 6. Set the initial conditions for v and x, and the boundary conditions on
276  // a beam-like mesh (see description above).
278  v.ProjectCoefficient(velo);
279  v.SetTrueVector();
281  x.ProjectCoefficient(deform);
282  x.SetTrueVector();
283 
284  Array<int> ess_bdr(fespace.GetMesh()->bdr_attributes.Max());
285  ess_bdr = 0;
286  ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed
287 
288  // 7. Initialize the hyperelastic operator, the GLVis visualization and print
289  // the initial energies.
290  HyperelasticOperator oper(fespace, ess_bdr, visc, mu, K);
291 
292  socketstream vis_v, vis_w;
293  if (visualization)
294  {
295  char vishost[] = "localhost";
296  int visport = 19916;
297  vis_v.open(vishost, visport);
298  vis_v.precision(8);
300  visualize(vis_v, mesh, &x, &v, "Velocity", true);
301  vis_w.open(vishost, visport);
302  if (vis_w)
303  {
304  oper.GetElasticEnergyDensity(x, w);
305  vis_w.precision(8);
306  visualize(vis_w, mesh, &x, &w, "Elastic energy density", true);
307  }
308  }
309 
310  double ee0 = oper.ElasticEnergy(x.GetTrueVector());
311  double ke0 = oper.KineticEnergy(v.GetTrueVector());
312  cout << "initial elastic energy (EE) = " << ee0 << endl;
313  cout << "initial kinetic energy (KE) = " << ke0 << endl;
314  cout << "initial total energy (TE) = " << (ee0 + ke0) << endl;
315 
316  double t = 0.0;
317  oper.SetTime(t);
318  ode_solver->Init(oper);
319 
320  // 8. Perform time-integration (looping over the time iterations, ti, with a
321  // time-step dt).
322  bool last_step = false;
323  for (int ti = 1; !last_step; ti++)
324  {
325  double dt_real = min(dt, t_final - t);
326 
327  ode_solver->Step(vx, t, dt_real);
328 
329  last_step = (t >= t_final - 1e-8*dt);
330 
331  if (last_step || (ti % vis_steps) == 0)
332  {
333  double ee = oper.ElasticEnergy(x.GetTrueVector());
334  double ke = oper.KineticEnergy(v.GetTrueVector());
335 
336  cout << "step " << ti << ", t = " << t << ", EE = " << ee << ", KE = "
337  << ke << ", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
338 
339  if (visualization)
340  {
342  visualize(vis_v, mesh, &x, &v);
343  if (vis_w)
344  {
345  oper.GetElasticEnergyDensity(x, w);
346  visualize(vis_w, mesh, &x, &w);
347  }
348  }
349  }
350  }
351 
352  // 9. Save the displaced mesh, the velocity and elastic energy.
353  {
355  GridFunction *nodes = &x;
356  int owns_nodes = 0;
357  mesh->SwapNodes(nodes, owns_nodes);
358  ofstream mesh_ofs("deformed.mesh");
359  mesh_ofs.precision(8);
360  mesh->Print(mesh_ofs);
361  mesh->SwapNodes(nodes, owns_nodes);
362  ofstream velo_ofs("velocity.sol");
363  velo_ofs.precision(8);
364  v.Save(velo_ofs);
365  ofstream ee_ofs("elastic_energy.sol");
366  ee_ofs.precision(8);
367  oper.GetElasticEnergyDensity(x, w);
368  w.Save(ee_ofs);
369  }
370 
371  // 10. Free the used memory.
372  delete ode_solver;
373  delete mesh;
374 
375  return 0;
376 }
377 
378 
379 void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes,
380  GridFunction *field, const char *field_name, bool init_vis)
381 {
382  if (!out)
383  {
384  return;
385  }
386 
387  GridFunction *nodes = deformed_nodes;
388  int owns_nodes = 0;
389 
390  mesh->SwapNodes(nodes, owns_nodes);
391 
392  out << "solution\n" << *mesh << *field;
393 
394  mesh->SwapNodes(nodes, owns_nodes);
395 
396  if (init_vis)
397  {
398  out << "window_size 800 800\n";
399  out << "window_title '" << field_name << "'\n";
400  if (mesh->SpaceDimension() == 2)
401  {
402  out << "view 0 0\n"; // view from top
403  out << "keys jl\n"; // turn off perspective and light
404  }
405  out << "keys cm\n"; // show colorbar and mesh
406  out << "autoscale value\n"; // update value-range; keep mesh-extents fixed
407  out << "pause\n";
408  }
409  out << flush;
410 }
411 
412 
413 ReducedSystemOperator::ReducedSystemOperator(
415  : Operator(M_->Height()), M(M_), S(S_), H(H_), Jacobian(NULL),
416  dt(0.0), v(NULL), x(NULL), w(height), z(height)
417 { }
418 
419 void ReducedSystemOperator::SetParameters(double dt_, const Vector *v_,
420  const Vector *x_)
421 {
422  dt = dt_; v = v_; x = x_;
423 }
424 
425 void ReducedSystemOperator::Mult(const Vector &k, Vector &y) const
426 {
427  // compute: y = H(x + dt*(v + dt*k)) + M*k + S*(v + dt*k)
428  add(*v, dt, k, w);
429  add(*x, dt, w, z);
430  H->Mult(z, y);
431  M->AddMult(k, y);
432  S->AddMult(w, y);
433 }
434 
435 Operator &ReducedSystemOperator::GetGradient(const Vector &k) const
436 {
437  delete Jacobian;
438  Jacobian = Add(1.0, M->SpMat(), dt, S->SpMat());
439  add(*v, dt, k, w);
440  add(*x, dt, w, z);
441  SparseMatrix *grad_H = dynamic_cast<SparseMatrix *>(&H->GetGradient(z));
442  Jacobian->Add(dt*dt, *grad_H);
443  return *Jacobian;
444 }
445 
446 ReducedSystemOperator::~ReducedSystemOperator()
447 {
448  delete Jacobian;
449 }
450 
451 
452 HyperelasticOperator::HyperelasticOperator(FiniteElementSpace &f,
453  Array<int> &ess_bdr, double visc,
454  double mu, double K)
455  : TimeDependentOperator(2*f.GetTrueVSize(), 0.0), fespace(f),
456  M(&fespace), S(&fespace), H(&fespace),
457  viscosity(visc), z(height/2)
458 {
459  const double rel_tol = 1e-8;
460  const int skip_zero_entries = 0;
461 
462  const double ref_density = 1.0; // density in the reference configuration
463  ConstantCoefficient rho0(ref_density);
464  M.AddDomainIntegrator(new VectorMassIntegrator(rho0));
465  M.Assemble(skip_zero_entries);
466  Array<int> ess_tdof_list;
467  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
468  SparseMatrix tmp;
469  M.FormSystemMatrix(ess_tdof_list, tmp);
470 
471  M_solver.iterative_mode = false;
472  M_solver.SetRelTol(rel_tol);
473  M_solver.SetAbsTol(0.0);
474  M_solver.SetMaxIter(30);
475  M_solver.SetPrintLevel(0);
476  M_solver.SetPreconditioner(M_prec);
477  M_solver.SetOperator(M.SpMat());
478 
479  model = new NeoHookeanModel(mu, K);
480  H.AddDomainIntegrator(new HyperelasticNLFIntegrator(model));
481  H.SetEssentialTrueDofs(ess_tdof_list);
482 
483  ConstantCoefficient visc_coeff(viscosity);
484  S.AddDomainIntegrator(new VectorDiffusionIntegrator(visc_coeff));
485  S.Assemble(skip_zero_entries);
486  S.FormSystemMatrix(ess_tdof_list, tmp);
487 
488  reduced_oper = new ReducedSystemOperator(&M, &S, &H);
489 
490 #ifndef MFEM_USE_SUITESPARSE
491  J_prec = new DSmoother(1);
492  MINRESSolver *J_minres = new MINRESSolver;
493  J_minres->SetRelTol(rel_tol);
494  J_minres->SetAbsTol(0.0);
495  J_minres->SetMaxIter(300);
496  J_minres->SetPrintLevel(-1);
497  J_minres->SetPreconditioner(*J_prec);
498  J_solver = J_minres;
499 #else
500  J_solver = new UMFPackSolver;
501  J_prec = NULL;
502 #endif
503 
504  newton_solver.iterative_mode = false;
505  newton_solver.SetSolver(*J_solver);
506  newton_solver.SetOperator(*reduced_oper);
507  newton_solver.SetPrintLevel(1); // print Newton iterations
508  newton_solver.SetRelTol(rel_tol);
509  newton_solver.SetAbsTol(0.0);
510  newton_solver.SetMaxIter(10);
511 }
512 
513 void HyperelasticOperator::Mult(const Vector &vx, Vector &dvx_dt) const
514 {
515  // Create views to the sub-vectors v, x of vx, and dv_dt, dx_dt of dvx_dt
516  int sc = height/2;
517  Vector v(vx.GetData() + 0, sc);
518  Vector x(vx.GetData() + sc, sc);
519  Vector dv_dt(dvx_dt.GetData() + 0, sc);
520  Vector dx_dt(dvx_dt.GetData() + sc, sc);
521 
522  H.Mult(x, z);
523  if (viscosity != 0.0)
524  {
525  S.AddMult(v, z);
526  }
527  z.Neg(); // z = -z
528  M_solver.Mult(z, dv_dt);
529 
530  dx_dt = v;
531 }
532 
533 void HyperelasticOperator::ImplicitSolve(const double dt,
534  const Vector &vx, Vector &dvx_dt)
535 {
536  int sc = height/2;
537  Vector v(vx.GetData() + 0, sc);
538  Vector x(vx.GetData() + sc, sc);
539  Vector dv_dt(dvx_dt.GetData() + 0, sc);
540  Vector dx_dt(dvx_dt.GetData() + sc, sc);
541 
542  // By eliminating kx from the coupled system:
543  // kv = -M^{-1}*[H(x + dt*kx) + S*(v + dt*kv)]
544  // kx = v + dt*kv
545  // we reduce it to a nonlinear equation for kv, represented by the
546  // reduced_oper. This equation is solved with the newton_solver
547  // object (using J_solver and J_prec internally).
548  reduced_oper->SetParameters(dt, &v, &x);
549  Vector zero; // empty vector is interpreted as zero r.h.s. by NewtonSolver
550  newton_solver.Mult(zero, dv_dt);
551  MFEM_VERIFY(newton_solver.GetConverged(), "Newton solver did not converge.");
552  add(v, dt, dv_dt, dx_dt);
553 }
554 
555 double HyperelasticOperator::ElasticEnergy(const Vector &x) const
556 {
557  return H.GetEnergy(x);
558 }
559 
560 double HyperelasticOperator::KineticEnergy(const Vector &v) const
561 {
562  return 0.5*M.InnerProduct(v, v);
563 }
564 
565 void HyperelasticOperator::GetElasticEnergyDensity(
566  const GridFunction &x, GridFunction &w) const
567 {
568  ElasticEnergyCoefficient w_coeff(*model, x);
569  w.ProjectCoefficient(w_coeff);
570 }
571 
572 HyperelasticOperator::~HyperelasticOperator()
573 {
574  delete J_solver;
575  delete J_prec;
576  delete reduced_oper;
577  delete model;
578 }
579 
580 
582  const IntegrationPoint &ip)
583 {
584  model.SetTransformation(T);
585  x.GetVectorGradient(T, J);
586  // return model.EvalW(J); // in reference configuration
587  return model.EvalW(J)/J.Det(); // in deformed configuration
588 }
589 
590 
591 void InitialDeformation(const Vector &x, Vector &y)
592 {
593  // set the initial configuration to be the same as the reference, stress
594  // free, configuration
595  y = x;
596 }
597 
598 void InitialVelocity(const Vector &x, Vector &v)
599 {
600  const int dim = x.Size();
601  const double s = 0.1/64.;
602 
603  v = 0.0;
604  v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
605  v(0) = -s*x(0)*x(0);
606 }
void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
Definition: ex10.cpp:379
double Eval(ElementTransformation &T, const IntegrationPoint &ip, double t)
Evaluate the coefficient in the element described by T at the point ip at time t. ...
Definition: coefficient.hpp:66
void Add(const int i, const int j, const double a)
Definition: sparsemat.cpp:2577
void InitialDeformation(const Vector &x, Vector &y)
Definition: ex10.cpp:591
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1387
Conjugate gradient method.
Definition: solvers.hpp:316
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
Data type for scaled Jacobi-type smoother of sparse matrix.
A class to handle Vectors in a block fashion.
Definition: blockvector.hpp:30
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:143
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
Base abstract class for first order time dependent operators.
Definition: operator.hpp:282
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Definition: mesh.cpp:7386
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
virtual void Step(Vector &x, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
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:652
MINRES method.
Definition: solvers.hpp:426
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:389
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:199
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:291
void InitialVelocity(const Vector &x, Vector &v)
Definition: ex10.cpp:598
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:839
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3484
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:1931
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:137
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:439
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
Data type sparse matrix.
Definition: sparsemat.hpp:41
constexpr char vishost[]
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:410
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
constexpr int visport
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
Definition: gridfunc.cpp:219
void SetMaxIter(int max_it)
Definition: solvers.hpp:100
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
Newton&#39;s method for solving F(x)=b for a given operator F.
Definition: solvers.hpp:464
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:543
int Dimension() const
Definition: mesh.hpp:911
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition: gridfunc.hpp:127
A general vector function coefficient.
int SpaceDimension() const
Definition: mesh.hpp:912
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
void SetAbsTol(double atol)
Definition: solvers.hpp:99
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:204
void SetRelTol(double rtol)
Definition: solvers.hpp:98
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
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)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:149
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Definition: gridfunc.cpp:1661
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:402
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
Class for integration point with weight.
Definition: intrules.hpp:25
double mu
Definition: ex25.cpp:139
int dim
Definition: ex24.cpp:53
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: gridfunc.cpp:2278
Abstract class for hyperelastic models.
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
RefCoord t[3]
Vector data type.
Definition: vector.hpp:60
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7343
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
RefCoord s[3]
Base class for solvers.
Definition: operator.hpp:648
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:66
The classical forward Euler method.
Definition: ode.hpp:116
Abstract operator.
Definition: operator.hpp:24
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
int main()
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:90
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:285
double f(const Vector &p)
void Neg()
(*this) = -(*this)
Definition: vector.cpp:283
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150