50#ifndef MFEM_USE_SUNDIALS 
   51#error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES 
   57class ReducedSystemOperator;
 
   84   ReducedSystemOperator *reduced_oper;
 
  104                        double visc, 
double mu, 
double K,
 
  105                        int kinsol_nls_type = -1, 
double kinsol_damping = 0.0,
 
  106                        int kinsol_aa_n = 0);
 
  109   virtual void Mult(
const Vector &vx, 
Vector &dvx_dt) 
const;
 
  113   virtual void ImplicitSolve(
const double dt, 
const Vector &x, 
Vector &k);
 
  139   virtual int SUNImplicitSetup(
const Vector &y, 
const Vector &fy,
 
  140                                int jok, 
int *jcur, 
double gamma);
 
  144   virtual int SUNImplicitSolve(
const Vector &
b, 
Vector &x, 
double tol);
 
  146   double ElasticEnergy(
const Vector &x) 
const;
 
  147   double KineticEnergy(
const Vector &v) 
const;
 
  150   virtual ~HyperelasticOperator();
 
  157class ReducedSystemOperator : 
public Operator 
  171   void SetParameters(
double dt_, 
const Vector *v_, 
const Vector *x_);
 
  179   virtual ~ReducedSystemOperator();
 
  194      : model(m), x(x_) { }
 
  196   virtual ~ElasticEnergyCoefficient() { }
 
  205               bool init_vis = 
false);
 
  208int main(
int argc, 
char *argv[])
 
  214   const char *mesh_file = 
"../../data/beam-quad.mesh";
 
  217   int ode_solver_type = 3;
 
  218   double t_final = 300.0;
 
  223   bool visualization = 
true;
 
  224   int nonlinear_solver_type = 0;
 
  226   double kinsol_damping = 0.0;
 
  227   int kinsol_aa_n = -1;
 
  230   const double reltol = 1e-1, abstol = 1e-1;
 
  234   const double cvode_eps_lin = 1e-4;
 
  236   const double arkode_eps_nonlin = 1e-6;
 
  239   args.
AddOption(&mesh_file, 
"-m", 
"--mesh",
 
  240                  "Mesh file to use.");
 
  241   args.
AddOption(&ref_levels, 
"-r", 
"--refine",
 
  242                  "Number of times to refine the mesh uniformly.");
 
  244                  "Order (degree) of the finite elements.");
 
  245   args.
AddOption(&ode_solver_type, 
"-s", 
"--ode-solver",
 
  247                  "1  - Backward Euler,\n\t" 
  248                  "2  - SDIRK2, L-stable\n\t" 
  249                  "3  - SDIRK3, L-stable\n\t" 
  250                  "4  - Implicit Midpoint,\n\t" 
  251                  "5  - SDIRK2, A-stable,\n\t" 
  252                  "6  - SDIRK3, A-stable,\n\t" 
  253                  "7  - Forward Euler,\n\t" 
  257                  "11 - CVODE implicit BDF, approximate Jacobian,\n\t" 
  258                  "12 - CVODE implicit BDF, specified Jacobian,\n\t" 
  259                  "13 - CVODE implicit ADAMS, approximate Jacobian,\n\t" 
  260                  "14 - CVODE implicit ADAMS, specified Jacobian,\n\t" 
  261                  "15 - ARKODE implicit, approximate Jacobian,\n\t" 
  262                  "16 - ARKODE implicit, specified Jacobian,\n\t" 
  263                  "17 - ARKODE explicit, 4th order.");
 
  264   args.
AddOption(&nonlinear_solver_type, 
"-nls", 
"--nonlinear-solver",
 
  265                  "Nonlinear system solver:\n\t" 
  266                  "0  - MFEM Newton method,\n\t" 
  267                  "1  - KINSOL Newton method,\n\t" 
  268                  "2  - KINSOL Newton method with globalization,\n\t" 
  269                  "3  - KINSOL fixed-point method (with or without AA),\n\t" 
  270                  "4  - KINSOL Picard method (with or without AA).");
 
  271   args.
