33 #ifndef MFEM_USE_PETSC
34 #error This example requires that MFEM is built with MFEM_USE_PETSC=YES
40 int main(
int argc,
char *argv[])
43 Mpi::Init(argc, argv);
44 int num_procs = Mpi::WorldSize();
45 int myid = Mpi::WorldRank();
49 const char *mesh_file =
"../../data/star.mesh";
51 bool visualization =
true;
52 int max_dofs = 100000;
53 bool use_petsc =
true;
54 const char *petscrc_file =
"";
55 bool use_nonoverlapping =
false;
58 args.
AddOption(&mesh_file,
"-m",
"--mesh",
61 "Finite element order (polynomial degree).");
62 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
64 "Enable or disable GLVis visualization.");
65 args.
AddOption(&max_dofs,
"-md",
"--max_dofs",
66 "Maximum number of dofs.");
67 args.
AddOption(&use_petsc,
"-usepetsc",
"--usepetsc",
"-no-petsc",
69 "Use or not PETSc to solve the linear system.");
70 args.
AddOption(&petscrc_file,
"-petscopts",
"--petscopts",
71 "PetscOptions file to use.");
72 args.
AddOption(&use_nonoverlapping,
"-nonoverlapping",
"--nonoverlapping",
73 "-no-nonoverlapping",
"--no-nonoverlapping",
74 "Use or not the block diagonal PETSc's matrix format "
75 "for non-overlapping domain decomposition.");
95 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
111 ParMesh pmesh(MPI_COMM_WORLD, *mesh);
115 "Boundary attributes required in the mesh.");
148 sout.
open(vishost, visport);
153 cout <<
"Unable to connect to GLVis server at "
154 << vishost <<
':' << visport << endl;
155 cout <<
"GLVis visualization disabled.\n";
157 visualization =
false;
186 for (
int it = 0; ; it++)
191 cout <<
"\nAMR iteration " << it << endl;
192 cout <<
"Number of unknowns: " << global_dofs << endl;
208 const int copy_interior = 1;
213 Operator::PETSC_MATIS : Operator::PETSC_MATAIJ);
216 MPI_Barrier(MPI_COMM_WORLD);
219 MPI_Barrier(MPI_COMM_WORLD);
221 if (myid == 0) { cout <<
"PETSc assembly timing : " << time << endl; }
229 MPI_Barrier(MPI_COMM_WORLD);
232 MPI_Barrier(MPI_COMM_WORLD);
234 if (myid == 0) { cout <<
"HYPRE assembly timing : " << time << endl; }
245 pcg.SetPrintLevel(3);
255 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
256 sout <<
"solution\n" << pmesh << x << flush;
259 if (global_dofs > max_dofs)
263 cout <<
"Reached the maximum number of dofs. Stop." << endl;
278 refiner.
Apply(pmesh);
283 cout <<
"Stopping criterion satisfied. Stop." << endl;
int Size() const
Return the logical size of the array.
Class for domain integration L(v) := (f, v)
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Conjugate gradient method.
MPI_Comm GetComm() const
MPI communicator.
A coefficient that is constant across space and time.
Wrapper for PETSc's matrix class.
virtual void Update(bool want_transform=true)
HYPRE_BigInt GlobalTrueVSize() const
bool Stop() const
Check if STOP action is requested, e.g. stopping criterion is satisfied.
Abstract parallel finite element space.
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
The BoomerAMG solver in hypre.
bool Apply(Mesh &mesh)
Perform the mesh operation.
bool Nonconforming() const
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void EnsureNCMesh(bool simplices_nonconforming=false)
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void SetPrintLevel(int print_level)
Mesh refinement operator using an error threshold.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void PrintUsage(std::ostream &out) const
Print the usage message.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
int SpaceDimension() const
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
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...
NURBSExtension * NURBSext
Optional NURBS mesh extension.
void MFEMInitializePetsc()
Convenience functions to initialize/finalize PETSc.
void PrintOptions(std::ostream &out) const
Print the options.
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Arbitrary order H1-conforming (continuous) finite elements.
Class for parallel grid function.
void SetTotalErrorFraction(double fraction)
Set the total fraction used in the computation of the threshold. The default value is 1/2...
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.
Arbitrary order "L2-conforming" discontinuous finite elements.
bool Good() const
Return true if the command line options were parsed successfully.