44 class ReducedSystemOperator;
71 ReducedSystemOperator *reduced_oper;
85 double visc,
double mu,
double K);
91 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
93 double ElasticEnergy(
Vector &x)
const;
94 double KineticEnergy(
Vector &v)
const;
97 virtual ~HyperelasticOperator();
104 class ReducedSystemOperator :
public Operator
118 void SetParameters(
double dt_,
const Vector *v_,
const Vector *x_);
126 virtual ~ReducedSystemOperator();
132 class ElasticEnergyCoefficient :
public Coefficient
141 : model(m), x(x_) { }
143 virtual ~ElasticEnergyCoefficient() { }
152 bool init_vis =
false);
155 int main(
int argc,
char *argv[])
158 const char *mesh_file =
"../data/beam-quad.mesh";
161 int ode_solver_type = 3;
162 double t_final = 300.0;
167 bool visualization =
true;
171 args.
AddOption(&mesh_file,
"-m",
"--mesh",
172 "Mesh file to use.");
173 args.
AddOption(&ref_levels,
"-r",
"--refine",
174 "Number of times to refine the mesh uniformly.");
176 "Order (degree) of the finite elements.");
177 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
178 "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
179 " 11 - Forward Euler, 12 - RK2,\n\t"
180 " 13 - RK3 SSP, 14 - RK4.");
181 args.
AddOption(&t_final,
"-tf",
"--t-final",
182 "Final time; start time is 0.");
183 args.
AddOption(&dt,
"-dt",
"--time-step",
185 args.
AddOption(&visc,
"-v",
"--viscosity",
186 "Viscosity coefficient.");
187 args.
AddOption(&mu,
"-mu",
"--shear-modulus",
188 "Shear modulus in the Neo-Hookean hyperelastic model.");
189 args.
AddOption(&K,
"-K",
"--bulk-modulus",
190 "Bulk modulus in the Neo-Hookean hyperelastic model.");
191 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
192 "--no-visualization",
193 "Enable or disable GLVis visualization.");
194 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
195 "Visualize every n-th timestep.");
206 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
213 switch (ode_solver_type)
221 case 12: ode_solver =
new RK2Solver(0.5);
break;
223 case 14: ode_solver =
new RK4Solver;
break;
229 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
236 for (
int lev = 0; lev < ref_levels; lev++)
251 cout <<
"Number of velocity/deformation unknowns: " << fe_size << endl;
254 fe_offset[1] = fe_size;
255 fe_offset[2] = 2*fe_size;
282 HyperelasticOperator oper(fespace, ess_bdr, visc, mu, K);
287 char vishost[] =
"localhost";
289 vis_v.
open(vishost, visport);
291 visualize(vis_v, mesh, &x, &v,
"Velocity",
true);
292 vis_w.
open(vishost, visport);
295 oper.GetElasticEnergyDensity(x, w);
297 visualize(vis_w, mesh, &x, &w,
"Elastic energy density",
true);
301 double ee0 = oper.ElasticEnergy(x);
302 double ke0 = oper.KineticEnergy(v);
303 cout <<
"initial elastic energy (EE) = " << ee0 << endl;
304 cout <<
"initial kinetic energy (KE) = " << ke0 << endl;
305 cout <<
"initial total energy (TE) = " << (ee0 + ke0) << endl;
309 ode_solver->
Init(oper);
313 bool last_step =
false;
314 for (
int ti = 1; !last_step; ti++)
316 double dt_real = min(dt, t_final - t);
318 ode_solver->
Step(vx, t, dt_real);
320 last_step = (t >= t_final - 1e-8*dt);
322 if (last_step || (ti % vis_steps) == 0)
324 double ee = oper.ElasticEnergy(x);
325 double ke = oper.KineticEnergy(v);
327 cout <<
"step " << ti <<
", t = " << t <<
", EE = " << ee <<
", KE = "
328 << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
335 oper.GetElasticEnergyDensity(x, w);
347 ofstream mesh_ofs(
"deformed.mesh");
348 mesh_ofs.precision(8);
349 mesh->
Print(mesh_ofs);
351 ofstream velo_ofs(
"velocity.sol");
352 velo_ofs.precision(8);
354 ofstream ee_ofs(
"elastic_energy.sol");
356 oper.GetElasticEnergyDensity(x, w);
369 GridFunction *field,
const char *field_name,
bool init_vis)
381 out <<
"solution\n" << *mesh << *field;
387 out <<
"window_size 800 800\n";
388 out <<
"window_title '" << field_name <<
"'\n";
395 out <<
"autoscale value\n";
402 ReducedSystemOperator::ReducedSystemOperator(
404 :
Operator(M_->Height()), M(M_), S(S_), H(H_), Jacobian(NULL),
405 dt(0.0), v(NULL), x(NULL), w(height), z(height)
408 void ReducedSystemOperator::SetParameters(
double dt_,
const Vector *v_,
411 dt = dt_; v = v_; x = x_;
424 Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const
427 Jacobian =
Add(1.0, M->SpMat(), dt, S->SpMat());
431 Jacobian->
Add(dt*dt, *grad_H);
435 ReducedSystemOperator::~ReducedSystemOperator()
445 M(&fespace), S(&fespace), H(&fespace),
446 viscosity(visc), z(height/2)
448 const double rel_tol = 1e-8;
449 const int skip_zero_entries = 0;
451 const double ref_density = 1.0;
454 M.Assemble(skip_zero_entries);
455 M.EliminateEssentialBC(ess_bdr);
456 M.Finalize(skip_zero_entries);
458 M_solver.iterative_mode =
false;
459 M_solver.SetRelTol(rel_tol);
460 M_solver.SetAbsTol(0.0);
461 M_solver.SetMaxIter(30);
462 M_solver.SetPrintLevel(0);
463 M_solver.SetPreconditioner(M_prec);
464 M_solver.SetOperator(M.SpMat());
468 H.SetEssentialBC(ess_bdr);
472 S.Assemble(skip_zero_entries);
473 S.EliminateEssentialBC(ess_bdr);
474 S.Finalize(skip_zero_entries);
476 reduced_oper =
new ReducedSystemOperator(&M, &S, &H);
478 #ifndef MFEM_USE_SUITESPARSE
493 newton_solver.SetSolver(*J_solver);
494 newton_solver.SetOperator(*reduced_oper);
495 newton_solver.SetPrintLevel(1);
496 newton_solver.SetRelTol(rel_tol);
497 newton_solver.SetAbsTol(0.0);
498 newton_solver.SetMaxIter(10);
511 if (viscosity != 0.0)
516 M_solver.Mult(z, dv_dt);
521 void HyperelasticOperator::ImplicitSolve(
const double dt,
536 reduced_oper->SetParameters(dt, &v, &x);
538 newton_solver.Mult(zero, dv_dt);
539 MFEM_VERIFY(newton_solver.GetConverged(),
"Newton solver did not converge.");
540 add(v, dt, dv_dt, dx_dt);
543 double HyperelasticOperator::ElasticEnergy(
Vector &x)
const
545 return H.GetEnergy(x);
548 double HyperelasticOperator::KineticEnergy(
Vector &v)
const
550 return 0.5*M.InnerProduct(v, v);
553 void HyperelasticOperator::GetElasticEnergyDensity(
556 ElasticEnergyCoefficient w_coeff(*model, x);
560 HyperelasticOperator::~HyperelasticOperator()
572 model.SetTransformation(T);
575 return model.EvalW(J)/J.Det();
589 const double s = 0.1/64.;
592 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)
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.
Data type for scaled Jacobi-type smoother of sparse matrix.
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)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Hyperelastic integrator for any given HyperelasticModel.
int main(int argc, char *argv[])
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.
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.
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.
Newton's method for solving F(x)=b for a given operator F.
void MakeRef(FiniteElementSpace *f, double *v)
Make the GridFunction reference external data on a new FiniteElementSpace.
void PrintUsage(std::ostream &out) const
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.
Implicit midpoint method. A-stable, not L-stable.
Class for integration point with weight.
void PrintOptions(std::ostream &out) const
void ProjectCoefficient(Coefficient &coeff)
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.
The classical forward Euler method.
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad)
Vector & GetBlock(int i)
Get the i-th vector in the block.
Arbitrary order "L2-conforming" discontinuous finite elements.
void Neg()
(*this) = -(*this)
virtual void Print(std::ostream &out=std::cout) const