40 #ifndef MFEM_USE_PETSC
41 #error This example requires that MFEM is built with MFEM_USE_PETSC=YES
47 class ReducedSystemOperator;
75 ReducedSystemOperator *reduced_oper;
92 double visc,
double mu,
double K,
bool use_petsc);
98 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
104 virtual ~HyperelasticOperator();
111 class ReducedSystemOperator :
public Operator
126 void SetParameters(
double dt_,
const Vector *v_,
const Vector *x_);
134 virtual ~ReducedSystemOperator();
140 class ElasticEnergyCoefficient :
public Coefficient
149 : model(m), x(x_) { }
151 virtual ~ElasticEnergyCoefficient() { }
160 bool init_vis =
false);
163 int main(
int argc,
char *argv[])
167 MPI_Init(&argc, &argv);
168 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
169 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
172 const char *mesh_file =
"../../data/beam-quad.mesh";
173 int ser_ref_levels = 2;
174 int par_ref_levels = 0;
176 int ode_solver_type = 3;
177 double t_final = 300.0;
182 bool visualization =
true;
184 bool use_petsc =
true;
185 const char *petscrc_file =
"";
188 args.
AddOption(&mesh_file,
"-m",
"--mesh",
189 "Mesh file to use.");
190 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
191 "Number of times to refine the mesh uniformly in serial.");
192 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
193 "Number of times to refine the mesh uniformly in parallel.");
195 "Order (degree) of the finite elements.");
196 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
197 "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
198 " 11 - Forward Euler, 12 - RK2,\n\t"
199 " 13 - RK3 SSP, 14 - RK4.");
200 args.
AddOption(&t_final,
"-tf",
"--t-final",
201 "Final time; start time is 0.");
202 args.
AddOption(&dt,
"-dt",
"--time-step",
204 args.
AddOption(&visc,
"-v",
"--viscosity",
205 "Viscosity coefficient.");
206 args.
AddOption(&mu,
"-mu",
"--shear-modulus",
207 "Shear modulus in the Neo-Hookean hyperelastic model.");
208 args.
AddOption(&K,
"-K",
"--bulk-modulus",
209 "Bulk modulus in the Neo-Hookean hyperelastic model.");
210 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
211 "--no-visualization",
212 "Enable or disable GLVis visualization.");
213 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
214 "Visualize every n-th timestep.");
215 args.
AddOption(&use_petsc,
"-usepetsc",
"--usepetsc",
"-no-petsc",
217 "Use or not PETSc to solve the nonlinear system.");
218 args.
AddOption(&petscrc_file,
"-petscopts",
"--petscopts",
219 "PetscOptions file to use.");
238 PetscInitialize(NULL,NULL,petscrc_file,NULL);
244 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
251 switch (ode_solver_type)
259 case 12: ode_solver =
new RK2Solver(0.5);
break;
261 case 14: ode_solver =
new RK4Solver;
break;
269 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
278 for (
int lev = 0; lev < ser_ref_levels; lev++)
288 for (
int lev = 0; lev < par_ref_levels; lev++)
305 cout <<
"Number of velocity/deformation unknowns: " << glob_size << endl;
310 true_offset[1] = true_size;
311 true_offset[2] = 2*true_size;
326 v_gf.ProjectCoefficient(velo);
339 HyperelasticOperator *oper =
new HyperelasticOperator(fespace, ess_bdr, visc,
345 char vishost[] =
"localhost";
347 vis_v.
open(vishost, visport);
349 visualize(vis_v, pmesh, &x_gf, &v_gf,
"Velocity",
true);
353 vis_w.
open(vishost, visport);
356 oper->GetElasticEnergyDensity(x_gf, w_gf);
358 visualize(vis_w, pmesh, &x_gf, &w_gf,
"Elastic energy density",
true);
362 double ee0 = oper->ElasticEnergy(x_gf);
363 double ke0 = oper->KineticEnergy(v_gf);
366 cout <<
"initial elastic energy (EE) = " << ee0 << endl;
367 cout <<
"initial kinetic energy (KE) = " << ke0 << endl;
368 cout <<
"initial total energy (TE) = " << (ee0 + ke0) << endl;
373 ode_solver->
Init(*oper);
377 bool last_step =
false;
378 for (
int ti = 1; !last_step; ti++)
380 double dt_real = min(dt, t_final - t);
382 ode_solver->
Step(vx, t, dt_real);
384 last_step = (t >= t_final - 1e-8*dt);
386 if (last_step || (ti % vis_steps) == 0)
391 double ee = oper->ElasticEnergy(x_gf);
392 double ke = oper->KineticEnergy(v_gf);
396 cout <<
"step " << ti <<
", t = " << t <<
", EE = " << ee
397 <<
", KE = " << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
405 oper->GetElasticEnergyDensity(x_gf, w_gf);
418 ostringstream mesh_name, velo_name, ee_name;
419 mesh_name <<
"deformed." << setfill(
'0') << setw(6) << myid;
420 velo_name <<
"velocity." << setfill(
'0') << setw(6) << myid;
421 ee_name <<
"elastic_energy." << setfill(
'0') << setw(6) << myid;
423 ofstream mesh_ofs(mesh_name.str().c_str());
424 mesh_ofs.precision(8);
425 pmesh->
Print(mesh_ofs);
427 ofstream velo_ofs(velo_name.str().c_str());
428 velo_ofs.precision(8);
430 ofstream ee_ofs(ee_name.str().c_str());
432 oper->GetElasticEnergyDensity(x_gf, w_gf);
442 if (use_petsc) { PetscFinalize(); }
463 out <<
"solution\n" << *mesh << *field;
469 out <<
"window_size 800 800\n";
470 out <<
"window_title '" << field_name <<
"'\n";
477 out <<
"autoscale value\n";
484 ReducedSystemOperator::ReducedSystemOperator(
486 :
Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
487 Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height)
490 void ReducedSystemOperator::SetParameters(
double dt_,
const Vector *v_,
493 dt = dt_; v = v_; x = x_;
502 M->TrueAddMult(k, y);
503 S->TrueAddMult(w, y);
506 Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const
512 localJ->
Add(dt*dt, H->GetLocalGradient(z));
515 Jacobian = M->ParallelAssemble(localJ);
520 ReducedSystemOperator::~ReducedSystemOperator()
528 double mu,
double K,
bool use_petsc)
530 M(&fespace), S(&fespace), H(&fespace),
531 viscosity(visc), M_solver(f.GetComm()),
532 newton_solver(f.GetComm()), pnewton_solver(NULL), z(height/2)
534 const double rel_tol = 1e-8;
535 const int skip_zero_entries = 0;
537 const double ref_density = 1.0;
540 M.Assemble(skip_zero_entries);
541 M.EliminateEssentialBC(ess_bdr);
542 M.Finalize(skip_zero_entries);
543 Mmat = M.ParallelAssemble();
545 M_solver.iterative_mode =
false;
546 M_solver.SetRelTol(rel_tol);
547 M_solver.SetAbsTol(0.0);
548 M_solver.SetMaxIter(30);
549 M_solver.SetPrintLevel(0);
550 M_prec.SetType(HypreSmoother::Jacobi);
551 M_solver.SetPreconditioner(M_prec);
552 M_solver.SetOperator(*Mmat);
556 H.SetEssentialBC(ess_bdr);
560 S.Assemble(skip_zero_entries);
561 S.EliminateEssentialBC(ess_bdr);
562 S.Finalize(skip_zero_entries);
564 reduced_oper =
new ReducedSystemOperator(&M, &S, &H);
568 J_hypreSmoother->
SetType(HypreSmoother::l1Jacobi);
570 J_prec = J_hypreSmoother;
581 newton_solver.SetSolver(*J_solver);
582 newton_solver.SetOperator(*reduced_oper);
583 newton_solver.SetPrintLevel(1);
584 newton_solver.SetRelTol(rel_tol);
585 newton_solver.SetAbsTol(0.0);
586 newton_solver.SetMaxIter(10);
597 pnewton_solver->SetRelTol(rel_tol);
598 pnewton_solver->SetAbsTol(0.0);
599 pnewton_solver->SetMaxIter(10);
613 if (viscosity != 0.0)
618 M_solver.Mult(z, dv_dt);
623 void HyperelasticOperator::ImplicitSolve(
const double dt,
638 reduced_oper->SetParameters(dt, &v, &x);
642 newton_solver.Mult(zero, dv_dt);
643 MFEM_VERIFY(newton_solver.GetConverged(),
644 "Newton solver did not converge.");
648 pnewton_solver->Mult(zero, dv_dt);
649 MFEM_VERIFY(pnewton_solver->GetConverged(),
650 "Newton solver did not converge.");
652 add(v, dt, dv_dt, dx_dt);
657 return H.GetEnergy(x);
662 double loc_energy = 0.5*M.InnerProduct(v, v);
664 MPI_Allreduce(&loc_energy, &energy, 1, MPI_DOUBLE, MPI_SUM,
669 void HyperelasticOperator::GetElasticEnergyDensity(
672 ElasticEnergyCoefficient w_coeff(*model, x);
676 HyperelasticOperator::~HyperelasticOperator()
683 delete pnewton_solver;
690 model.SetTransformation(T);
693 return model.EvalW(J)/J.Det();
707 const double s = 0.1/64.;
710 v(dim-1) = s*x(0)*x(0)*(8.0-x(0));
void visualize(ostream &out, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
void SetPrintLevel(int plev)
double Eval(ElementTransformation &T, const IntegrationPoint &ip, double 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.
Subclass constant coefficient.
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Base abstract class for 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]...
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)
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 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.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
void SetPrintLevel(int print_lvl)
Mesh * GetMesh() const
Returns the mesh.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
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
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
Abstract class for PETSc's nonlinear solvers.
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 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 Distribute(const Vector *tv)
Implicit midpoint method. A-stable, not L-stable.
Class for integration point with weight.
void PrintOptions(std::ostream &out) const
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
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.
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 GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad)
Vector & GetBlock(int i)
Get the i-th vector in the block.
Class for parallel meshes.
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Arbitrary order "L2-conforming" discontinuous finite elements.
void Neg()
(*this) = -(*this)