46class ReducedSystemOperator;
75 ReducedSystemOperator *reduced_oper;
92 virtual void Mult(
const Vector &vx,
Vector &dvx_dt)
const;
102 virtual ~HyperelasticOperator();
109class ReducedSystemOperator :
public Operator
133 virtual ~ReducedSystemOperator();
148 : model(m), x(x_) { }
150 virtual ~ElasticEnergyCoefficient() { }
160 bool init_vis =
false);
163int main(
int argc,
char *argv[])
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;
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 real_t ee0 = oper.ElasticEnergy(x_gf);
362 real_t 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 real_t 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 real_t ee = oper.ElasticEnergy(x_gf);
390 real_t 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";
480ReducedSystemOperator::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_)
488void ReducedSystemOperator::SetParameters(
real_t dt_,
const Vector *v_,
491 dt = dt_; v = v_; x = x_;
494void ReducedSystemOperator::Mult(
const Vector &k,
Vector &y)
const
500 M->TrueAddMult(k, y);
501 S->TrueAddMult(w, y);
505Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const
511 localJ->
Add(dt*dt, H->GetLocalGradient(z));
512 Jacobian = M->ParallelAssemble(localJ);
519ReducedSystemOperator::~ReducedSystemOperator()
529 M(&fespace), S(&fespace), H(&fespace),
530 viscosity(visc), M_solver(
f.GetComm()), newton_solver(
f.GetComm()),
533#if defined(MFEM_USE_DOUBLE)
534 const real_t rel_tol = 1e-8;
535 const real_t newton_abs_tol = 0.0;
536#elif defined(MFEM_USE_SINGLE)
537 const real_t rel_tol = 1e-3;
538 const real_t newton_abs_tol = 1e-4;
540#error "Only single and double precision are supported!"
544 const int skip_zero_entries = 0;
546 const real_t ref_density = 1.0;
549 M.Assemble(skip_zero_entries);
550 M.Finalize(skip_zero_entries);
551 Mmat = M.ParallelAssemble();
552 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
556 M_solver.iterative_mode =
false;
557 M_solver.SetRelTol(rel_tol);
558 M_solver.SetAbsTol(0.0);
559 M_solver.SetMaxIter(30);
560 M_solver.SetPrintLevel(0);
562 M_solver.SetPreconditioner(M_prec);
563 M_solver.SetOperator(*Mmat);
567 H.SetEssentialTrueDofs(ess_tdof_list);
571 S.Assemble(skip_zero_entries);
572 S.Finalize(skip_zero_entries);
574 reduced_oper =
new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
579 J_prec = J_hypreSmoother;
590 newton_solver.SetSolver(*J_solver);
591 newton_solver.SetOperator(*reduced_oper);
592 newton_solver.SetPrintLevel(1);
593 newton_solver.SetRelTol(rel_tol);
594 newton_solver.SetAbsTol(newton_abs_tol);
595 newton_solver.SetAdaptiveLinRtol(2, 0.5, 0.9);
596 newton_solver.SetMaxIter(10);
599void HyperelasticOperator::Mult(
const Vector &vx,
Vector &dvx_dt)
const
609 if (viscosity != 0.0)
615 M_solver.
Mult(z, dv_dt);
620void HyperelasticOperator::ImplicitSolve(
const real_t dt,
635 reduced_oper->SetParameters(dt, &v, &x);
637 newton_solver.
Mult(zero, dv_dt);
638 MFEM_VERIFY(newton_solver.
GetConverged(),
"Newton solver did not converge.");
639 add(v, dt, dv_dt, dx_dt);
649 real_t energy = 0.5*M.ParInnerProduct(v, v);
653void HyperelasticOperator::GetElasticEnergyDensity(
656 ElasticEnergyCoefficient w_coeff(*model, x);
660HyperelasticOperator::~HyperelasticOperator()
693 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.
virtual void Mult(const Vector &b, Vector &x) const
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.
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.
virtual void SetPreconditioner(Solver &pr)
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.
virtual void Mult(const Vector &b, Vector &x) const
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].
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
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 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)
void add(const Vector &v1, const Vector &v2, Vector &v)
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.