45#error This example requires that MFEM is built with MFEM_USE_PETSC=YES
51class ReducedSystemOperator;
80 ReducedSystemOperator *reduced_oper;
100 bool use_petsc,
bool petsc_use_jfnk);
103 virtual void Mult(
const Vector &vx,
Vector &dvx_dt)
const;
113 virtual ~HyperelasticOperator();
120class ReducedSystemOperator :
public Operator
144 virtual ~ReducedSystemOperator();
155 PreconditionerFactory(
const ReducedSystemOperator& op_,
const string& name_)
158 virtual ~PreconditionerFactory() {}
172 : model(m), x(x_) { }
174 virtual ~ElasticEnergyCoefficient() { }
183 bool init_vis =
false);
186int main(
int argc,
char *argv[])
194 const char *mesh_file =
"../../data/beam-quad.mesh";
195 int ser_ref_levels = 2;
196 int par_ref_levels = 0;
198 int ode_solver_type = 3;
204 bool visualization =
true;
206 bool use_petsc =
true;
207 const char *petscrc_file =
"";
208 bool petsc_use_jfnk =
false;
209 const char *device_config =
"cpu";
212 args.
AddOption(&mesh_file,
"-m",
"--mesh",
213 "Mesh file to use.");
214 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
215 "Number of times to refine the mesh uniformly in serial.");
216 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
217 "Number of times to refine the mesh uniformly in parallel.");
219 "Order (degree) of the finite elements.");
220 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
221 "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
222 " 11 - Forward Euler, 12 - RK2,\n\t"
223 " 13 - RK3 SSP, 14 - RK4.");
224 args.
AddOption(&t_final,
"-tf",
"--t-final",
225 "Final time; start time is 0.");
226 args.
AddOption(&dt,
"-dt",
"--time-step",
228 args.
AddOption(&visc,
"-v",
"--viscosity",
229 "Viscosity coefficient.");
231 "Shear modulus in the Neo-Hookean hyperelastic model.");
232 args.
AddOption(&K,
"-K",
"--bulk-modulus",
233 "Bulk modulus in the Neo-Hookean hyperelastic model.");
234 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
235 "--no-visualization",
236 "Enable or disable GLVis visualization.");
237 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
238 "Visualize every n-th timestep.");
239 args.
AddOption(&use_petsc,
"-usepetsc",
"--usepetsc",
"-no-petsc",
241 "Use or not PETSc to solve the nonlinear system.");
242 args.
AddOption(&petscrc_file,
"-petscopts",
"--petscopts",
243 "PetscOptions file to use.");
244 args.
AddOption(&petsc_use_jfnk,
"-jfnk",
"--jfnk",
"-no-jfnk",
246 "Use JFNK with user-defined preconditioner factory.");
247 args.
AddOption(&device_config,
"-d",
"--device",
248 "Device configuration string, see Device::Configure().");
265 Device device(device_config);
266 if (myid == 0) { device.
Print(); }
277 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
284 switch (ode_solver_type)
292 case 12: ode_solver =
new RK2Solver(0.5);
break;
294 case 14: ode_solver =
new RK4Solver;
break;
302 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
310 for (
int lev = 0; lev < ser_ref_levels; lev++)
320 for (
int lev = 0; lev < par_ref_levels; lev++)
337 cout <<
"Number of velocity/deformation unknowns: " << glob_size << endl;
342 true_offset[1] = true_size;
343 true_offset[2] = 2*true_size;
347 v_gf.
MakeTRef(&fespace, vx, true_offset[0]);
348 x_gf.
MakeTRef(&fespace, vx, true_offset[1]);
374 HyperelasticOperator *oper =
new HyperelasticOperator(fespace, ess_bdr, visc,
385 visualize(vis_v, pmesh, &x_gf, &v_gf,
"Velocity",
true);
392 oper->GetElasticEnergyDensity(x_gf, w_gf);
394 visualize(vis_w, pmesh, &x_gf, &w_gf,
"Elastic energy density",
true);
398 real_t ee0 = oper->ElasticEnergy(x_gf);
399 real_t ke0 = oper->KineticEnergy(v_gf);
402 cout <<
"initial elastic energy (EE) = " << ee0 << endl;
403 cout <<
"initial kinetic energy (KE) = " << ke0 << endl;
404 cout <<
"initial total energy (TE) = " << (ee0 + ke0) << endl;
409 ode_solver->
Init(*oper);
413 bool last_step =
false;
414 for (
int ti = 1; !last_step; ti++)
416 real_t dt_real = min(dt, t_final - t);
418 ode_solver->
Step(vx, t, dt_real);
420 last_step = (t >= t_final - 1e-8*dt);
422 if (last_step || (ti % vis_steps) == 0)
426 real_t ee = oper->ElasticEnergy(x_gf);
427 real_t ke = oper->KineticEnergy(v_gf);
431 cout <<
"step " << ti <<
", t = " << t <<
", EE = " << ee
432 <<
", KE = " << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
440 oper->GetElasticEnergyDensity(x_gf, w_gf);
454 ostringstream mesh_name, velo_name, ee_name;
455 mesh_name <<
"deformed." << setfill(
'0') << setw(6) << myid;
456 velo_name <<
"velocity." << setfill(
'0') << setw(6) << myid;
457 ee_name <<
"elastic_energy." << setfill(
'0') << setw(6) << myid;
459 ofstream mesh_ofs(mesh_name.str().c_str());
460 mesh_ofs.precision(8);
461 pmesh->
Print(mesh_ofs);
463 ofstream velo_ofs(velo_name.str().c_str());
464 velo_ofs.precision(8);
466 ofstream ee_ofs(ee_name.str().c_str());
468 oper->GetElasticEnergyDensity(x_gf, w_gf);
497 os <<
"solution\n" << *mesh << *field;
503 os <<
"window_size 800 800\n";
504 os <<
"window_title '" << field_name <<
"'\n";
511 os <<
"autoscale value\n";
518ReducedSystemOperator::ReducedSystemOperator(
521 :
Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
522 Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height),
523 ess_tdof_list(ess_tdof_list_)
526void ReducedSystemOperator::SetParameters(
real_t dt_,
const Vector *v_,
529 dt = dt_; v = v_; x = x_;
532void ReducedSystemOperator::Mult(
const Vector &k,
Vector &y)
const
538 M->TrueAddMult(k, y);
539 S->TrueAddMult(w, y);
543Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const
549 localJ->
Add(dt*dt, H->GetLocalGradient(z));
552 Jacobian = M->ParallelAssemble(localJ);
559ReducedSystemOperator::~ReducedSystemOperator()
568 bool use_petsc_factory)
570 M(&fespace), S(&fespace), H(&fespace),
571 viscosity(visc), M_solver(
f.GetComm()),
572 newton_solver(
f.GetComm()), pnewton_solver(NULL), z(height/2)
574#if defined(MFEM_USE_DOUBLE)
575 const real_t rel_tol = 1e-8;
576 const real_t newton_abs_tol = 0.0;
577#elif defined(MFEM_USE_SINGLE)
578 const real_t rel_tol = 1e-3;
579 const real_t newton_abs_tol = 1e-4;
581 const int skip_zero_entries = 0;
583 const real_t ref_density = 1.0;
586 M.Assemble(skip_zero_entries);
587 M.Finalize(skip_zero_entries);
588 Mmat = M.ParallelAssemble();
589 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
593 M_solver.iterative_mode =
false;
594 M_solver.SetRelTol(rel_tol);
595 M_solver.SetAbsTol(0.0);
596 M_solver.SetMaxIter(30);
597 M_solver.SetPrintLevel(0);
599 M_solver.SetPreconditioner(M_prec);
600 M_solver.SetOperator(*Mmat);
604 H.SetEssentialTrueDofs(ess_tdof_list);
608 S.Assemble(skip_zero_entries);
609 S.Finalize(skip_zero_entries);
611 reduced_oper =
new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
617 J_prec = J_hypreSmoother;
630 newton_solver.SetSolver(*J_solver);
631 newton_solver.SetOperator(*reduced_oper);
632 newton_solver.SetPrintLevel(1);
633 newton_solver.SetRelTol(rel_tol);
634 newton_solver.SetAbsTol(newton_abs_tol);
635 newton_solver.SetMaxIter(10);
648 if (use_petsc_factory)
650 J_factory =
new PreconditionerFactory(*reduced_oper,
"JFNK preconditioner");
651 pnewton_solver->SetPreconditionerFactory(J_factory);
653 pnewton_solver->SetPrintLevel(1);
654 pnewton_solver->SetRelTol(rel_tol);
655 pnewton_solver->SetAbsTol(newton_abs_tol);
656 pnewton_solver->SetMaxIter(10);
660void HyperelasticOperator::Mult(
const Vector &vx,
Vector &dvx_dt)
const
670 if (viscosity != 0.0)
676 M_solver.
Mult(z, dv_dt);
681void HyperelasticOperator::ImplicitSolve(
const real_t dt,
696 reduced_oper->SetParameters(dt, &v, &x);
700 newton_solver.
Mult(zero, dv_dt);
702 "Newton solver did not converge.");
706 pnewton_solver->
Mult(zero, dv_dt);
708 "Newton solver did not converge.");
710 add(v, dt, dv_dt, dx_dt);
720 real_t energy = 0.5*M.ParInnerProduct(v, v);
724void HyperelasticOperator::GetElasticEnergyDensity(
727 ElasticEnergyCoefficient w_coeff(*model, x);
731HyperelasticOperator::~HyperelasticOperator()
739 delete pnewton_solver;
779 v(
dim-1) = s*x(0)*x(0)*(8.0-x(0));
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.
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
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.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
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.
void SetRelTol(real_t rtol)
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)
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.
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.
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].
Pointer to an Operator of a specified type.
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
int height
Dimension of the output / number of rows in the matrix.
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
Abstract class for PETSc's nonlinear solvers.
void Mult(const Vector &b, Vector &x) const override
Application of the solver.
Wrapper for PETSc's matrix class.
Abstract class for PETSc's preconditioners.
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.
void Add(const int i, const int j, const real_t val)
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)
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
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)
std::array< int, NCMesh::MaxFaceNodes > nodes