44 class ReducedSystemOperator;
72 ReducedSystemOperator *reduced_oper;
86 double visc,
double mu,
double K);
92 virtual void ImplicitSolve(
const double dt,
const Vector &x,
Vector &k);
98 virtual ~HyperelasticOperator();
105 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[])
161 MPI_Init(&argc, &argv);
162 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
163 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
166 const char *mesh_file =
"../data/beam-quad.mesh";
167 int ser_ref_levels = 2;
168 int par_ref_levels = 0;
170 int ode_solver_type = 3;
171 double t_final = 300.0;
176 bool visualization =
true;
180 args.
AddOption(&mesh_file,
"-m",
"--mesh",
181 "Mesh file to use.");
182 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
183 "Number of times to refine the mesh uniformly in serial.");
184 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
185 "Number of times to refine the mesh uniformly in parallel.");
187 "Order (degree) of the finite elements.");
188 args.
AddOption(&ode_solver_type,
"-s",
"--ode-solver",
189 "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3,\n\t"
190 " 11 - Forward Euler, 12 - RK2,\n\t"
191 " 13 - RK3 SSP, 14 - RK4.");
192 args.
AddOption(&t_final,
"-tf",
"--t-final",
193 "Final time; start time is 0.");
194 args.
AddOption(&dt,
"-dt",
"--time-step",
196 args.
AddOption(&visc,
"-v",
"--viscosity",
197 "Viscosity coefficient.");
198 args.
AddOption(&mu,
"-mu",
"--shear-modulus",
199 "Shear modulus in the Neo-Hookean hyperelastic model.");
200 args.
AddOption(&K,
"-K",
"--bulk-modulus",
201 "Bulk modulus in the Neo-Hookean hyperelastic model.");
202 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
203 "--no-visualization",
204 "Enable or disable GLVis visualization.");
205 args.
AddOption(&vis_steps,
"-vs",
"--visualization-steps",
206 "Visualize every n-th timestep.");
225 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
232 switch (ode_solver_type)
240 case 12: ode_solver =
new RK2Solver(0.5);
break;
242 case 14: ode_solver =
new RK4Solver;
break;
250 cout <<
"Unknown ODE solver type: " << ode_solver_type <<
'\n';
259 for (
int lev = 0; lev < ser_ref_levels; lev++)
269 for (
int lev = 0; lev < par_ref_levels; lev++)
286 cout <<
"Number of velocity/deformation unknowns: " << glob_size << endl;
291 true_offset[1] = true_size;
292 true_offset[2] = 2*true_size;
307 v_gf.ProjectCoefficient(velo);
320 HyperelasticOperator oper(fespace, ess_bdr, visc, mu, K);
325 char vishost[] =
"localhost";
327 vis_v.
open(vishost, visport);
329 visualize(vis_v, pmesh, &x_gf, &v_gf,
"Velocity",
true);
333 vis_w.
open(vishost, visport);
336 oper.GetElasticEnergyDensity(x_gf, w_gf);
338 visualize(vis_w, pmesh, &x_gf, &w_gf,
"Elastic energy density",
true);
342 double ee0 = oper.ElasticEnergy(x_gf);
343 double ke0 = oper.KineticEnergy(v_gf);
346 cout <<
"initial elastic energy (EE) = " << ee0 << endl;
347 cout <<
"initial kinetic energy (KE) = " << ke0 << endl;
348 cout <<
"initial total energy (TE) = " << (ee0 + ke0) << endl;
353 ode_solver->
Init(oper);
357 bool last_step =
false;
358 for (
int ti = 1; !last_step; ti++)
360 double dt_real = min(dt, t_final - t);
362 ode_solver->
Step(vx, t, dt_real);
364 last_step = (t >= t_final - 1e-8*dt);
366 if (last_step || (ti % vis_steps) == 0)
371 double ee = oper.ElasticEnergy(x_gf);
372 double ke = oper.KineticEnergy(v_gf);
376 cout <<
"step " << ti <<
", t = " << t <<
", EE = " << ee
377 <<
", KE = " << ke <<
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
385 oper.GetElasticEnergyDensity(x_gf, w_gf);
398 ostringstream mesh_name, velo_name, ee_name;
399 mesh_name <<
"deformed." << setfill(
'0') << setw(6) << myid;
400 velo_name <<
"velocity." << setfill(
'0') << setw(6) << myid;
401 ee_name <<
"elastic_energy." << setfill(
'0') << setw(6) << myid;
403 ofstream mesh_ofs(mesh_name.str().c_str());
404 mesh_ofs.precision(8);
405 pmesh->
Print(mesh_ofs);
407 ofstream velo_ofs(velo_name.str().c_str());
408 velo_ofs.precision(8);
410 ofstream ee_ofs(ee_name.str().c_str());
412 oper.GetElasticEnergyDensity(x_gf, w_gf);
439 out <<
"solution\n" << *mesh << *field;
445 out <<
"window_size 800 800\n";
446 out <<
"window_title '" << field_name <<
"'\n";
453 out <<
"autoscale value\n";
460 ReducedSystemOperator::ReducedSystemOperator(
462 :
Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
463 Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height)
466 void ReducedSystemOperator::SetParameters(
double dt_,
const Vector *v_,
469 dt = dt_; v = v_; x = x_;
478 M->TrueAddMult(k, y);
479 S->TrueAddMult(w, y);
482 Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
const
488 localJ->
Add(dt*dt, H->GetLocalGradient(z));
489 Jacobian = M->ParallelAssemble(localJ);
494 ReducedSystemOperator::~ReducedSystemOperator()
504 M(&fespace), S(&fespace), H(&fespace),
505 viscosity(visc), M_solver(f.GetComm()), newton_solver(f.GetComm()),
508 const double rel_tol = 1e-8;
509 const int skip_zero_entries = 0;
511 const double ref_density = 1.0;
514 M.Assemble(skip_zero_entries);
515 M.EliminateEssentialBC(ess_bdr);
516 M.Finalize(skip_zero_entries);
517 Mmat = M.ParallelAssemble();
519 M_solver.iterative_mode =
false;
520 M_solver.SetRelTol(rel_tol);
521 M_solver.SetAbsTol(0.0);
522 M_solver.SetMaxIter(30);
523 M_solver.SetPrintLevel(0);
524 M_prec.SetType(HypreSmoother::Jacobi);
525 M_solver.SetPreconditioner(M_prec);
526 M_solver.SetOperator(*Mmat);
530 H.SetEssentialBC(ess_bdr);
534 S.Assemble(skip_zero_entries);
535 S.EliminateEssentialBC(ess_bdr);
536 S.Finalize(skip_zero_entries);
538 reduced_oper =
new ReducedSystemOperator(&M, &S, &H);
541 J_hypreSmoother->
SetType(HypreSmoother::l1Jacobi);
543 J_prec = J_hypreSmoother;
554 newton_solver.SetSolver(*J_solver);
555 newton_solver.SetOperator(*reduced_oper);
556 newton_solver.SetPrintLevel(1);
557 newton_solver.SetRelTol(rel_tol);
558 newton_solver.SetAbsTol(0.0);
559 newton_solver.SetMaxIter(10);
572 if (viscosity != 0.0)
577 M_solver.Mult(z, dv_dt);
582 void HyperelasticOperator::ImplicitSolve(
const double dt,
597 reduced_oper->SetParameters(dt, &v, &x);
599 newton_solver.Mult(zero, dv_dt);
600 MFEM_VERIFY(newton_solver.GetConverged(),
"Newton solver did not converge.");
601 add(v, dt, dv_dt, dx_dt);
606 return H.GetEnergy(x);
611 double loc_energy = 0.5*M.InnerProduct(v, v);
613 MPI_Allreduce(&loc_energy, &energy, 1, MPI_DOUBLE, MPI_SUM,
618 void HyperelasticOperator::GetElasticEnergyDensity(
621 ElasticEnergyCoefficient w_coeff(*model, x);
625 HyperelasticOperator::~HyperelasticOperator()
638 model.SetTransformation(T);
641 return model.EvalW(J)/J.Det();
655 const double s = 0.1/64.;
658 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.
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]...
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.
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 SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
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.
Parallel smoothers in hypre.
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.
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
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
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