36#error This example requires that MFEM is built with MFEM_USE_PETSC=YES
49int main(
int argc,
char *argv[])
58 bool verbose = (myid == 0);
61 const char *mesh_file =
"../../data/star.mesh";
62 int ser_ref_levels = -1;
63 int par_ref_levels = 2;
65 bool par_format =
false;
66 bool visualization = 1;
67 bool use_petsc =
true;
68 bool use_nonoverlapping =
false;
69 bool local_bdr_spec =
false;
70 const char *petscrc_file =
"";
71 const char *device_config =
"cpu";
74 args.
AddOption(&mesh_file,
"-m",
"--mesh",
76 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
77 "Number of times to refine the mesh uniformly in serial.");
78 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
79 "Number of times to refine the mesh uniformly in parallel.");
81 "Finite element order (polynomial degree).");
82 args.
AddOption(&par_format,
"-pf",
"--parallel-format",
"-sf",
84 "Format to use when saving the results for VisIt.");
85 args.
AddOption(&device_config,
"-d",
"--device",
86 "Device configuration string, see Device::Configure().");
87 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
89 "Enable or disable GLVis visualization.");
90 args.
AddOption(&use_petsc,
"-usepetsc",
"--usepetsc",
"-no-petsc",
92 "Use or not PETSc to solve the linear system.");
93 args.
AddOption(&petscrc_file,
"-petscopts",
"--petscopts",
94 "PetscOptions file to use.");
95 args.
AddOption(&use_nonoverlapping,
"-nonoverlapping",
"--nonoverlapping",
96 "-no-nonoverlapping",
"--no-nonoverlapping",
97 "Use or not the block diagonal PETSc's matrix format "
98 "for non-overlapping domain decomposition.");
99 args.
AddOption(&local_bdr_spec,
"-local-bdr",
"--local-bdr",
"-no-local-bdr",
101 "Specify boundary dofs in local (Vdofs) ordering.");
118 Device device(device_config);
119 if (myid == 0) { device.
Print(); }
127 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
135 if (ser_ref_levels < 0)
137 ser_ref_levels = (int)floor(log(10000./mesh->
GetNE())/log(2.)/
dim);
139 for (
int l = 0; l < ser_ref_levels; l++)
151 for (
int l = 0; l < par_ref_levels; l++)
170 std::cout <<
"***********************************************************\n";
171 std::cout <<
"dim(R) = " << dimR <<
"\n";
172 std::cout <<
"dim(W) = " << dimW <<
"\n";
173 std::cout <<
"dim(R+W) = " << dimR + dimW <<
"\n";
174 std::cout <<
"***********************************************************\n";
183 block_offsets[0] = 0;
184 block_offsets[1] = R_space->
GetVSize();
185 block_offsets[2] = W_space->
GetVSize();
189 block_trueOffsets[0] = 0;
190 block_trueOffsets[1] = R_space->
TrueVSize();
191 block_trueOffsets[2] = W_space->
TrueVSize();
207 BlockVector x(block_offsets, mt), rhs(block_offsets, mt);
208 BlockVector trueX(block_trueOffsets, mt), trueRhs(block_trueOffsets, mt);
249 if (!use_petsc) { Mh.Get(M); }
251 Mh.SetOperatorOwner(
false);
269 if (!use_petsc) { BT = B->Transpose(); }
270 else { pBT = pB->Transpose(); };
323 invS->iterative_mode =
false;
331 if (use_nonoverlapping)
362 MFEM_WARNING(
"Missing boundary dofs. This may cause solver failures.");
396 solver.
Mult(trueRhs, trueX);
402 <<
" iterations with a residual norm of "
407 std::cout <<
"MINRES did not converge in "
409 <<
" iterations. Residual norm is "
412 std::cout <<
"MINRES solver took " << chrono.
RealTime() <<
"s. \n";
418 std::string solvertype;
420 if (use_nonoverlapping)
429 solvertype =
"MINRES";
437 solver->
Mult(trueRhs, trueX);
442 std::cout << solvertype <<
" converged in "
444 <<
" iterations with a residual norm of "
449 std::cout << solvertype <<
" did not converge in "
451 <<
" iterations. Residual norm is "
454 std::cout << solvertype <<
" solver took "
466 u->MakeRef(R_space, x.GetBlock(0), 0);
467 p->MakeRef(W_space, x.GetBlock(1), 0);
468 u->Distribute(&(trueX.GetBlock(0)));
469 p->Distribute(&(trueX.GetBlock(1)));
471 int order_quad = max(2, 2*order+1);
478 real_t err_u =
u->ComputeL2Error(ucoeff, irs);
480 real_t err_p =
p->ComputeL2Error(pcoeff, irs);
485 std::cout <<
"|| u_h - u_ex || / || u_ex || = " << err_u / norm_u <<
"\n";
486 std::cout <<
"|| p_h - p_ex || / || p_ex || = " << err_p / norm_p <<
"\n";
492 ostringstream mesh_name, u_name, p_name;
493 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
494 u_name <<
"sol_u." << setfill(
'0') << setw(6) << myid;
495 p_name <<
"sol_p." << setfill(
'0') << setw(6) << myid;
497 ofstream mesh_ofs(mesh_name.str().c_str());
498 mesh_ofs.precision(8);
499 pmesh->
Print(mesh_ofs);
501 ofstream u_ofs(u_name.str().c_str());
505 ofstream p_ofs(p_name.str().c_str());
525 u_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
527 u_sock <<
"solution\n" << *pmesh << *
u <<
"window_title 'Velocity'"
533 p_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
535 p_sock <<
"solution\n" << *pmesh << *
p <<
"window_title 'Pressure'"
49int main(
int argc,
char *argv[]) {
…}
583 u(0) = - exp(xi)*sin(yi)*cos(zi);
584 u(1) = - exp(xi)*cos(yi)*cos(zi);
588 u(2) = exp(xi)*sin(yi)*sin(zi);
604 return exp(xi)*sin(yi)*cos(zi);
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
A class to handle Block diagonal preconditioners in a matrix-free implementation.
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
A class to handle Block systems in a matrix-free implementation.
void SetBlock(int iRow, int iCol, Operator *op, real_t c=1.0)
Add a block op in the block-entry (iblock, jblock).
A class to handle Vectors in a block fashion.
Vector & GetBlock(int i)
Get the i-th vector in the block.
A coefficient that is constant across space and time.
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
static MemoryType GetMemoryType()
(DEPRECATED) Equivalent to GetDeviceMemoryType().
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...
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
A general function coefficient.
The BoomerAMG solver in hypre.
Jacobi preconditioner in hypre.
Wrapper for hypre's ParCSR matrix class.
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
HYPRE_BigInt * GetRowStarts() const
Return the parallel row partitioning array.
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
HypreParMatrix * Transpose() const
Returns the transpose of *this.
Wrapper for hypre's parallel vector class.
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(),...
Class for an integration rule - an Array of IntegrationPoint.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
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)
Arbitrary order "L2-conforming" discontinuous finite elements.
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the MINRES method.
void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
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 GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
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.
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
Type
Enumeration defining IDs for some classes derived from Operator.
@ PETSC_MATIS
ID for class PetscParMatrix, MATIS format.
@ Hypre_ParCSR
ID for class HypreParMatrix.
@ PETSC_MATAIJ
ID for class PetscParMatrix, MATAIJ format.
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
void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const override
Determine the boundary degrees of freedom.
int TrueVSize() const
Obsolete, kept for backward compatibility.
Class for parallel grid function.
Class for parallel meshes.
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Auxiliary class for BDDC customization.
void SetNatBdrDofs(const Array< int > *natdofs, bool loc=false)
Specify dofs on the natural boundary.
Abstract class for PETSc's linear solvers.
void SetOperator(const Operator &op) override
void SetPreconditioner(Solver &precond)
void Mult(const Vector &b, Vector &x) const override
Application of the solver.
Wrapper for PETSc's matrix class.
MPI_Comm GetComm() const
Get the associated MPI communicator.
Abstract class for PETSc's preconditioners.
void SetAbsTol(real_t tol)
void SetPrintLevel(int plev)
void SetMaxIter(int max_iter)
void SetRelTol(real_t tol)
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
void Start()
Start the stopwatch. The elapsed time is not cleared.
void Stop()
Stop the stopwatch.
void Clear()
Clear the elapsed time on the stopwatch and restart it if it's running.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
A general vector function coefficient.
void SyncAliasMemory(const Vector &v) const
Update the alias memory location of the vector to match v.
int Size() const
Returns the size of the vector.
Data collection with VisIt I/O routines.
void Save() override
Save the collection and a VisIt root file.
void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
real_t u(const Vector &xvec)
real_t ComputeGlobalLpNorm(real_t p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
MemoryType
Memory types supported by MFEM.
std::function< real_t(const Vector &)> f(real_t mass_coeff)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
real_t p(const Vector &x, real_t t)
real_t pFun_ex(const Vector &x)
void uFun_ex(const Vector &x, Vector &u)
real_t gFun(const Vector &x)
real_t f_natural(const Vector &x)
void fFun(const Vector &x, Vector &f)