66 static double a_ = 0.2;
77 double alpha,
double beta,
double gamma,
80 int main(
int argc,
char *argv[])
83 int ser_ref_levels = 2;
88 bool visualization =
true;
93 double rbc_a_val = 1.0;
94 double rbc_b_val = 1.0;
97 args.
AddOption(&h1,
"-h1",
"--continuous",
"-dg",
"--discontinuous",
98 "Select continuous \"H1\" or discontinuous \"DG\" basis.");
100 "Finite element order (polynomial degree) or -1 for"
101 " isoparametric space.");
103 "One of the two DG penalty parameters, typically +1/-1."
104 " See the documentation of class DGDiffusionIntegrator.");
106 "One of the two DG penalty parameters, should be positive."
107 " Negative values are replaced with (order+1)^2.");
108 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
109 "Number of times to refine the mesh uniformly in serial.");
110 args.
AddOption(&mat_val,
"-mat",
"--material-value",
111 "Constant value for material coefficient "
112 "in the Laplace operator.");
113 args.
AddOption(&dbc_val,
"-dbc",
"--dirichlet-value",
114 "Constant value for Dirichlet Boundary Condition.");
115 args.
AddOption(&nbc_val,
"-nbc",
"--neumann-value",
116 "Constant value for Neumann Boundary Condition.");
117 args.
AddOption(&rbc_a_val,
"-rbc-a",
"--robin-a-value",
118 "Constant 'a' value for Robin Boundary Condition: "
119 "du/dn + a * u = b.");
120 args.
AddOption(&rbc_b_val,
"-rbc-b",
"--robin-b-value",
121 "Constant 'b' value for Robin Boundary Condition: "
122 "du/dn + a * u = b.");
124 "Radius of holes in the mesh.");
125 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
126 "--no-visualization",
127 "Enable or disable GLVis visualization.");
134 if (kappa < 0 && !h1)
136 kappa = (order+1)*(order+1);
142 mfem::out <<
"Hole radius too small, resetting to 0.01.\n";
147 mfem::out <<
"Hole radius too large, resetting to 0.49.\n";
163 mfem::out <<
"Number of finite element unknowns: " << size << endl;
173 nbc_bdr = 0; nbc_bdr[0] = 1;
174 rbc_bdr = 0; rbc_bdr[1] = 1;
175 dbc_bdr = 0; dbc_bdr[2] = 1;
270 #ifndef MFEM_USE_SUITESPARSE
278 PCG(*A, M, B, X, 1, 500, 1e-12, 0.0);
282 GMRES(*A, M, B, X, 1, 500, 10, 1e-12, 0.0);
289 umf_solver.
Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
291 umf_solver.
Mult(B, X);
300 <<
"Verifying boundary conditions" << endl
301 <<
"=============================" << endl;
305 double error, avg =
IntegrateBC(u, dbc_bdr, 0.0, 1.0, dbc_val, error);
307 bool hom_dbc = (dbc_val == 0.0);
308 error /= hom_dbc ? 1.0 : fabs(dbc_val);
309 mfem::out <<
"Average of solution on Gamma_dbc:\t"
311 << (hom_dbc ?
"absolute" :
"relative")
312 <<
" error " << error << endl;
317 double error, avg =
IntegrateBC(u, nbc_bdr, 1.0, 0.0, nbc_val, error);
319 bool hom_nbc = (nbc_val == 0.0);
320 error /= hom_nbc ? 1.0 : fabs(nbc_val);
321 mfem::out <<
"Average of n.Grad(u) on Gamma_nbc:\t"
323 << (hom_nbc ?
"absolute" :
"relative")
324 <<
" error " << error << endl;
333 double error, avg =
IntegrateBC(u, nbc0_bdr, 1.0, 0.0, 0.0, error);
336 mfem::out <<
"Average of n.Grad(u) on Gamma_nbc0:\t"
338 << (hom_nbc ?
"absolute" :
"relative")
339 <<
" error " << error << endl;
345 double avg =
IntegrateBC(u, rbc_bdr, 1.0, rbc_a_val, rbc_b_val, error);
347 bool hom_rbc = (rbc_b_val == 0.0);
348 error /= hom_rbc ? 1.0 : fabs(rbc_b_val);
349 mfem::out <<
"Average of n.Grad(u)+a*u on Gamma_rbc:\t"
351 << (hom_rbc ?
"absolute" :
"relative")
352 <<
" error " << error << endl;
358 ofstream mesh_ofs(
"refined.mesh");
359 mesh_ofs.precision(8);
360 mesh->
Print(mesh_ofs);
361 ofstream sol_ofs(
"sol.gf");
362 sol_ofs.precision(8);
369 string title_str = h1 ?
"H1" :
"DG";
373 sol_sock.precision(8);
374 sol_sock <<
"solution\n" << *mesh << u
375 <<
"window_title '" << title_str <<
" Solution'"
376 <<
" keys 'mmc'" << flush;
386 void quad_trans(
double u,
double v,
double &x,
double &y,
bool log =
false)
390 double d = 4.0 * a * (M_SQRT2 - 2.0 *
a) * (1.0 - 2.0 * v);
392 double v0 = (1.0 + M_SQRT2) * (M_SQRT2 * a - 2.0 * v) *
393 ((4.0 - 3 * M_SQRT2) * a +
394 (8.0 * (M_SQRT2 - 1.0) * a - 2.0) * v) / d;
396 double r = 2.0 * ((M_SQRT2 - 1.0) * a * a * (1.0 - 4.0 *v) +
397 2.0 * (1.0 + M_SQRT2 *
398 (1.0 + 2.0 * (2.0 * a - M_SQRT2 - 1.0) *
a)) * v * v
401 double t = asin(v / r) * u / v;
405 << u <<
" " << v <<
" " << r <<
" " << v0 <<
" " << t
416 if (u[1] > 0.5 - tol || u[1] < -0.5 + tol)
421 if (u[0] > 1.0 - tol || u[0] < -1.0 + tol || fabs(u[0]) < tol)
429 if (u[1] > fabs(u[0] - 0.5))
435 if (u[1] < -fabs(u[0] - 0.5))
442 if (u[0] - 0.5 > fabs(u[1]))
448 if (u[0] - 0.5 < -fabs(u[1]))
458 if (u[1] > fabs(u[0] + 0.5))
464 if (u[1] < -fabs(u[0] + 0.5))
471 if (u[0] + 0.5 > fabs(u[1]))
477 if (u[0] + 0.5 < -fabs(u[1]))
490 Mesh * mesh =
new Mesh(2, 29, 16, 24, 2);
494 for (
int i=0; i<2; i++)
497 vi[0] = o + 0; vi[1] = o + 3; vi[2] = o + 4; vi[3] = o + 1;
500 vi[0] = o + 1; vi[1] = o + 4; vi[2] = o + 5; vi[3] = o + 2;
503 vi[0] = o + 5; vi[1] = o + 8; vi[2] = o + 9; vi[3] = o + 2;
506 vi[0] = o + 8; vi[1] = o + 12; vi[2] = o + 15; vi[3] = o + 9;
509 vi[0] = o + 11; vi[1] = o + 14; vi[2] = o + 15; vi[3] = o + 12;
512 vi[0] = o + 10; vi[1] = o + 13; vi[2] = o + 14; vi[3] = o + 11;
515 vi[0] = o + 6; vi[1] = o + 13; vi[2] = o + 10; vi[3] = o + 7;
518 vi[0] = o + 0; vi[1] = o + 6; vi[2] = o + 7; vi[3] = o + 3;
532 for (
int i=0; i<2; i++)
535 vi[0] = o + 7; vi[1] = o + 3; mesh->
AddBdrSegment(vi, 3 + i);
536 vi[0] = o + 10; vi[1] = o + 7; mesh->
AddBdrSegment(vi, 3 + i);
537 vi[0] = o + 11; vi[1] = o + 10; mesh->
AddBdrSegment(vi, 3 + i);
538 vi[0] = o + 12; vi[1] = o + 11; mesh->
AddBdrSegment(vi, 3 + i);
539 vi[0] = o + 8; vi[1] = o + 12; mesh->
AddBdrSegment(vi, 3 + i);
540 vi[0] = o + 5; vi[1] = o + 8; mesh->
AddBdrSegment(vi, 3 + i);
541 vi[0] = o + 4; vi[1] = o + 5; mesh->
AddBdrSegment(vi, 3 + i);
542 vi[0] = o + 3; vi[1] = o + 4; mesh->
AddBdrSegment(vi, 3 + i);
546 double a = a_ / M_SQRT2;
548 d[0] = -1.0; d[1] = -0.5; mesh->
AddVertex(d);
549 d[0] = -1.0; d[1] = 0.0; mesh->
AddVertex(d);
550 d[0] = -1.0; d[1] = 0.5; mesh->
AddVertex(d);
553 d[0] = -0.5 -
a; d[1] = 0.0; mesh->
AddVertex(d);
556 d[0] = -0.5; d[1] = -0.5; mesh->
AddVertex(d);
559 d[0] = -0.5; d[1] = 0.5; mesh->
AddVertex(d);
562 d[0] = -0.5 +
a; d[1] = 0.0; mesh->
AddVertex(d);
565 d[0] = 0.0; d[1] = -0.5; mesh->
AddVertex(d);
566 d[0] = 0.0; d[1] = 0.0; mesh->
AddVertex(d);
567 d[0] = 0.0; d[1] = 0.5; mesh->
AddVertex(d);
570 d[0] = 0.5 -
a; d[1] = 0.0; mesh->
AddVertex(d);
573 d[0] = 0.5; d[1] = -0.5; mesh->
AddVertex(d);
576 d[0] = 0.5; d[1] = 0.5; mesh->
AddVertex(d);
579 d[0] = 0.5 +
a; d[1] = 0.0; mesh->
AddVertex(d);
582 d[0] = 1.0; d[1] = -0.5; mesh->
AddVertex(d);
583 d[0] = 1.0; d[1] = 0.0; mesh->
AddVertex(d);
584 d[0] = 1.0; d[1] = 0.5; mesh->
AddVertex(d);
593 for (
int i = 0; i < v2v.Size() - 3; i++)
598 v2v[v2v.Size() - 3] = 0;
599 v2v[v2v.Size() - 2] = 1;
600 v2v[v2v.Size() - 1] = 2;
603 for (
int i = 0; i < mesh->
GetNE(); i++)
608 for (
int j = 0; j < nv; j++)
614 for (
int i = 0; i < mesh->
GetNBE(); i++)
619 for (
int j = 0; j < nv; j++)
629 for (
int l = 0; l < ref; l++)
640 double alpha,
double beta,
double gamma,
647 const bool a_is_zero = alpha == 0.0;
648 const bool b_is_zero = beta == 0.0;
651 MFEM_ASSERT(fes.
GetVDim() == 1,
"");
653 Vector shape, loc_dofs, w_nor;
656 for (
int i = 0; i < mesh.GetNBE(); i++)
658 if (bdr[mesh.GetBdrAttribute(i)-1] == 0) {
continue; }
661 if (FTr ==
nullptr) {
continue; }
664 MFEM_ASSERT(fe.
GetMapType() == FiniteElement::VALUE,
"");
665 const int int_order = 2*fe.
GetOrder() + 3;
693 val += alpha * dshape.
InnerProduct(w_nor, loc_dofs) / face_weight;
698 val += beta * (shape * loc_dofs);
702 nrm += ip.
weight * face_weight;
705 avg += val * ip.
weight * face_weight;
709 error += (val*val) * ip.
weight * face_weight;
714 if (std::abs(nrm) > 0.0)
722 error = (error >= 0.0) ? sqrt(error) : -sqrt(-error);
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for all finite elements.
int Size() const
Return the logical size of the array.
void trans(const Vector &u, Vector &x)
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
A coefficient that is constant across space and time.
virtual void GetVertices(Array< int > &v) const =0
Returns element's vertices.
void n4Vec(const Vector &x, Vector &n)
void SetSize(int s)
Resize the vector to size s.
int GetNBE() const
Returns number of boundary elements.
double Norml2() const
Returns the l2 norm of the vector.
Pointer to an Operator of a specified type.
Scalar coefficient defined as the product of two scalar coefficients or a scalar and a scalar coeffic...
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
void quad_trans(double u, double v, double &x, double &y, bool log=false)
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Data type dense matrix using column-major storage.
Class for boundary integration L(v) := (g, v)
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 ...
int GetNE() const
Returns number of elements.
int GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
void CalcOrtho(const DenseMatrix &J, Vector &n)
void RemoveInternalBoundaries()
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
int AddVertex(double x, double y=0.0, double z=0.0)
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Mesh * GetMesh() const
Returns the mesh.
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
int AddBdrSegment(int v1, int v2, int attr=1)
const Element * GetElement(int i) const
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)
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
void CalcPhysDShape(ElementTransformation &Trans, DenseMatrix &dshape) const
Evaluate the gradients of all shape functions of a scalar finite element in physical space at the poi...
int GetVDim() const
Returns vector dimension.
FiniteElementSpace * FESpace()
void PrintUsage(std::ostream &out) const
Print the usage message.
virtual DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element 'elem'.
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
double IntegrateBC(const GridFunction &sol, const Array< int > &bdr_marker, double alpha, double beta, double gamma, double &error)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
double Control[UMFPACK_CONTROL]
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
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 GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
Class for integration point with weight.
virtual void Print(std::ostream &os=mfem::out) const
Mesh * GenerateSerialMesh(int ref)
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
void PrintOptions(std::ostream &out) const
Print the options.
virtual int GetNVertices() const =0
Arbitrary order H1-conforming (continuous) finite elements.
double u(const Vector &xvec)
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Abstract data type element.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
const Element * GetBdrElement(int i) const
int GMRES(const Operator &A, Vector &x, const Vector &b, Solver &M, int &max_iter, int m, double &tol, double atol, int printit)
GMRES method. (tolerances are squared)
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
double sigma(const Vector &x)
bool Good() const
Return true if the command line options were parsed successfully.