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);
310 <<
"Verifying boundary conditions" << endl
311 <<
"=============================" << endl;
315 double err, avg =
IntegrateBC(u, dbc_bdr, 0.0, 1.0, dbc_val, err);
317 bool hom_dbc = (dbc_val == 0.0);
318 err /= hom_dbc ? 1.0 : fabs(dbc_val);
319 mfem::out <<
"Average of solution on Gamma_dbc:\t"
321 << (hom_dbc ?
"absolute" :
"relative")
322 <<
" error " << err << endl;
327 double err, avg =
IntegrateBC(u, nbc_bdr, 1.0, 0.0, nbc_val, err);
329 bool hom_nbc = (nbc_val == 0.0);
330 err /= hom_nbc ? 1.0 : fabs(nbc_val);
331 mfem::out <<
"Average of n.Grad(u) on Gamma_nbc:\t"
333 << (hom_nbc ?
"absolute" :
"relative")
334 <<
" error " << err << endl;
346 mfem::out <<
"Average of n.Grad(u) on Gamma_nbc0:\t"
348 << (hom_nbc ?
"absolute" :
"relative")
349 <<
" error " << err << endl;
354 double err, avg =
IntegrateBC(u, rbc_bdr, 1.0, rbc_a_val, rbc_b_val, err);
356 bool hom_rbc = (rbc_b_val == 0.0);
357 err /= hom_rbc ? 1.0 : fabs(rbc_b_val);
358 mfem::out <<
"Average of n.Grad(u)+a*u on Gamma_rbc:\t"
360 << (hom_rbc ?
"absolute" :
"relative")
361 <<
" error " << err << endl;
367 ofstream mesh_ofs(
"refined.mesh");
368 mesh_ofs.precision(8);
369 mesh->
Print(mesh_ofs);
370 ofstream sol_ofs(
"sol.gf");
371 sol_ofs.precision(8);
378 string title_str = h1 ?
"H1" :
"DG";
382 sol_sock.precision(8);
383 sol_sock <<
"solution\n" << *mesh << u
384 <<
"window_title '" << title_str <<
" Solution'"
385 <<
" keys 'mmc'" << flush;
395 void quad_trans(
double u,
double v,
double &x,
double &y,
bool log =
false)
399 double d = 4.0 * a * (M_SQRT2 - 2.0 *
a) * (1.0 - 2.0 * v);
401 double v0 = (1.0 + M_SQRT2) * (M_SQRT2 * a - 2.0 * v) *
402 ((4.0 - 3 * M_SQRT2) * a +
403 (8.0 * (M_SQRT2 - 1.0) * a - 2.0) * v) / d;
405 double r = 2.0 * ((M_SQRT2 - 1.0) * a * a * (1.0 - 4.0 *v) +
406 2.0 * (1.0 + M_SQRT2 *
407 (1.0 + 2.0 * (2.0 * a - M_SQRT2 - 1.0) *
a)) * v * v
410 double t = asin(v / r) * u / v;
414 << u <<
" " << v <<
" " << r <<
" " << v0 <<
" " << t
425 if (u[1] > 0.5 - tol || u[1] < -0.5 + tol)
430 if (u[0] > 1.0 - tol || u[0] < -1.0 + tol || fabs(u[0]) < tol)
438 if (u[1] > fabs(u[0] - 0.5))
444 if (u[1] < -fabs(u[0] - 0.5))
451 if (u[0] - 0.5 > fabs(u[1]))
457 if (u[0] - 0.5 < -fabs(u[1]))
467 if (u[1] > fabs(u[0] + 0.5))
473 if (u[1] < -fabs(u[0] + 0.5))
480 if (u[0] + 0.5 > fabs(u[1]))
486 if (u[0] + 0.5 < -fabs(u[1]))
499 Mesh * mesh =
new Mesh(2, 29, 16, 24, 2);
503 for (
int i=0; i<2; i++)
506 vi[0] = o + 0; vi[1] = o + 3; vi[2] = o + 4; vi[3] = o + 1;
509 vi[0] = o + 1; vi[1] = o + 4; vi[2] = o + 5; vi[3] = o + 2;
512 vi[0] = o + 5; vi[1] = o + 8; vi[2] = o + 9; vi[3] = o + 2;
515 vi[0] = o + 8; vi[1] = o + 12; vi[2] = o + 15; vi[3] = o + 9;
518 vi[0] = o + 11; vi[1] = o + 14; vi[2] = o + 15; vi[3] = o + 12;
521 vi[0] = o + 10; vi[1] = o + 13; vi[2] = o + 14; vi[3] = o + 11;
524 vi[0] = o + 6; vi[1] = o + 13; vi[2] = o + 10; vi[3] = o + 7;
527 vi[0] = o + 0; vi[1] = o + 6; vi[2] = o + 7; vi[3] = o + 3;
541 for (
int i=0; i<2; i++)
544 vi[0] = o + 7; vi[1] = o + 3; mesh->
AddBdrSegment(vi, 3 + i);
545 vi[0] = o + 10; vi[1] = o + 7; mesh->
AddBdrSegment(vi, 3 + i);
546 vi[0] = o + 11; vi[1] = o + 10; mesh->
AddBdrSegment(vi, 3 + i);
547 vi[0] = o + 12; vi[1] = o + 11; mesh->
AddBdrSegment(vi, 3 + i);
548 vi[0] = o + 8; vi[1] = o + 12; mesh->
AddBdrSegment(vi, 3 + i);
549 vi[0] = o + 5; vi[1] = o + 8; mesh->
AddBdrSegment(vi, 3 + i);
550 vi[0] = o + 4; vi[1] = o + 5; mesh->
AddBdrSegment(vi, 3 + i);
551 vi[0] = o + 3; vi[1] = o + 4; mesh->
AddBdrSegment(vi, 3 + i);
555 double a = a_ / M_SQRT2;
557 d[0] = -1.0; d[1] = -0.5; mesh->
AddVertex(d);
558 d[0] = -1.0; d[1] = 0.0; mesh->
AddVertex(d);
559 d[0] = -1.0; d[1] = 0.5; mesh->
AddVertex(d);
562 d[0] = -0.5 -
a; d[1] = 0.0; mesh->
AddVertex(d);
565 d[0] = -0.5; d[1] = -0.5; mesh->
AddVertex(d);
568 d[0] = -0.5; d[1] = 0.5; mesh->
AddVertex(d);
571 d[0] = -0.5 +
a; d[1] = 0.0; mesh->
AddVertex(d);
574 d[0] = 0.0; d[1] = -0.5; mesh->
AddVertex(d);
575 d[0] = 0.0; d[1] = 0.0; mesh->
AddVertex(d);
576 d[0] = 0.0; d[1] = 0.5; mesh->
AddVertex(d);
579 d[0] = 0.5 -
a; d[1] = 0.0; mesh->
AddVertex(d);
582 d[0] = 0.5; d[1] = -0.5; mesh->
AddVertex(d);
585 d[0] = 0.5; d[1] = 0.5; mesh->
AddVertex(d);
588 d[0] = 0.5 +
a; d[1] = 0.0; mesh->
AddVertex(d);
591 d[0] = 1.0; d[1] = -0.5; mesh->
AddVertex(d);
592 d[0] = 1.0; d[1] = 0.0; mesh->
AddVertex(d);
593 d[0] = 1.0; d[1] = 0.5; mesh->
AddVertex(d);
602 for (
int i = 0; i < v2v.Size() - 3; i++)
607 v2v[v2v.Size() - 3] = 0;
608 v2v[v2v.Size() - 2] = 1;
609 v2v[v2v.Size() - 1] = 2;
612 for (
int i = 0; i < mesh->
GetNE(); i++)
617 for (
int j = 0; j < nv; j++)
623 for (
int i = 0; i < mesh->
GetNBE(); i++)
628 for (
int j = 0; j < nv; j++)
638 for (
int l = 0; l < ref; l++)
649 double alpha,
double beta,
double gamma,
656 const bool a_is_zero = alpha == 0.0;
657 const bool b_is_zero = beta == 0.0;
660 MFEM_ASSERT(fes.
GetVDim() == 1,
"");
662 Vector shape, loc_dofs, w_nor;
665 for (
int i = 0; i < mesh.GetNBE(); i++)
667 if (bdr[mesh.GetBdrAttribute(i)-1] == 0) {
continue; }
670 if (FTr ==
nullptr) {
continue; }
673 MFEM_ASSERT(fe.
GetMapType() == FiniteElement::VALUE,
"");
674 const int int_order = 2*fe.
GetOrder() + 3;
702 val += alpha * dshape.
InnerProduct(w_nor, loc_dofs) / face_weight;
707 val += beta * (shape * loc_dofs);
711 nrm += ip.
weight * face_weight;
714 avg += val * ip.
weight * face_weight;
718 err += (val*val) * ip.
weight * face_weight;
723 if (std::abs(nrm) > 0.0)
731 err = (err >= 0.0) ? sqrt(err) : -sqrt(-err);
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.
virtual void Print(std::ostream &out=mfem::out) const
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...
int main(int argc, char *argv[])
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.
virtual void GetElementDofs(int i, Array< int > &dofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
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.
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
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.
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
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...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
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.
Mesh * GenerateSerialMesh(int ref)
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.
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)
double IntegrateBC(const GridFunction &sol, const Array< int > &bdr_marker, double alpha, double beta, double gamma, double &err)
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.