MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
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
43using namespace std;
44using namespace mfem;
45
46class 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. */
58class HyperelasticOperator : public TimeDependentOperator
59{
60protected:
61 FiniteElementSpace &fespace;
62
63 BilinearForm M, S;
65 real_t 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
85public:
86 HyperelasticOperator(FiniteElementSpace &f, Array<int> &ess_bdr,
87 real_t visc, real_t mu, real_t 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 real_t dt, const Vector &x, Vector &k);
94
95 real_t ElasticEnergy(const Vector &x) const;
96 real_t 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. */
106class ReducedSystemOperator : public Operator
107{
108private:
109 BilinearForm *M, *S;
110 NonlinearForm *H;
111 mutable SparseMatrix *Jacobian;
112 real_t dt;
113 const Vector *v, *x;
114 mutable Vector w, z;
115
116public:
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(real_t 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. */
134class ElasticEnergyCoefficient : public Coefficient
135{
136private:
137 HyperelasticModel &model;
138 const GridFunction &x;
139 DenseMatrix J;
140
141public:
142 ElasticEnergyCoefficient(HyperelasticModel &m, const GridFunction &x_)
143 : model(m), x(x_) { }
144 virtual real_t Eval(ElementTransformation &T, const IntegrationPoint &ip);
145 virtual ~ElasticEnergyCoefficient() { }
146};
147
148void InitialDeformation(const Vector &x, Vector &y);
149
150void InitialVelocity(const Vector &x, Vector &v);
151
152void visualize(ostream &os, Mesh *mesh, GridFunction *deformed_nodes,
153 GridFunction *field, const char *field_name = NULL,
154 bool init_vis = false);
155
156
157int 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 real_t t_final = 300.0;
165 real_t dt = 3.0;
166 real_t visc = 1e-2;
167 real_t mu = 0.25;
168 real_t 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 cout << "GLVis visualization paused."
309 << " Press space (in the GLVis window) to resume it.\n";
310 }
311
312 real_t ee0 = oper.ElasticEnergy(x.GetTrueVector());
313 real_t ke0 = oper.KineticEnergy(v.GetTrueVector());
314 cout << "initial elastic energy (EE) = " << ee0 << endl;
315 cout << "initial kinetic energy (KE) = " << ke0 << endl;
316 cout << "initial total energy (TE) = " << (ee0 + ke0) << endl;
317
318 real_t t = 0.0;
319 oper.SetTime(t);
320 ode_solver->Init(oper);
321
322 // 8. Perform time-integration (looping over the time iterations, ti, with a
323 // time-step dt).
324 bool last_step = false;
325 for (int ti = 1; !last_step; ti++)
326 {
327 real_t dt_real = min(dt, t_final - t);
328
329 ode_solver->Step(vx, t, dt_real);
330
331 last_step = (t >= t_final - 1e-8*dt);
332
333 if (last_step || (ti % vis_steps) == 0)
334 {
335 real_t ee = oper.ElasticEnergy(x.GetTrueVector());
336 real_t ke = oper.KineticEnergy(v.GetTrueVector());
337
338 cout << "step " << ti << ", t = " << t << ", EE = " << ee << ", KE = "
339 << ke << ", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
340
341 if (visualization)
342 {
344 visualize(vis_v, mesh, &x, &v);
345 if (vis_w)
346 {
347 oper.GetElasticEnergyDensity(x, w);
348 visualize(vis_w, mesh, &x, &w);
349 }
350 }
351 }
352 }
353
354 // 9. Save the displaced mesh, the velocity and elastic energy.
355 {
357 GridFunction *nodes = &x;
358 int owns_nodes = 0;
359 mesh->SwapNodes(nodes, owns_nodes);
360 ofstream mesh_ofs("deformed.mesh");
361 mesh_ofs.precision(8);
362 mesh->Print(mesh_ofs);
363 mesh->SwapNodes(nodes, owns_nodes);
364 ofstream velo_ofs("velocity.sol");
365 velo_ofs.precision(8);
366 v.Save(velo_ofs);
367 ofstream ee_ofs("elastic_energy.sol");
368 ee_ofs.precision(8);
369 oper.GetElasticEnergyDensity(x, w);
370 w.Save(ee_ofs);
371 }
372
373 // 10. Free the used memory.
374 delete ode_solver;
375 delete mesh;
376
377 return 0;
378}
379
380
381void visualize(ostream &os, Mesh *mesh, GridFunction *deformed_nodes,
382 GridFunction *field, const char *field_name, bool init_vis)
383{
384 if (!os)
385 {
386 return;
387 }
388
389 GridFunction *nodes = deformed_nodes;
390 int owns_nodes = 0;
391
392 mesh->SwapNodes(nodes, owns_nodes);
393
394 os << "solution\n" << *mesh << *field;
395
396 mesh->SwapNodes(nodes, owns_nodes);
397
398 if (init_vis)
399 {
400 os << "window_size 800 800\n";
401 os << "window_title '" << field_name << "'\n";
402 if (mesh->SpaceDimension() == 2)
403 {
404 os << "view 0 0\n"; // view from top
405 os << "keys jl\n"; // turn off perspective and light
406 }
407 os << "keys cm\n"; // show colorbar and mesh
408 // update value-range; keep mesh-extents fixed
409 os << "autoscale value\n";
410 os << "pause\n";
411 }
412 os << flush;
413}
414
415
416ReducedSystemOperator::ReducedSystemOperator(
418 : Operator(M_->Height()), M(M_), S(S_), H(H_), Jacobian(NULL),
419 dt(0.0), v(NULL), x(NULL), w(height), z(height)
420{ }
421
422void ReducedSystemOperator::SetParameters(real_t dt_, const Vector *v_,
423 const Vector *x_)
424{
425 dt = dt_; v = v_; x = x_;
426}
427
428void ReducedSystemOperator::Mult(const Vector &k, Vector &y) const
429{
430 // compute: y = H(x + dt*(v + dt*k)) + M*k + S*(v + dt*k)
431 add(*v, dt, k, w);
432 add(*x, dt, w, z);
433 H->Mult(z, y);
434 M->AddMult(k, y);
435 S->AddMult(w, y);
436}
437
438Operator &ReducedSystemOperator::GetGradient(const Vector &k) const
439{
440 delete Jacobian;
441 Jacobian = Add(1.0, M->SpMat(), dt, S->SpMat());
442 add(*v, dt, k, w);
443 add(*x, dt, w, z);
444 SparseMatrix *grad_H = dynamic_cast<SparseMatrix *>(&H->GetGradient(z));
445 Jacobian->Add(dt*dt, *grad_H);
446 return *Jacobian;
447}
448
449ReducedSystemOperator::~ReducedSystemOperator()
450{
451 delete Jacobian;
452}
453
454
455HyperelasticOperator::HyperelasticOperator(FiniteElementSpace &f,
456 Array<int> &ess_bdr, real_t visc,
457 real_t mu, real_t K)
458 : TimeDependentOperator(2*f.GetTrueVSize(), (real_t) 0.0), fespace(f),
459 M(&fespace), S(&fespace), H(&fespace),
460 viscosity(visc), z(height/2)
461{
462#if defined(MFEM_USE_DOUBLE)
463 const real_t rel_tol = 1e-8;
464 const real_t newton_abs_tol = 0.0;
465#elif defined(MFEM_USE_SINGLE)
466 const real_t rel_tol = 1e-3;
467 const real_t newton_abs_tol = 1e-4;
468#else
469#error "Only single and double precision are supported!"
470 const real_t rel_tol = real_t(1);
471 const real_t newton_abs_tol = real_t(0);
472#endif
473 const int skip_zero_entries = 0;
474
475 const real_t ref_density = 1.0; // density in the reference configuration
476 ConstantCoefficient rho0(ref_density);
477 M.AddDomainIntegrator(new VectorMassIntegrator(rho0));
478 M.Assemble(skip_zero_entries);
479 Array<int> ess_tdof_list;
480 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
481 SparseMatrix tmp;
482 M.FormSystemMatrix(ess_tdof_list, tmp);
483
484 M_solver.iterative_mode = false;
485 M_solver.SetRelTol(rel_tol);
486 M_solver.SetAbsTol(0.0);
487 M_solver.SetMaxIter(30);
488 M_solver.SetPrintLevel(0);
489 M_solver.SetPreconditioner(M_prec);
490 M_solver.SetOperator(M.SpMat());
491
492 model = new NeoHookeanModel(mu, K);
493 H.AddDomainIntegrator(new HyperelasticNLFIntegrator(model));
494 H.SetEssentialTrueDofs(ess_tdof_list);
495
496 ConstantCoefficient visc_coeff(viscosity);
497 S.AddDomainIntegrator(new VectorDiffusionIntegrator(visc_coeff));
498 S.Assemble(skip_zero_entries);
499 S.FormSystemMatrix(ess_tdof_list, tmp);
500
501 reduced_oper = new ReducedSystemOperator(&M, &S, &H);
502
503#ifndef MFEM_USE_SUITESPARSE
504 J_prec = new DSmoother(1);
505 MINRESSolver *J_minres = new MINRESSolver;
506 J_minres->SetRelTol(rel_tol);
507 J_minres->SetAbsTol(0.0);
508 J_minres->SetMaxIter(300);
509 J_minres->SetPrintLevel(-1);
510 J_minres->SetPreconditioner(*J_prec);
511 J_solver = J_minres;
512#else
513 J_solver = new UMFPackSolver;
514 J_prec = NULL;
515#endif
516
517 newton_solver.iterative_mode = false;
518 newton_solver.SetSolver(*J_solver);
519 newton_solver.SetOperator(*reduced_oper);
520 newton_solver.SetPrintLevel(1); // print Newton iterations
521 newton_solver.SetRelTol(rel_tol);
522 newton_solver.SetAbsTol(newton_abs_tol);
523 newton_solver.SetMaxIter(10);
524}
525
526void HyperelasticOperator::Mult(const Vector &vx, Vector &dvx_dt) const
527{
528 // Create views to the sub-vectors v, x of vx, and dv_dt, dx_dt of dvx_dt
529 int sc = height/2;
530 Vector v(vx.GetData() + 0, sc);
531 Vector x(vx.GetData() + sc, sc);
532 Vector dv_dt(dvx_dt.GetData() + 0, sc);
533 Vector dx_dt(dvx_dt.GetData() + sc, sc);
534
535 H.Mult(x, z);
536 if (viscosity != 0.0)
537 {
538 S.AddMult(v, z);
539 }
540 z.Neg(); // z = -z
541 M_solver.Mult(z, dv_dt);
542
543 dx_dt = v;
544}
545
546void HyperelasticOperator::ImplicitSolve(const real_t dt,
547 const Vector &vx, Vector &dvx_dt)
548{
549 int sc = height/2;
550 Vector v(vx.GetData() + 0, sc);
551 Vector x(vx.GetData() + sc, sc);
552 Vector dv_dt(dvx_dt.GetData() + 0, sc);
553 Vector dx_dt(dvx_dt.GetData() + sc, sc);
554
555 // By eliminating kx from the coupled system:
556 // kv = -M^{-1}*[H(x + dt*kx) + S*(v + dt*kv)]
557 // kx = v + dt*kv
558 // we reduce it to a nonlinear equation for kv, represented by the
559 // reduced_oper. This equation is solved with the newton_solver
560 // object (using J_solver and J_prec internally).
561 reduced_oper->SetParameters(dt, &v, &x);
562 Vector zero; // empty vector is interpreted as zero r.h.s. by NewtonSolver
563 newton_solver.Mult(zero, dv_dt);
564 MFEM_VERIFY(newton_solver.GetConverged(), "Newton solver did not converge.");
565 add(v, dt, dv_dt, dx_dt);
566}
567
568real_t HyperelasticOperator::ElasticEnergy(const Vector &x) const
569{
570 return H.GetEnergy(x);
571}
572
573real_t HyperelasticOperator::KineticEnergy(const Vector &v) const
574{
575 return 0.5*M.InnerProduct(v, v);
576}
577
578void HyperelasticOperator::GetElasticEnergyDensity(
579 const GridFunction &x, GridFunction &w) const
580{
581 ElasticEnergyCoefficient w_coeff(*model, x);
582 w.ProjectCoefficient(w_coeff);
583}
584
585HyperelasticOperator::~HyperelasticOperator()
586{
587 delete J_solver;
588 delete J_prec;
589 delete reduced_oper;
590 delete model;
591}
592
593
594real_t ElasticEnergyCoefficient::Eval(ElementTransformation &T,
595 const IntegrationPoint &ip)
596{
597 model.SetTransformation(T);
598 x.GetVectorGradient(T, J);
599 // return model.EvalW(J); // in reference configuration
600 return model.EvalW(J)/J.Det(); // in deformed configuration
601}
602
603
605{
606 // set the initial configuration to be the same as the reference, stress
607 // free, configuration
608 y = x;
609}
610
611void InitialVelocity(const Vector &x, Vector &v)
612{
613 const int dim = x.Size();
614 const real_t s = 0.1/64.;
615
616 v = 0.0;
617 v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
618 v(0) = -s*x(0)*x(0);
619}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
Backward Euler ODE solver. L-stable.
Definition ode.hpp:412
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
real_t InnerProduct(const Vector &x, const Vector &y) const
Compute .
virtual void AddMult(const Vector &x, Vector &y, const real_t a=1.0) const
Add the matrix vector multiple to a vector: .
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix: .
A class to handle Vectors in a block fashion.
Vector & GetBlock(int i)
Get the i-th vector in the block.
Conjugate gradient method.
Definition solvers.hpp:513
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:718
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
A coefficient that is constant across space and time.
Data type for scaled Jacobi-type smoother of sparse matrix.
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
real_t Det() const
Definition densemat.cpp:535
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:716
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
The classical forward Euler method.
Definition ode.hpp:118
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Compute the vector gradient with respect to the physical element variable.
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition gridfunc.hpp:144
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
void MakeTRef(FiniteElementSpace *f, real_t *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
Definition gridfunc.cpp:221
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Definition gridfunc.hpp:150
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Definition gridfunc.hpp:130
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
Abstract class for hyperelastic models.
virtual real_t EvalW(const DenseMatrix &Jpt) const =0
Evaluate the strain energy density function, W = W(Jpt).
void SetTransformation(ElementTransformation &Ttr_)
Implicit midpoint method. A-stable, not L-stable.
Definition ode.hpp:425
Class for integration point with weight.
Definition intrules.hpp:35
void SetRelTol(real_t rtol)
Definition solvers.hpp:209
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:71
void SetMaxIter(int max_it)
Definition solvers.hpp:211
bool GetConverged() const
Returns true if the last call to Mult() converged successfully.
Definition solvers.hpp:262
void SetAbsTol(real_t atol)
Definition solvers.hpp:210
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition fe_coll.hpp:330
MINRES method.
Definition solvers.hpp:628
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.hpp:640
Mesh data type.
Definition mesh.hpp:56
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Definition mesh.hpp:2288
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
int SpaceDimension() const
Dimension of the physical space containing the mesh.
Definition mesh.hpp:1163
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:8973
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Swap the internal node GridFunction pointer and ownership flag members with the given ones.
Definition mesh.cpp:9022
Newton's method for solving F(x)=b for a given operator F.
Definition solvers.hpp:667
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Definition solvers.cpp:1821
virtual void Mult(const Vector &x, Vector &y) const
Evaluate the action of the NonlinearForm.
virtual real_t GetEnergy(const Vector &x) const
Compute the energy corresponding to the state x.
virtual Operator & GetGradient(const Vector &x) const
Compute the gradient Operator of the NonlinearForm corresponding to the state x.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition ode.hpp:24
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:18
virtual void Step(Vector &x, real_t &t, real_t &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Abstract operator.
Definition operator.hpp:25
int height
Dimension of the output / number of rows in the matrix.
Definition operator.hpp:27
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition ode.hpp:151
The classical explicit forth-order Runge-Kutta method, RK4.
Definition ode.hpp:164
Base class for solvers.
Definition operator.hpp:683
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition operator.hpp:686
Data type sparse matrix.
Definition sparsemat.hpp:51
void Add(const int i, const int j, const real_t val)
Base abstract class for first order time dependent operators.
Definition operator.hpp:317
virtual void SetTime(const real_t t_)
Set the current time.
Definition operator.hpp:360
Direct sparse solver using UMFPACK.
Definition solvers.hpp:1096
A general vector function coefficient.
Vector data type.
Definition vector.hpp:80
void Neg()
(*this) = -(*this)
Definition vector.cpp:300
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:227
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
void visualize(ostream &os, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
Definition ex10.cpp:381
void InitialDeformation(const Vector &x, Vector &y)
Definition ex10.cpp:604
void InitialVelocity(const Vector &x, Vector &v)
Definition ex10.cpp:611
int dim
Definition ex24.cpp:53
real_t mu
Definition ex25.cpp:140
int main()
const int visport
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition vector.cpp:316
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
void Add(const DenseMatrix &A, const DenseMatrix &B, real_t alpha, DenseMatrix &C)
C = A + alpha*B.
const char vishost[]
RefCoord t[3]
RefCoord s[3]