35 int main(
int argc,
char *argv[])
38 Mpi::Init(argc, argv);
39 int num_procs = Mpi::WorldSize();
40 int myid = Mpi::WorldRank();
48 bool always_snap =
false;
49 bool visualization = 1;
50 const char *device_config =
"cpu";
53 args.
AddOption(&elem_type,
"-e",
"--elem",
54 "Type of elements to use: 0 - triangles, 1 - quads.");
56 "Finite element order (polynomial degree).");
57 args.
AddOption(&ref_levels,
"-r",
"--refine",
58 "Number of times to refine the mesh uniformly.");
59 args.
AddOption(&amr,
"-amr",
"--refine-locally",
60 "Additional local (non-conforming) refinement:"
61 " 1 = refine around north pole, 2 = refine randomly.");
62 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
64 "Enable or disable GLVis visualization.");
65 args.
AddOption(&always_snap,
"-snap",
"--always-snap",
"-no-snap",
67 "If true, snap nodes to the sphere initially and after each refinement "
68 "otherwise, snap only after the last refinement");
69 args.
AddOption(&device_config,
"-d",
"--device",
70 "Device configuration string, see Device::Configure().");
87 Device device(device_config);
88 if (myid == 0) { device.
Print(); }
94 int Nvert = 8, Nelem = 6;
100 Mesh *mesh =
new Mesh(2, Nvert, Nelem, 0, 3);
104 const double tri_v[6][3] =
106 { 1, 0, 0}, { 0, 1, 0}, {-1, 0, 0},
107 { 0, -1, 0}, { 0, 0, 1}, { 0, 0, -1}
109 const int tri_e[8][3] =
111 {0, 1, 4}, {1, 2, 4}, {2, 3, 4}, {3, 0, 4},
112 {1, 0, 5}, {2, 1, 5}, {3, 2, 5}, {0, 3, 5}
115 for (
int j = 0; j < Nvert; j++)
119 for (
int j = 0; j < Nelem; j++)
121 int attribute = j + 1;
128 const double quad_v[8][3] =
130 {-1, -1, -1}, {+1, -1, -1}, {+1, +1, -1}, {-1, +1, -1},
131 {-1, -1, +1}, {+1, -1, +1}, {+1, +1, +1}, {-1, +1, +1}
133 const int quad_e[6][4] =
135 {3, 2, 1, 0}, {0, 1, 5, 4}, {1, 2, 6, 5},
136 {2, 3, 7, 6}, {3, 0, 4, 7}, {4, 5, 6, 7}
139 for (
int j = 0; j < Nvert; j++)
143 for (
int j = 0; j < Nelem; j++)
145 int attribute = j + 1;
146 mesh->
AddQuad(quad_e[j], attribute);
158 for (
int l = 0; l <= ref_levels; l++)
174 Vertex target(0.0, 0.0, 1.0);
175 for (
int l = 0; l < 3; l++)
183 for (
int l = 0; l < 2; l++)
193 int par_ref_levels = 2;
194 for (
int l = 0; l < par_ref_levels; l++)
204 if (!always_snap || par_ref_levels < 1)
212 Vertex target(0.0, 0.0, 1.0);
213 for (
int l = 0; l < 2; l++)
221 for (
int l = 0; l < 2; l++)
234 cout <<
"Number of unknowns: " << size << endl;
287 cout <<
"\nL2 norm of error: " << error << endl;
293 ostringstream mesh_name, sol_name;
294 mesh_name <<
"sphere_refined." << setfill(
'0') << setw(6) << myid;
295 sol_name <<
"sol." << setfill(
'0') << setw(6) << myid;
297 ofstream mesh_ofs(mesh_name.str().c_str());
298 mesh_ofs.precision(8);
299 pmesh->
Print(mesh_ofs);
301 ofstream sol_ofs(sol_name.str().c_str());
302 sol_ofs.precision(8);
312 sol_sock <<
"parallel " << num_procs <<
" " << myid <<
"\n";
313 sol_sock.precision(8);
314 sol_sock <<
"solution\n" << *pmesh << x << flush;
328 double l2 = x(0)*x(0) + x(1)*x(1) + x(2)*x(2);
334 double l2 = x(0)*x(0) + x(1)*x(1) + x(2)*x(2);
335 return 7*x(0)*x(1)/l2;
Class for domain integration L(v) := (f, v)
int GetNDofs() const
Returns number of degrees of freedom.
Class for grid function - Vector with associated FE space.
int DofToVDof(int dof, int vd, int ndofs=-1) const
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
A coefficient that is constant across space and time.
double Norml2() const
Returns the l2 norm of the vector.
HYPRE_BigInt GlobalTrueVSize() const
int AddTriangle(int v1, int v2, int v3, int attr=1)
virtual void Save(std::ostream &out) const
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Abstract parallel finite element space.
void SetPrintLevel(int print_lvl)
The BoomerAMG solver in hypre.
int AddVertex(double x, double y=0.0, double z=0.0)
bool Nonconforming() const
void FinalizeTriMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a triangular Mesh.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
double analytic_solution(const Vector &x)
void RandomRefinement(double prob, bool aniso=false, int nonconforming=-1, int nc_limit=0)
Refine each element with given probability. Uses GeneralRefinement.
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
void SetNodalFESpace(FiniteElementSpace *nfes)
FiniteElementSpace * FESpace()
void PrintUsage(std::ostream &out) const
Print the usage message.
void SetMaxIter(int max_iter)
double analytic_rhs(const Vector &x)
int SpaceDimension() const
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
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 RefineAtVertex(const Vertex &vert, double eps=0.0, int nonconforming=-1)
Refine elements sharing the specified vertex. Uses GeneralRefinement.
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
void PrintOptions(std::ostream &out) const
Print the options.
void SnapNodes(Mesh &mesh)
Abstract class for hypre's solvers and preconditioners.
A general function coefficient.
void GetNodes(Vector &node_coord) const
Arbitrary order H1-conforming (continuous) finite elements.
void Print(std::ostream &out=mfem::out) const override
Class for parallel grid function.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Wrapper for hypre's ParCSR matrix class.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
Class for parallel meshes.
bool Good() const
Return true if the command line options were parsed successfully.