45 class BackwardEulerOperator;
73 BackwardEulerOperator *backward_euler_oper;
85 double visc,
double mu,
double K);
90 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
96 virtual ~HyperelasticOperator();
103 class BackwardEulerOperator :
public Operator
116 void SetParameters(
double dt_,
const Vector *v_,
const Vector *x_);
119 virtual ~BackwardEulerOperator();
124 class ElasticEnergyCoefficient :
public Coefficient
133 : model(m), x(x_) { }
135 virtual ~ElasticEnergyCoefficient() { }
144 bool init_vis =
false);
147 int main(
int argc,
char *argv[])
151 MPI_Init(&argc, &argv);
152 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
153 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
156 const char *mesh_file =
"../data/beam-quad.mesh";
157 int ser_ref_levels = 2;
158 int par_ref_levels = 0;
160 int ode_solver_type = 3;
161 double t_final = 300.0;
166 bool visualization =
true;
170 args.
AddOption(&mesh_file,
"-m",
"--mesh",
171 "Mesh file to use.");
172 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
173 "Number of times to refine the mesh uniformly in serial.");
174 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
175 "Number of times to refine the mesh uniformly in parallel.");
177 "Order (degree) of the finite elements.");
178 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
179 "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
180 "\t 11 - Forward Euler, 12 - RK2, 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.");
214 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
221 switch (ode_solver_type)
229 case 12: ode_solver =
new RK2Solver(0.5);
break;
231 case 14: ode_solver =
new RK4Solver;
break;
239 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
248 for (
int lev = 0; lev < ser_ref_levels; lev++)
258 for (
int lev = 0; lev < par_ref_levels; lev++)
275 cout <<
"Number of velocity/deformation unknowns: " << glob_size << endl;
280 true_offset[1] = true_size;
281 true_offset[2] = 2*true_size;
296 v_gf.ProjectCoefficient(velo);
309 HyperelasticOperator oper(fespace, ess_bdr, visc, mu, K);
314 char vishost[] =
"localhost";
316 vis_v.
open(vishost, visport);
318 visualize(vis_v, pmesh, &x_gf, &v_gf,
"Velocity",
true);
322 vis_w.
open(vishost, visport);
325 oper.GetElasticEnergyDensity(x_gf, w_gf);
327 visualize(vis_w, pmesh, &x_gf, &w_gf,
"Elastic energy density",
true);
331 double ee0 = oper.ElasticEnergy(x_gf);
332 double ke0 = oper.KineticEnergy(v_gf);
335 cout <<
"initial elastic energy (EE) = " << ee0 << endl;
336 cout <<
"initial kinetic energy (KE) = " << ke0 << endl;
337 cout <<
"initial total energy (TE) = " << (ee0 + ke0) << endl;
342 ode_solver->
Init(oper);
345 bool last_step =
false;
346 for (
int ti = 1; !last_step; ti++)
348 if (t + dt >= t_final - dt/2)
353 ode_solver->
Step(vx, t, dt);
355 if (last_step || (ti % vis_steps) == 0)
360 double ee = oper.ElasticEnergy(x_gf);
361 double ke = oper.KineticEnergy(v_gf);
364 cout <<
"step " << ti <<
", t = " << t <<
", EE = " << ee
365 <<
", KE = " << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
372 oper.GetElasticEnergyDensity(x_gf, w_gf);
385 ostringstream mesh_name, velo_name, ee_name;
386 mesh_name <<
"deformed." << setfill(
'0') << setw(6) << myid;
387 velo_name <<
"velocity." << setfill(
'0') << setw(6) << myid;
388 ee_name <<
"elastic_energy." << setfill(
'0') << setw(6) << myid;
390 ofstream mesh_ofs(mesh_name.str().c_str());
391 mesh_ofs.precision(8);
392 pmesh->
Print(mesh_ofs);
394 ofstream velo_ofs(velo_name.str().c_str());
395 velo_ofs.precision(8);
397 ofstream ee_ofs(ee_name.str().c_str());
399 oper.GetElasticEnergyDensity(x_gf, w_gf);
426 out <<
"solution\n" << *mesh << *field;
432 out <<
"window_size 800 800\n";
433 out <<
"window_title '" << field_name <<
"'\n";
440 out <<
"autoscale value\n";
446 BackwardEulerOperator::BackwardEulerOperator(
448 :
Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
449 Jacobian(NULL), v(NULL), x(NULL), dt(0.0), w(height), z(height)
452 void BackwardEulerOperator::SetParameters(
double dt_,
const Vector *v_,
455 dt = dt_; v = v_; x = x_;
464 M->TrueAddMult(k, y);
465 S->TrueAddMult(w, y);
468 Operator &BackwardEulerOperator::GetGradient(
const Vector &k)
const
474 localJ->
Add(dt*dt, H->GetLocalGradient(z));
475 Jacobian = M->ParallelAssemble(localJ);
480 BackwardEulerOperator::~BackwardEulerOperator()
490 M(&fespace), S(&fespace), H(&fespace), M_solver(f.GetComm()),
491 newton_solver(f.GetComm()), z(height/2)
493 const double rel_tol = 1e-8;
494 const int skip_zero_entries = 0;
496 const double ref_density = 1.0;
499 M.Assemble(skip_zero_entries);
500 M.EliminateEssentialBC(ess_bdr);
501 M.Finalize(skip_zero_entries);
502 Mmat = M.ParallelAssemble();
504 M_solver.iterative_mode =
false;
505 M_solver.SetRelTol(rel_tol);
506 M_solver.SetAbsTol(0.0);
507 M_solver.SetMaxIter(30);
508 M_solver.SetPrintLevel(0);
509 M_prec.SetType(HypreSmoother::Jacobi);
510 M_solver.SetPreconditioner(M_prec);
511 M_solver.SetOperator(*Mmat);
515 H.SetEssentialBC(ess_bdr);
520 S.Assemble(skip_zero_entries);
521 S.EliminateEssentialBC(ess_bdr);
522 S.Finalize(skip_zero_entries);
524 backward_euler_oper =
new BackwardEulerOperator(&M, &S, &H);
527 J_hypreSmoother->
SetType(HypreSmoother::l1Jacobi);
528 J_prec = J_hypreSmoother;
539 newton_solver.SetSolver(*J_solver);
540 newton_solver.SetOperator(*backward_euler_oper);
541 newton_solver.SetPrintLevel(1);
542 newton_solver.SetRelTol(rel_tol);
543 newton_solver.SetAbsTol(0.0);
544 newton_solver.SetMaxIter(10);
557 if (viscosity != 0.0)
562 M_solver.Mult(z, dv_dt);
567 void HyperelasticOperator::ImplicitSolve(
const double dt,
582 backward_euler_oper->SetParameters(dt, &v, &x);
584 newton_solver.Mult(zero, dv_dt);
585 add(v, dt, dv_dt, dx_dt);
587 MFEM_VERIFY(newton_solver.GetConverged(),
"Newton Solver did not converge.");
592 return H.GetEnergy(x);
597 double loc_energy = 0.5*M.InnerProduct(v, v);
599 MPI_Allreduce(&loc_energy, &energy, 1, MPI_DOUBLE, MPI_SUM,
604 void HyperelasticOperator::GetElasticEnergyDensity(
607 ElasticEnergyCoefficient w_coeff(*model, x);
611 HyperelasticOperator::~HyperelasticOperator()
614 delete backward_euler_oper;
624 model.SetTransformation(T);
627 return model.EvalW(J)/J.Det();
641 const double s = 0.1/64.;
644 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.
Subclass constant coefficient.
Base abstract class for time dependent operators: (x,t) -> f(x,t)
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
HYPRE_Int GlobalTrueVSize()
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)
virtual void Save(std::ostream &out) const
Abstract parallel finite element space.
void ProjectCoefficient(Coefficient &coeff)
bool iterative_mode
If true, use the second argument of Mult as an initial guess.
virtual void Init(TimeDependentOperator &_f)
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)
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)
Parallel smoothers in hypre.
void PrintUsage(std::ostream &out) const
void GetTrueDofs(Vector &tv) const
Returns the true dofs in a Vector.
int SpaceDimension() const
The classical explicit forth-order Runge-Kutta method, RK4.
void SetAbsTol(double atol)
Array< int > bdr_attributes
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.
void Distribute(const Vector *tv)
Implicit midpoint method. A-stable, not L-stable.
Class for integration point with weight.
void PrintOptions(std::ostream &out) const
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.
Class for parallel grid function.
The classical forward Euler method.
Wrapper for hypre's ParCSR matrix class.
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad)
Vector & GetBlock(int i)
Get the i-th vector in the block.
Class for parallel meshes.
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Arbitrary order "L2-conforming" discontinuous finite elements.
void Neg()
(*this) = -(*this)
virtual void Print(std::ostream &out=std::cout) const