35int main(
int argc,
char *argv[])
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 real_t 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 real_t 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;
265 a->FormLinearSystem(empty_tdof_list, x, *
b, A, X, B);
278 a->RecoverFEMSolution(X, *
b, x);
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 real_t l2 = x(0)*x(0) + x(1)*x(1) + x(2)*x(2);
334 real_t l2 = x(0)*x(0) + x(1)*x(1) + x(2)*x(2);
335 return 7*x(0)*x(1)/l2;
342 for (
int i = 0; i < nodes.
FESpace()->GetNDofs(); i++)
A coefficient that is constant across space and time.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Class for domain integration .
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
int DofToVDof(int dof, int vd, int ndofs=-1) const
Compute a single vdof corresponding to the index dof and the vector index vd.
A general function coefficient.
Class for grid function - Vector with associated FE space.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
FiniteElementSpace * FESpace()
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
Arbitrary order H1-conforming (continuous) finite elements.
The BoomerAMG solver in hypre.
void SetPrintLevel(int print_lvl)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
void SetMaxIter(int max_iter)
Wrapper for hypre's ParCSR matrix class.
Abstract class for hypre's solvers and preconditioners.
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Adds a quadrilateral to the mesh given by 4 vertices v1 through v4.
bool Nonconforming() const
int AddTriangle(int v1, int v2, int v3, int attr=1)
Adds a triangle to the mesh given by 3 vertices v1 through v3.
int AddVertex(real_t x, real_t y=0.0, real_t z=0.0)
int Dimension() const
Dimension of the reference space used within the elements.
void RandomRefinement(real_t prob, bool aniso=false, int nonconforming=-1, int nc_limit=0)
Refine each element with given probability. Uses GeneralRefinement.
void RefineAtVertex(const Vertex &vert, real_t eps=0.0, int nonconforming=-1)
Refine elements sharing the specified vertex. Uses GeneralRefinement.
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
int SpaceDimension() const
Dimension of the physical space containing the mesh.
void GetNodes(Vector &node_coord) const
void FinalizeTriMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a triangular Mesh.
virtual void SetNodalFESpace(FiniteElementSpace *nfes)
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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...
bool Good() const
Return true if the command line options were parsed successfully.
Abstract parallel finite element space.
HYPRE_BigInt GlobalTrueVSize() const
Class for parallel grid function.
void Save(std::ostream &out) const override
real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const override
Class for parallel meshes.
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
real_t Norml2() const
Returns the l2 norm of the vector.
void SnapNodes(Mesh &mesh)
real_t analytic_rhs(const Vector &x)
real_t analytic_solution(const Vector &x)