30#error This example requires that MFEM is built with MFEM_USE_PETSC=YES
41int main(
int argc,
char *argv[])
50 const char *mesh_file =
"../../data/star.mesh";
51 int ser_ref_levels = -1;
52 int par_ref_levels = 2;
55 bool static_cond =
false;
56 bool hybridization =
false;
57 bool visualization = 1;
58 bool use_petsc =
true;
59 const char *petscrc_file =
"";
60 bool use_nonoverlapping =
false;
63 args.
AddOption(&mesh_file,
"-m",
"--mesh",
65 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
66 "Number of times to refine the mesh uniformly in serial.");
67 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
68 "Number of times to refine the mesh uniformly in parallel.");
70 "Finite element order (polynomial degree).");
71 args.
AddOption(&set_bc,
"-bc",
"--impose-bc",
"-no-bc",
"--dont-impose-bc",
72 "Impose or not essential boundary conditions.");
73 args.
AddOption(&
freq,
"-f",
"--frequency",
"Set the frequency for the exact"
75 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
76 "--no-static-condensation",
"Enable static condensation.");
77 args.
AddOption(&hybridization,
"-hb",
"--hybridization",
"-no-hb",
78 "--no-hybridization",
"Enable hybridization.");
79 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
81 "Enable or disable GLVis visualization.");
82 args.
AddOption(&use_petsc,
"-usepetsc",
"--usepetsc",
"-no-petsc",
84 "Use or not PETSc to solve the linear system.");
85 args.
AddOption(&petscrc_file,
"-petscopts",
"--petscopts",
86 "PetscOptions file to use.");
87 args.
AddOption(&use_nonoverlapping,
"-nonoverlapping",
"--nonoverlapping",
88 "-no-nonoverlapping",
"--no-nonoverlapping",
89 "Use or not the block diagonal PETSc's matrix format "
90 "for non-overlapping domain decomposition.");
111 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
120 if (ser_ref_levels < 0)
122 ser_ref_levels = (int)floor(log(1000./mesh->
GetNE())/log(2.)/
dim);
124 for (
int l = 0; l < ser_ref_levels; l++)
136 for (
int l = 0; l < par_ref_levels; l++)
149 cout <<
"Number of finite element unknowns: " << size << endl;
160 ess_bdr = set_bc ? 1 : 0;
200 a->EnableStaticCondensation();
202 else if (hybridization)
225 cout <<
"Size of linear system: " << glob_size << endl;
237 (
a->StaticCondensationIsEnabled() ?
a->SCParFESpace() : fespace);
238 if (
dim == 2) { prec =
new HypreAMS(A, prec_fespace); }
239 else { prec =
new HypreADS(A, prec_fespace); }
249 a->SetOperatorType(use_nonoverlapping ?
255 cout <<
"Size of linear system: " << A.
M() << endl;
259 if (use_nonoverlapping)
262 (
a->StaticCondensationIsEnabled() ?
a->SCParFESpace() :
263 (hfes ? NULL : fespace));
294 cout <<
"\n|| F_h - F ||_{L^2} = " <<
err <<
'\n' << endl;
301 ostringstream mesh_name, sol_name;
302 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
303 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
305 ofstream mesh_ofs(mesh_name.str().c_str());
306 mesh_ofs.precision(8);
307 pmesh->
Print(mesh_ofs);
309 ofstream sol_ofs(sol_name.str().c_str());
310 sol_ofs.precision(8);
320 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
321 sol_sock.precision(8);
322 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.
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.
for Raviart-Thomas elements
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
The Auxiliary-space Divergence Solver in hypre.
The Auxiliary-space Maxwell Solver in hypre.
The BoomerAMG solver in hypre.
Wrapper for hypre's ParCSR matrix class.
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Abstract class for hypre's solvers and preconditioners.
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.
int SpaceDimension() const
Dimension of the physical space containing the mesh.
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).
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
Form a constrained linear system using a matrix-free approach.
@ PETSC_MATIS
ID for class PetscParMatrix, MATIS format.
@ PETSC_MATAIJ
ID for class PetscParMatrix, MATAIJ format.
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
Class for parallel grid function.
void Save(std::ostream &out) const override
real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) 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
Auxiliary class for BDDC customization.
void SetEssBdrDofs(const Array< int > *essdofs, bool loc=false)
Specify dofs on the essential boundary.
void SetSpace(ParFiniteElementSpace *fe)
Wrapper for PETSc's matrix class.
PetscInt M() const
Returns the global number of rows.
Abstract class for PETSc's preconditioners.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
A general vector function coefficient.
int Size() const
Returns the size of the vector.
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
std::function< real_t(const Vector &)> f(real_t mass_coeff)
real_t p(const Vector &x, real_t t)
void f_exact(const Vector &, Vector &)
void F_exact(const Vector &, Vector &)