48 #ifndef MFEM_USE_SUNDIALS
49 #error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES
55 class ReducedSystemOperator;
84 ReducedSystemOperator *reduced_oper;
103 enum NonlinearSolverType
110 double visc,
double mu,
double K,
111 NonlinearSolverType nls_type);
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();
163 class ReducedSystemOperator :
public Operator
179 void SetParameters(
double dt_,
const Vector *v_,
const Vector *x_);
187 virtual ~ReducedSystemOperator();
193 class ElasticEnergyCoefficient :
public Coefficient
202 : model(m), x(x_) { }
204 virtual ~ElasticEnergyCoefficient() { }
213 bool init_vis =
false);
216 int main(
int argc,
char *argv[])
220 MPI_Init(&argc, &argv);
221 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
222 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
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.");
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.");
310 if (ode_solver_type < 1 || ode_solver_type > 17)
314 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
323 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
327 map<string,HyperelasticOperator::NonlinearSolverType> nls_map;
328 nls_map[
"newton"] = HyperelasticOperator::NEWTON;
329 nls_map[
"kinsol"] = HyperelasticOperator::KINSOL;
330 if (nls_map.find(nls) == nls_map.end())
334 cout <<
"Unknown type of nonlinear solver: " << nls << endl;
344 for (
int lev = 0; lev < ser_ref_levels; lev++)
354 for (
int lev = 0; lev < par_ref_levels; lev++)
371 cout <<
"Number of velocity/deformation unknowns: " << glob_size << endl;
376 true_offset[1] = true_size;
377 true_offset[2] = 2*true_size;
381 v_gf.
MakeTRef(&fespace, vx, true_offset[0]);
382 x_gf.
MakeTRef(&fespace, vx, true_offset[1]);
408 HyperelasticOperator oper(fespace, ess_bdr, visc, mu, K, nls_map[nls]);
413 char vishost[] =
"localhost";
415 vis_v.
open(vishost, visport);
417 visualize(vis_v, pmesh, &x_gf, &v_gf,
"Velocity",
true);
421 vis_w.
open(vishost, visport);
424 oper.GetElasticEnergyDensity(x_gf, w_gf);
426 visualize(vis_w, pmesh, &x_gf, &w_gf,
"Elastic energy density",
true);
430 double ee0 = oper.ElasticEnergy(x_gf);
431 double ke0 = oper.KineticEnergy(v_gf);
434 cout <<
"initial elastic energy (EE) = " << ee0 << endl;
435 cout <<
"initial kinetic energy (KE) = " << ke0 << endl;
436 cout <<
"initial total energy (TE) = " << (ee0 + ke0) << endl;
448 switch (ode_solver_type)
460 case 8: ode_solver =
new RK2Solver(0.5);
break;
462 case 10: ode_solver =
new RK4Solver;
break;
469 CVodeSetEpsLin(cvode->
GetMem(), cvode_eps_lin);
471 if (ode_solver_type == 11)
475 ode_solver = cvode;
break;
482 CVodeSetEpsLin(cvode->
GetMem(), cvode_eps_lin);
484 if (ode_solver_type == 13)
488 ode_solver = cvode;
break;
492 arkode =
new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::IMPLICIT);
495 ARKStepSetNonlinConvCoef(arkode->
GetMem(), arkode_eps_nonlin);
497 if (ode_solver_type == 15)
501 ode_solver = arkode;
break;
504 arkode =
new ARKStepSolver(MPI_COMM_WORLD, ARKStepSolver::EXPLICIT);
508 ode_solver = arkode;
break;
512 if (ode_solver_type < 11) { ode_solver->
Init(oper); }
516 bool last_step =
false;
517 for (
int ti = 1; !last_step; ti++)
519 double dt_real = min(dt, t_final - t);
521 ode_solver->
Step(vx, t, dt_real);
523 last_step = (t >= t_final - 1e-8*dt);
525 if (last_step || (ti % vis_steps) == 0)
529 double ee = oper.ElasticEnergy(x_gf);
530 double ke = oper.KineticEnergy(v_gf);
534 cout <<
"step " << ti <<
", t = " << t <<
", EE = " << ee
535 <<
", KE = " << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
538 else if (arkode) { arkode->
PrintInfo(); }
546 oper.GetElasticEnergyDensity(x_gf, w_gf);
560 ostringstream mesh_name, velo_name, ee_name;
561 mesh_name <<
"deformed." << setfill(
'0') << setw(6) << myid;
562 velo_name <<
"velocity." << setfill(
'0') << setw(6) << myid;
563 ee_name <<
"elastic_energy." << setfill(
'0') << setw(6) << myid;
565 ofstream mesh_ofs(mesh_name.str().c_str());
566 mesh_ofs.precision(8);
567 pmesh->
Print(mesh_ofs);
569 ofstream velo_ofs(velo_name.str().c_str());
570 velo_ofs.precision(8);
572 ofstream ee_ofs(ee_name.str().c_str());
574 oper.GetElasticEnergyDensity(x_gf, w_gf);
601 out <<
"solution\n" << *mesh << *field;
607 out <<
"window_size 800 800\n";
608 out <<
"window_title '" << field_name <<
"'\n";
615 out <<
"autoscale value\n";
622 ReducedSystemOperator::ReducedSystemOperator(
625 :
Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
626 Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height),
627 ess_tdof_list(ess_tdof_list_)
630 void ReducedSystemOperator::SetParameters(
double dt_,
const Vector *v_,
633 dt = dt_; v = v_; x = x_;
642 M->TrueAddMult(k, y);
643 S->TrueAddMult(w, y);
647 Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const
653 localJ->
Add(dt*dt, H->GetLocalGradient(z));
654 Jacobian = M->ParallelAssemble(localJ);
661 ReducedSystemOperator::~ReducedSystemOperator()
670 NonlinearSolverType nls_type)
672 M(&fespace), S(&fespace), H(&fespace),
673 viscosity(visc), M_solver(f.GetComm()), z(height/2),
674 local_grad_H(NULL), Jacobian(NULL)
676 const double rel_tol = 1e-8;
677 const int skip_zero_entries = 0;
679 const double ref_density = 1.0;
682 M.Assemble(skip_zero_entries);
683 M.Finalize(skip_zero_entries);
684 Mmat = M.ParallelAssemble();
685 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
689 M_solver.iterative_mode =
false;
690 M_solver.SetRelTol(rel_tol);
691 M_solver.SetAbsTol(0.0);
692 M_solver.SetMaxIter(30);
693 M_solver.SetPrintLevel(0);
694 M_prec.SetType(HypreSmoother::Jacobi);
695 M_solver.SetPreconditioner(M_prec);
696 M_solver.SetOperator(*Mmat);
700 H.SetEssentialTrueDofs(ess_tdof_list);
704 S.Assemble(skip_zero_entries);
705 S.Finalize(skip_zero_entries);
707 reduced_oper =
new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
710 J_hypreSmoother->
SetType(HypreSmoother::l1Jacobi);
712 J_prec = J_hypreSmoother;
722 if (nls_type == KINSOL)
725 newton_solver = kinsolver;
754 if (viscosity != 0.0)
760 M_solver.Mult(z, dv_dt);
765 void HyperelasticOperator::ImplicitSolve(
const double dt,
780 reduced_oper->SetParameters(dt, &v, &x);
782 newton_solver->
Mult(zero, dv_dt);
784 "Nonlinear solver did not converge.");
786 if (fespace.GetMyRank() == 0)
789 <<
", final norm = " << newton_solver->
GetFinalNorm() <<
'\n';
792 add(v, dt, dv_dt, dx_dt);
795 int HyperelasticOperator::SUNImplicitSetup(
const Vector &y,
796 const Vector &fy,
int jok,
int *jcur,
799 int sc = y.
Size() / 2;
803 if (Jacobian) {
delete Jacobian; }
805 local_grad_H = &H.GetLocalGradient(x);
806 localJ->
Add(gamma*gamma, *local_grad_H);
807 Jacobian = M.ParallelAssemble(localJ);
813 J_solver->SetOperator(*Jacobian);
825 int HyperelasticOperator::SUNImplicitSolve(
const Vector &
b,
Vector &x,
828 int sc = b.
Size() / 2;
840 lb_x.Distribute(b_x);
841 local_grad_H->Mult(lb_x, lrhs);
842 lrhs.ParallelAssemble(rhs);
844 M.TrueAddMult(b_v, rhs);
845 rhs.SetSubVector(ess_tdof_list, 0.0);
847 J_solver->iterative_mode =
false;
848 J_solver->Mult(rhs, x_v);
850 add(b_x, saved_gamma, x_v, x_x);
855 double HyperelasticOperator::ElasticEnergy(
const ParGridFunction &x)
const
857 return H.GetEnergy(x);
860 double HyperelasticOperator::KineticEnergy(
const ParGridFunction &v)
const
862 double loc_energy = 0.5*M.InnerProduct(v, v);
864 MPI_Allreduce(&loc_energy, &energy, 1, MPI_DOUBLE, MPI_SUM,
869 void HyperelasticOperator::GetElasticEnergyDensity(
872 ElasticEnergyCoefficient w_coeff(*model, x);
876 HyperelasticOperator::~HyperelasticOperator()
879 delete newton_solver;
891 model.SetTransformation(T);
894 return model.EvalW(J)/J.Det();
908 const double s = 0.1/64.;
911 v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
double Eval(ElementTransformation &T, const IntegrationPoint &ip, double t)
Evaluate the coefficient in the element described by T at the point ip at time t. ...
void Add(const int i, const int j, const double a)
void InitialDeformation(const Vector &x, Vector &y)
Conjugate gradient method.
Class for grid function - Vector with associated FE space.
int GetNumIterations() const
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
A class to handle Vectors in a block fashion.
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
Subclass constant coefficient.
Base abstract class for first order time dependent operators.
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
virtual void Step(Vector &x, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
void * GetMem() const
Access the SUNDIALS memory structure.
Data type dense matrix using column-major storage.
int Size() const
Returns the size of the vector.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
virtual void Save(std::ostream &out) const
void SetMaxStep(double dt_max)
Set the maximum time step.
Abstract parallel finite element space.
virtual void ProjectCoefficient(Coefficient &coeff)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
int main(int argc, char *argv[])
Backward Euler ODE solver. L-stable.
double * GetData() const
Return a pointer to the beginning of the Vector data.
void add(const Vector &v1, const Vector &v2, Vector &v)
void SetMaxSetupCalls(int max_calls)
Set maximum number of nonlinear iterations without a Jacobian update.
double GetFinalNorm() const
void InitialVelocity(const Vector &x, Vector &v)
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
void SetPrintLevel(int print_lvl)
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Mesh * GetMesh() const
Returns the mesh.
Interface to the CVODE library – linear multi-step methods.
Interface to the KINSOL library – nonlinear solver methods.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void MakeTRef(FiniteElementSpace *f, double *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
void SetMaxIter(int max_it)
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
HYPRE_Int GlobalTrueVSize() const
void SetMaxStep(double dt_max)
Set the maximum time step.
Newton's method for solving F(x)=b for a given operator F.
virtual void Print(std::ostream &out=mfem::out) const
Parallel smoothers in hypre.
void PrintUsage(std::ostream &out) const
static int Mult(const N_Vector u, N_Vector fu, void *user_data)
Wrapper to compute the nonlinear residual .
int SpaceDimension() const
The classical explicit forth-order Runge-Kutta method, RK4.
void SetAbsTol(double atol)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void SetRelTol(double rtol)
virtual void SetOperator(const Operator &op)
Set the nonlinear Operator of the system and initialize KINSOL.
Base class Coefficient that may optionally depend on time.
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)
Third-order, strong stability preserving (SSP) Runge-Kutta method.
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Implicit midpoint method. A-stable, not L-stable.
Class for integration point with weight.
void PrintInfo() const
Print various ARKStep statistics.
void PrintOptions(std::ostream &out) const
Abstract class for hyperelastic models.
int open(const char hostname[], int port)
void GetNodes(Vector &node_coord) const
Arbitrary order H1-conforming (continuous) finite elements.
void PrintInfo() const
Print various CVODE statistics.
Class for parallel grid function.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
The classical forward Euler method.
Wrapper for hypre's ParCSR matrix class.
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
Class for parallel meshes.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
void UseSundialsLinearSolver()
Attach a SUNDIALS GMRES linear solver to ARKode.
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
Arbitrary order "L2-conforming" discontinuous finite elements.