39 return (0.25 * (2.0 + x[0]) - x[2]) * (x[2] + 0.25 * (2.0 + x[0]));
45 du[0] = 0.125 * (2.0 + x[0]) * x[1] * x[1];
46 du[1] = -0.125 * (2.0 + x[0]) * x[0] * x[1];
64 int main(
int argc,
char *argv[])
71 bool static_cond =
false;
72 bool visualization =
true;
75 args.
AddOption(&mesh_type,
"-mt",
"--mesh-type",
76 "Mesh type: 3 - Triangular, 4 - Quadrilateral.");
77 args.
AddOption(&mesh_order,
"-mo",
"--mesh-order",
78 "Geometric order of the curved mesh.");
79 args.
AddOption(&ref_levels,
"-r",
"--refine",
80 "Number of times to refine the mesh uniformly in serial.");
82 "Finite element order (polynomial degree).");
83 args.
AddOption(&static_cond,
"-sc",
"--static-condensation",
"-no-sc",
84 "--no-static-condensation",
"Enable static condensation.");
85 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
87 "Enable or disable GLVis visualization.");
97 for (
int l = 0; l < ref_levels; l++)
110 cout <<
"Number of finite element unknowns: " 145 a.AddDomainIntegrator(integ);
151 if (static_cond) {
a.EnableStaticCondensation(); }
156 a.FormLinearSystem(ess_tdof_list, x,
b, A, X, B);
158 cout <<
"Size of linear system: " << A->
Height() << endl;
163 PCG(*A, M, B, X, 1, 200, 1e-12, 0.0);
166 a.RecoverFEMSolution(X,
b, x);
172 cout <<
"|u - u_h|_2 = " << error << endl;
181 cout <<
"|f - f_h|_2 = " << flux_err << endl;
185 ofstream mesh_ofs(
"refined.mesh");
186 mesh_ofs.precision(8);
187 mesh->
Print(mesh_ofs);
188 ofstream sol_ofs(
"sol.gf");
189 sol_ofs.precision(8);
198 sol_sock.precision(8);
199 sol_sock <<
"solution\n" << *mesh << x
200 <<
"window_title 'Solution'\n" << flush;
203 flux_sock.precision(8);
204 flux_sock <<
"solution\n" << *mesh << flux
206 <<
"window_geometry 402 0 400 350\n" 207 <<
"window_title 'Flux'\n" << flush;
224 mesh =
new Mesh(2, 12, 16, 8, 3);
267 mesh =
new Mesh(2, 8, 4, 8, 3);
294 MFEM_ABORT(
"Unrecognized mesh type " << type <<
"!");
309 if (fabs(x[1] + 1.0) < tol)
311 theta = 0.25 * M_PI * (x[0] - 2.0);
313 else if (fabs(x[0] - 1.0) < tol)
315 theta = 0.25 * M_PI * x[1];
317 else if (fabs(x[1] - 1.0) < tol)
319 theta = 0.25 * M_PI * (2.0 - x[0]);
321 else if (fabs(x[0] + 1.0) < tol)
323 theta = 0.25 * M_PI * (4.0 - x[1]);
327 cerr <<
"side not recognized " 328 << x[0] <<
" " << x[1] <<
" " << x[2] << endl;
333 r[2] = 0.25 * (2.0 * x[2] - 1.0) * (r[0] + 2.0);
340 double a = 17.0 - 2.0 * x[0] * (1.0 + x[0]);
341 s(0,0) = 0.5 + x[0] * x[0] * (8.0 /
a - 0.5);
342 s(0,1) = x[0] * x[1] * (8.0 /
a - 0.5);
345 s(1,1) = 0.5 * x[0] * x[0] + 8.0 * x[1] * x[1] /
a;
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Class for domain integration L(v) := (f, v)
void duExact(const Vector &x, Vector &du)
Class for grid function - Vector with associated FE space.
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
A coefficient that is constant across space and time.
void SetSize(int s)
Resize the vector to size s.
Pointer to an Operator of a specified type.
void sigmaFunc(const Vector &x, DenseMatrix &s)
int AddTriangle(int v1, int v2, int v3, int attr=1)
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Data type dense matrix using column-major storage.
void Transform(void(*f)(const Vector &, Vector &))
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Data type for Gauss-Seidel smoother of sparse matrix.
int AddVertex(double x, double y=0.0, double z=0.0)
A matrix coefficient with an optional scalar coefficient multiplier q. The matrix function can either...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
int AddBdrSegment(int v1, int v2, int attr=1)
int main(int argc, char *argv[])
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
A general vector function coefficient.
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a 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...
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
void trans(const Vector &x, Vector &r)
double uExact(const Vector &x)
int Size() const
Return the logical size of the array.
A general function coefficient.
virtual void Print(std::ostream &os=mfem::out) const
virtual void ComputeFlux(BilinearFormIntegrator &blfi, GridFunction &flux, bool wcoef=true, int subdomain=-1)
Arbitrary order H1-conforming (continuous) finite elements.
void fluxExact(const Vector &x, Vector &f)
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
void ParseCheck(std::ostream &out=mfem::out)
double f(const Vector &p)
double sigma(const Vector &x)