49 int main(
int argc,
char *argv[])
53 MPI_Init(&argc, &argv);
54 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
55 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
58 const char *mesh_file =
"../data/star.mesh";
61 const char *device_config =
"cpu";
62 bool visualization =
true;
65 args.
AddOption(&mesh_file,
"-m",
"--mesh",
68 "Finite element order (polynomial degree).");
69 args.
AddOption(&pa,
"-pa",
"--partial-assembly",
"-no-pa",
70 "--no-partial-assembly",
"Enable Partial Assembly.");
71 args.
AddOption(&device_config,
"-d",
"--device",
72 "Device configuration string, see Device::Configure().");
73 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
75 "Enable or disable GLVis visualization.");
93 Device device(device_config);
94 if (myid == 0) { device.
Print(); }
99 Mesh *mesh =
new Mesh(mesh_file, 1, 1);
115 ParMesh pmesh(MPI_COMM_WORLD, *mesh);
119 "Boundary attributes required in the mesh.");
157 sout.
open(vishost, visport);
162 cout <<
"Unable to connect to GLVis server at "
163 << vishost <<
':' << visport << endl;
164 cout <<
"GLVis visualization disabled.\n";
166 visualization =
false;
195 const int max_dofs = 100000;
196 for (
int it = 0; ; it++)
201 cout <<
"\nAMR iteration " << it << endl;
202 cout <<
"Number of unknowns: " << global_dofs << endl;
222 const int copy_interior = 1;
256 sout <<
"parallel " << num_procs <<
" " << myid <<
"\n";
257 sout <<
"solution\n" << pmesh << x << flush;
260 if (global_dofs > max_dofs)
264 cout <<
"Reached the maximum number of dofs. Stop." << endl;
273 refiner.
Apply(pmesh);
278 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.
A coefficient that is constant across space and time.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
virtual void Update(bool want_transform=true)
Pointer to an Operator of a specified type.
bool Stop() const
Check if STOP action is requested, e.g. stopping criterion is satisfied.
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Abstract parallel finite element space.
int main(int argc, char *argv[])
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 SetPrintLevel(int print_lvl)
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).
Jacobi smoothing for a given bilinear form (no matrix necessary).
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void SetPrintLevel(int print_level)
Mesh refinement operator using an error threshold.
void SetMaxIter(int max_it)
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.
HYPRE_Int GlobalTrueVSize() const
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 SetRelTol(double rtol)
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 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 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.
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...
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Class for parallel meshes.
Arbitrary order "L2-conforming" discontinuous finite elements.
bool Good() const
Return true if the command line options were parsed successfully.