80int main(
int argc,
char *argv[])
88 int ser_ref_levels = 2;
89 int par_ref_levels = 1;
94 bool visualization =
true;
103 args.
AddOption(&h1,
"-h1",
"--continuous",
"-dg",
"--discontinuous",
104 "Select continuous \"H1\" or discontinuous \"DG\" basis.");
106 "Finite element order (polynomial degree) or -1 for"
107 " isoparametric space.");
109 "One of the two DG penalty parameters, typically +1/-1."
110 " See the documentation of class DGDiffusionIntegrator.");
112 "One of the two DG penalty parameters, should be positive."
113 " Negative values are replaced with (order+1)^2.");
114 args.
AddOption(&ser_ref_levels,
"-rs",
"--refine-serial",
115 "Number of times to refine the mesh uniformly in serial.");
116 args.
AddOption(&par_ref_levels,
"-rp",
"--refine-parallel",
117 "Number of times to refine the mesh uniformly in parallel.");
118 args.
AddOption(&mat_val,
"-mat",
"--material-value",
119 "Constant value for material coefficient "
120 "in the Laplace operator.");
121 args.
AddOption(&dbc_val,
"-dbc",
"--dirichlet-value",
122 "Constant value for Dirichlet Boundary Condition.");
123 args.
AddOption(&nbc_val,
"-nbc",
"--neumann-value",
124 "Constant value for Neumann Boundary Condition.");
125 args.
AddOption(&rbc_a_val,
"-rbc-a",
"--robin-a-value",
126 "Constant 'a' value for Robin Boundary Condition: "
127 "du/dn + a * u = b.");
128 args.
AddOption(&rbc_b_val,
"-rbc-b",
"--robin-b-value",
129 "Constant 'b' value for Robin Boundary Condition: "
130 "du/dn + a * u = b.");
132 "Radius of holes in the mesh.");
133 args.
AddOption(&visualization,
"-vis",
"--visualization",
"-no-vis",
134 "--no-visualization",
135 "Enable or disable GLVis visualization.");
142 if (
kappa < 0 && !h1)
144 kappa = (order+1)*(order+1);
150 mfem::out <<
"Hole radius too small, resetting to 0.01.\n";
155 mfem::out <<
"Hole radius too large, resetting to 0.49.\n";
166 ParMesh pmesh(MPI_COMM_WORLD, *mesh);
168 for (
int l = 0; l < par_ref_levels; l++)
181 mfem::out <<
"Number of finite element unknowns: " << size << endl;
191 nbc_bdr = 0; nbc_bdr[0] = 1;
192 rbc_bdr = 0; rbc_bdr[1] = 1;
193 dbc_bdr = 0; dbc_bdr[2] = 1;
258 u.ProjectBdrCoefficient(dbcCoef, dbc_bdr);
286 a.FormLinearSystem(ess_tdof_list,
u,
b, A, X, B);
290 if (h1 ||
sigma == -1.0)
320 <<
"Verifying boundary conditions" << endl
321 <<
"=============================" << endl;
327 bool hom_dbc = (dbc_val == 0.0);
328 error /= hom_dbc ? 1.0 : fabs(dbc_val);
329 mfem::out <<
"Average of solution on Gamma_dbc:\t"
331 << (hom_dbc ?
"absolute" :
"relative")
332 <<
" error " << error << endl;
339 bool hom_nbc = (nbc_val == 0.0);
340 error /= hom_nbc ? 1.0 : fabs(nbc_val);
341 mfem::out <<
"Average of n.Grad(u) on Gamma_nbc:\t"
343 << (hom_nbc ?
"absolute" :
"relative")
344 <<
" error " << error << endl;
356 mfem::out <<
"Average of n.Grad(u) on Gamma_nbc0:\t"
358 << (hom_nbc ?
"absolute" :
"relative")
359 <<
" error " << error << endl;
367 bool hom_rbc = (rbc_b_val == 0.0);
368 error /= hom_rbc ? 1.0 : fabs(rbc_b_val);
369 mfem::out <<
"Average of n.Grad(u)+a*u on Gamma_rbc:\t"
371 << (hom_rbc ?
"absolute" :
"relative")
372 <<
" error " << error << endl;
378 ostringstream mesh_name, sol_name;
379 mesh_name <<
"mesh." << setfill(
'0') << setw(6) <<
Mpi::WorldRank();
380 sol_name <<
"sol." << setfill(
'0') << setw(6) <<
Mpi::WorldRank();
382 ofstream mesh_ofs(mesh_name.str().c_str());
383 mesh_ofs.precision(8);
384 pmesh.
Print(mesh_ofs);
386 ofstream sol_ofs(sol_name.str().c_str());
387 sol_ofs.precision(8);
394 string title_str = h1 ?
"H1" :
"DG";
400 sol_sock.precision(8);
401 sol_sock <<
"solution\n" << pmesh <<
u
402 <<
"window_title '" << title_str <<
" Solution'"
403 <<
" keys 'mmc'" << flush;
416 real_t d = 4.0 *
a * (M_SQRT2 - 2.0 *
a) * (1.0 - 2.0 * v);
418 real_t v0 = (1.0 + M_SQRT2) * (M_SQRT2 *
a - 2.0 * v) *
419 ((4.0 - 3 * M_SQRT2) *
a +
420 (8.0 * (M_SQRT2 - 1.0) *
a - 2.0) * v) / d;
422 real_t r = 2.0 * ((M_SQRT2 - 1.0) *
a *
a * (1.0 - 4.0 *v) +
423 2.0 * (1.0 + M_SQRT2 *
424 (1.0 + 2.0 * (2.0 *
a - M_SQRT2 - 1.0) *
a)) * v * v
431 <<
u <<
" " << v <<
" " << r <<
" " << v0 <<
" " <<
t
442 if (
u[1] > 0.5 - tol ||
u[1] < -0.5 + tol)
447 if (
u[0] > 1.0 - tol ||
u[0] < -1.0 + tol || fabs(
u[0]) < tol)
455 if (
u[1] > fabs(
u[0] - 0.5))
461 if (
u[1] < -fabs(
u[0] - 0.5))
468 if (
u[0] - 0.5 > fabs(
u[1]))
474 if (
u[0] - 0.5 < -fabs(
u[1]))
484 if (
u[1] > fabs(
u[0] + 0.5))
490 if (
u[1] < -fabs(
u[0] + 0.5))
497 if (
u[0] + 0.5 > fabs(
u[1]))
503 if (
u[0] + 0.5 < -fabs(
u[1]))
516 Mesh * mesh =
new Mesh(2, 29, 16, 24, 2);
520 for (
int i=0; i<2; i++)
523 vi[0] = o + 0; vi[1] = o + 3; vi[2] = o + 4; vi[3] = o + 1;
526 vi[0] = o + 1; vi[1] = o + 4; vi[2] = o + 5; vi[3] = o + 2;
529 vi[0] = o + 5; vi[1] = o + 8; vi[2] = o + 9; vi[3] = o + 2;
532 vi[0] = o + 8; vi[1] = o + 12; vi[2] = o + 15; vi[3] = o + 9;
535 vi[0] = o + 11; vi[1] = o + 14; vi[2] = o + 15; vi[3] = o + 12;
538 vi[0] = o + 10; vi[1] = o + 13; vi[2] = o + 14; vi[3] = o + 11;
541 vi[0] = o + 6; vi[1] = o + 13; vi[2] = o + 10; vi[3] = o + 7;
544 vi[0] = o + 0; vi[1] = o + 6; vi[2] = o + 7; vi[3] = o + 3;
558 for (
int i=0; i<2; i++)
561 vi[0] = o + 7; vi[1] = o + 3; mesh->
AddBdrSegment(vi, 3 + i);
562 vi[0] = o + 10; vi[1] = o + 7; mesh->
AddBdrSegment(vi, 3 + i);
563 vi[0] = o + 11; vi[1] = o + 10; mesh->
AddBdrSegment(vi, 3 + i);
564 vi[0] = o + 12; vi[1] = o + 11; mesh->
AddBdrSegment(vi, 3 + i);
565 vi[0] = o + 8; vi[1] = o + 12; mesh->
AddBdrSegment(vi, 3 + i);
566 vi[0] = o + 5; vi[1] = o + 8; mesh->
AddBdrSegment(vi, 3 + i);
567 vi[0] = o + 4; vi[1] = o + 5; mesh->
AddBdrSegment(vi, 3 + i);
568 vi[0] = o + 3; vi[1] = o + 4; mesh->
AddBdrSegment(vi, 3 + i);
574 d[0] = -1.0; d[1] = -0.5; mesh->
AddVertex(d);
575 d[0] = -1.0; d[1] = 0.0; mesh->
AddVertex(d);
576 d[0] = -1.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] = 0.0; d[1] = -0.5; mesh->
AddVertex(d);
592 d[0] = 0.0; d[1] = 0.0; mesh->
AddVertex(d);
593 d[0] = 0.0; d[1] = 0.5; mesh->
AddVertex(d);
596 d[0] = 0.5 -
a; d[1] = 0.0; mesh->
AddVertex(d);
599 d[0] = 0.5; d[1] = -0.5; mesh->
AddVertex(d);
602 d[0] = 0.5; d[1] = 0.5; mesh->
AddVertex(d);
605 d[0] = 0.5 +
a; d[1] = 0.0; mesh->
AddVertex(d);
608 d[0] = 1.0; d[1] = -0.5; mesh->
AddVertex(d);
609 d[0] = 1.0; d[1] = 0.0; mesh->
AddVertex(d);
610 d[0] = 1.0; d[1] = 0.5; mesh->
AddVertex(d);
619 for (
int i = 0; i < v2v.
Size() - 3; i++)
624 v2v[v2v.
Size() - 3] = 0;
625 v2v[v2v.
Size() - 2] = 1;
626 v2v[v2v.
Size() - 1] = 2;
629 for (
int i = 0; i < mesh->
GetNE(); i++)
634 for (
int j = 0; j < nv; j++)
640 for (
int i = 0; i < mesh->
GetNBE(); i++)
645 for (
int j = 0; j < nv; j++)
655 for (
int l = 0; l < ref; l++)
670 real_t &nrm = loc_vals[0];
671 real_t &avg = loc_vals[1];
672 real_t &error = loc_vals[2];
678 const bool a_is_zero =
alpha == 0.0;
679 const bool b_is_zero =
beta == 0.0;
682 MFEM_ASSERT(fes.
GetVDim() == 1,
"");
684 Vector shape, loc_dofs, w_nor;
687 for (
int i = 0; i < mesh.
GetNBE(); i++)
692 if (FTr ==
nullptr) {
continue; }
696 const int int_order = 2*fe.
GetOrder() + 3;
729 val +=
beta * (shape * loc_dofs);
733 nrm += ip.
weight * face_weight;
736 avg += val * ip.
weight * face_weight;
740 error += (val*val) * ip.
weight * face_weight;
748 real_t glb_nrm = glb_vals[0];
749 real_t glb_avg = glb_vals[1];
750 glb_err = glb_vals[2];
753 if (std::abs(glb_nrm) > 0.0)
761 glb_err = (glb_err >= 0.0) ? sqrt(glb_err) : -sqrt(-glb_err);
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
int Size() const
Return the logical size of the array.
Class for boundary integration .
A coefficient that is constant across space and time.
Data type dense matrix using column-major storage.
real_t InnerProduct(const real_t *x, const real_t *y) const
Compute y^t A x.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Abstract data type element.
virtual void GetVertices(Array< int > &v) const =0
Get the indices defining the vertices.
virtual int GetNVertices() const =0
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
int GetVDim() const
Returns vector dimension.
Abstract class for all finite elements.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
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 GetMapType() const
Returns the FiniteElement::MapType of the element describing how reference functions are mapped to ph...
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...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
void SetKDim(int dim)
Set the number of iteration to perform between restarts, default is 50.
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the GMRES method.
Arbitrary order H1-conforming (continuous) finite elements.
The BoomerAMG solver in hypre.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
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)
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(),...
Class for integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
int GetNPoints() const
Returns the number of the points in the integration rule.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
void SetRelTol(real_t rtol)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetMaxIter(int max_it)
void SetAbsTol(real_t atol)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
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.
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
const Element * GetElement(int i) const
Return pointer to the i'th element object.
void FinalizeTopology(bool generate_bdr=true)
Finalize the construction of the secondary topology (connectivity) data of a Mesh.
int AddVertex(real_t x, real_t y=0.0, real_t z=0.0)
int GetNE() const
Returns number of elements.
int Dimension() const
Dimension of the reference space used within the elements.
const Element * GetBdrElement(int i) const
Return pointer to the i'th boundary element object.
int AddBdrSegment(int v1, int v2, int attr=1)
FaceElementTransformations * GetBdrFaceTransformations(int BdrElemNo)
Builds the transformation defining the given boundary face.
int GetNV() const
Returns number of vertices. Vertices are only at the corners of elements, where you would expect them...
int GetNBE() const
Returns number of boundary elements.
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.
void Transform(void(*f)(const Vector &, Vector &))
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
void RemoveInternalBoundaries()
void RemoveUnusedVertices()
Remove unused vertices and rebuild mesh connectivity.
static bool Root()
Return true if the rank in MPI_COMM_WORLD is zero.
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).
Pointer to an Operator of a specified type.
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear sys...
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.
void Disable()
Disable output.
Abstract parallel finite element space.
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalTrueVSize() const
void GetElementDofs(int i, Array< int > &dofs, DofTransformation &doftrans) const override
The same as GetElementDofs(), but with a user-allocated DofTransformation object. doftrans must be al...
ParMesh * GetParMesh() const
const FiniteElement * GetFE(int i) const override
Class for parallel grid function.
ParFiniteElementSpace * ParFESpace() const
Class for parallel meshes.
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Scalar coefficient defined as the product of two scalar coefficients or a scalar and a scalar coeffic...
real_t Norml2() const
Returns the l2 norm of the vector.
void SetSize(int s)
Resize the vector to size s.
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
real_t sigma(const Vector &x)
void n4Vec(const Vector &x, Vector &n)
real_t IntegrateBC(const ParGridFunction &sol, const Array< int > &bdr_marker, real_t alpha, real_t beta, real_t gamma, real_t &error)
void trans(const Vector &u, Vector &x)
void quad_trans(real_t u, real_t v, real_t &x, real_t &y, bool log=false)
Mesh * GenerateSerialMesh(int ref)
void CalcOrtho(const DenseMatrix &J, Vector &n)
real_t u(const Vector &xvec)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
L2_FECollection DG_FECollection
Declare an alternative name for L2_FECollection = DG_FECollection.
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Helper struct to convert a C++ type to an MPI type.