46 class ReducedSystemOperator;
73 ReducedSystemOperator *reduced_oper;
87 double visc,
double mu,
double K);
93 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
95 double ElasticEnergy(
const Vector &x)
const;
96 double KineticEnergy(
const Vector &v)
const;
99 virtual ~HyperelasticOperator();
106 class ReducedSystemOperator :
public Operator 120 void SetParameters(
double dt_,
const Vector *v_,
const Vector *x_);
128 virtual ~ReducedSystemOperator();
134 class ElasticEnergyCoefficient :
public Coefficient 143 : model(m), x(x_) { }
145 virtual ~ElasticEnergyCoefficient() { }
154 bool init_vis =
false);
157 int main(
int argc,
char *argv[])
160 const char *mesh_file =
"../data/beam-quad.mesh";
163 int ode_solver_type = 3;
164 double t_final = 300.0;
169 bool visualization =
true;
173 args.
AddOption(&mesh_file,
"-m",
"--mesh",
174 "Mesh file to use.");
175 args.
AddOption(&ref_levels,
"-r",
"--refine",
176 "Number of times to refine the mesh uniformly.");
178 "Order (degree) of the finite elements.");
179 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
180 "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t" 181 " 11 - Forward Euler, 12 - RK2,\n\t" 182 " 13 - RK3 SSP, 14 - RK4." 183 " 22 - Implicit Midpoint Method,\n\t" 184 " 23 - SDIRK23 (A-stable), 24 - SDIRK34");
185 args.
AddOption(&t_final,
"-tf",
"--t-final",
186 "Final time; start time is 0.");
187 args.
AddOption(&dt,
"-dt",
"--time-step",
189 args.
AddOption(&visc,
"-v",
"--viscosity",
190 "Viscosity coefficient.");
192 "Shear modulus in the Neo-Hookean hyperelastic model.");
193 args.
AddOption(&K,
"-K",
"--bulk-modulus",
194 "Bulk modulus in the Neo-Hookean hyperelastic model.");
195 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
196 "--no-visualization",
197 "Enable or disable GLVis visualization.");
198 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
199 "Visualize every n-th timestep.");
210 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
217 switch (ode_solver_type)
225 case 12: ode_solver =
new RK2Solver(0.5);
break;
227 case 14: ode_solver =
new RK4Solver;
break;
234 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
242 for (
int lev = 0; lev < ref_levels; lev++)
257 cout <<
"Number of velocity/deformation unknowns: " << fe_size << endl;
260 fe_offset[1] = fe_size;
261 fe_offset[2] = 2*fe_size;
290 HyperelasticOperator oper(fespace, ess_bdr, visc,
mu, K);
300 visualize(vis_v, mesh, &x, &v,
"Velocity",
true);
304 oper.GetElasticEnergyDensity(x, w);
306 visualize(vis_w, mesh, &x, &w,
"Elastic energy density",
true);
308 cout <<
"GLVis visualization paused." 309 <<
" Press space (in the GLVis window) to resume it.\n";
314 cout <<
"initial elastic energy (EE) = " << ee0 << endl;
315 cout <<
"initial kinetic energy (KE) = " << ke0 << endl;
316 cout <<
"initial total energy (TE) = " << (ee0 + ke0) << endl;
320 ode_solver->
Init(oper);
324 bool last_step =
false;
325 for (
int ti = 1; !last_step; ti++)
327 double dt_real = min(dt, t_final -
t);
329 ode_solver->
Step(vx,
t, dt_real);
331 last_step = (
t >= t_final - 1e-8*dt);
333 if (last_step || (ti % vis_steps) == 0)
338 cout <<
"step " << ti <<
", t = " <<
t <<
", EE = " << ee <<
", KE = " 339 << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
347 oper.GetElasticEnergyDensity(x, w);
360 ofstream mesh_ofs(
"deformed.mesh");
361 mesh_ofs.precision(8);
362 mesh->
Print(mesh_ofs);
364 ofstream velo_ofs(
"velocity.sol");
365 velo_ofs.precision(8);
367 ofstream ee_ofs(
"elastic_energy.sol");
369 oper.GetElasticEnergyDensity(x, w);
382 GridFunction *field,
const char *field_name,
bool init_vis)
394 os <<
"solution\n" << *mesh << *field;
400 os <<
"window_size 800 800\n";
401 os <<
"window_title '" << field_name <<
"'\n";
409 os <<
"autoscale value\n";
416 ReducedSystemOperator::ReducedSystemOperator(
418 :
Operator(M_->Height()), M(M_), S(S_), H(H_), Jacobian(NULL),
419 dt(0.0), v(NULL), x(NULL), w(height), z(height)
422 void ReducedSystemOperator::SetParameters(
double dt_,
const Vector *v_,
425 dt = dt_; v = v_; x = x_;
438 Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const 441 Jacobian =
Add(1.0, M->SpMat(), dt, S->SpMat());
445 Jacobian->
Add(dt*dt, *grad_H);
449 ReducedSystemOperator::~ReducedSystemOperator()
459 M(&fespace), S(&fespace), H(&fespace),
460 viscosity(visc), z(height/2)
462 const double rel_tol = 1e-8;
463 const int skip_zero_entries = 0;
465 const double ref_density = 1.0;
468 M.Assemble(skip_zero_entries);
470 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
472 M.FormSystemMatrix(ess_tdof_list, tmp);
474 M_solver.iterative_mode =
false;
475 M_solver.SetRelTol(rel_tol);
476 M_solver.SetAbsTol(0.0);
477 M_solver.SetMaxIter(30);
478 M_solver.SetPrintLevel(0);
479 M_solver.SetPreconditioner(M_prec);
480 M_solver.SetOperator(M.SpMat());
484 H.SetEssentialTrueDofs(ess_tdof_list);
488 S.Assemble(skip_zero_entries);
489 S.FormSystemMatrix(ess_tdof_list, tmp);
491 reduced_oper =
new ReducedSystemOperator(&M, &S, &H);
493 #ifndef MFEM_USE_SUITESPARSE 508 newton_solver.SetSolver(*J_solver);
509 newton_solver.SetOperator(*reduced_oper);
510 newton_solver.SetPrintLevel(1);
511 newton_solver.SetRelTol(rel_tol);
512 newton_solver.SetAbsTol(0.0);
513 newton_solver.SetMaxIter(10);
526 if (viscosity != 0.0)
531 M_solver.Mult(z, dv_dt);
536 void HyperelasticOperator::ImplicitSolve(
const double dt,
551 reduced_oper->SetParameters(dt, &v, &x);
553 newton_solver.Mult(zero, dv_dt);
554 MFEM_VERIFY(newton_solver.GetConverged(),
"Newton solver did not converge.");
555 add(v, dt, dv_dt, dx_dt);
558 double HyperelasticOperator::ElasticEnergy(
const Vector &x)
const 560 return H.GetEnergy(x);
563 double HyperelasticOperator::KineticEnergy(
const Vector &v)
const 565 return 0.5*M.InnerProduct(v, v);
568 void HyperelasticOperator::GetElasticEnergyDensity(
571 ElasticEnergyCoefficient w_coeff(*model, x);
575 HyperelasticOperator::~HyperelasticOperator()
587 model.SetTransformation(T);
590 return model.EvalW(J)/J.Det();
604 const double s = 0.1/64.;
607 v(
dim-1) =
s*x(0)*x(0)*(8.0-x(0));
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.
Data type for scaled Jacobi-type smoother of sparse matrix.
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 Dimension() const
Dimension of the reference space used within the elements.
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.
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)
Direct sparse solver using UMFPACK.
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)
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.
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.
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
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)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
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.
int main(int argc, char *argv[])
int SpaceDimension() const
Dimension of the physical space containing the mesh.
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
Implicit midpoint method. A-stable, not L-stable.
Class for integration point with weight.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Abstract class for hyperelastic models.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
virtual void Print(std::ostream &os=mfem::out) const
Arbitrary order H1-conforming (continuous) finite elements.
void GetNodes(Vector &node_coord) const
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
The classical forward Euler method.
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Vector & GetBlock(int i)
Get the i-th vector in the block.
Arbitrary order "L2-conforming" discontinuous finite elements.
double f(const Vector &p)
void Neg()
(*this) = -(*this)