AddOption(&kinsol_damping, 
"-damp", 
"--kinsol-damping",
 
  272                  "Picard or Fixed-Point damping parameter (only valid with KINSOL): " 
  274   args.
AddOption(&kinsol_aa_n, 
"-aan", 
"--anderson-subspace",
 
  275                  "Anderson Acceleration subspace size (only valid with KINSOL)");
 
  276   args.
AddOption(&t_final, 
"-tf", 
"--t-final",
 
  277                  "Final time; start time is 0.");
 
  278   args.
AddOption(&dt, 
"-dt", 
"--time-step",
 
  280   args.
AddOption(&visc, 
"-v", 
"--viscosity",
 
  281                  "Viscosity coefficient.");
 
  283                  "Shear modulus in the Neo-Hookean hyperelastic model.");
 
  284   args.
AddOption(&K, 
"-K", 
"--bulk-modulus",
 
  285                  "Bulk modulus in the Neo-Hookean hyperelastic model.");
 
  286   args.
AddOption(&visualization, 
"-vis", 
"--visualization", 
"-no-vis",
 
  287                  "--no-visualization",
 
  288                  "Enable or disable GLVis visualization.");
 
  289   args.
AddOption(&vis_steps, 
"-vs", 
"--visualization-steps",
 
  290                  "Visualize every n-th timestep.");
 
  300   if (ode_solver_type < 1 || ode_solver_type > 17)
 
  302      cout << 
"Unknown ODE solver type: " << ode_solver_type << 
'\n';
 
  307   if (nonlinear_solver_type < 0 || nonlinear_solver_type > 4)
 
  309      cout << 
"Unknown nonlinear solver type: " << nonlinear_solver_type << 
"\n";
 
  312   if (kinsol_damping > 0.0 &&
 
  313       !(nonlinear_solver_type == 3 || nonlinear_solver_type == 4))
 
  315      cout << 
"Only KINSOL fixed-point and Picard methods can use damping\n";
 
  318   if (kinsol_aa_n > 0 &&
 
  319       !(nonlinear_solver_type == 3 || nonlinear_solver_type == 4))
 
  321      cout << 
"Only KINSOL fixed-point and Picard methods can use AA\n";
 
  328   Mesh *mesh = 
new Mesh(mesh_file, 1, 1);
 
  334   for (
int lev = 0; lev < ref_levels; lev++)
 
  349   cout << 
"Number of velocity/deformation unknowns: " << fe_size << endl;
 
  352   fe_offset[1] = fe_size;
 
  353   fe_offset[2] = 2*fe_size;
 
  382   std::unique_ptr<HyperelasticOperator> oper;
 
  383   if (nonlinear_solver_type == 0)
 
  384      oper = std::make_unique<HyperelasticOperator>(fespace, ess_bdr, visc, 
mu,
 
  388      switch (nonlinear_solver_type)
 
  391            oper = std::make_unique<HyperelasticOperator>(fespace, ess_bdr,
 
  392                                                          visc, 
mu, K, KIN_NONE);
 
  395            oper = std::make_unique<HyperelasticOperator>(fespace, ess_bdr,
 
  396                                                          visc, 
mu, K, KIN_LINESEARCH);
 
  399            oper = std::make_unique<HyperelasticOperator>(fespace, ess_bdr,
 
  400                                                          visc, 
mu, K, KIN_FP, kinsol_damping, kinsol_aa_n);
 
  403            oper = std::make_unique<HyperelasticOperator>(fespace, ess_bdr,
 
  404                                                          visc, 
mu, K, KIN_PICARD, kinsol_damping, kinsol_aa_n);
 
  417      visualize(vis_v, mesh, &x, &v, 
"Velocity", 
true);
 
  421         oper->GetElasticEnergyDensity(x, w);
 
  423         visualize(vis_w, mesh, &x, &w, 
"Elastic energy density", 
true);
 
  429   cout << 
"initial elastic energy (EE) = " << ee0 << endl;
 
  430   cout << 
"initial kinetic energy (KE) = " << ke0 << endl;
 
  431   cout << 
"initial   total energy (TE) = " << (ee0 + ke0) << endl;
 
  442   switch (ode_solver_type)
 
  454      case 8:  ode_solver = 
new RK2Solver(0.5); 
break; 
 
  456      case 10: ode_solver = 
new RK4Solver; 
break;
 
  463         CVodeSetEpsLin(cvode->
GetMem(), cvode_eps_lin);
 
  465         if (ode_solver_type == 11)
 
  469         ode_solver = cvode; 
break;
 
  476         CVodeSetEpsLin(cvode->
GetMem(), cvode_eps_lin);
 
  478         if (ode_solver_type == 13)
 
  482         ode_solver = cvode; 
break;
 
  489#if MFEM_SUNDIALS_VERSION < 70100 
  490         ARKStepSetNonlinConvCoef(arkode->
GetMem(), arkode_eps_nonlin);
 
  492         ARKodeSetNonlinConvCoef(arkode->
GetMem(), arkode_eps_nonlin);
 
  495         if (ode_solver_type == 15)
 
  499         ode_solver = arkode; 
break;
 
  506         ode_solver = arkode; 
break;
 
  510   if (ode_solver_type < 11) { ode_solver->
Init(*oper); }
 
  514   bool last_step = 
false;
 
  515   for (
int ti = 1; !last_step; ti++)
 
  517      double dt_real = min(dt, t_final - t);
 
  519      ode_solver->
Step(vx, t, dt_real);
 
  521      last_step = (t >= t_final - 1e-8*dt);
 
  523      if (last_step || (ti % vis_steps) == 0)
 
  528         cout << 
"step " << ti << 
", t = " << t << 
", EE = " << ee << 
", KE = " 
  529              << ke << 
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
 
  532         else if (arkode) { arkode->
PrintInfo(); }
 
  540               oper->GetElasticEnergyDensity(x, w);
 
  553      ofstream mesh_ofs(
"deformed.mesh");
 
  554      mesh_ofs.precision(8);
 
  555      mesh->
Print(mesh_ofs);
 
  557      ofstream velo_ofs(
"velocity.sol");
 
  558      velo_ofs.precision(8);
 
  560      ofstream ee_ofs(
"elastic_energy.sol");
 
  562      oper->GetElasticEnergyDensity(x, w);
 
 
  575               GridFunction *field, 
const char *field_name, 
bool init_vis)
 
  587   os << 
"solution\n" << *mesh << *field;
 
  593      os << 
"window_size 800 800\n";
 
  594      os << 
"window_title '" << field_name << 
"'\n";
 
  601      os << 
"autoscale value\n"; 
 
 
  608ReducedSystemOperator::ReducedSystemOperator(
 
  610   : 
Operator(M_->Height()), M(M_), S(S_), H(H_), Jacobian(NULL),
 
  611     dt(0.0), v(NULL), x(NULL), w(height), z(height)
 
  614void ReducedSystemOperator::SetParameters(
double dt_, 
const Vector *v_,
 
  617   dt = dt_;  v = v_;  x = x_;
 
  620void ReducedSystemOperator::Mult(
const Vector &k, 
Vector &y)
 const 
  630Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
 const 
  637   Jacobian->
Add(dt*dt, *grad_H);
 
  641ReducedSystemOperator::~ReducedSystemOperator()
 
  651                                           double kinsol_damping,
 
  654     M(&fespace), S(&fespace), H(&fespace),
 
  655     viscosity(visc), z(height/2),
 
  656     grad_H(NULL), Jacobian(NULL)
 
  658   const double rel_tol = 1e-8;
 
  659   const int skip_zero_entries = 0;
 
  661   const double ref_density = 1.0; 
 
  664   M.Assemble(skip_zero_entries);
 
  666   fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
 
  668   M.FormSystemMatrix(ess_tdof_list, tmp);
 
  670   M_solver.iterative_mode = 
false;
 
  671   M_solver.SetRelTol(rel_tol);
 
  672   M_solver.SetAbsTol(0.0);
 
  673   M_solver.SetMaxIter(30);
 
  674   M_solver.SetPrintLevel(0);
 
  675   M_solver.SetPreconditioner(M_prec);
 
  676   M_solver.SetOperator(M.SpMat());
 
  680   H.SetEssentialTrueDofs(ess_tdof_list);
 
  684   S.Assemble(skip_zero_entries);
 
  685   S.FormSystemMatrix(ess_tdof_list, tmp);
 
  687   reduced_oper = 
new ReducedSystemOperator(&M, &S, &H);
 
  689#ifndef MFEM_USE_SUITESPARSE 
  703   if (kinsol_nls_type > 0)
 
  706      if (kinsol_nls_type != KIN_PICARD)
 
  715      newton_solver = kinsolver;
 
  717      newton_solver->SetMaxIter(200);
 
  718      newton_solver->SetRelTol(rel_tol);
 
  719      newton_solver->SetPrintLevel(0);
 
  721      if (kinsol_damping > 0.0)
 
  729      newton_solver->SetOperator(*reduced_oper);
 
  730      newton_solver->SetMaxIter(10);
 
  731      newton_solver->SetRelTol(rel_tol);
 
  732      newton_solver->SetPrintLevel(-1);
 
  734   newton_solver->SetSolver(*J_solver);
 
  735   newton_solver->iterative_mode = 
false;
 
  738void HyperelasticOperator::Mult(
const Vector &vx, 
Vector &dvx_dt)
 const 
  748   if (viscosity != 0.0)
 
  753   M_solver.
Mult(z, dv_dt);
 
  758void HyperelasticOperator::ImplicitSolve(
const double dt,
 
  773   reduced_oper->SetParameters(dt, &v, &x);
 
  775   newton_solver->
Mult(zero, dv_dt);
 
  777               "Nonlinear solver did not converge.");
 
  780        << 
", final norm = " << newton_solver->
GetFinalNorm() << 
'\n';
 
  782   add(v, dt, dv_dt, dx_dt);
 
  785int HyperelasticOperator::SUNImplicitSetup(
const Vector &y,
 
  786                                           const Vector &fy, 
int jok, 
int *jcur,
 
  789   int sc = y.
Size() / 2;
 
  793   if (Jacobian) { 
delete Jacobian; }
 
  796   Jacobian->
Add(gamma * gamma, *grad_H);
 
  811int HyperelasticOperator::SUNImplicitSolve(
const Vector &
b, 
Vector &x,
 
  814   int sc = 
b.Size() / 2;
 
  815   Vector b_v(
b.GetData() +  0, sc);
 
  816   Vector b_x(
b.GetData() + sc, sc);
 
  822   grad_H->
Mult(b_x, rhs);
 
  827   J_solver->
Mult(rhs, x_v);
 
  829   add(b_x, saved_gamma, x_v, x_x);
 
  834double HyperelasticOperator::ElasticEnergy(
const Vector &x)
 const 
  839double HyperelasticOperator::KineticEnergy(
const Vector &v)
 const 
  844void HyperelasticOperator::GetElasticEnergyDensity(
 
  847   ElasticEnergyCoefficient w_coeff(*model, x);
 
  851HyperelasticOperator::~HyperelasticOperator()
 
  854   delete newton_solver;
 
  882   const double s = 0.1/64.;
 
  885   v(
dim-1) = s*x(0)*x(0)*(8.0-x(0));
 
 
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
 
void SetMaxStep(double dt_max)
Set the maximum time step.
 
void PrintInfo() const
Print various ARKStep statistics.
 
@ IMPLICIT
Implicit RK method.
 
@ EXPLICIT
Explicit RK method.
 
void UseSundialsLinearSolver()
Attach a SUNDIALS GMRES linear solver to ARKode.
 
void Init(TimeDependentOperator &f_) override
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults.
 
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.
 
Backward Euler ODE solver. L-stable.
 
A class to handle Vectors in a block fashion.
 
Vector & GetBlock(int i)
Get the i-th vector in the block.
 
Conjugate gradient method.
 
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
 
Interface to the CVODE library – linear multi-step methods.
 
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
 
void Init(TimeDependentOperator &f_) override
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
 
void SetMaxStep(double dt_max)
Set the maximum time step.
 
void PrintInfo() const
Print various CVODE statistics.
 
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
 
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.
 
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
 
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
 
Mesh * GetMesh() const
Returns the mesh.
 
The classical forward Euler method.
 
Class for grid function - Vector with associated FE space.
 
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.
 
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.
 
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
 
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
 
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.
 
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.
 
Class for integration point with weight.
 
real_t GetFinalNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult().
 
void SetRelTol(real_t rtol)
 
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
 
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
 
void SetMaxIter(int max_it)
 
bool GetConverged() const
Returns true if the last call to Mult() converged successfully.
 
void SetAbsTol(real_t atol)
 
Interface to the KINSOL library – nonlinear solver methods.
 
void SetJFNK(bool use_jfnk)
Set the Jacobian Free Newton Krylov flag. The default is false.
 
void SetOperator(const Operator &op) override
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 SetDamping(double damping)
 
void EnableAndersonAcc(int n, int orth=KIN_ORTH_MGS, int delay=0, double damping=1.0)
Enable Anderson Acceleration for KIN_FP or KIN_PICARD.
 
void SetLSMaxIter(int m)
Set the maximum number of linear solver iterations.
 
Arbitrary order "L2-conforming" discontinuous finite elements.
 
void SetPreconditioner(Solver &pr) override
This should be called before SetOperator.
 
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
 
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
 
int Dimension() const
Dimension of the reference space used within the elements.
 
int SpaceDimension() const
Dimension of the physical space containing the mesh.
 
void GetNodes(Vector &node_coord) const
 
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
 
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Swap the internal node GridFunction pointer and ownership flag members with the given ones.
 
Newton's method for solving F(x)=b for a given operator F.
 
void Mult(const Vector &b, Vector &x) const override
Solve the nonlinear system with right-hand side b.
 
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
 
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
 
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].
 
int height
Dimension of the output / number of rows in the matrix.
 
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...
 
bool Good() const
Return true if the command line options were parsed successfully.
 
Third-order, strong stability preserving (SSP) Runge-Kutta method.
 
The classical explicit forth-order Runge-Kutta method, RK4.
 
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
 
virtual void SetOperator(const Operator &op)=0
Set/update the solver for the given operator.
 
void Add(const int i, const int j, const real_t val)
 
void Mult(const Vector &x, Vector &y) const override
Matrix vector multiplication.
 
void * GetMem() const
Access the SUNDIALS memory structure.
 
Base abstract class for first order time dependent operators.
 
Direct sparse solver using UMFPACK.
 
A general vector function coefficient.
 
void Neg()
(*this) = -(*this)
 
int Size() const
Returns the size of the vector.
 
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
 
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
 
void add(const Vector &v1, const Vector &v2, Vector &v)
 
std::function< real_t(const Vector &)> f(real_t mass_coeff)
 
void Add(const DenseMatrix &A, const DenseMatrix &B, real_t alpha, DenseMatrix &C)
C = A + alpha*B.
 
std::array< int, NCMesh::MaxFaceNodes > nodes
 
void visualize(ostream &os, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
 
void InitialDeformation(const Vector &x, Vector &y)
 
void InitialVelocity(const Vector &x, Vector &v)