46 class ReducedSystemOperator;
75 ReducedSystemOperator *reduced_oper;
89 double visc,
double mu,
double K);
95 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
102 virtual ~HyperelasticOperator();
109 class ReducedSystemOperator :
public Operator 125 void SetParameters(
double dt_,
const Vector *v_,
const Vector *x_);
133 virtual ~ReducedSystemOperator();
139 class ElasticEnergyCoefficient :
public Coefficient 148 : model(m), x(x_) { }
150 virtual ~ElasticEnergyCoefficient() { }
160 bool init_vis =
false);
163 int main(
int argc,
char *argv[])
166 Mpi::Init(argc, argv);
167 int myid = Mpi::WorldRank();
171 const char *mesh_file =
"../data/beam-quad.mesh";
172 int ser_ref_levels = 2;
173 int par_ref_levels = 0;
175 int ode_solver_type = 3;
176 double t_final = 300.0;
181 bool adaptive_lin_rtol =
true;
182 bool visualization =
true;
186 args.
AddOption(&mesh_file,
"-m",
"--mesh",
187 "Mesh file to use.");
188 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
189 "Number of times to refine the mesh uniformly in serial.");
190 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
191 "Number of times to refine the mesh uniformly in parallel.");
193 "Order (degree) of the finite elements.");
194 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
195 "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t" 196 " 11 - Forward Euler, 12 - RK2,\n\t" 197 " 13 - RK3 SSP, 14 - RK4." 198 " 22 - Implicit Midpoint Method,\n\t" 199 " 23 - SDIRK23 (A-stable), 24 - SDIRK34");
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.");
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(&adaptive_lin_rtol,
"-alrtol",
"--adaptive-lin-rtol",
211 "-no-alrtol",
"--no-adaptive-lin-rtol",
212 "Enable or disable adaptive linear solver rtol.");
213 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
214 "--no-visualization",
215 "Enable or disable GLVis visualization.");
216 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
217 "Visualize every n-th timestep.");
235 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
242 switch (ode_solver_type)
250 case 12: ode_solver =
new RK2Solver(0.5);
break;
252 case 14: ode_solver =
new RK4Solver;
break;
261 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
270 for (
int lev = 0; lev < ser_ref_levels; lev++)
280 for (
int lev = 0; lev < par_ref_levels; lev++)
297 cout <<
"Number of velocity/deformation unknowns: " << glob_size << endl;
302 true_offset[1] = true_size;
303 true_offset[2] = 2*true_size;
307 v_gf.
MakeTRef(&fespace, vx, true_offset[0]);
308 x_gf.
MakeTRef(&fespace, vx, true_offset[1]);
334 HyperelasticOperator oper(fespace, ess_bdr, visc,
mu, K);
343 visualize(vis_v, pmesh, &x_gf, &v_gf,
"Velocity",
true);
350 oper.GetElasticEnergyDensity(x_gf, w_gf);
352 visualize(vis_w, pmesh, &x_gf, &w_gf,
"Elastic energy density",
true);
356 cout <<
"GLVis visualization paused." 357 <<
" Press space (in the GLVis window) to resume it.\n";
361 double ee0 = oper.ElasticEnergy(x_gf);
362 double ke0 = oper.KineticEnergy(v_gf);
365 cout <<
"initial elastic energy (EE) = " << ee0 << endl;
366 cout <<
"initial kinetic energy (KE) = " << ke0 << endl;
367 cout <<
"initial total energy (TE) = " << (ee0 + ke0) << endl;
372 ode_solver->
Init(oper);
376 bool last_step =
false;
377 for (
int ti = 1; !last_step; ti++)
379 double dt_real = min(dt, t_final -
t);
381 ode_solver->
Step(vx,
t, dt_real);
383 last_step = (
t >= t_final - 1e-8*dt);
385 if (last_step || (ti % vis_steps) == 0)
389 double ee = oper.ElasticEnergy(x_gf);
390 double ke = oper.KineticEnergy(v_gf);
394 cout <<
"step " << ti <<
", t = " <<
t <<
", EE = " << ee
395 <<
", KE = " << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
403 oper.GetElasticEnergyDensity(x_gf, w_gf);
417 ostringstream mesh_name, velo_name, ee_name;
418 mesh_name <<
"deformed." << setfill(
'0') << setw(6) << myid;
419 velo_name <<
"velocity." << setfill(
'0') << setw(6) << myid;
420 ee_name <<
"elastic_energy." << setfill(
'0') << setw(6) << myid;
422 ofstream mesh_ofs(mesh_name.str().c_str());
423 mesh_ofs.precision(8);
424 pmesh->
Print(mesh_ofs);
426 ofstream velo_ofs(velo_name.str().c_str());
427 velo_ofs.precision(8);
429 ofstream ee_ofs(ee_name.str().c_str());
431 oper.GetElasticEnergyDensity(x_gf, w_gf);
458 os <<
"solution\n" << *mesh << *field;
464 os <<
"window_size 800 800\n";
465 os <<
"window_title '" << field_name <<
"'\n";
473 os <<
"autoscale value\n";
480 ReducedSystemOperator::ReducedSystemOperator(
483 :
Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
484 Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height),
485 ess_tdof_list(ess_tdof_list_)
488 void ReducedSystemOperator::SetParameters(
double dt_,
const Vector *v_,
491 dt = dt_; v = v_; x = x_;
500 M->TrueAddMult(k, y);
501 S->TrueAddMult(w, y);
505 Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const 511 localJ->
Add(dt*dt, H->GetLocalGradient(z));
512 Jacobian = M->ParallelAssemble(localJ);
519 ReducedSystemOperator::~ReducedSystemOperator()
529 M(&fespace), S(&fespace), H(&fespace),
530 viscosity(visc), M_solver(
f.GetComm()), newton_solver(
f.GetComm()),
533 const double rel_tol = 1e-8;
534 const int skip_zero_entries = 0;
536 const double ref_density = 1.0;
539 M.Assemble(skip_zero_entries);
540 M.Finalize(skip_zero_entries);
541 Mmat = M.ParallelAssemble();
542 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
546 M_solver.iterative_mode =
false;
547 M_solver.SetRelTol(rel_tol);
548 M_solver.SetAbsTol(0.0);
549 M_solver.SetMaxIter(30);
550 M_solver.SetPrintLevel(0);
551 M_prec.SetType(HypreSmoother::Jacobi);
552 M_solver.SetPreconditioner(M_prec);
553 M_solver.SetOperator(*Mmat);
557 H.SetEssentialTrueDofs(ess_tdof_list);
561 S.Assemble(skip_zero_entries);
562 S.Finalize(skip_zero_entries);
564 reduced_oper =
new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
567 J_hypreSmoother->
SetType(HypreSmoother::l1Jacobi);
569 J_prec = J_hypreSmoother;
580 newton_solver.SetSolver(*J_solver);
581 newton_solver.SetOperator(*reduced_oper);
582 newton_solver.SetPrintLevel(1);
583 newton_solver.SetRelTol(rel_tol);
584 newton_solver.SetAbsTol(0.0);
585 newton_solver.SetAdaptiveLinRtol(2, 0.5, 0.9);
586 newton_solver.SetMaxIter(10);
599 if (viscosity != 0.0)
605 M_solver.Mult(z, dv_dt);
610 void HyperelasticOperator::ImplicitSolve(
const double dt,
625 reduced_oper->SetParameters(dt, &v, &x);
627 newton_solver.Mult(zero, dv_dt);
628 MFEM_VERIFY(newton_solver.GetConverged(),
"Newton solver did not converge.");
629 add(v, dt, dv_dt, dx_dt);
632 double HyperelasticOperator::ElasticEnergy(
const ParGridFunction &x)
const 634 return H.GetEnergy(x);
637 double HyperelasticOperator::KineticEnergy(
const ParGridFunction &v)
const 639 double loc_energy = 0.5*M.InnerProduct(v, v);
641 MPI_Allreduce(&loc_energy, &energy, 1, MPI_DOUBLE, MPI_SUM,
646 void HyperelasticOperator::GetElasticEnergyDensity(
649 ElasticEnergyCoefficient w_coeff(*model, x);
653 HyperelasticOperator::~HyperelasticOperator()
666 model.SetTransformation(T);
669 return model.EvalW(J)/J.Det();
683 const double s = 0.1/64.;
686 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.
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.
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)
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.
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)
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
virtual void Save(std::ostream &out) const
Implicit midpoint method. A-stable, not L-stable.
Class for integration point with weight.
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)