50#ifndef MFEM_USE_SUNDIALS 
   51#error This example requires that MFEM is built with MFEM_USE_SUNDIALS=YES 
   57class ReducedSystemOperator;
 
   86   ReducedSystemOperator *reduced_oper;
 
  106                        double visc, 
double mu, 
double K,
 
  107                        int kinsol_nls_type = -1, 
double kinsol_damping = 0.0,
 
  108                        int kinsol_aa_n = 0);
 
  111   virtual void Mult(
const Vector &vx, 
Vector &dvx_dt) 
const;
 
  115   virtual void ImplicitSolve(
const double dt, 
const Vector &x, 
Vector &k);
 
  141   virtual int SUNImplicitSetup(
const Vector &y, 
const Vector &fy,
 
  142                                int jok, 
int *jcur, 
double gamma);
 
  146   virtual int SUNImplicitSolve(
const Vector &
b, 
Vector &x, 
double tol);
 
  153   virtual ~HyperelasticOperator();
 
  160class ReducedSystemOperator : 
public Operator 
  176   void SetParameters(
double dt_, 
const Vector *v_, 
const Vector *x_);
 
  184   virtual ~ReducedSystemOperator();
 
  199      : model(m), x(x_) { }
 
  201   virtual ~ElasticEnergyCoefficient() { }
 
  210               bool init_vis = 
false);
 
  213int main(
int argc, 
char *argv[])
 
  222   const char *mesh_file = 
"../../data/beam-quad.mesh";
 
  223   int ser_ref_levels = 2;
 
  224   int par_ref_levels = 0;
 
  226   int ode_solver_type = 3;
 
  227   double t_final = 300.0;
 
  232   bool visualization = 
true;
 
  233   int nonlinear_solver_type = 0;
 
  235   double kinsol_damping = 0.0;
 
  236   int kinsol_aa_n = -1;
 
  239   const double reltol = 1e-1, abstol = 1e-1;
 
  243   const double cvode_eps_lin = 1e-4;
 
  245   const double arkode_eps_nonlin = 1e-6;
 
  248   args.
AddOption(&mesh_file, 
"-m", 
"--mesh",
 
  249                  "Mesh file to use.");
 
  250   args.
AddOption(&ser_ref_levels, 
"-rs", 
"--refine-serial",
 
  251                  "Number of times to refine the mesh uniformly in serial.");
 
  252   args.
AddOption(&par_ref_levels, 
"-rp", 
"--refine-parallel",
 
  253                  "Number of times to refine the mesh uniformly in parallel.");
 
  255                  "Order (degree) of the finite elements.");
 
  256   args.
AddOption(&ode_solver_type, 
"-s", 
"--ode-solver",
 
  258                  "1  - Backward Euler,\n\t" 
  259                  "2  - SDIRK2, L-stable\n\t" 
  260                  "3  - SDIRK3, L-stable\n\t" 
  261                  "4  - Implicit Midpoint,\n\t" 
  262                  "5  - SDIRK2, A-stable,\n\t" 
  263                  "6  - SDIRK3, A-stable,\n\t" 
  264                  "7  - Forward Euler,\n\t" 
  268                  "11 - CVODE implicit BDF, approximate Jacobian,\n\t" 
  269                  "12 - CVODE implicit BDF, specified Jacobian,\n\t" 
  270                  "13 - CVODE implicit ADAMS, approximate Jacobian,\n\t" 
  271                  "14 - CVODE implicit ADAMS, specified Jacobian,\n\t" 
  272                  "15 - ARKODE implicit, approximate Jacobian,\n\t" 
  273                  "16 - ARKODE implicit, specified Jacobian,\n\t" 
  274                  "17 - ARKODE explicit, 4th order.");
 
  275   args.
