50 ___ ___ ________ ________ ____ __.___________ 51 / | \\_____ \ \_____ \ | |/ _|\_ _____/ 52 / ~ \/ | \ / | \| < | __)_ 53 \ Y / | \/ | \ | \ | \ 54 \___|_ /\_______ /\_______ /____|__ \/_______ / 60 int main(
int argc,
char *argv[])
62 Mpi::Init(argc, argv);
63 int num_procs = Mpi::WorldSize();
64 int myid = Mpi::WorldRank();
67 const char *device_config =
"cpu";
68 int diagpc_type = ElasticityDiagonalPreconditioner::Type::Diagonal;
69 int serial_refinement_levels = 0;
70 bool visualization =
true;
71 bool paraview =
false;
80 "Finite element order (polynomial degree).");
81 args.
AddOption(&device_config,
"-d",
"--device",
82 "Device configuration string, see Device::Configure().");
83 args.
AddOption(&diagpc_type,
"-pc",
"--pctype",
84 "Select diagonal preconditioner type" 85 " (0:Diagonal, 1:BlockDiagonal).");
86 args.
AddOption(&serial_refinement_levels,
"-rs",
"--ref-serial",
87 "Number of uniform refinements on the serial mesh.");
88 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
90 "Enable or disable GLVis visualization.");
91 args.
AddOption(¶view,
"-pv",
"--paraview",
"-no-pv",
93 "Enable or disable ParaView DataCollection output.");
96 Device device(device_config);
103 Mesh::MakeCartesian3D(8, 2, 2, Element::HEXAHEDRON, 8.0, 1.0, 1.0);
106 MFEM_ABORT(
"This example only works in 3D.");
110 for (
int l = 0; l < serial_refinement_levels; l++)
112 mesh.UniformRefinement();
115 ParMesh pmesh(MPI_COMM_WORLD, mesh);
153 displaced_attr[2] = 1;
168 static_cast<ElasticityDiagonalPreconditioner::Type>(diagpc_type));
184 newton.
Mult(zero, U);
193 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
194 sol_sock.precision(8);
195 sol_sock <<
"solution\n" << pmesh << U_gf << flush;
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Conjugate gradient method.
void SetMaterial(const material_type &material)
Set the material type.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
void SetDataFormat(VTKFormat fmt)
Helper class for ParaView visualization data.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
int material(Vector &x, Vector &xmin, Vector &xmax)
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetMaxIter(int max_it)
Newton's method for solving F(x)=b for a given operator F.
void SetHighOrderOutput(bool high_order_output_)
void SetTime(double t)
Set physical time (for time-dependent simulations)
virtual void Mult(const Vector &b, Vector &x) const
Solve the nonlinear system with right-hand side b.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void SetRelTol(double rtol)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
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 Distribute(const Vector *tv)
const Array< int > & GetPrescribedDisplacementTDofs()
Return the T-vector degrees of freedom that have been marked as displaced.
void display_banner(ostream &os)
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
int Size() const
Return the logical size of the array.
int main(int argc, char *argv[])
void SetLevelsOfDetail(int levels_of_detail_)
void SetPrescribedDisplacement(const Array< int > attr)
Set the attributes which mark the degrees of freedom that have a fixed displacement.
void SetEssentialAttributes(const Array< int > attr)
Set the essential attributes which mark degrees of freedom for the solving process.
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
ElasticityDiagonalPreconditioner acts as a matrix-free preconditioner for ElasticityOperator.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
virtual void Save() override
Class for parallel grid function.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
ParFiniteElementSpace h1_fes_
H1 finite element space.
Class for parallel meshes.
void ParseCheck(std::ostream &out=mfem::out)
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.