44 #ifndef MFEM_USE_PETSC 45 #error This example requires that MFEM is built with MFEM_USE_PETSC=YES 51 class ReducedSystemOperator;
80 ReducedSystemOperator *reduced_oper;
99 double visc,
double mu,
double K,
100 bool use_petsc,
bool petsc_use_jfnk);
106 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
113 virtual ~HyperelasticOperator();
120 class ReducedSystemOperator :
public Operator 136 void SetParameters(
double dt_,
const Vector *v_,
const Vector *x_);
144 virtual ~ReducedSystemOperator();
155 PreconditionerFactory(
const ReducedSystemOperator& op_,
const string& name_)
158 virtual ~PreconditionerFactory() {}
163 class ElasticEnergyCoefficient :
public Coefficient 172 : model(m), x(x_) { }
174 virtual ~ElasticEnergyCoefficient() { }
183 bool init_vis =
false);
186 int main(
int argc,
char *argv[])
189 Mpi::Init(argc, argv);
190 int myid = Mpi::WorldRank();
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;
199 double t_final = 300.0;
204 bool visualization =
true;
206 bool use_petsc =
true;
207 const char *petscrc_file =
"";
208 bool petsc_use_jfnk =
false;
211 args.
AddOption(&mesh_file,
"-m",
"--mesh",
212 "Mesh file to use.");
213 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
214 "Number of times to refine the mesh uniformly in serial.");
215 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
216 "Number of times to refine the mesh uniformly in parallel.");
218 "Order (degree) of the finite elements.");
219 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
220 "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t" 221 " 11 - Forward Euler, 12 - RK2,\n\t" 222 " 13 - RK3 SSP, 14 - RK4.");
223 args.
AddOption(&t_final,
"-tf",
"--t-final",
224 "Final time; start time is 0.");
225 args.
AddOption(&dt,
"-dt",
"--time-step",
227 args.
AddOption(&visc,
"-v",
"--viscosity",
228 "Viscosity coefficient.");
230 "Shear modulus in the Neo-Hookean hyperelastic model.");
231 args.
AddOption(&K,
"-K",
"--bulk-modulus",
232 "Bulk modulus in the Neo-Hookean hyperelastic model.");
233 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
234 "--no-visualization",
235 "Enable or disable GLVis visualization.");
236 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
237 "Visualize every n-th timestep.");
238 args.
AddOption(&use_petsc,
"-usepetsc",
"--usepetsc",
"-no-petsc",
240 "Use or not PETSc to solve the nonlinear system.");
241 args.
AddOption(&petscrc_file,
"-petscopts",
"--petscopts",
242 "PetscOptions file to use.");
243 args.
AddOption(&petsc_use_jfnk,
"-jfnk",
"--jfnk",
"-no-jfnk",
245 "Use JFNK with user-defined preconditioner factory.");
269 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
276 switch (ode_solver_type)
284 case 12: ode_solver =
new RK2Solver(0.5);
break;
286 case 14: ode_solver =
new RK4Solver;
break;
294 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
302 for (
int lev = 0; lev < ser_ref_levels; lev++)
312 for (
int lev = 0; lev < par_ref_levels; lev++)
329 cout <<
"Number of velocity/deformation unknowns: " << glob_size << endl;
334 true_offset[1] = true_size;
335 true_offset[2] = 2*true_size;
339 v_gf.
MakeTRef(&fespace, vx, true_offset[0]);
340 x_gf.
MakeTRef(&fespace, vx, true_offset[1]);
366 HyperelasticOperator *oper =
new HyperelasticOperator(fespace, ess_bdr, visc,
377 visualize(vis_v, pmesh, &x_gf, &v_gf,
"Velocity",
true);
384 oper->GetElasticEnergyDensity(x_gf, w_gf);
386 visualize(vis_w, pmesh, &x_gf, &w_gf,
"Elastic energy density",
true);
390 double ee0 = oper->ElasticEnergy(x_gf);
391 double ke0 = oper->KineticEnergy(v_gf);
394 cout <<
"initial elastic energy (EE) = " << ee0 << endl;
395 cout <<
"initial kinetic energy (KE) = " << ke0 << endl;
396 cout <<
"initial total energy (TE) = " << (ee0 + ke0) << endl;
401 ode_solver->
Init(*oper);
405 bool last_step =
false;
406 for (
int ti = 1; !last_step; ti++)
408 double dt_real = min(dt, t_final -
t);
410 ode_solver->
Step(vx,
t, dt_real);
412 last_step = (
t >= t_final - 1e-8*dt);
414 if (last_step || (ti % vis_steps) == 0)
418 double ee = oper->ElasticEnergy(x_gf);
419 double ke = oper->KineticEnergy(v_gf);
423 cout <<
"step " << ti <<
", t = " <<
t <<
", EE = " << ee
424 <<
", KE = " << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
432 oper->GetElasticEnergyDensity(x_gf, w_gf);
446 ostringstream mesh_name, velo_name, ee_name;
447 mesh_name <<
"deformed." << setfill(
'0') << setw(6) << myid;
448 velo_name <<
"velocity." << setfill(
'0') << setw(6) << myid;
449 ee_name <<
"elastic_energy." << setfill(
'0') << setw(6) << myid;
451 ofstream mesh_ofs(mesh_name.str().c_str());
452 mesh_ofs.precision(8);
453 pmesh->
Print(mesh_ofs);
455 ofstream velo_ofs(velo_name.str().c_str());
456 velo_ofs.precision(8);
458 ofstream ee_ofs(ee_name.str().c_str());
460 oper->GetElasticEnergyDensity(x_gf, w_gf);
489 os <<
"solution\n" << *mesh << *field;
495 os <<
"window_size 800 800\n";
496 os <<
"window_title '" << field_name <<
"'\n";
503 os <<
"autoscale value\n";
510 ReducedSystemOperator::ReducedSystemOperator(
513 :
Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
514 Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height),
515 ess_tdof_list(ess_tdof_list_)
518 void ReducedSystemOperator::SetParameters(
double dt_,
const Vector *v_,
521 dt = dt_; v = v_; x = x_;
530 M->TrueAddMult(k, y);
531 S->TrueAddMult(w, y);
535 Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const 541 localJ->
Add(dt*dt, H->GetLocalGradient(z));
544 Jacobian = M->ParallelAssemble(localJ);
551 ReducedSystemOperator::~ReducedSystemOperator()
559 double mu,
double K,
bool use_petsc,
560 bool use_petsc_factory)
562 M(&fespace), S(&fespace), H(&fespace),
563 viscosity(visc), M_solver(
f.GetComm()),
564 newton_solver(
f.GetComm()), pnewton_solver(NULL), z(height/2)
566 const double rel_tol = 1e-8;
567 const int skip_zero_entries = 0;
569 const double ref_density = 1.0;
572 M.Assemble(skip_zero_entries);
573 M.Finalize(skip_zero_entries);
574 Mmat = M.ParallelAssemble();
575 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
579 M_solver.iterative_mode =
false;
580 M_solver.SetRelTol(rel_tol);
581 M_solver.SetAbsTol(0.0);
582 M_solver.SetMaxIter(30);
583 M_solver.SetPrintLevel(0);
584 M_prec.SetType(HypreSmoother::Jacobi);
585 M_solver.SetPreconditioner(M_prec);
586 M_solver.SetOperator(*Mmat);
590 H.SetEssentialTrueDofs(ess_tdof_list);
594 S.Assemble(skip_zero_entries);
595 S.Finalize(skip_zero_entries);
597 reduced_oper =
new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
601 J_hypreSmoother->
SetType(HypreSmoother::l1Jacobi);
603 J_prec = J_hypreSmoother;
616 newton_solver.SetSolver(*J_solver);
617 newton_solver.SetOperator(*reduced_oper);
618 newton_solver.SetPrintLevel(1);
619 newton_solver.SetRelTol(rel_tol);
620 newton_solver.SetAbsTol(0.0);
621 newton_solver.SetMaxIter(10);
634 if (use_petsc_factory)
636 J_factory =
new PreconditionerFactory(*reduced_oper,
"JFNK preconditioner");
637 pnewton_solver->SetPreconditionerFactory(J_factory);
639 pnewton_solver->SetPrintLevel(1);
640 pnewton_solver->SetRelTol(rel_tol);
641 pnewton_solver->SetAbsTol(0.0);
642 pnewton_solver->SetMaxIter(10);
656 if (viscosity != 0.0)
662 M_solver.Mult(z, dv_dt);
667 void HyperelasticOperator::ImplicitSolve(
const double dt,
682 reduced_oper->SetParameters(dt, &v, &x);
686 newton_solver.Mult(zero, dv_dt);
687 MFEM_VERIFY(newton_solver.GetConverged(),
688 "Newton solver did not converge.");
692 pnewton_solver->Mult(zero, dv_dt);
693 MFEM_VERIFY(pnewton_solver->GetConverged(),
694 "Newton solver did not converge.");
696 add(v, dt, dv_dt, dx_dt);
699 double HyperelasticOperator::ElasticEnergy(
const ParGridFunction &x)
const 701 return H.GetEnergy(x);
704 double HyperelasticOperator::KineticEnergy(
const ParGridFunction &v)
const 706 double loc_energy = 0.5*M.InnerProduct(v, v);
708 MPI_Allreduce(&loc_energy, &energy, 1, MPI_DOUBLE, MPI_SUM,
713 void HyperelasticOperator::GetElasticEnergyDensity(
716 ElasticEnergyCoefficient w_coeff(*model, x);
720 HyperelasticOperator::~HyperelasticOperator()
728 delete pnewton_solver;
748 model.SetTransformation(T);
751 return model.EvalW(J)/J.Det();
765 const double s = 0.1/64.;
768 v(
dim-1) =
s*x(0)*x(0)*(8.0-x(0));
void InitialDeformation(const Vector &x, Vector &y)
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. ...
Conjugate gradient method.
Class for grid function - Vector with associated FE space.
Abstract class for PETSc's preconditioners.
A class to handle Vectors in a block fashion.
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
A coefficient that is constant across space and time.
void PrintOptions(std::ostream &out) const
Print the options.
int TrueVSize() const
Obsolete, kept for backward compatibility.
int Dimension() const
Dimension of the reference space used within the elements.
Wrapper for PETSc's matrix class.
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Compute the vector gradient with respect to the physical element variable.
void PrintUsage(std::ostream &out) const
Print the usage message.
Base abstract class for first order time dependent operators.
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Swap the internal node GridFunction pointer and ownership flag members with the given ones...
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Pointer to an Operator of a specified type.
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]...
int Size() const
Returns the size of the vector.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Data type dense matrix using column-major storage.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
bool Good() const
Return true if the command line options were parsed successfully.
Abstract parallel finite element space.
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
int main(int argc, char *argv[])
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Backward Euler ODE solver. L-stable.
void add(const Vector &v1, const Vector &v2, Vector &v)
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
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.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
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)
Newton's method for solving F(x)=b for a given operator F.
Parallel smoothers in hypre.
HYPRE_BigInt GlobalTrueVSize() const
void Add(const int i, const int j, const double val)
Abstract class for PETSc's nonlinear solvers.
double * GetData() const
Return a pointer to the beginning of the Vector data.
A general vector function coefficient.
Mesh * GetMesh() const
Returns the mesh.
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)
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
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...
void visualize(ostream &os, ParMesh *mesh, ParGridFunction *deformed_nodes, ParGridFunction *field, const char *field_name=NULL, bool init_vis=false)
Third-order, strong stability preserving (SSP) Runge-Kutta method.
void InitialVelocity(const Vector &x, Vector &v)
int SpaceDimension() const
Dimension of the physical space containing the mesh.
virtual void Save(std::ostream &out) const
Implicit midpoint method. A-stable, not L-stable.
Class for integration point with weight.
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
Abstract class for hyperelastic models.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
Arbitrary order H1-conforming (continuous) finite elements.
void GetNodes(Vector &node_coord) const
void Print(std::ostream &out=mfem::out) const override
Class for parallel grid function.
The classical forward Euler method.
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Arbitrary order "L2-conforming" discontinuous finite elements.
double f(const Vector &p)