AddOption(&nonlinear_solver_type, 
"-nls", 
"--nonlinear-solver",
 
  276                  "Nonlinear system solver:\n\t" 
  277                  "0  - MFEM Newton method,\n\t" 
  278                  "1  - KINSOL Newton method,\n\t" 
  279                  "2  - KINSOL Newton method with globalization,\n\t" 
  280                  "3  - KINSOL fixed-point method (with or without AA),\n\t" 
  281                  "4  - KINSOL Picard method (with or without AA).");
 
  282   args.
AddOption(&kinsol_damping, 
"-damp", 
"--kinsol-damping",
 
  283                  "Picard or Fixed-Point damping parameter (only valid with KINSOL): " 
  285   args.
AddOption(&kinsol_aa_n, 
"-aan", 
"--anderson-subspace",
 
  286                  "Anderson Acceleration subspace size (only valid with KINSOL)");
 
  287   args.
AddOption(&t_final, 
"-tf", 
"--t-final",
 
  288                  "Final time; start time is 0.");
 
  289   args.
AddOption(&dt, 
"-dt", 
"--time-step",
 
  291   args.
AddOption(&visc, 
"-v", 
"--viscosity",
 
  292                  "Viscosity coefficient.");
 
  294                  "Shear modulus in the Neo-Hookean hyperelastic model.");
 
  295   args.
AddOption(&K, 
"-K", 
"--bulk-modulus",
 
  296                  "Bulk modulus in the Neo-Hookean hyperelastic model.");
 
  297   args.
AddOption(&visualization, 
"-vis", 
"--visualization", 
"-no-vis",
 
  298                  "--no-visualization",
 
  299                  "Enable or disable GLVis visualization.");
 
  300   args.
AddOption(&vis_steps, 
"-vs", 
"--visualization-steps",
 
  301                  "Visualize every n-th timestep.");
 
  317   if (ode_solver_type < 1 || ode_solver_type > 17)
 
  321         cout << 
"Unknown ODE solver type: " << ode_solver_type << 
'\n';
 
  327   if (nonlinear_solver_type < 0 || nonlinear_solver_type > 4)
 
  331         cout << 
"Unknown nonlinear solver type: " << nonlinear_solver_type
 
  336   if (kinsol_damping > 0.0 &&
 
  337       !(nonlinear_solver_type == 3 || nonlinear_solver_type == 4))
 
  341         cout << 
"Only KINSOL fixed-point and Picard methods can use damping\n";
 
  345   if (kinsol_aa_n > 0 &&
 
  346       !(nonlinear_solver_type == 3 || nonlinear_solver_type == 4))
 
  350         cout << 
"Only KINSOL fixed-point and Picard methods can use AA\n";
 
  358   Mesh *mesh = 
new Mesh(mesh_file, 1, 1);
 
  364   for (
int lev = 0; lev < ser_ref_levels; lev++)
 
  374   for (
int lev = 0; lev < par_ref_levels; lev++)
 
  391      cout << 
"Number of velocity/deformation unknowns: " << glob_size << endl;
 
  396   true_offset[1] = true_size;
 
  397   true_offset[2] = 2*true_size;
 
  401   v_gf.
MakeTRef(&fespace, vx, true_offset[0]);
 
  402   x_gf.
MakeTRef(&fespace, vx, true_offset[1]);
 
  428   std::unique_ptr<HyperelasticOperator> oper;
 
  429   if (nonlinear_solver_type == 0)
 
  430      oper = std::make_unique<HyperelasticOperator>(fespace, ess_bdr, visc, 
mu,
 
  434      switch (nonlinear_solver_type)
 
  437            oper = std::make_unique<HyperelasticOperator>(fespace, ess_bdr,
 
  438                                                          visc, 
mu, K, KIN_NONE);
 
  441            oper = std::make_unique<HyperelasticOperator>(fespace, ess_bdr,
 
  442                                                          visc, 
mu, K, KIN_LINESEARCH);
 
  445            oper = std::make_unique<HyperelasticOperator>(fespace, ess_bdr,
 
  446                                                          visc, 
mu, K, KIN_FP, kinsol_damping, kinsol_aa_n);
 
  449            oper = std::make_unique<HyperelasticOperator>(fespace, ess_bdr,
 
  450                                                          visc, 
mu, K, KIN_PICARD, kinsol_damping, kinsol_aa_n);
 
  453            cout << 
"Unknown type of nonlinear solver: " 
  454                 << nonlinear_solver_type << endl;
 
  466      visualize(vis_v, pmesh, &x_gf, &v_gf, 
"Velocity", 
true);
 
  473         oper->GetElasticEnergyDensity(x_gf, w_gf);
 
  475         visualize(vis_w, pmesh, &x_gf, &w_gf, 
"Elastic energy density", 
true);
 
  479   double ee0 = oper->ElasticEnergy(x_gf);
 
  480   double ke0 = oper->KineticEnergy(v_gf);
 
  483      cout << 
"initial elastic energy (EE) = " << ee0 << endl;
 
  484      cout << 
"initial kinetic energy (KE) = " << ke0 << endl;
 
  485      cout << 
"initial   total energy (TE) = " << (ee0 + ke0) << endl;
 
  497   switch (ode_solver_type)
 
  509      case 8:  ode_solver = 
new RK2Solver(0.5); 
break; 
 
  511      case 10: ode_solver = 
new RK4Solver; 
break;
 
  518         CVodeSetEpsLin(cvode->
GetMem(), cvode_eps_lin);
 
  520         if (ode_solver_type == 11)
 
  524         ode_solver = cvode; 
break;
 
  531         CVodeSetEpsLin(cvode->
GetMem(), cvode_eps_lin);
 
  533         if (ode_solver_type == 13)
 
  537         ode_solver = cvode; 
break;
 
  544#if MFEM_SUNDIALS_VERSION < 70100 
  545         ARKStepSetNonlinConvCoef(arkode->
GetMem(), arkode_eps_nonlin);
 
  547         ARKodeSetNonlinConvCoef(arkode->
GetMem(), arkode_eps_nonlin);
 
  550         if (ode_solver_type == 15)
 
  554         ode_solver = arkode; 
break;
 
  561         ode_solver = arkode; 
break;
 
  565   if (ode_solver_type < 11) { ode_solver->
Init(*oper); }
 
  569   bool last_step = 
false;
 
  570   for (
int ti = 1; !last_step; ti++)
 
  572      double dt_real = min(dt, t_final - t);
 
  574      ode_solver->
Step(vx, t, dt_real);
 
  576      last_step = (t >= t_final - 1e-8*dt);
 
  578      if (last_step || (ti % vis_steps) == 0)
 
  582         double ee = oper->ElasticEnergy(x_gf);
 
  583         double ke = oper->KineticEnergy(v_gf);
 
  587            cout << 
"step " << ti << 
", t = " << t << 
", EE = " << ee
 
  588                 << 
", KE = " << ke << 
", ΔTE = " << (ee+ke)-(ee0+ke0) << endl;
 
  591            else if (arkode) { arkode->
PrintInfo(); }
 
  599               oper->GetElasticEnergyDensity(x_gf, w_gf);
 
  613      ostringstream mesh_name, velo_name, ee_name;
 
  614      mesh_name << 
"deformed." << setfill(
'0') << setw(6) << myid;
 
  615      velo_name << 
"velocity." << setfill(
'0') << setw(6) << myid;
 
  616      ee_name << 
"elastic_energy." << setfill(
'0') << setw(6) << myid;
 
  618      ofstream mesh_ofs(mesh_name.str().c_str());
 
  619      mesh_ofs.precision(8);
 
  620      pmesh->
Print(mesh_ofs);
 
  622      ofstream velo_ofs(velo_name.str().c_str());
 
  623      velo_ofs.precision(8);
 
  625      ofstream ee_ofs(ee_name.str().c_str());
 
  627      oper->GetElasticEnergyDensity(x_gf, w_gf);
 
 
  652   os << 
"solution\n" << *mesh << *field;
 
  658      os << 
"window_size 800 800\n";
 
  659      os << 
"window_title '" << field_name << 
"'\n";
 
  666      os << 
"autoscale value\n"; 
 
 
  673ReducedSystemOperator::ReducedSystemOperator(
 
  676   : 
Operator(M_->ParFESpace()->TrueVSize()), M(M_), S(S_), H(H_),
 
  677     Jacobian(NULL), dt(0.0), v(NULL), x(NULL), w(height), z(height),
 
  678     ess_tdof_list(ess_tdof_list_)
 
  681void ReducedSystemOperator::SetParameters(
double dt_, 
const Vector *v_,
 
  684   dt = dt_;  v = v_;  x = x_;
 
  687void ReducedSystemOperator::Mult(
const Vector &k, 
Vector &y)
 const 
  693   M->TrueAddMult(k, y);
 
  694   S->TrueAddMult(w, y);
 
  698Operator &ReducedSystemOperator::GetGradient(
const Vector &k)
 const 
  704   localJ->
Add(dt*dt, H->GetLocalGradient(z));
 
  705   Jacobian = M->ParallelAssemble(localJ);
 
  712ReducedSystemOperator::~ReducedSystemOperator()
 
  722                                           double kinsol_damping,
 
  726     M(&fespace), S(&fespace), H(&fespace),
 
  727     viscosity(visc), M_solver(
f.GetComm()), z(height/2),
 
  728     local_grad_H(NULL), Jacobian(NULL)
 
  730   const double rel_tol = 1e-8;
 
  731   const int skip_zero_entries = 0;
 
  733   const double ref_density = 1.0; 
 
  736   M.Assemble(skip_zero_entries);
 
  737   M.Finalize(skip_zero_entries);
 
  738   Mmat = M.ParallelAssemble();
 
  739   fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
 
  743   M_solver.iterative_mode = 
false;
 
  744   M_solver.SetRelTol(rel_tol);
 
  745   M_solver.SetAbsTol(0.0);
 
  746   M_solver.SetMaxIter(30);
 
  747   M_solver.SetPrintLevel(0);
 
  749   M_solver.SetPreconditioner(M_prec);
 
  750   M_solver.SetOperator(*Mmat);
 
  754   H.SetEssentialTrueDofs(ess_tdof_list);
 
  758   S.Assemble(skip_zero_entries);
 
  759   S.Finalize(skip_zero_entries);
 
  761   reduced_oper = 
new ReducedSystemOperator(&M, &S, &H, ess_tdof_list);
 
  766   J_prec = J_hypreSmoother;
 
  776   if (kinsol_nls_type > 0)
 
  779      if (kinsol_nls_type != KIN_PICARD)
 
  788      newton_solver = kinsolver;
 
  790      newton_solver->SetMaxIter(200);
 
  791      newton_solver->SetRelTol(rel_tol);
 
  792      newton_solver->SetPrintLevel(0);
 
  794      if (kinsol_damping > 0.0)
 
  802      newton_solver->SetOperator(*reduced_oper);
 
  803      newton_solver->SetMaxIter(10);
 
  804      newton_solver->SetRelTol(rel_tol);
 
  805      newton_solver->SetPrintLevel(-1);
 
  807   newton_solver->SetSolver(*J_solver);
 
  808   newton_solver->iterative_mode = 
false;
 
  811void HyperelasticOperator::Mult(
const Vector &vx, 
Vector &dvx_dt)
 const 
  821   if (viscosity != 0.0)
 
  827   M_solver.
Mult(z, dv_dt);
 
  832void HyperelasticOperator::ImplicitSolve(
const double dt,
 
  847   reduced_oper->SetParameters(dt, &v, &x);
 
  849   newton_solver->
Mult(zero, dv_dt);
 
  851               "Nonlinear solver did not converge.");
 
  853   if (fespace.GetMyRank() == 0)
 
  856           << 
", final norm = " << newton_solver->
GetFinalNorm() << 
'\n';
 
  859   add(v, dt, dv_dt, dx_dt);
 
  862int HyperelasticOperator::SUNImplicitSetup(
const Vector &y,
 
  863                                           const Vector &fy, 
int jok, 
int *jcur,
 
  866   int sc = y.
Size() / 2;
 
  870   if (Jacobian) { 
delete Jacobian; }
 
  872   local_grad_H = &H.GetLocalGradient(x);
 
  873   localJ->
Add(gamma*gamma, *local_grad_H);
 
  874   Jacobian = M.ParallelAssemble(localJ);
 
  892int HyperelasticOperator::SUNImplicitSolve(
const Vector &
b, 
Vector &x,
 
  895   int sc = 
b.Size() / 2;
 
  897   Vector b_v(
b.GetData() +  0, sc);
 
  898   Vector b_x(
b.GetData() + sc, sc);
 
  907   lb_x.Distribute(b_x);
 
  908   local_grad_H->
Mult(lb_x, lrhs);
 
  909   lrhs.ParallelAssemble(rhs);
 
  911   M.TrueAddMult(b_v, rhs);
 
  912   rhs.SetSubVector(ess_tdof_list, 0.0);
 
  915   J_solver->
Mult(rhs, x_v);
 
  917   add(b_x, saved_gamma, x_v, x_x);
 
  922double HyperelasticOperator::ElasticEnergy(
const ParGridFunction &x)
 const 
  927double HyperelasticOperator::KineticEnergy(
const ParGridFunction &v)
 const 
  929   double energy = 0.5*M.ParInnerProduct(v, v);
 
  933void HyperelasticOperator::GetElasticEnergyDensity(
 
  936   ElasticEnergyCoefficient w_coeff(*model, x);
 
  940HyperelasticOperator::~HyperelasticOperator()
 
  943   delete newton_solver;
 
  972   const double s = 0.1/64.;
 
  975   v(
dim-1) = s*x(0)*x(0)*(8.0-x(0));
 
 
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
 
void SetMaxStep(double dt_max)
Set the maximum time step.
 
void PrintInfo() const
Print various ARKStep statistics.
 
@ IMPLICIT
Implicit RK method.
 
@ EXPLICIT
Explicit RK method.
 
void UseSundialsLinearSolver()
Attach a SUNDIALS GMRES linear solver to ARKode.
 
void Init(TimeDependentOperator &f_) override
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults.
 
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
 
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
 
Backward Euler ODE solver. L-stable.
 
A class to handle Vectors in a block fashion.
 
Conjugate gradient method.
 
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
 
Interface to the CVODE library – linear multi-step methods.
 
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
 
void Init(TimeDependentOperator &f_) override
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
 
void SetMaxStep(double dt_max)
Set the maximum time step.
 
void PrintInfo() const
Print various CVODE statistics.
 
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
 
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
 
A coefficient that is constant across space and time.
 
Data type dense matrix using column-major storage.
 
Mesh * GetMesh() const
Returns the mesh.
 
The classical forward Euler method.
 
Class for grid function - Vector with associated FE space.
 
void GetVectorGradient(ElementTransformation &tr, DenseMatrix &grad) const
Compute the vector gradient with respect to the physical element variable.
 
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
 
void MakeTRef(FiniteElementSpace *f, real_t *tv)
Associate a new FiniteElementSpace and new true-dof data with the GridFunction.
 
void SetFromTrueVector()
Shortcut for calling SetFromTrueDofs() with GetTrueVector() as argument.
 
Arbitrary order H1-conforming (continuous) finite elements.
 
Abstract class for hyperelastic models.
 
virtual real_t EvalW(const DenseMatrix &Jpt) const =0
Evaluate the strain energy density function, W = W(Jpt).
 
void SetTransformation(ElementTransformation &Ttr_)
 
Wrapper for hypre's ParCSR matrix class.
 
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
 
Parallel smoothers in hypre.
 
void SetPositiveDiagonal(bool pos=true)
After computing l1-norms, replace them with their absolute values.
 
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
 
@ l1Jacobi
l1-scaled Jacobi
 
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
 
Implicit midpoint method. A-stable, not L-stable.
 
Class for integration point with weight.
 
real_t GetFinalNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult().
 
void SetRelTol(real_t rtol)
 
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
 
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
 
void SetMaxIter(int max_it)
 
bool GetConverged() const
Returns true if the last call to Mult() converged successfully.
 
void SetAbsTol(real_t atol)
 
Interface to the KINSOL library – nonlinear solver methods.
 
void SetJFNK(bool use_jfnk)
Set the Jacobian Free Newton Krylov flag. The default is false.
 
void SetOperator(const Operator &op) override
Set the nonlinear Operator of the system and initialize KINSOL.
 
void SetMaxSetupCalls(int max_calls)
Set maximum number of nonlinear iterations without a Jacobian update.
 
void SetDamping(double damping)
 
void EnableAndersonAcc(int n, int orth=KIN_ORTH_MGS, int delay=0, double damping=1.0)
Enable Anderson Acceleration for KIN_FP or KIN_PICARD.
 
void SetLSMaxIter(int m)
Set the maximum number of linear solver iterations.
 
Arbitrary order "L2-conforming" discontinuous finite elements.
 
void SetPreconditioner(Solver &pr) override
This should be called before SetOperator.
 
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
 
int Dimension() const
Dimension of the reference space used within the elements.
 
int SpaceDimension() const
Dimension of the physical space containing the mesh.
 
void GetNodes(Vector &node_coord) const
 
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
 
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Swap the internal node GridFunction pointer and ownership flag members with the given ones.
 
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
 
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
 
Newton's method for solving F(x)=b for a given operator F.
 
void Mult(const Vector &b, Vector &x) const override
Solve the nonlinear system with right-hand side b.
 
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
 
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
 
virtual void Step(Vector &x, real_t &t, real_t &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
 
int height
Dimension of the output / number of rows in the matrix.
 
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
 
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
 
void PrintUsage(std::ostream &out) const
Print the usage message.
 
void PrintOptions(std::ostream &out) const
Print the options.
 
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...
 
bool Good() const
Return true if the command line options were parsed successfully.
 
Abstract parallel finite element space.
 
HYPRE_BigInt GlobalTrueVSize() const
 
int TrueVSize() const
Obsolete, kept for backward compatibility.
 
Class for parallel grid function.
 
void Save(std::ostream &out) const override
 
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
 
Class for parallel meshes.
 
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
 
Third-order, strong stability preserving (SSP) Runge-Kutta method.
 
The classical explicit forth-order Runge-Kutta method, RK4.
 
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
 
virtual void SetOperator(const Operator &op)=0
Set/update the solver for the given operator.
 
void Add(const int i, const int j, const real_t val)
 
void Mult(const Vector &x, Vector &y) const override
Matrix vector multiplication.
 
void * GetMem() const
Access the SUNDIALS memory structure.
 
Base abstract class for first order time dependent operators.
 
A general vector function coefficient.
 
void Neg()
(*this) = -(*this)
 
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
 
int Size() const
Returns the size of the vector.
 
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
 
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
 
void add(const Vector &v1, const Vector &v2, Vector &v)
 
std::function< real_t(const Vector &)> f(real_t mass_coeff)
 
void Add(const DenseMatrix &A, const DenseMatrix &B, real_t alpha, DenseMatrix &C)
C = A + alpha*B.
 
std::array< int, NCMesh::MaxFaceNodes > nodes
 
void visualize(ostream &os, ParMesh *mesh, ParGridFunction *deformed_nodes, ParGridFunction *field, const char *field_name=NULL, bool init_vis=false)
 
void InitialDeformation(const Vector &x, Vector &y)
 
void InitialVelocity(const Vector &x, Vector &v)