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.");
191 args.
AddOption(&mu,
"-mu",
"--shear-modulus",
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);
297 vis_v.
open(vishost, visport);
300 visualize(vis_v, mesh, &x, &v,
"Velocity",
true);
301 vis_w.
open(vishost, visport);
304 oper.GetElasticEnergyDensity(x, w);
306 visualize(vis_w, mesh, &x, &w,
"Elastic energy density",
true);
312 cout <<
"initial elastic energy (EE) = " << ee0 << endl;
313 cout <<
"initial kinetic energy (KE) = " << ke0 << endl;
314 cout <<
"initial total energy (TE) = " << (ee0 + ke0) << endl;
318 ode_solver->
Init(oper);
322 bool last_step =
false;
323 for (
int ti = 1; !last_step; ti++)
325 double dt_real = min(dt, t_final - t);
327 ode_solver->
Step(vx, t, dt_real);
329 last_step = (t >= t_final - 1e-8*dt);
331 if (last_step || (ti % vis_steps) == 0)
336 cout <<
"step " << ti <<
", t = " << t <<
", EE = " << ee <<
", KE = "
337 << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
345 oper.GetElasticEnergyDensity(x, w);
358 ofstream mesh_ofs(
"deformed.mesh");
359 mesh_ofs.precision(8);
360 mesh->
Print(mesh_ofs);
362 ofstream velo_ofs(
"velocity.sol");
363 velo_ofs.precision(8);
365 ofstream ee_ofs(
"elastic_energy.sol");
367 oper.GetElasticEnergyDensity(x, w);
380 GridFunction *field,
const char *field_name,
bool init_vis)
392 out <<
"solution\n" << *mesh << *field;
398 out <<
"window_size 800 800\n";
399 out <<
"window_title '" << field_name <<
"'\n";
406 out <<
"autoscale value\n";
413 ReducedSystemOperator::ReducedSystemOperator(
415 :
Operator(M_->Height()), M(M_), S(S_), H(H_), Jacobian(NULL),
416 dt(0.0), v(NULL), x(NULL), w(height), z(height)
419 void ReducedSystemOperator::SetParameters(
double dt_,
const Vector *v_,
422 dt = dt_; v = v_; x = x_;
435 Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const
438 Jacobian =
Add(1.0, M->SpMat(), dt, S->SpMat());
442 Jacobian->
Add(dt*dt, *grad_H);
446 ReducedSystemOperator::~ReducedSystemOperator()
456 M(&fespace), S(&fespace), H(&fespace),
457 viscosity(visc), z(height/2)
459 const double rel_tol = 1e-8;
460 const int skip_zero_entries = 0;
462 const double ref_density = 1.0;
465 M.Assemble(skip_zero_entries);
467 fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
469 M.FormSystemMatrix(ess_tdof_list, tmp);
471 M_solver.iterative_mode =
false;
472 M_solver.SetRelTol(rel_tol);
473 M_solver.SetAbsTol(0.0);
474 M_solver.SetMaxIter(30);
475 M_solver.SetPrintLevel(0);
476 M_solver.SetPreconditioner(M_prec);
477 M_solver.SetOperator(M.SpMat());
481 H.SetEssentialTrueDofs(ess_tdof_list);
485 S.Assemble(skip_zero_entries);
486 S.FormSystemMatrix(ess_tdof_list, tmp);
488 reduced_oper =
new ReducedSystemOperator(&M, &S, &H);
490 #ifndef MFEM_USE_SUITESPARSE
505 newton_solver.SetSolver(*J_solver);
506 newton_solver.SetOperator(*reduced_oper);
507 newton_solver.SetPrintLevel(1);
508 newton_solver.SetRelTol(rel_tol);
509 newton_solver.SetAbsTol(0.0);
510 newton_solver.SetMaxIter(10);
523 if (viscosity != 0.0)
528 M_solver.Mult(z, dv_dt);
533 void HyperelasticOperator::ImplicitSolve(
const double dt,
548 reduced_oper->SetParameters(dt, &v, &x);
550 newton_solver.Mult(zero, dv_dt);
551 MFEM_VERIFY(newton_solver.GetConverged(),
"Newton solver did not converge.");
552 add(v, dt, dv_dt, dx_dt);
555 double HyperelasticOperator::ElasticEnergy(
const Vector &x)
const
557 return H.GetEnergy(x);
560 double HyperelasticOperator::KineticEnergy(
const Vector &v)
const
562 return 0.5*M.InnerProduct(v, v);
565 void HyperelasticOperator::GetElasticEnergyDensity(
568 ElasticEnergyCoefficient w_coeff(*model, x);
572 HyperelasticOperator::~HyperelasticOperator()
584 model.SetTransformation(T);
587 return model.EvalW(J)/J.Det();
601 const double s = 0.1/64.;
604 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)
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 Add(const int i, const int j, const double a)
void InitialDeformation(const Vector &x, Vector &y)
virtual void Print(std::ostream &out=mfem::out) const
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.
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)
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)
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)
Direct sparse solver using UMFPACK.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
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.
void SetPrintLevel(int print_lvl)
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.
Newton's method for solving F(x)=b for a given operator F.
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
void PrintUsage(std::ostream &out) const
Print the usage message.
const Vector & GetTrueVector() const
Read only access to the (optional) internal true-dof Vector.
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)
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.
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Implicit midpoint method. A-stable, not L-stable.
Class for integration point with weight.
void PrintOptions(std::ostream &out) const
Print the options.
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'.
void GetNodes(Vector &node_coord) const
Arbitrary order H1-conforming (continuous) finite elements.
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.
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)
bool Good() const
Return true if the command line options were parsed successfully.