35int main(
int argc,
char *argv[])
42 bool always_snap =
false;
43 bool visualization = 1;
46 args.
AddOption(&elem_type,
"-e",
"--elem",
47 "Type of elements to use: 0 - triangles, 1 - quads.");
49 "Finite element order (polynomial degree).");
50 args.
AddOption(&ref_levels,
"-r",
"--refine",
51 "Number of times to refine the mesh uniformly.");
52 args.
AddOption(&amr,
"-amr",
"--refine-locally",
53 "Additional local (non-conforming) refinement:"
54 " 1 = refine around north pole, 2 = refine randomly.");
55 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
57 "Enable or disable GLVis visualization.");
58 args.
AddOption(&always_snap,
"-snap",
"--always-snap",
"-no-snap",
60 "If true, snap nodes to the sphere initially and after each refinement "
61 "otherwise, snap only after the last refinement");
74 int Nvert = 8, Nelem = 6;
80 Mesh *mesh =
new Mesh(2, Nvert, Nelem, 0, 3);
86 { 1, 0, 0}, { 0, 1, 0}, {-1, 0, 0},
87 { 0, -1, 0}, { 0, 0, 1}, { 0, 0, -1}
89 const int tri_e[8][3] =
91 {0, 1, 4}, {1, 2, 4}, {2, 3, 4}, {3, 0, 4},
92 {1, 0, 5}, {2, 1, 5}, {3, 2, 5}, {0, 3, 5}
95 for (
int j = 0; j < Nvert; j++)
99 for (
int j = 0; j < Nelem; j++)
101 int attribute = j + 1;
108 const real_t quad_v[8][3] =
110 {-1, -1, -1}, {+1, -1, -1}, {+1, +1, -1}, {-1, +1, -1},
111 {-1, -1, +1}, {+1, -1, +1}, {+1, +1, +1}, {-1, +1, +1}
113 const int quad_e[6][4] =
115 {3, 2, 1, 0}, {0, 1, 5, 4}, {1, 2, 6, 5},
116 {2, 3, 7, 6}, {3, 0, 4, 7}, {4, 5, 6, 7}
119 for (
int j = 0; j < Nvert; j++)
123 for (
int j = 0; j < Nelem; j++)
125 int attribute = j + 1;
126 mesh->
AddQuad(quad_e[j], attribute);
137 for (
int l = 0; l <= ref_levels; l++)
145 if (always_snap || l == ref_levels)
153 Vertex target(0.0, 0.0, 1.0);
154 for (
int l = 0; l < 5; l++)
162 for (
int l = 0; l < 4; l++)
172 cout <<
"Number of unknowns: " << fespace->
GetTrueVSize() << endl;
201 a->FormLinearSystem(empty_tdof_list, x, *
b, A, X, B);
203#ifndef MFEM_USE_SUITESPARSE
207 PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
211 umf_solver.
Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
213 umf_solver.
Mult(B, X);
217 a->RecoverFEMSolution(X, *
b, x);
220 cout<<
"\nL2 norm of error: " << x.
ComputeL2Error(sol_coef) << endl;
225 ofstream mesh_ofs(
"sphere_refined.mesh");
226 mesh_ofs.precision(8);
227 mesh->
Print(mesh_ofs);
228 ofstream sol_ofs(
"sol.gf");
229 sol_ofs.precision(8);
239 sol_sock.precision(8);
240 sol_sock <<
"solution\n" << *mesh << x << flush;
254 real_t l2 = x(0)*x(0) + x(1)*x(1) + x(2)*x(2);
260 real_t l2 = x(0)*x(0) + x(1)*x(1) + x(2)*x(2);
261 return 7*x(0)*x(1)/l2;
268 for (
int i = 0; i < nodes.
FESpace()->GetNDofs(); i++)
A coefficient that is constant across space and time.
Class for domain integration .
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
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.
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
virtual real_t ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
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.
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.
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
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 *)
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.
Direct sparse solver using UMFPACK.
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
real_t Control[UMFPACK_CONTROL]
virtual void Mult(const Vector &b, Vector &x) const
Direct solution of the linear system using UMFPACK.
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)
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)