MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex10p.cpp
Go to the documentation of this file.
1// MFEM Example 10 - Parallel Version
2// SUNDIALS Modification
3//
4// Compile with: make ex10p
5//
6// Sample runs:
7// mpirun -np 4 ex10p -m ../../data/beam-quad.mesh -rp 1 -o 2 -s 12 -dt 0.15 -vs 10
8// mpirun -np 4 ex10p -m ../../data/beam-tri.mesh -rp 1 -o 2 -s 16 -dt 0.25 -vs 10
9// mpirun -np 4 ex10p -m ../../data/beam-hex.mesh -rp 0 -o 2 -s 12 -dt 0.15 -vs 10
10// mpirun -np 4 ex10p -m ../../data/beam-tri.mesh -rp 1 -o 2 -s 2 -dt 3 -nls kinsol
11// mpirun -np 4 ex10p -m ../../data/beam-quad.mesh -rp 1 -o 2 -s 2 -dt 3 -nls kinsol
12// mpirun -np 4 ex10p -m ../../data/beam-hex.mesh -rs 1 -o 2 -s 2 -dt 3 -nls kinsol
13// mpirun -np 4 ex10p -m ../../data/beam-quad.mesh -rp 1 -o 2 -s 14 -dt 0.15 -vs 10
14// mpirun -np 4 ex10p -m ../../data/beam-tri.mesh -rp 1 -o 2 -s 17 -dt 5e-3 -vs 60
15// mpirun -np 4 ex10p -m ../../data/beam-hex.mesh -rp 0 -o 2 -s 14 -dt 0.15 -vs 10
16// mpirun -np 4 ex10p -m ../../data/beam-quad-amr.mesh -rp 1 -o 2 -s 12 -dt 0.15 -vs 10
17//
18// Description: This examples solves a time dependent nonlinear elasticity
19// problem of the form dv/dt = H(x) + S v, dx/dt = v, where H is a
20// hyperelastic model and S is a viscosity operator of Laplacian
21// type. The geometry of the domain is assumed to be as follows:
22//
23// +---------------------+
24// boundary --->| |
25// attribute 1 | |
26// (fixed) +---------------------+
27//
28// The example demonstrates the use of nonlinear operators (the
29// class HyperelasticOperator defining H(x)), as well as their
30// implicit time integration using a Newton method for solving an
31// associated reduced backward-Euler type nonlinear equation
32// (class ReducedSystemOperator). Each Newton step requires the
33// inversion of a Jacobian matrix, which is done through a
34// (preconditioned) inner solver. Note that implementing the
35// method HyperelasticOperator::ImplicitSolve is the only
36// requirement for high-order implicit (SDIRK) time integration.
37//
38// We recommend viewing examples 2 and 9 before viewing this
39// example.
40
41#include "mfem.hpp"
42#include <memory>
43#include <iostream>
44#include <fstream>
45#include <string>
46#include <map>
47
48#ifndef MFEM_USE_SUNDIALS
49#error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES
50#endif
51
52using namespace std;
53using namespace mfem;
54
55class ReducedSystemOperator;
56
57/** After spatial discretization, the hyperelastic model can be written as a
58 * system of ODEs:
59 * dv/dt = -M^{-1}*(H(x) + S*v)
60 * dx/dt = v,
61 * where x is the vector representing the deformation, v is the velocity field,
62 * M is the mass matrix, S is the viscosity matrix, and H(x) is the nonlinear
63 * hyperelastic operator.
64 *
65 * Class HyperelasticOperator represents the right-hand side of the above
66 * system of ODEs. */
67class HyperelasticOperator : public TimeDependentOperator
68{
69protected:
70 ParFiniteElementSpace &fespace;
71 Array<int> ess_tdof_list;
72
73 ParBilinearForm M, S;
75 double viscosity;
76 HyperelasticModel *model;
77
78 HypreParMatrix *Mmat; // Mass matrix from ParallelAssemble()
79 CGSolver M_solver; // Krylov solver for inverting the mass matrix M
80 HypreSmoother M_prec; // Preconditioner for the mass matrix M
81
82 /** Nonlinear operator defining the reduced backward Euler equation for the
83 velocity. Used in the implementation of method ImplicitSolve. */
84 ReducedSystemOperator *reduced_oper;
85
86 /// Newton solver for the reduced backward Euler equation
87 NewtonSolver *newton_solver;
88
89 /// Solver for the Jacobian solve in the Newton method
90 Solver *J_solver;
91 /// Preconditioner for the Jacobian solve in the Newton method
92 Solver *J_prec;
93
94 mutable Vector z; // auxiliary vector
95
96 const SparseMatrix *local_grad_H;
97 HypreParMatrix *Jacobian;
98
99 double saved_gamma; // saved gamma value from implicit setup
100
101public:
102 /// Solver type to use in the ImplicitSolve() method, used by SDIRK methods.
103 enum NonlinearSolverType
104 {
105 NEWTON = 0, ///< Use MFEM's plain NewtonSolver
106 KINSOL = 1 ///< Use SUNDIALS' KINSOL (through MFEM's class KINSolver)
107 };
108
109 HyperelasticOperator(ParFiniteElementSpace &f, Array<int> &ess_bdr,
110 double visc, double mu, double K,
111 NonlinearSolverType nls_type);
112
113 /// Compute the right-hand side of the ODE system.
114 virtual void Mult(const Vector &vx, Vector &dvx_dt) const;
115
116 /** Solve the Backward-Euler equation: k = f(x + dt*k, t), for the unknown k.
117 This is the only requirement for high-order SDIRK implicit integration.*/
118 virtual void ImplicitSolve(const double dt, const Vector &x, Vector &k);
119
120
121 /// Custom Jacobian system solver for the SUNDIALS time integrators.
122 /** For the ODE system represented by HyperelasticOperator
123
124 M dv/dt = -(H(x) + S*v)
125 dx/dt = v,
126
127 this class facilitates the solution of linear systems of the form
128
129 (M + γS) yv + γJ yx = M bv, J=(dH/dx)(x)
130 - γ yv + yx = bx
131
132 for given bv, bx, x, and γ = GetTimeStep(). */
133
134 /** Linear solve applicable to the SUNDIALS format.
135 Solves (Mass - dt J) y = Mass b, where in our case:
136 Mass = | M 0 | J = | -S -grad_H | y = | v_hat | b = | b_v |
137 | 0 I | | I 0 | | x_hat | | b_x |
138 The result replaces the rhs b.
139 We substitute x_hat = b_x + dt v_hat and solve
140 (M + dt S + dt^2 grad_H) v_hat = M b_v - dt grad_H b_x. */
141
142 /** Setup the linear system. This method is used by the implicit
143 SUNDIALS solvers. */
144 virtual int SUNImplicitSetup(const Vector &y, const Vector &fy,
145 int jok, int *jcur, double gamma);
146
147 /** Solve the linear system. This method is used by the implicit
148 SUNDIALS solvers. */
149 virtual int SUNImplicitSolve(const Vector &b, Vector &x, double tol);
150
151 double ElasticEnergy(const ParGridFunction &x) const;
152 double KineticEnergy(const ParGridFunction &v) const;
153 void GetElasticEnergyDensity(const ParGridFunction &x,
154 ParGridFunction &w) const;
155
156 virtual ~HyperelasticOperator();
157};
158
159/** Nonlinear operator of the form:
160 k --> (M + dt*S)*k + H(x + dt*v + dt^2*k) + S*v,
161 where M and S are given BilinearForms, H is a given NonlinearForm, v and x
162 are given vectors, and dt is a scalar. */
163class ReducedSystemOperator : public Operator
164{
165private:
166 ParBilinearForm *M, *S;
168 mutable HypreParMatrix *Jacobian;
169 double dt;
170 const Vector *v, *x;
171 mutable Vector w, z;
172 const Array<int> &ess_tdof_list;
173
174public:
175 ReducedSystemOperator(ParBilinearForm *M_, ParBilinearForm *S_,
176 ParNonlinearForm *H_, const Array<int> &ess_tdof_list);
177
178 /// Set current dt, v, x values - needed to compute action and Jacobian.
179 void SetParameters(double dt_, const Vector *v_, const Vector *x_);
180
181 /// Compute y = H(x + dt (v + dt k)) + M k + S (v + dt k).
182 virtual void Mult(const Vector &k, Vector &y) const;
183
184 /// Compute J = M + dt S + dt^2 grad_H(x + dt (v + dt k)).
185 virtual Operator &GetGradient(const Vector &k) const;
186
187 virtual ~ReducedSystemOperator();
188};
189
190
191/** Function representing the elastic energy density for the given hyperelastic
192 model+deformation. Used in HyperelasticOperator::GetElasticEnergyDensity. */
193class ElasticEnergyCoefficient : public Coefficient
194{
195private:
196 HyperelasticModel &model;
197 const ParGridFunction &x;
198 DenseMatrix J;
199
200public:
201 ElasticEnergyCoefficient(HyperelasticModel &m, const ParGridFunction &x_)
202 : model(m), x(x_) { }
203 virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip);
204 virtual ~ElasticEnergyCoefficient() { }
205};
206
207void InitialDeformation(const Vector &x, Vector &y);
208
209void InitialVelocity(const Vector &x, Vector &v);
210
211void visualize(ostream &os, ParMesh *mesh, ParGridFunction *deformed_nodes,
212 ParGridFunction *field, const char *field_name = NULL,
213 bool init_vis = false);
214
215
216int main(int argc, char *argv[])
217{
218 // 1. Initialize MPI, HYPRE, and SUNDIALS.
219 Mpi::Init(argc, argv);
220 int myid = Mpi::WorldRank();
221 Hypre::Init();
223
224 // 2. Parse command-line options.
225 const char *mesh_file = "../../data/beam-quad.mesh";
226 int ser_ref_levels = 2;
227 int par_ref_levels = 0;
228 int order = 2;
229 int ode_solver_type = 3;
230 double t_final = 300.0;
231 double dt = 3.0;
232 double visc = 1e-2;
233 double mu = 0.25;
234 double K = 5.0;
235 bool visualization = true;
236 const char *nls = "newton";
237 int vis_steps = 1;
238
239 // Relative and absolute tolerances for CVODE and ARKODE.
240 const double reltol = 1e-1, abstol = 1e-1;
241 // Since this example uses the loose tolerances defined above, it is
242 // necessary to lower the linear solver tolerance for CVODE which is relative
243 // to the above tolerances.
244 const double cvode_eps_lin = 1e-4;
245 // Similarly, the nonlinear tolerance for ARKODE needs to be tightened.
246 const double arkode_eps_nonlin = 1e-6;
247
248 OptionsParser args(argc, argv);
249 args.AddOption(&mesh_file, "-m", "--mesh",
250 "Mesh file to use.");
251 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
252 "Number of times to refine the mesh uniformly in serial.");
253 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
254 "Number of times to refine the mesh uniformly in parallel.");
255 args.AddOption(&order, "-o", "--order",
256 "Order (degree) of the finite elements.");
257 args.AddOption(&ode_solver_type, "-s", "--ode-solver",
258 "ODE solver:\n\t"
259 "1 - Backward Euler,\n\t"
260 "2 - SDIRK2, L-stable\n\t"
261 "3 - SDIRK3, L-stable\n\t"
262 "4 - Implicit Midpoint,\n\t"
263 "5 - SDIRK2, A-stable,\n\t"
264 "6 - SDIRK3, A-stable,\n\t"
265 "7 - Forward Euler,\n\t"
266 "8 - RK2,\n\t"
267 "9 - RK3 SSP,\n\t"
268 "10 - RK4,\n\t"
269 "11 - CVODE implicit BDF, approximate Jacobian,\n\t"
270 "12 - CVODE implicit BDF, specified Jacobian,\n\t"
271 "13 - CVODE implicit ADAMS, approximate Jacobian,\n\t"
272 "14 - CVODE implicit ADAMS, specified Jacobian,\n\t"
273 "15 - ARKODE implicit, approximate Jacobian,\n\t"
274 "16 - ARKODE implicit, specified Jacobian,\n\t"
275 "17 - ARKODE explicit, 4th order.");
276 args.AddOption(&nls, "-nls", "--nonlinear-solver",
277 "Nonlinear systems solver: "
278 "\"newton\" (plain Newton) or \"kinsol\" (KINSOL).");
279 args.AddOption(&t_final, "-tf", "--t-final",
280 "Final time; start time is 0.");
281 args.AddOption(&dt, "-dt", "--time-step",
282 "Time step.");
283 args.AddOption(&visc, "-v", "--viscosity",
284 "Viscosity coefficient.");
285 args.AddOption(&mu, "-mu", "--shear-modulus",
286 "Shear modulus in the Neo-Hookean hyperelastic model.");
287 args.AddOption(&K, "-K", "--bulk-modulus",
288 "Bulk modulus in the Neo-Hookean hyperelastic model.");
289 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
290 "--no-visualization",
291 "Enable or disable GLVis visualization.");
292 args.AddOption(&vis_steps, "-vs", "--visualization-steps",
293 "Visualize every n-th timestep.");
294 args.Parse();
295 if (!args.Good())
296 {
297 if (myid == 0)
298 {
299 args.PrintUsage(cout);
300 }
301 return 1;
302 }
303 if (myid == 0)
304 {
305 args.PrintOptions(cout);
306 }
307
308 // check for valid ODE solver option
309 if (ode_solver_type < 1 || ode_solver_type > 17)
310 {
311 if (myid == 0)
312 {
313 cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
314 }
315 return 1;
316 }
317
318 // 3. Read the serial mesh from the given mesh file on all processors. We can
319 // handle triangular, quadrilateral, tetrahedral and hexahedral meshes
320 // with the same code.
321 Mesh *mesh = new Mesh(mesh_file, 1, 1);
322 int dim = mesh->Dimension();
323
324 // 4. Nonlinear solver
325 map<string,HyperelasticOperator::NonlinearSolverType> nls_map;
326 nls_map["newton"] = HyperelasticOperator::NEWTON;
327 nls_map["kinsol"] = HyperelasticOperator::KINSOL;
328 if (nls_map.find(nls) == nls_map.end())
329 {
330 if (myid == 0)
331 {
332 cout << "Unknown type of nonlinear solver: " << nls << endl;
333 }
334 delete mesh;
335 return 4;
336 }
337
338 // 5. Refine the mesh in serial to increase the resolution. In this example
339 // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
340 // a command-line parameter.
341 for (int lev = 0; lev < ser_ref_levels; lev++)
342 {
343 mesh->UniformRefinement();
344 }
345
346 // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
347 // this mesh further in parallel to increase the resolution. Once the
348 // parallel mesh is defined, the serial mesh can be deleted.
349 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
350 delete mesh;
351 for (int lev = 0; lev < par_ref_levels; lev++)
352 {
353 pmesh->UniformRefinement();
354 }
355
356 // 7. Define the parallel vector finite element spaces representing the mesh
357 // deformation x_gf, the velocity v_gf, and the initial configuration,
358 // x_ref. Define also the elastic energy density, w_gf, which is in a
359 // discontinuous higher-order space. Since x and v are integrated in time
360 // as a system, we group them together in block vector vx, on the unique
361 // parallel degrees of freedom, with offsets given by array true_offset.
362 H1_FECollection fe_coll(order, dim);
363 ParFiniteElementSpace fespace(pmesh, &fe_coll, dim);
364
365 HYPRE_BigInt glob_size = fespace.GlobalTrueVSize();
366 if (myid == 0)
367 {
368 cout << "Number of velocity/deformation unknowns: " << glob_size << endl;
369 }
370 int true_size = fespace.TrueVSize();
371 Array<int> true_offset(3);
372 true_offset[0] = 0;
373 true_offset[1] = true_size;
374 true_offset[2] = 2*true_size;
375
376 BlockVector vx(true_offset);
377 ParGridFunction v_gf, x_gf;
378 v_gf.MakeTRef(&fespace, vx, true_offset[0]);
379 x_gf.MakeTRef(&fespace, vx, true_offset[1]);
380
381 ParGridFunction x_ref(&fespace);
382 pmesh->GetNodes(x_ref);
383
384 L2_FECollection w_fec(order + 1, dim);
385 ParFiniteElementSpace w_fespace(pmesh, &w_fec);
386 ParGridFunction w_gf(&w_fespace);
387
388 // 8. Set the initial conditions for v_gf, x_gf and vx, and define the
389 // boundary conditions on a beam-like mesh (see description above).
391 v_gf.ProjectCoefficient(velo);
392 v_gf.SetTrueVector();
394 x_gf.ProjectCoefficient(deform);
395 x_gf.SetTrueVector();
396
398
399 Array<int> ess_bdr(fespace.GetMesh()->bdr_attributes.Max());
400 ess_bdr = 0;
401 ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed
402
403 // 9. Initialize the hyperelastic operator, the GLVis visualization and print
404 // the initial energies.
405 HyperelasticOperator oper(fespace, ess_bdr, visc, mu, K, nls_map[nls]);
406
407 socketstream vis_v, vis_w;
408 if (visualization)
409 {
410 char vishost[] = "localhost";
411 int visport = 19916;
412 vis_v.open(vishost, visport);
413 vis_v.precision(8);
414 visualize(vis_v, pmesh, &x_gf, &v_gf, "Velocity", true);
415 // Make sure all ranks have sent their 'v' solution before initiating
416 // another set of GLVis connections (one from each rank):
417 MPI_Barrier(pmesh->GetComm());
418 vis_w.open(vishost, visport);
419 if (vis_w)
420 {
421 oper.GetElasticEnergyDensity(x_gf, w_gf);
422 vis_w.precision(8);
423 visualize(vis_w, pmesh, &x_gf, &w_gf, "Elastic energy density", true);
424 }
425 }
426
427 double ee0 = oper.ElasticEnergy(x_gf);
428 double ke0 = oper.KineticEnergy(v_gf);
429 if (myid == 0)
430 {
431 cout << "initial elastic energy (EE) = " << ee0 << endl;
432 cout << "initial kinetic energy (KE) = " << ke0 << endl;
433 cout << "initial total energy (TE) = " << (ee0 + ke0) << endl;
434 }
435
436 // 10. Define the ODE solver used for time integration. Several implicit
437 // singly diagonal implicit Runge-Kutta (SDIRK) methods, as well as
438 // explicit Runge-Kutta methods are available.
439 double t = 0.0;
440 oper.SetTime(t);
441
442 ODESolver *ode_solver = NULL;
443 CVODESolver *cvode = NULL;
444 ARKStepSolver *arkode = NULL;
445 switch (ode_solver_type)
446 {
447 // Implicit L-stable methods
448 case 1: ode_solver = new BackwardEulerSolver; break;
449 case 2: ode_solver = new SDIRK23Solver(2); break;
450 case 3: ode_solver = new SDIRK33Solver; break;
451 // Implicit A-stable methods (not L-stable)
452 case 4: ode_solver = new ImplicitMidpointSolver; break;
453 case 5: ode_solver = new SDIRK23Solver; break;
454 case 6: ode_solver = new SDIRK34Solver; break;
455 // Explicit methods
456 case 7: ode_solver = new ForwardEulerSolver; break;
457 case 8: ode_solver = new RK2Solver(0.5); break; // midpoint method
458 case 9: ode_solver = new RK3SSPSolver; break;
459 case 10: ode_solver = new RK4Solver; break;
460 // CVODE BDF
461 case 11:
462 case 12:
463 cvode = new CVODESolver(MPI_COMM_WORLD, CV_BDF);
464 cvode->Init(oper);
465 cvode->SetSStolerances(reltol, abstol);
466 CVodeSetEpsLin(cvode->GetMem(), cvode_eps_lin);
467 cvode->SetMaxStep(dt);
468 if (ode_solver_type == 11)
469 {
471 }
472 ode_solver = cvode; break;
473 // CVODE Adams
474 case 13:
475 case 14:
476 cvode = new CVODESolver(MPI_COMM_WORLD, CV_ADAMS);
477 cvode->Init(oper);
478 cvode->SetSStolerances(reltol, abstol);
479 CVodeSetEpsLin(cvode->GetMem(), cvode_eps_lin);
480 cvode->SetMaxStep(dt);
481 if (ode_solver_type == 13)
482 {
484 }
485 ode_solver = cvode; break;
486 // ARKStep Implicit methods
487 case 15:
488 case 16:
489 arkode = new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::IMPLICIT);
490 arkode->Init(oper);
491 arkode->SetSStolerances(reltol, abstol);
492 ARKStepSetNonlinConvCoef(arkode->GetMem(), arkode_eps_nonlin);
493 arkode->SetMaxStep(dt);
494 if (ode_solver_type == 15)
495 {
496 arkode->UseSundialsLinearSolver();
497 }
498 ode_solver = arkode; break;
499 // ARKStep Explicit methods
500 case 17:
501 arkode = new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::EXPLICIT);
502 arkode->Init(oper);
503 arkode->SetSStolerances(reltol, abstol);
504 arkode->SetMaxStep(dt);
505 ode_solver = arkode; break;
506 }
507
508 // Initialize MFEM integrators, SUNDIALS integrators are initialized above
509 if (ode_solver_type < 11) { ode_solver->Init(oper); }
510
511 // 11. Perform time-integration
512 // (looping over the time iterations, ti, with a time-step dt).
513 bool last_step = false;
514 for (int ti = 1; !last_step; ti++)
515 {
516 double dt_real = min(dt, t_final - t);
517
518 ode_solver->Step(vx, t, dt_real);
519
520 last_step = (t >= t_final - 1e-8*dt);
521
522 if (last_step || (ti % vis_steps) == 0)
523 {
525
526 double ee = oper.ElasticEnergy(x_gf);
527 double ke = oper.KineticEnergy(v_gf);
528
529 if (myid == 0)
530 {
531 cout << "step " << ti << ", t = " << t << ", EE = " << ee
532 << ", KE = " << ke << ", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
533
534 if (cvode) { cvode->PrintInfo(); }
535 else if (arkode) { arkode->PrintInfo(); }
536 }
537
538 if (visualization)
539 {
540 visualize(vis_v, pmesh, &x_gf, &v_gf);
541 if (vis_w)
542 {
543 oper.GetElasticEnergyDensity(x_gf, w_gf);
544 visualize(vis_w, pmesh, &x_gf, &w_gf);
545 }
546 }
547 }
548 }
549
550 // 12. Save the displaced mesh, the velocity and elastic energy.
551 {
553 GridFunction *nodes = &x_gf;
554 int owns_nodes = 0;
555 pmesh->SwapNodes(nodes, owns_nodes);
556
557 ostringstream mesh_name, velo_name, ee_name;
558 mesh_name << "deformed." << setfill('0') << setw(6) << myid;
559 velo_name << "velocity." << setfill('0') << setw(6) << myid;
560 ee_name << "elastic_energy." << setfill('0') << setw(6) << myid;
561
562 ofstream mesh_ofs(mesh_name.str().c_str());
563 mesh_ofs.precision(8);
564 pmesh->Print(mesh_ofs);
565 pmesh->SwapNodes(nodes, owns_nodes);
566 ofstream velo_ofs(velo_name.str().c_str());
567 velo_ofs.precision(8);
568 v_gf.Save(velo_ofs);
569 ofstream ee_ofs(ee_name.str().c_str());
570 ee_ofs.precision(8);
571 oper.GetElasticEnergyDensity(x_gf, w_gf);
572 w_gf.Save(ee_ofs);
573 }
574
575 // 13. Free the used memory.
576 delete ode_solver;
577 delete pmesh;
578
579 return 0;
580}
581
582void visualize(ostream &os, ParMesh *mesh, ParGridFunction *deformed_nodes,
583 ParGridFunction *field, const char *field_name, bool init_vis)
584{
585 if (!os)
586 {
587 return;
588 }
589
590 GridFunction *nodes = deformed_nodes;
591 int owns_nodes = 0;
592
593 mesh->SwapNodes(nodes, owns_nodes);
594
595 os << "parallel " << mesh->GetNRanks() << " " << mesh->GetMyRank() << "\n";
596 os << "solution\n" << *mesh << *field;
597
598 mesh->SwapNodes(nodes, owns_nodes);
599
600 if (init_vis)
601 {
602 os << "window_size 800 800\n";
603 os << "window_title '" << field_name << "'\n";
604 if (mesh->SpaceDimension() == 2)
605 {
606 os << "view 0 0\n"; // view from top
607 os << "keys jl\n"; // turn off perspective and light
608 }
609 os << "keys cm\n"; // show colorbar and mesh
610 os << "autoscale value\n"; // update value-range; keep mesh-extents fixed
611 os << "pause\n";
612 }
613 os << flush;
614}
615
616
617ReducedSystemOperator::ReducedSystemOperator(
619 const Array<int> &ess_tdof_list_)
620 : Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
621 Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height),
622 ess_tdof_list(ess_tdof_list_)
623{ }
624
625void ReducedSystemOperator::SetParameters(double dt_, const Vector *v_,
626 const Vector *x_)
627{
628 dt = dt_; v = v_; x = x_;
629}
630
631void ReducedSystemOperator::Mult(const Vector &k, Vector &y) const
632{
633 // compute: y = H(x + dt*(v + dt*k)) + M*k + S*(v + dt*k)
634 add(*v, dt, k, w);
635 add(*x, dt, w, z);
636 H->Mult(z, y);
637 M->TrueAddMult(k, y);
638 S->TrueAddMult(w, y);
639 y.SetSubVector(ess_tdof_list, 0.0);
640}
641
642Operator &ReducedSystemOperator::GetGradient(const Vector &k) const
643{
644 delete Jacobian;
645 SparseMatrix *localJ = Add(1.0, M->SpMat(), dt, S->SpMat());
646 add(*v, dt, k, w);
647 add(*x, dt, w, z);
648 localJ->Add(dt*dt, H->GetLocalGradient(z));
649 Jacobian = M->ParallelAssemble(localJ);
650 delete localJ;
651 HypreParMatrix *Je = Jacobian->EliminateRowsCols(ess_tdof_list);
652 delete Je;
653 return *Jacobian;
654}
655
656ReducedSystemOperator::~ReducedSystemOperator()
657{
658 delete Jacobian;
659}
660
661
662HyperelasticOperator::HyperelasticOperator(ParFiniteElementSpace &f,
663 Array<int> &ess_bdr, double visc,
664 double mu, double K,
665 NonlinearSolverType nls_type)
666 : TimeDependentOperator(2*f.TrueVSize(), 0.0), fespace(f),
667 M(&fespace), S(&fespace), H(&fespace),
668 viscosity(visc), M_solver(f.GetComm()), z(height/2),
669 local_grad_H(NULL), Jacobian(NULL)
670{
671 const double rel_tol = 1e-8;
672 const int skip_zero_entries = 0;
673
674 const double ref_density = 1.0; // density in the reference configuration
675 ConstantCoefficient rho0(ref_density);
676 M.AddDomainIntegrator(new VectorMassIntegrator(rho0));
677 M.Assemble(skip_zero_entries);
678 M.Finalize(skip_zero_entries);
679 Mmat = M.ParallelAssemble();
680 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
681 HypreParMatrix *Me = Mmat->EliminateRowsCols(ess_tdof_list);
682 delete Me;
683
684 M_solver.iterative_mode = false;
685 M_solver.SetRelTol(rel_tol);
686 M_solver.SetAbsTol(0.0);
687 M_solver.SetMaxIter(30);
688 M_solver.SetPrintLevel(0);
689 M_prec.SetType(HypreSmoother::Jacobi);
690 M_solver.SetPreconditioner(M_prec);
691 M_solver.SetOperator(*Mmat);
692
693 model = new NeoHookeanModel(mu, K);
694 H.AddDomainIntegrator(new HyperelasticNLFIntegrator(model));
695 H.SetEssentialTrueDofs(ess_tdof_list);
696
697 ConstantCoefficient visc_coeff(viscosity);
698 S.AddDomainIntegrator(new VectorDiffusionIntegrator(visc_coeff));
699 S.Assemble(skip_zero_entries);
700 S.Finalize(skip_zero_entries);
701
702 reduced_oper = new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
703
704 HypreSmoother *J_hypreSmoother = new HypreSmoother;
705 J_hypreSmoother->SetType(HypreSmoother::l1Jacobi);
706 J_hypreSmoother->SetPositiveDiagonal(true);
707 J_prec = J_hypreSmoother;
708
709 MINRESSolver *J_minres = new MINRESSolver(f.GetComm());
710 J_minres->SetRelTol(rel_tol);
711 J_minres->SetAbsTol(0.0);
712 J_minres->SetMaxIter(300);
713 J_minres->SetPrintLevel(-1);
714 J_minres->SetPreconditioner(*J_prec);
715 J_solver = J_minres;
716
717 if (nls_type == KINSOL)
718 {
719 KINSolver *kinsolver = new KINSolver(f.GetComm(), KIN_LINESEARCH, true);
720 kinsolver->SetJFNK(true);
721 kinsolver->SetLSMaxIter(100);
722 newton_solver = kinsolver;
723 newton_solver->SetOperator(*reduced_oper);
724 newton_solver->SetMaxIter(200);
725 newton_solver->SetRelTol(rel_tol);
726 newton_solver->SetPrintLevel(1);
727 kinsolver->SetMaxSetupCalls(4);
728 }
729 else
730 {
731 newton_solver = new NewtonSolver(f.GetComm());
732 newton_solver->SetOperator(*reduced_oper);
733 newton_solver->SetMaxIter(10);
734 newton_solver->SetRelTol(rel_tol);
735 newton_solver->SetPrintLevel(-1);
736 }
737 newton_solver->SetSolver(*J_solver);
738 newton_solver->iterative_mode = false;
739}
740
741void HyperelasticOperator::Mult(const Vector &vx, Vector &dvx_dt) const
742{
743 // Create views to the sub-vectors v, x of vx, and dv_dt, dx_dt of dvx_dt
744 int sc = height/2;
745 Vector v(vx.GetData() + 0, sc);
746 Vector x(vx.GetData() + sc, sc);
747 Vector dv_dt(dvx_dt.GetData() + 0, sc);
748 Vector dx_dt(dvx_dt.GetData() + sc, sc);
749
750 H.Mult(x, z);
751 if (viscosity != 0.0)
752 {
753 S.TrueAddMult(v, z);
754 z.SetSubVector(ess_tdof_list, 0.0);
755 }
756 z.Neg(); // z = -z
757 M_solver.Mult(z, dv_dt);
758
759 dx_dt = v;
760}
761
762void HyperelasticOperator::ImplicitSolve(const double dt,
763 const Vector &vx, Vector &dvx_dt)
764{
765 int sc = height/2;
766 Vector v(vx.GetData() + 0, sc);
767 Vector x(vx.GetData() + sc, sc);
768 Vector dv_dt(dvx_dt.GetData() + 0, sc);
769 Vector dx_dt(dvx_dt.GetData() + sc, sc);
770
771 // By eliminating kx from the coupled system:
772 // kv = -M^{-1}*[H(x + dt*kx) + S*(v + dt*kv)]
773 // kx = v + dt*kv
774 // we reduce it to a nonlinear equation for kv, represented by the
775 // reduced_oper. This equation is solved with the newton_solver
776 // object (using J_solver and J_prec internally).
777 reduced_oper->SetParameters(dt, &v, &x);
778 Vector zero; // empty vector is interpreted as zero r.h.s. by NewtonSolver
779 newton_solver->Mult(zero, dv_dt);
780 MFEM_VERIFY(newton_solver->GetConverged(),
781 "Nonlinear solver did not converge.");
782#ifdef MFEM_DEBUG
783 if (fespace.GetMyRank() == 0)
784 {
785 cout << " num nonlin sol iters = " << newton_solver->GetNumIterations()
786 << ", final norm = " << newton_solver->GetFinalNorm() << '\n';
787 }
788#endif
789 add(v, dt, dv_dt, dx_dt);
790}
791
792int HyperelasticOperator::SUNImplicitSetup(const Vector &y,
793 const Vector &fy, int jok, int *jcur,
794 double gamma)
795{
796 int sc = y.Size() / 2;
797 const Vector x(y.GetData() + sc, sc);
798
799 // J = M + dt*(S + dt*grad(H))
800 if (Jacobian) { delete Jacobian; }
801 SparseMatrix *localJ = Add(1.0, M.SpMat(), gamma, S.SpMat());
802 local_grad_H = &H.GetLocalGradient(x);
803 localJ->Add(gamma*gamma, *local_grad_H);
804 Jacobian = M.ParallelAssemble(localJ);
805 delete localJ;
806 HypreParMatrix *Je = Jacobian->EliminateRowsCols(ess_tdof_list);
807 delete Je;
808
809 // Set Jacobian solve operator
810 J_solver->SetOperator(*Jacobian);
811
812 // Indicate that the Jacobian was updated
813 *jcur = 1;
814
815 // Save gamma for use in solve
816 saved_gamma = gamma;
817
818 // Return success
819 return 0;
820}
821
822int HyperelasticOperator::SUNImplicitSolve(const Vector &b, Vector &x,
823 double tol)
824{
825 int sc = b.Size() / 2;
826 ParFiniteElementSpace *fes = H.ParFESpace();
827 Vector b_v(b.GetData() + 0, sc);
828 Vector b_x(b.GetData() + sc, sc);
829 Vector x_v(x.GetData() + 0, sc);
830 Vector x_x(x.GetData() + sc, sc);
831 Vector rhs(sc);
832
833 // We can assume that b_v and b_x have zeros at essential tdofs.
834
835 // rhs = M b_v - dt*grad(H) b_x
836 ParGridFunction lb_x(fes), lrhs(fes);
837 lb_x.Distribute(b_x);
838 local_grad_H->Mult(lb_x, lrhs);
839 lrhs.ParallelAssemble(rhs);
840 rhs *= -saved_gamma;
841 M.TrueAddMult(b_v, rhs);
842 rhs.SetSubVector(ess_tdof_list, 0.0);
843
844 J_solver->iterative_mode = false;
845 J_solver->Mult(rhs, x_v);
846
847 add(b_x, saved_gamma, x_v, x_x);
848
849 return 0;
850}
851
852double HyperelasticOperator::ElasticEnergy(const ParGridFunction &x) const
853{
854 return H.GetEnergy(x);
855}
856
857double HyperelasticOperator::KineticEnergy(const ParGridFunction &v) const
858{
859 double energy = 0.5*M.ParInnerProduct(v, v);
860 return energy;
861}
862
863void HyperelasticOperator::GetElasticEnergyDensity(
864 const ParGridFunction &x, ParGridFunction &w) const
865{
866 ElasticEnergyCoefficient w_coeff(*model, x);
867 w.ProjectCoefficient(w_coeff);
868}
869
870HyperelasticOperator::~HyperelasticOperator()
871{
872 delete Jacobian;
873 delete newton_solver;
874 delete J_solver;
875 delete J_prec;
876 delete reduced_oper;
877 delete model;
878 delete Mmat;
879}
880
881
882double ElasticEnergyCoefficient::Eval(ElementTransformation &T,
883 const IntegrationPoint &ip)
884{
885 model.SetTransformation(T);
886 x.GetVectorGradient(T, J);
887 // return model.EvalW(J); // in reference configuration
888 return model.EvalW(J)/J.Det(); // in deformed configuration
889}
890
891
893{
894 // set the initial configuration to be the same as the reference, stress
895 // free, configuration
896 y = x;
897}
898
899void InitialVelocity(const Vector &x, Vector &v)
900{
901 const int dim = x.Size();
902 const double s = 0.1/64.;
903
904 v = 0.0;
905 v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
906 v(0) = -s*x(0)*x(0);
907}
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
Definition sundials.hpp:673
void SetMaxStep(double dt_max)
Set the maximum time step.
void PrintInfo() const
Print various ARKStep statistics.
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults.
@ IMPLICIT
Implicit RK method.
Definition sundials.hpp:679
@ EXPLICIT
Explicit RK method.
Definition sundials.hpp:678
void UseSundialsLinearSolver()
Attach a SUNDIALS GMRES linear solver to ARKode.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
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
const SparseMatrix & SpMat() const
Returns a const reference to the sparse matrix: .
A class to handle Vectors in a block fashion.
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
Interface to the CVODE library – linear multi-step methods.
Definition sundials.hpp:386
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition sundials.cpp:875
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition sundials.cpp:893
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
Definition sundials.cpp:709
void PrintInfo() const
Print various CVODE statistics.
Definition sundials.cpp:919
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
Definition sundials.cpp:855
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 dense matrix using column-major storage.
Definition densemat.hpp:24
real_t Det() const
Definition densemat.cpp:535
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
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
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_)
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Definition hypre.cpp:2356
Parallel smoothers in hypre.
Definition hypre.hpp:1020
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
Definition hypre.hpp:1136
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition hypre.cpp:3550
@ l1Jacobi
l1-scaled Jacobi
Definition hypre.hpp:1080
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.hpp:74
Implicit midpoint method. A-stable, not L-stable.
Definition ode.hpp:425
Class for integration point with weight.
Definition intrules.hpp:35
real_t GetFinalNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult().
Definition solvers.hpp:275
void SetRelTol(real_t rtol)
Definition solvers.hpp:209
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
Definition solvers.hpp:260
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
Interface to the KINSOL library – nonlinear solver methods.
Definition sundials.hpp:847
void SetJFNK(bool use_jfnk)
Set the Jacobian Free Newton Krylov flag. The default is false.
Definition sundials.hpp:959
virtual void SetOperator(const Operator &op)
Set the nonlinear Operator of the system and initialize KINSOL.
void SetMaxSetupCalls(int max_calls)
Set maximum number of nonlinear iterations without a Jacobian update.
void SetLSMaxIter(int m)
Set the maximum number of linear solver iterations.
Definition sundials.hpp:963
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
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
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
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.
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
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
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.
Class for parallel bilinear form.
Abstract parallel finite element space.
Definition pfespace.hpp:29
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
int TrueVSize() const
Obsolete, kept for backward compatibility.
Definition pfespace.hpp:436
Class for parallel grid function.
Definition pgridfunc.hpp:33
void Save(std::ostream &out) const override
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Class for parallel meshes.
Definition pmesh.hpp:34
MPI_Comm GetComm() const
Definition pmesh.hpp:402
int GetMyRank() const
Definition pmesh.hpp:404
int GetNRanks() const
Definition pmesh.hpp:403
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4801
Parallel non-linear operator on the true dofs.
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
virtual void SetOperator(const Operator &op)=0
Set/update the solver for the given operator.
Data type sparse matrix.
Definition sparsemat.hpp:51
void Add(const int i, const int j, const real_t val)
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void * GetMem() const
Access the SUNDIALS memory structure.
Definition sundials.hpp:373
static void Init()
Definition sundials.cpp:164
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
A general vector function coefficient.
Vector data type.
Definition vector.hpp:80
void Neg()
(*this) = -(*this)
Definition vector.cpp:300
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
Definition vector.cpp:604
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'.
int dim
Definition ex24.cpp:53
real_t mu
Definition ex25.cpp:140
int main()
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
const int visport
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition vector.cpp:316
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]
void visualize(ostream &os, ParMesh *mesh, ParGridFunction *deformed_nodes, ParGridFunction *field, const char *field_name=NULL, bool init_vis=false)
Definition ex10p.cpp:582
void InitialDeformation(const Vector &x, Vector &y)
Definition ex10p.cpp:892
void InitialVelocity(const Vector &x, Vector &v)
Definition ex10p.cpp:899