35 #ifndef MFEM_USE_PETSC 36 #error This example requires that MFEM is built with MFEM_USE_PETSC=YES 49 int main(
int argc,
char *argv[])
54 Mpi::Init(argc, argv);
55 int num_procs = Mpi::WorldSize();
56 int myid = Mpi::WorldRank();
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);
220 !use_petsc ? Operator::Hypre_ParCSR :
221 (use_nonoverlapping ? Operator::PETSC_MATIS : Operator::PETSC_MATAIJ);
228 if (!use_petsc) { Mh.Get(M); }
230 Mh.SetOperatorOwner(
false);
248 if (!use_petsc) { BT = B->Transpose(); }
249 else { pBT = pB->Transpose(); };
270 use_nonoverlapping ? Operator::PETSC_MATIS :
271 Operator::PETSC_MATAIJ);
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);
452 for (
int i=0; i < Geometry::NumGeom; ++i)
457 double err_u =
u->ComputeL2Error(ucoeff, irs);
459 double 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());
494 DataCollection::SERIAL_FORMAT :
495 DataCollection::PARALLEL_FORMAT);
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);
void uFun_ex(const Vector &x, Vector &u)
void SetPreconditioner(Solver &precond)
void SetPrintLevel(int plev)
Class for domain integration L(v) := (f, v)
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Class for an integration rule - an Array of IntegrationPoint.
Abstract class for PETSc's preconditioners.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
A class to handle Vectors in a block fashion.
virtual void SetOperator(const Operator &op)
A coefficient that is constant across space and time.
void PrintOptions(std::ostream &out) const
Print the options.
int TrueVSize() const
Obsolete, kept for backward compatibility.
Wrapper for PETSc's matrix class.
Abstract class for PETSc's linear solvers.
void PrintUsage(std::ostream &out) const
Print the usage message.
double pFun_ex(const Vector &x)
Pointer to an Operator of a specified type.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
int Size() const
Returns the size of the vector.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
bool Good() const
Return true if the command line options were parsed successfully.
Abstract parallel finite element space.
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
void SetRelTol(double tol)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void SetMaxIter(int max_iter)
MPI_Comm GetComm() const
Get the associated MPI communicator.
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
double GetFinalNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult()...
void Stop()
Stop the stopwatch.
double f_natural(const Vector &x)
The BoomerAMG solver in hypre.
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const
Determine the boundary degrees of freedom.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
int main(int argc, char *argv[])
Jacobi preconditioner in hypre.
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Auxiliary class for BDDC customization.
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
HYPRE_BigInt GlobalTrueVSize() const
void Start()
Start the stopwatch. The elapsed time is not cleared.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
A general vector function coefficient.
Wrapper for hypre's parallel vector class.
Type
Enumeration defining IDs for some classes derived from Operator.
void SetAbsTol(double atol)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
double p(const Vector &x, double t)
void SetRelTol(double rtol)
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
virtual void Save() override
Save the collection and a VisIt root file.
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 SetSize(int nsize)
Change the logical size of the array, keep existing entries.
bool GetConverged() const
Returns true if the last call to Mult() converged successfully.
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
double ComputeGlobalLpNorm(double p, Coefficient &coeff, ParMesh &pmesh, const IntegrationRule *irs[])
Compute the global Lp norm of a function f. .
HypreParMatrix * Transpose() const
Returns the transpose of *this.
int GetNE() const
Returns number of elements.
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
double gFun(const Vector &x)
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
void SetAbsTol(double tol)
void fFun(const Vector &x, Vector &f)
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Abstract class for hypre's solvers and preconditioners.
int Size() const
Return the logical size of the array.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
A general function coefficient.
void Print(std::ostream &out=mfem::out) const override
double u(const Vector &xvec)
Class for parallel grid function.
void SetNatBdrDofs(const Array< int > *natdofs, bool loc=false)
Specify dofs on the natural boundary.
Wrapper for hypre's ParCSR matrix class.
A class to handle Block systems in a matrix-free implementation.
Class for parallel meshes.
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
void Clear()
Clear the elapsed time on the stopwatch and restart it if it's running.
void SetBlock(int iRow, int iCol, Operator *op, double c=1.0)
Add a block op in the block-entry (iblock, jblock).
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
Vector & GetBlock(int i)
Get the i-th vector in the block.
Arbitrary order "L2-conforming" discontinuous finite elements.
double f(const Vector &p)
virtual void RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
HYPRE_BigInt * GetRowStarts() const
Return the parallel row partitioning array.