48#ifndef MFEM_USE_SUNDIALS
49#error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES
55class ReducedSystemOperator;
84 ReducedSystemOperator *reduced_oper;
103 enum NonlinearSolverType
110 double visc,
double mu,
double K,
111 NonlinearSolverType nls_type);
114 virtual void Mult(
const Vector &vx,
Vector &dvx_dt)
const;
118 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
144 virtual int SUNImplicitSetup(
const Vector &y,
const Vector &fy,
145 int jok,
int *jcur,
double gamma);
149 virtual int SUNImplicitSolve(
const Vector &
b,
Vector &x,
double tol);
156 virtual ~HyperelasticOperator();
163class ReducedSystemOperator :
public Operator
179 void SetParameters(
double dt_,
const Vector *v_,
const Vector *x_);
187 virtual ~ReducedSystemOperator();
202 : model(m), x(x_) { }
204 virtual ~ElasticEnergyCoefficient() { }
213 bool init_vis =
false);
216int main(
int argc,
char *argv[])
225 const char *mesh_file =
"../../data/beam-quad.mesh";
226 int ser_ref_levels = 2;
227 int par_ref_levels = 0;
229 int ode_solver_type = 3;
230 double t_final = 300.0;
235 bool visualization =
true;
236 const char *nls =
"newton";
240 const double reltol = 1e-1, abstol = 1e-1;
244 const double cvode_eps_lin = 1e-4;
246 const double arkode_eps_nonlin = 1e-6;
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.");
256 "Order (degree) of the finite elements.");
257 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
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"
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",
283 args.
AddOption(&visc,
"-v",
"--viscosity",
284 "Viscosity coefficient.");
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.");
309 if (ode_solver_type < 1 || ode_solver_type > 17)
313 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
321 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
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())
332 cout <<
"Unknown type of nonlinear solver: " << nls << endl;
341 for (
int lev = 0; lev < ser_ref_levels; lev++)
351 for (
int lev = 0; lev < par_ref_levels; lev++)
368 cout <<
"Number of velocity/deformation unknowns: " << glob_size << endl;
373 true_offset[1] = true_size;
374 true_offset[2] = 2*true_size;
378 v_gf.
MakeTRef(&fespace, vx, true_offset[0]);
379 x_gf.
MakeTRef(&fespace, vx, true_offset[1]);
405 HyperelasticOperator oper(fespace, ess_bdr, visc,
mu, K, nls_map[nls]);
414 visualize(vis_v, pmesh, &x_gf, &v_gf,
"Velocity",
true);
421 oper.GetElasticEnergyDensity(x_gf, w_gf);
423 visualize(vis_w, pmesh, &x_gf, &w_gf,
"Elastic energy density",
true);
427 double ee0 = oper.ElasticEnergy(x_gf);
428 double ke0 = oper.KineticEnergy(v_gf);
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;
445 switch (ode_solver_type)
457 case 8: ode_solver =
new RK2Solver(0.5);
break;
459 case 10: ode_solver =
new RK4Solver;
break;
466 CVodeSetEpsLin(cvode->
GetMem(), cvode_eps_lin);
468 if (ode_solver_type == 11)
472 ode_solver = cvode;
break;
479 CVodeSetEpsLin(cvode->
GetMem(), cvode_eps_lin);
481 if (ode_solver_type == 13)
485 ode_solver = cvode;
break;
492 ARKStepSetNonlinConvCoef(arkode->
GetMem(), arkode_eps_nonlin);
494 if (ode_solver_type == 15)
498 ode_solver = arkode;
break;
505 ode_solver = arkode;
break;
509 if (ode_solver_type < 11) { ode_solver->
Init(oper); }
513 bool last_step =
false;
514 for (
int ti = 1; !last_step; ti++)
516 double dt_real = min(dt, t_final -
t);
518 ode_solver->
Step(vx,
t, dt_real);
520 last_step = (
t >= t_final - 1e-8*dt);
522 if (last_step || (ti % vis_steps) == 0)
526 double ee = oper.ElasticEnergy(x_gf);
527 double ke = oper.KineticEnergy(v_gf);
531 cout <<
"step " << ti <<
", t = " <<
t <<
", EE = " << ee
532 <<
", KE = " << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
535 else if (arkode) { arkode->
PrintInfo(); }
543 oper.GetElasticEnergyDensity(x_gf, w_gf);
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;
562 ofstream mesh_ofs(mesh_name.str().c_str());
563 mesh_ofs.precision(8);
564 pmesh->
Print(mesh_ofs);
566 ofstream velo_ofs(velo_name.str().c_str());
567 velo_ofs.precision(8);
569 ofstream ee_ofs(ee_name.str().c_str());
571 oper.GetElasticEnergyDensity(x_gf, w_gf);
596 os <<
"solution\n" << *mesh << *field;
602 os <<
"window_size 800 800\n";
603 os <<
"window_title '" << field_name <<
"'\n";
610 os <<
"autoscale value\n";
617ReducedSystemOperator::ReducedSystemOperator(
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_)
625void ReducedSystemOperator::SetParameters(
double dt_,
const Vector *v_,
628 dt = dt_; v = v_; x = x_;
631void ReducedSystemOperator::Mult(
const Vector &k,
Vector &y)
const
637 M->TrueAddMult(k, y);
638 S->TrueAddMult(w, y);
642Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const
648 localJ->
Add(dt*dt, H->GetLocalGradient(z));
649 Jacobian = M->ParallelAssemble(localJ);
656ReducedSystemOperator::~ReducedSystemOperator()
665 NonlinearSolverType nls_type)
667 M(&fespace), S(&fespace), H(&fespace),
668 viscosity(visc), M_solver(
f.GetComm()), z(height/2),
669 local_grad_H(NULL), Jacobian(NULL)
671 const double rel_tol = 1e-8;
672 const int skip_zero_entries = 0;
674 const double ref_density = 1.0;
677 M.Assemble(skip_zero_entries);
678 M.Finalize(skip_zero_entries);
679 Mmat = M.ParallelAssemble();
680 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
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);
690 M_solver.SetPreconditioner(M_prec);
691 M_solver.SetOperator(*Mmat);
695 H.SetEssentialTrueDofs(ess_tdof_list);
699 S.Assemble(skip_zero_entries);
700 S.Finalize(skip_zero_entries);
702 reduced_oper =
new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
707 J_prec = J_hypreSmoother;
717 if (nls_type == KINSOL)
722 newton_solver = kinsolver;
724 newton_solver->SetMaxIter(200);
725 newton_solver->SetRelTol(rel_tol);
726 newton_solver->SetPrintLevel(1);
732 newton_solver->SetOperator(*reduced_oper);
733 newton_solver->SetMaxIter(10);
734 newton_solver->SetRelTol(rel_tol);
735 newton_solver->SetPrintLevel(-1);
737 newton_solver->SetSolver(*J_solver);
738 newton_solver->iterative_mode =
false;
741void HyperelasticOperator::Mult(
const Vector &vx,
Vector &dvx_dt)
const
751 if (viscosity != 0.0)
757 M_solver.
Mult(z, dv_dt);
762void HyperelasticOperator::ImplicitSolve(
const double dt,
777 reduced_oper->SetParameters(dt, &v, &x);
779 newton_solver->
Mult(zero, dv_dt);
781 "Nonlinear solver did not converge.");
783 if (fespace.GetMyRank() == 0)
786 <<
", final norm = " << newton_solver->
GetFinalNorm() <<
'\n';
789 add(v, dt, dv_dt, dx_dt);
792int HyperelasticOperator::SUNImplicitSetup(
const Vector &y,
793 const Vector &fy,
int jok,
int *jcur,
796 int sc = y.
Size() / 2;
800 if (Jacobian) {
delete Jacobian; }
802 local_grad_H = &H.GetLocalGradient(x);
803 localJ->
Add(gamma*gamma, *local_grad_H);
804 Jacobian = M.ParallelAssemble(localJ);
822int HyperelasticOperator::SUNImplicitSolve(
const Vector &
b,
Vector &x,
825 int sc =
b.Size() / 2;
827 Vector b_v(
b.GetData() + 0, sc);
828 Vector b_x(
b.GetData() + sc, sc);
837 lb_x.Distribute(b_x);
838 local_grad_H->
Mult(lb_x, lrhs);
839 lrhs.ParallelAssemble(rhs);
841 M.TrueAddMult(b_v, rhs);
842 rhs.SetSubVector(ess_tdof_list, 0.0);
845 J_solver->
Mult(rhs, x_v);
847 add(b_x, saved_gamma, x_v, x_x);
852double HyperelasticOperator::ElasticEnergy(
const ParGridFunction &x)
const
857double HyperelasticOperator::KineticEnergy(
const ParGridFunction &v)
const
859 double energy = 0.5*M.ParInnerProduct(v, v);
863void HyperelasticOperator::GetElasticEnergyDensity(
866 ElasticEnergyCoefficient w_coeff(*model, x);
870HyperelasticOperator::~HyperelasticOperator()
873 delete newton_solver;
902 const double s = 0.1/64.;
905 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.
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults.
@ IMPLICIT
Implicit RK method.
@ EXPLICIT
Explicit RK method.
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.
Backward Euler ODE solver. L-stable.
A class to handle Vectors in a block fashion.
Conjugate gradient method.
virtual void Mult(const Vector &b, Vector &x) const
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 SetMaxStep(double dt_max)
Set the maximum time step.
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
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 dense matrix using column-major storage.
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.
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.
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_)
Wrapper for hypre's ParCSR matrix class.
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
Parallel smoothers in hypre.
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
@ l1Jacobi
l1-scaled Jacobi
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
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.
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.
Arbitrary order "L2-conforming" discontinuous finite elements.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
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.
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.
virtual void Mult(const Vector &b, Vector &x) const
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.
Abstract parallel finite element space.
HYPRE_BigInt GlobalTrueVSize() const
int TrueVSize() const
Obsolete, kept for backward compatibility.
Class for parallel grid function.
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.
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
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)
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication.
void * GetMem() const
Access the SUNDIALS memory structure.
Base abstract class for first order time dependent operators.
virtual void SetTime(const real_t t_)
Set the current time.
A general vector function coefficient.
void Neg()
(*this) = -(*this)
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
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.
void visualize(ostream &os, ParMesh *mesh, ParGridFunction *deformed_nodes, ParGridFunction *field, const char *field_name=NULL, bool init_vis=false)
void InitialDeformation(const Vector &x, Vector &y)
void InitialVelocity(const Vector &x, Vector &v)