37#error This example requires that MFEM is built with MFEM_USE_PETSC=YES
56 if (!it || it > 5) {
return; }
60 a->RecoverFEMSolution(X, *b, x);
66 MPI_Comm_size(mesh->
GetComm(),&num_procs);
67 MPI_Comm_rank(mesh->
GetComm(),&myid);
69 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
70 sol_sock.precision(8);
71 sol_sock <<
"solution\n" << *mesh << x
72 <<
"window_title 'Iteration no " << it <<
"'" << flush;
76int main(
int argc,
char *argv[])
85 const char *mesh_file =
"../../data/star.mesh";
87 bool static_cond =
false;
89 bool visualization =
false;
90 const char *device_config =
"cpu";
91 bool use_petsc =
true;
92 const char *petscrc_file =
"";
93 bool petscmonitor =
false;
94 bool forcewrap =
false;
98 args.
AddOption(&mesh_file,
"-m",
"--mesh",
101 "Finite element order (polynomial degree) or -1 for"
102 " isoparametric space.");
103 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
104 "--no-static-condensation",
"Enable static condensation.");
105 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
106 "--no-partial-assembly",
"Enable Partial Assembly.");
107 args.
AddOption(&device_config,
"-d",
"--device",
108 "Device configuration string, see Device::Configure().");
109 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
110 "--no-visualization",
111 "Enable or disable GLVis visualization.");
112 args.
AddOption(&use_petsc,
"-usepetsc",
"--usepetsc",
"-no-petsc",
114 "Use or not PETSc to solve the linear system.");
115 args.
AddOption(&petscrc_file,
"-petscopts",
"--petscopts",
116 "PetscOptions file to use.");
117 args.
AddOption(&petscmonitor,
"-petscmonitor",
"--petscmonitor",
118 "-no-petscmonitor",
"--no-petscmonitor",
119 "Enable or disable GLVis visualization of residual.");
120 args.
AddOption(&forcewrap,
"-forcewrap",
"--forcewrap",
121 "-noforce-wrap",
"--noforce-wrap",
122 "Force matrix-free.");
123 args.
AddOption(&useh2,
"-useh2",
"--useh2",
"-no-h2",
125 "Use or not the H2 matrix solver.");
142 Device device(device_config);
143 if (myid == 0) { device.
Print(); }
151 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
160 (int)floor(log(10000./mesh->
GetNE())/log(2.)/
dim);
161 for (
int l = 0; l < ref_levels; l++)
173 int par_ref_levels = 2;
174 for (
int l = 0; l < par_ref_levels; l++)
196 cout <<
"Using isoparametric FEs: " << fec->
Name() << endl;
208 cout <<
"Number of finite element unknowns: " << size << endl;
241 if (pa) {
a->SetAssemblyLevel(AssemblyLevel::PARTIAL); }
248 if (static_cond) {
a->EnableStaticCondensation(); }
253 a->FormLinearSystem(ess_tdof_list, x, *
b, A, X, B);
292 bool wrap = forcewrap ? true : (pa ? true : !strlen(petscrc_file));
326 UserMonitor mymon(
a,
b);
327 if (visualization && petscmonitor)
344 ostringstream mesh_name, sol_name;
345 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
346 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
348 ofstream mesh_ofs(mesh_name.str().c_str());
349 mesh_ofs.precision(8);
350 pmesh->
Print(mesh_ofs);
352 ofstream sol_ofs(sol_name.str().c_str());
353 sol_ofs.precision(8);
363 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
364 sol_sock.precision(8);
365 sol_sock <<
"solution\n" << *pmesh << x << flush;
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
int Size() const
Return the logical size of the array.
Conjugate gradient method.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
A coefficient that is constant across space and time.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Class for domain integration .
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
virtual const char * Name() const
Arbitrary order H1-conforming (continuous) finite elements.
The BoomerAMG solver in hypre.
Wrapper for hypre's ParCSR matrix class.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
void SetRelTol(real_t rtol)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetMaxIter(int max_it)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
void GetNodes(Vector &node_coord) const
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of 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).
Pointer to an Operator of a specified type.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
Operator * Ptr() const
Access the underlying Operator pointer.
Jacobi smoothing for a given bilinear form (no matrix necessary).
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear sys...
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.
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalTrueVSize() const
ParMesh * GetParMesh() const
Class for parallel grid function.
void Save(std::ostream &out) const override
Class for parallel meshes.
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
virtual void SetOperator(const Operator &op)
void SetPreconditioner(Solver &precond)
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
Abstract class for monitoring PETSc's solvers.
void SetAbsTol(real_t tol)
void SetPrintLevel(int plev)
void SetMaxIter(int max_iter)
void SetRelTol(real_t tol)
void SetMonitor(PetscSolverMonitor *ctx)
Sets user-defined monitoring routine.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void Randomize(int seed=0)
Set random values in the vector.
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.