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";
63 bool par_format =
false;
64 bool visualization = 1;
65 bool use_petsc =
true;
66 bool use_nonoverlapping =
false;
67 bool local_bdr_spec =
false;
68 const char *petscrc_file =
"";
71 args.
AddOption(&mesh_file,
"-m",
"--mesh",
74 "Finite element order (polynomial degree).");
75 args.
AddOption(&par_format,
"-pf",
"--parallel-format",
"-sf",
77 "Format to use when saving the results for VisIt.");
78 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
80 "Enable or disable GLVis visualization.");
81 args.
AddOption(&use_petsc,
"-usepetsc",
"--usepetsc",
"-no-petsc",
83 "Use or not PETSc to solve the linear system.");
84 args.
AddOption(&petscrc_file,
"-petscopts",
"--petscopts",
85 "PetscOptions file to use.");
86 args.
AddOption(&use_nonoverlapping,
"-nonoverlapping",
"--nonoverlapping",
87 "-no-nonoverlapping",
"--no-nonoverlapping",
88 "Use or not the block diagonal PETSc's matrix format "
89 "for non-overlapping domain decomposition.");
90 args.
AddOption(&local_bdr_spec,
"-local-bdr",
"--local-bdr",
"-no-local-bdr",
92 "Specify boundary dofs in local (Vdofs) ordering.");
112 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
121 (int)floor(log(10000./mesh->
GetNE())/log(2.)/
dim);
122 for (
int l = 0; l < ref_levels; l++)
134 int par_ref_levels = 2;
135 for (
int l = 0; l < par_ref_levels; l++)
154 std::cout <<
"***********************************************************\n";
155 std::cout <<
"dim(R) = " << dimR <<
"\n";
156 std::cout <<
"dim(W) = " << dimW <<
"\n";
157 std::cout <<
"dim(R+W) = " << dimR + dimW <<
"\n";
158 std::cout <<
"***********************************************************\n";
167 block_offsets[0] = 0;
168 block_offsets[1] = R_space->
GetVSize();
169 block_offsets[2] = W_space->
GetVSize();
173 block_trueOffsets[0] = 0;
174 block_trueOffsets[1] = R_space->
TrueVSize();
175 block_trueOffsets[2] = W_space->
TrueVSize();
191 BlockVector trueX(block_trueOffsets), trueRhs(block_trueOffsets);
228 if (!use_petsc) { Mh.Get(M); }
230 Mh.SetOperatorOwner(
false);
248 if (!use_petsc) { BT = B->Transpose(); }
249 else { pBT = pB->Transpose(); };
302 invS->iterative_mode =
false;
310 if (use_nonoverlapping)
341 MFEM_WARNING(
"Missing boundary dofs. This may cause solver failures.");
375 solver.
Mult(trueRhs, trueX);
381 <<
" iterations with a residual norm of "
386 std::cout <<
"MINRES did not converge in "
388 <<
" iterations. Residual norm is "
391 std::cout <<
"MINRES solver took " << chrono.
RealTime() <<
"s. \n";
397 std::string solvertype;
399 if (use_nonoverlapping)
408 solvertype =
"MINRES";
416 solver->
Mult(trueRhs, trueX);
421 std::cout << solvertype <<
" converged in "
423 <<
" iterations with a residual norm of "
428 std::cout << solvertype <<
" did not converge in "
430 <<
" iterations. Residual norm is "
433 std::cout << solvertype <<
" solver took "
445 u->MakeRef(R_space, x.GetBlock(0), 0);
446 p->MakeRef(W_space, x.GetBlock(1), 0);
447 u->Distribute(&(trueX.GetBlock(0)));
448 p->Distribute(&(trueX.GetBlock(1)));
450 int order_quad = max(2, 2*order+1);
457 real_t err_u =
u->ComputeL2Error(ucoeff, irs);
459 real_t err_p =
p->ComputeL2Error(pcoeff, irs);
464 std::cout <<
"|| u_h - u_ex || / || u_ex || = " << err_u / norm_u <<
"\n";
465 std::cout <<
"|| p_h - p_ex || / || p_ex || = " << err_p / norm_p <<
"\n";
471 ostringstream mesh_name, u_name, p_name;
472 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
473 u_name <<
"sol_u." << setfill(
'0') << setw(6) << myid;
474 p_name <<
"sol_p." << setfill(
'0') << setw(6) << myid;
476 ofstream mesh_ofs(mesh_name.str().c_str());
477 mesh_ofs.precision(8);
478 pmesh->
Print(mesh_ofs);
480 ofstream u_ofs(u_name.str().c_str());
484 ofstream p_ofs(p_name.str().c_str());
504 u_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
506 u_sock <<
"solution\n" << *pmesh << *
u <<
"window_title 'Velocity'"
512 p_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
514 p_sock <<
"solution\n" << *pmesh << *
p <<
"window_title 'Pressure'"
562 u(0) = - exp(xi)*sin(yi)*cos(zi);
563 u(1) = - exp(xi)*cos(yi)*cos(zi);
567 u(2) = exp(xi)*sin(yi)*sin(zi);
583 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.
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.
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the MINRES method.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
virtual void SetPreconditioner(Solver &pr)
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.
virtual void SetOperator(const Operator &op)
void SetPreconditioner(Solver &precond)
virtual void Mult(const Vector &b, Vector &x) const
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.
int Size() const
Returns the size of the vector.
Data collection with VisIt I/O routines.
virtual void Save() override
Save the collection and a VisIt root file.
virtual 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)
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)