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 num_procs = Mpi::WorldSize();
191 int myid = Mpi::WorldRank();
195 const char *mesh_file =
"../../data/beam-quad.mesh";
196 int ser_ref_levels = 2;
197 int par_ref_levels = 0;
199 int ode_solver_type = 3;
200 double t_final = 300.0;
205 bool visualization =
true;
207 bool use_petsc =
true;
208 const char *petscrc_file =
"";
209 bool petsc_use_jfnk =
false;
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.");
230 args.
AddOption(&mu,
"-mu",
"--shear-modulus",
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.");
270 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
277 switch (ode_solver_type)
285 case 12: ode_solver =
new RK2Solver(0.5);
break;
287 case 14: ode_solver =
new RK4Solver;
break;
295 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
303 for (
int lev = 0; lev < ser_ref_levels; lev++)
313 for (
int lev = 0; lev < par_ref_levels; lev++)
330 cout <<
"Number of velocity/deformation unknowns: " << glob_size << endl;
335 true_offset[1] = true_size;
336 true_offset[2] = 2*true_size;
340 v_gf.
MakeTRef(&fespace, vx, true_offset[0]);
341 x_gf.
MakeTRef(&fespace, vx, true_offset[1]);
367 HyperelasticOperator *oper =
new HyperelasticOperator(fespace, ess_bdr, visc,
376 vis_v.
open(vishost, visport);
378 visualize(vis_v, pmesh, &x_gf, &v_gf,
"Velocity",
true);
382 vis_w.
open(vishost, visport);
385 oper->GetElasticEnergyDensity(x_gf, w_gf);
387 visualize(vis_w, pmesh, &x_gf, &w_gf,
"Elastic energy density",
true);
391 double ee0 = oper->ElasticEnergy(x_gf);
392 double ke0 = oper->KineticEnergy(v_gf);
395 cout <<
"initial elastic energy (EE) = " << ee0 << endl;
396 cout <<
"initial kinetic energy (KE) = " << ke0 << endl;
397 cout <<
"initial total energy (TE) = " << (ee0 + ke0) << endl;
402 ode_solver->
Init(*oper);
406 bool last_step =
false;
407 for (
int ti = 1; !last_step; ti++)
409 double dt_real = min(dt, t_final - t);
411 ode_solver->
Step(vx, t, dt_real);
413 last_step = (t >= t_final - 1e-8*dt);
415 if (last_step || (ti % vis_steps) == 0)
419 double ee = oper->ElasticEnergy(x_gf);
420 double ke = oper->KineticEnergy(v_gf);
424 cout <<
"step " << ti <<
", t = " << t <<
", EE = " << ee
425 <<
", KE = " << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
433 oper->GetElasticEnergyDensity(x_gf, w_gf);
447 ostringstream mesh_name, velo_name, ee_name;
448 mesh_name <<
"deformed." << setfill(
'0') << setw(6) << myid;
449 velo_name <<
"velocity." << setfill(
'0') << setw(6) << myid;
450 ee_name <<
"elastic_energy." << setfill(
'0') << setw(6) << myid;
452 ofstream mesh_ofs(mesh_name.str().c_str());
453 mesh_ofs.precision(8);
454 pmesh->
Print(mesh_ofs);
456 ofstream velo_ofs(velo_name.str().c_str());
457 velo_ofs.precision(8);
459 ofstream ee_ofs(ee_name.str().c_str());
461 oper->GetElasticEnergyDensity(x_gf, w_gf);
490 os <<
"solution\n" << *mesh << *field;
496 os <<
"window_size 800 800\n";
497 os <<
"window_title '" << field_name <<
"'\n";
504 os <<
"autoscale value\n";
511 ReducedSystemOperator::ReducedSystemOperator(
514 :
Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
515 Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height),
516 ess_tdof_list(ess_tdof_list_)
519 void ReducedSystemOperator::SetParameters(
double dt_,
const Vector *v_,
522 dt = dt_; v = v_; x = x_;
531 M->TrueAddMult(k, y);
532 S->TrueAddMult(w, y);
536 Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const
542 localJ->
Add(dt*dt, H->GetLocalGradient(z));
545 Jacobian = M->ParallelAssemble(localJ);
552 ReducedSystemOperator::~ReducedSystemOperator()
560 double mu,
double K,
bool use_petsc,
561 bool use_petsc_factory)
563 M(&fespace), S(&fespace), H(&fespace),
564 viscosity(visc), M_solver(f.GetComm()),
565 newton_solver(f.GetComm()), pnewton_solver(NULL), z(height/2)
567 const double rel_tol = 1e-8;
568 const int skip_zero_entries = 0;
570 const double ref_density = 1.0;
573 M.Assemble(skip_zero_entries);
574 M.Finalize(skip_zero_entries);
575 Mmat = M.ParallelAssemble();
576 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
580 M_solver.iterative_mode =
false;
581 M_solver.SetRelTol(rel_tol);
582 M_solver.SetAbsTol(0.0);
583 M_solver.SetMaxIter(30);
584 M_solver.SetPrintLevel(0);
585 M_prec.SetType(HypreSmoother::Jacobi);
586 M_solver.SetPreconditioner(M_prec);
587 M_solver.SetOperator(*Mmat);
591 H.SetEssentialTrueDofs(ess_tdof_list);
595 S.Assemble(skip_zero_entries);
596 S.Finalize(skip_zero_entries);
598 reduced_oper =
new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
602 J_hypreSmoother->
SetType(HypreSmoother::l1Jacobi);
604 J_prec = J_hypreSmoother;
617 newton_solver.SetSolver(*J_solver);
618 newton_solver.SetOperator(*reduced_oper);
619 newton_solver.SetPrintLevel(1);
620 newton_solver.SetRelTol(rel_tol);
621 newton_solver.SetAbsTol(0.0);
622 newton_solver.SetMaxIter(10);
635 if (use_petsc_factory)
637 J_factory =
new PreconditionerFactory(*reduced_oper,
"JFNK preconditioner");
638 pnewton_solver->SetPreconditionerFactory(J_factory);
640 pnewton_solver->SetPrintLevel(1);
641 pnewton_solver->SetRelTol(rel_tol);
642 pnewton_solver->SetAbsTol(0.0);
643 pnewton_solver->SetMaxIter(10);
657 if (viscosity != 0.0)
663 M_solver.Mult(z, dv_dt);
668 void HyperelasticOperator::ImplicitSolve(
const double dt,
683 reduced_oper->SetParameters(dt, &v, &x);
687 newton_solver.Mult(zero, dv_dt);
688 MFEM_VERIFY(newton_solver.GetConverged(),
689 "Newton solver did not converge.");
693 pnewton_solver->Mult(zero, dv_dt);
694 MFEM_VERIFY(pnewton_solver->GetConverged(),
695 "Newton solver did not converge.");
697 add(v, dt, dv_dt, dx_dt);
700 double HyperelasticOperator::ElasticEnergy(
const ParGridFunction &x)
const
702 return H.GetEnergy(x);
705 double HyperelasticOperator::KineticEnergy(
const ParGridFunction &v)
const
707 double loc_energy = 0.5*M.InnerProduct(v, v);
709 MPI_Allreduce(&loc_energy, &energy, 1, MPI_DOUBLE, MPI_SUM,
714 void HyperelasticOperator::GetElasticEnergyDensity(
717 ElasticEnergyCoefficient w_coeff(*model, x);
721 HyperelasticOperator::~HyperelasticOperator()
729 delete pnewton_solver;
749 model.SetTransformation(T);
752 return model.EvalW(J)/J.Det();
766 const double s = 0.1/64.;
769 v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
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 InitialDeformation(const Vector &x, Vector &y)
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.
Wrapper for PETSc's matrix class.
int TrueVSize() const
Obsolete, kept for backward compatibility.
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)
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]...
HYPRE_BigInt GlobalTrueVSize() const
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
Abstract parallel finite element space.
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.
double * GetData() const
Return a pointer to the beginning of the Vector data.
void add(const Vector &v1, const Vector &v2, Vector &v)
void InitialVelocity(const Vector &x, 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 ...
Mesh * GetMesh() const
Returns the mesh.
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.
void visualize(ostream &os, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
Newton's method for solving F(x)=b for a given operator F.
Parallel smoothers in hypre.
void PrintUsage(std::ostream &out) const
Print the usage message.
void Add(const int i, const int j, const double val)
Abstract class for PETSc's nonlinear solvers.
A general vector function coefficient.
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)
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...
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 MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
void PrintOptions(std::ostream &out) const
Print the options.
Abstract class for hyperelastic models.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
void GetNodes(Vector &node_coord) const
Arbitrary order H1-conforming (continuous) finite elements.
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)
bool Good() const
Return true if the command line options were parsed successfully.