45 GeneralResidualMonitor(
const std::string& prefix_,
int print_lvl)
48 print_level = print_lvl;
51 virtual void MonitorResidual(
int it,
double norm,
const Vector &r,
bool final);
54 const std::string prefix;
59 void GeneralResidualMonitor::MonitorResidual(
int it,
double norm,
60 const Vector &r,
bool final)
62 if (print_level == 1 || (print_level == 3 && (
final || it == 0)))
64 mfem::out << prefix <<
" iteration " << setw(2) << it
65 <<
" : ||r|| = " << norm;
68 mfem::out <<
", ||r||/||r_0|| = " << norm/norm0;
93 class JacobianPreconditioner :
public Solver 120 virtual void SetOperator(
const Operator &op);
122 virtual ~JacobianPreconditioner();
129 class RubberOperator :
public Operator 143 GeneralResidualMonitor newton_monitor;
147 GeneralResidualMonitor j_monitor;
160 Array<int> &block_trueOffsets,
double rel_tol,
double abs_tol,
168 void Solve(
Vector &xp)
const;
170 virtual ~RubberOperator();
176 bool init_vis =
false);
183 int main(
int argc,
char *argv[])
186 const char *mesh_file =
"../data/beam-tet.mesh";
189 bool visualization =
true;
190 double newton_rel_tol = 1e-4;
191 double newton_abs_tol = 1e-6;
192 int newton_iter = 500;
196 args.
AddOption(&mesh_file,
"-m",
"--mesh",
197 "Mesh file to use.");
198 args.
AddOption(&ref_levels,
"-r",
"--refine",
199 "Number of times to refine the mesh uniformly.");
201 "Order (degree) of the finite elements.");
202 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
203 "--no-visualization",
204 "Enable or disable GLVis visualization.");
205 args.
AddOption(&newton_rel_tol,
"-rel",
"--relative-tolerance",
206 "Relative tolerance for the Newton solve.");
207 args.
AddOption(&newton_abs_tol,
"-abs",
"--absolute-tolerance",
208 "Absolute tolerance for the Newton solve.");
209 args.
AddOption(&newton_iter,
"-it",
"--newton-iterations",
210 "Maximum iterations for the Newton solve.");
212 "Shear modulus for the neo-Hookean material.");
223 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
229 for (
int lev = 0; lev < ref_levels; lev++)
247 spaces[0] = &R_space;
248 spaces[1] = &W_space;
264 ess_bdr[0] = &ess_bdr_u;
265 ess_bdr[1] = &ess_bdr_p;
268 std::cout <<
"***********************************************************\n";
269 std::cout <<
"dim(u) = " << R_size <<
"\n";
270 std::cout <<
"dim(p) = " << W_size <<
"\n";
271 std::cout <<
"dim(u+p) = " << R_size + W_size <<
"\n";
272 std::cout <<
"***********************************************************\n";
276 block_trueOffsets[0] = 0;
304 RubberOperator oper(spaces, ess_bdr, block_trueOffsets,
305 newton_rel_tol, newton_abs_tol, newton_iter, c_mu);
323 visualize(vis_u, mesh, &x_gf, &x_def,
"Deformation",
true);
326 visualize(vis_p, mesh, &x_gf, &p_gf,
"Pressure",
true);
335 ofstream mesh_ofs(
"deformed.mesh");
336 mesh_ofs.precision(8);
337 mesh->
Print(mesh_ofs);
339 ofstream pressure_ofs(
"pressure.sol");
340 pressure_ofs.precision(8);
341 p_gf.
Save(pressure_ofs);
343 ofstream deformation_ofs(
"deformation.sol");
344 deformation_ofs.precision(8);
345 x_def.
Save(deformation_ofs);
358 :
Solver(offsets[2]), block_trueOffsets(offsets), pressure_mass(&mass)
368 mass_prec = mass_prec_gs;
379 mass_pcg = mass_pcg_iter;
391 block_trueOffsets[1]-block_trueOffsets[0]);
393 block_trueOffsets[2]-block_trueOffsets[1]);
396 block_trueOffsets[1]-block_trueOffsets[0]);
398 block_trueOffsets[2]-block_trueOffsets[1]);
400 Vector temp(block_trueOffsets[1]-block_trueOffsets[0]);
401 Vector temp2(block_trueOffsets[1]-block_trueOffsets[0]);
404 mass_pcg->Mult(pres_in, pres_out);
407 jacobian->GetBlock(0,1).Mult(pres_out, temp);
410 stiff_pcg->Mult(temp2, disp_out);
413 void JacobianPreconditioner::SetOperator(
const Operator &op)
418 if (stiff_prec == NULL)
422 stiff_prec = stiff_prec_gs;
432 stiff_pcg = stiff_pcg_iter;
440 JacobianPreconditioner::~JacobianPreconditioner()
456 :
Operator(fes[0]->GetTrueVSize() + fes[1]->GetTrueVSize()),
457 newton_solver(), newton_monitor(
"Newton", 1),
458 j_monitor(
" GMRES", 3),
mu(c_mu), block_trueOffsets(offsets)
472 Hform->SetEssentialBC(ess_bdr, rhs);
483 a->FormSystemMatrix(p_ess_tdofs, op);
484 pressure_mass =
a->LoseMat();
488 JacobianPreconditioner *jac_prec =
489 new JacobianPreconditioner(fes, *pressure_mass, block_trueOffsets);
505 newton_solver.SetSolver(*j_solver);
506 newton_solver.SetOperator(*
this);
507 newton_solver.SetPrintLevel(-1);
508 newton_solver.SetMonitor(newton_monitor);
509 newton_solver.SetRelTol(rel_tol);
510 newton_solver.SetAbsTol(abs_tol);
511 newton_solver.SetMaxIter(iter);
515 void RubberOperator::Solve(
Vector &xp)
const 518 newton_solver.Mult(zero, xp);
519 MFEM_VERIFY(newton_solver.GetConverged(),
520 "Newton Solver did not converge.");
535 RubberOperator::~RubberOperator()
538 delete pressure_mass;
546 GridFunction *field,
const char *field_name,
bool init_vis)
558 os <<
"solution\n" << *mesh << *field;
564 os <<
"window_size 800 800\n";
565 os <<
"window_title '" << field_name <<
"'\n";
574 os <<
"autoscale value\n";
590 y[1] = x[1] + 0.25*x[0];
Conjugate gradient method.
void ReferenceConfiguration(const Vector &x, Vector &y)
Class for grid function - Vector with associated FE space.
virtual Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
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.
void PrintOptions(std::ostream &out) const
Print the options.
void PrintUsage(std::ostream &out) const
Print the usage message.
void SwapNodes(GridFunction *&nodes, int &own_nodes_)
Swap the internal node GridFunction pointer and ownership flag members with the given ones...
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Pointer to an Operator of a specified type.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
bool Good() const
Return true if the command line options were parsed successfully.
void visualize(ostream &os, Mesh *mesh, GridFunction *deformed_nodes, GridFunction *field, const char *field_name=NULL, bool init_vis=false)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void InitialDeformation(const Vector &x, Vector &y)
Data type for Gauss-Seidel smoother of sparse matrix.
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
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)
Newton's method for solving F(x)=b for a given operator F.
double * GetData() const
Return a pointer to the beginning of the Vector data.
A general vector function coefficient.
Abstract base class for an iterative solver monitor.
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Mesh * GetMesh() const
Returns the mesh.
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...
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void subtract(const Vector &x, const Vector &y, Vector &z)
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...
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
int SpaceDimension() const
void SetMonitor(IterativeSolverMonitor &m)
Set the iterative solver monitor.
void Copy(Array ©) const
Create a copy of the internal array to the provided copy.
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void Print(std::ostream &os=mfem::out) const
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Arbitrary order H1-conforming (continuous) finite elements.
int main(int argc, char *argv[])
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
A class to handle Block systems in a matrix-free implementation.
Vector & GetBlock(int i)
Get the i-th vector in the block.