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[])
55 MPI_Init(&argc, &argv);
56 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
57 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
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.");
113 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
122 (int)floor(log(10000./mesh->
GetNE())/log(2.)/
dim);
123 for (
int l = 0; l < ref_levels; l++)
135 int par_ref_levels = 2;
136 for (
int l = 0; l < par_ref_levels; l++)
155 std::cout <<
"***********************************************************\n";
156 std::cout <<
"dim(R) = " << dimR <<
"\n";
157 std::cout <<
"dim(W) = " << dimW <<
"\n";
158 std::cout <<
"dim(R+W) = " << dimR + dimW <<
"\n";
159 std::cout <<
"***********************************************************\n";
168 block_offsets[0] = 0;
169 block_offsets[1] = R_space->
GetVSize();
170 block_offsets[2] = W_space->
GetVSize();
174 block_trueOffsets[0] = 0;
175 block_trueOffsets[1] = R_space->
TrueVSize();
176 block_trueOffsets[2] = W_space->
TrueVSize();
192 BlockVector trueX(block_trueOffsets), trueRhs(block_trueOffsets);
221 !use_petsc ? Operator::Hypre_ParCSR :
222 (use_nonoverlapping ? Operator::PETSC_MATIS : Operator::PETSC_MATAIJ);
229 if (!use_petsc) { Mh.Get(M); }
231 Mh.SetOperatorOwner(
false);
249 if (!use_petsc) { BT = B->Transpose(); }
250 else { pBT = pB->Transpose(); };
271 use_nonoverlapping ? Operator::PETSC_MATIS :
272 Operator::PETSC_MATAIJ);
303 invS->iterative_mode =
false;
311 if (use_nonoverlapping)
342 MFEM_WARNING(
"Missing boundary dofs. This may cause solver failures.");
376 solver.
Mult(trueRhs, trueX);
382 <<
" iterations with a residual norm of "
387 std::cout <<
"MINRES did not converge in "
389 <<
" iterations. Residual norm is "
392 std::cout <<
"MINRES solver took " << chrono.
RealTime() <<
"s. \n";
398 std::string solvertype;
400 if (use_nonoverlapping)
409 solvertype =
"MINRES";
417 solver->
Mult(trueRhs, trueX);
422 std::cout << solvertype <<
" converged in "
424 <<
" iterations with a residual norm of "
429 std::cout << solvertype <<
" did not converge in "
431 <<
" iterations. Residual norm is "
434 std::cout << solvertype <<
" solver took "
446 u->
MakeRef(R_space, x.GetBlock(0), 0);
447 p->
MakeRef(W_space, x.GetBlock(1), 0);
451 int order_quad = max(2, 2*order+1);
453 for (
int i=0; i < Geometry::NumGeom; ++i)
465 std::cout <<
"|| u_h - u_ex || / || u_ex || = " << err_u / norm_u <<
"\n";
466 std::cout <<
"|| p_h - p_ex || / || p_ex || = " << err_p / norm_p <<
"\n";
472 ostringstream mesh_name, u_name, p_name;
473 mesh_name <<
"mesh." << setfill(
'0') << setw(6) << myid;
474 u_name <<
"sol_u." << setfill(
'0') << setw(6) << myid;
475 p_name <<
"sol_p." << setfill(
'0') << setw(6) << myid;
477 ofstream mesh_ofs(mesh_name.str().c_str());
478 mesh_ofs.precision(8);
479 pmesh->
Print(mesh_ofs);
481 ofstream u_ofs(u_name.str().c_str());
485 ofstream p_ofs(p_name.str().c_str());
495 DataCollection::SERIAL_FORMAT :
496 DataCollection::PARALLEL_FORMAT);
505 u_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
507 u_sock <<
"solution\n" << *pmesh << *u <<
"window_title 'Velocity'"
513 p_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
515 p_sock <<
"solution\n" << *pmesh << *p <<
"window_title 'Pressure'"
565 u(0) = - exp(xi)*sin(yi)*cos(zi);
566 u(1) = - exp(xi)*cos(yi)*cos(zi);
570 u(2) = exp(xi)*sin(yi)*sin(zi);
586 return exp(xi)*sin(yi)*cos(zi);
int Size() const
Return the logical size of the array.
void SetPreconditioner(Solver &precond)
void SetPrintLevel(int plev)
Class for domain integration L(v) := (f, v)
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
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.
int GetNumIterations() const
Abstract class for PETSc's preconditioners.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void fFun(const Vector &x, Vector &f)
A class to handle Vectors in a block fashion.
virtual void SetOperator(const Operator &op)
virtual void GetEssentialVDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_dofs, int component=-1) const
Determine the boundary degrees of freedom.
A coefficient that is constant across space and time.
Wrapper for PETSc's matrix class.
Abstract class for PETSc's linear solvers.
int TrueVSize() const
Obsolete, kept for backward compatibility.
Pointer to an Operator of a specified type.
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
Abstract parallel finite element space.
void SetRelTol(double tol)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void SetMaxIter(int max_iter)
int main(int argc, char *argv[])
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
double GetFinalNorm() const
void Stop()
Stop the stopwatch.
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
The BoomerAMG solver in hypre.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
void SetPrintLevel(int print_lvl)
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. ...
Jacobi preconditioner in hypre.
double f_natural(const Vector &x)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Auxiliary class for BDDC customization.
MPI_Comm GetComm() const
Get the associated MPI communicator.
Data collection with VisIt I/O routines.
void SetMaxIter(int max_it)
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
HYPRE_Int GlobalTrueVSize() const
virtual void Print(std::ostream &out=mfem::out) const
virtual void SetFormat(int fmt)
Set the desired output mesh and data format.
double pFun_ex(const Vector &x)
HypreParMatrix * Transpose() const
Returns the transpose of *this.
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
void PrintUsage(std::ostream &out) const
Print the usage message.
void Start()
Clear the elapsed time and start the stopwatch.
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...
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.
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. .
void Distribute(const Vector *tv)
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
void uFun_ex(const Vector &x, Vector &u)
double gFun(const Vector &x)
void SetOperatorOwner(bool own=true)
Set the ownership flag for the held Operator.
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
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 PrintOptions(std::ostream &out) const
Print the options.
Abstract class for hypre's solvers and preconditioners.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
A general function coefficient.
void Get(OpType *&A) const
Return the Operator pointer statically cast to a given OpType.
HYPRE_Int * GetRowStarts() const
Return the parallel row partitioning array.
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.
virtual void Mult(const Vector &b, Vector &x) const
Application of the solver.
Class for parallel meshes.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
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)
bool Good() const
Return true if the command line options were parsed successfully.