53 const double solver_tol = 1e-12;
54 const int solver_max_it = 200;
62 #endif // MFEM_USE_MPI
68 : attributes(attributes_)
69 , flux_integrator(&di_)
71 , flux_space(&flux_fespace_)
72 , own_flux_fespace(false)
84 : attributes(attributes_)
85 , flux_integrator(&di_)
87 , flux_space(flux_fespace_)
88 , own_flux_fespace(true)
106 compute_element_coefficient = [](
Mesh* mesh,
const int e)
111 compute_face_coefficient = [](
Mesh* mesh,
const int f,
112 const bool shared_face)
119 return dynamic_cast<ParMesh*
>(mesh)->GetSharedFaceTransformations(f);
121 #endif // MFEM_USE_MPI
127 double diameter = 0.0;
137 for (
int i = 0; i < nip; i++)
140 auto fip1 = vtx_intrule->IntPoint(i);
141 FT->Transform(fip1, p1);
143 for (
int j = 0; j < nip; j++)
145 auto fip2 = vtx_intrule->IntPoint(j);
146 FT->Transform(fip2, p2);
148 diameter = std::max<double>(diameter, p2.DistanceTo(p1));
151 return diameter/(2.0*order);
155 void KellyErrorEstimator::ComputeEstimates()
169 flux_space->
Update(
false);
171 auto xfes = solution->
FESpace();
172 MFEM_ASSERT(xfes->GetVDim() == 1,
173 "Estimation for vector-valued problems not implemented yet.");
176 this->error_estimates.
SetSize(xfes->GetNE());
177 this->error_estimates = 0.0;
190 if (attributes.
Size())
195 Array<int> xdofs, fdofs;
197 for (
int e = 0; e < xfes->GetNE(); e++)
199 auto attr = xfes->GetAttribute(e);
205 xfes->GetElementVDofs(e, xdofs);
208 ElementTransformation* Transf = xfes->GetElementTransformation(e);
210 *flux_space->
GetFE(e), el_f,
true);
217 for (
int f = 0;
f < mesh->GetNumFaces();
f++)
219 auto FT = mesh->GetFaceElementTransformations(
f);
221 auto &int_rule =
IntRules.
Get(FT->FaceGeom, 2 * xfes->GetFaceOrder(
f));
224 if (mesh->FaceIsInterior(
f))
226 int Inf1, Inf2, NCFace;
227 mesh->GetFaceInfos(
f, &Inf1, &Inf2, &NCFace);
234 bool isNCSlave = FT->Elem2No >= 0 && NCFace >= 0;
235 bool isConforming = FT->Elem2No >= 0 && NCFace == -1;
236 if ((FT->Elem1No < FT->Elem2No && isConforming) || isNCSlave)
238 if (attributes.
Size() &&
239 (attributes.
FindSorted(FT->Elem1->Attribute) == -1
240 || attributes.
FindSorted(FT->Elem2->Attribute) == -1))
250 for (
int i = 0; i < nip; i++)
253 auto &fip = int_rule.IntPoint(i);
255 FT->Loc1.Transform(fip, ip);
257 Vector val(flux_space->
GetVDim());
261 Vector normal(mesh->SpaceDimension());
262 FT->Face->SetIntPoint(&fip);
263 if (mesh->Dimension() == mesh->SpaceDimension())
269 Vector ref_normal(mesh->Dimension());
270 FT->Loc1.Transf.SetIntPoint(&fip);
271 CalcOrtho(FT->Loc1.Transf.Jacobian(), ref_normal);
272 auto &e1 = FT->GetElement1Transformation();
273 e1.AdjugateJacobian().MultTranspose(ref_normal, normal);
274 normal /= e1.Weight();
276 jumps(i) = val * normal * fip.weight * FT->Face->Weight();
281 for (
int i = 0; i < nip; i++)
284 auto &fip = int_rule.IntPoint(i);
286 FT->Loc2.Transform(fip, ip);
288 Vector val(flux_space->
GetVDim());
292 Vector normal(mesh->SpaceDimension());
293 FT->Face->SetIntPoint(&fip);
294 if (mesh->Dimension() == mesh->SpaceDimension())
300 Vector ref_normal(mesh->Dimension());
301 FT->Loc1.Transf.SetIntPoint(&fip);
302 CalcOrtho(FT->Loc1.Transf.Jacobian(), ref_normal);
303 auto &e1 = FT->GetElement1Transformation();
304 e1.AdjugateJacobian().MultTranspose(ref_normal, normal);
305 normal /= e1.Weight();
308 jumps(i) -= val * normal * fip.weight * FT->Face->Weight();
312 for (
int i = 0; i < nip; i++)
314 jumps(i) *= jumps(i);
316 auto h_k_face = compute_face_coefficient(mesh,
f,
false);
317 double jump_integral = h_k_face*jumps.Sum();
323 error_estimates(FT->Elem1No) += jump_integral;
324 error_estimates(FT->Elem2No) += jump_integral;
333 #endif // MFEM_USE_MPI
336 for (
int e = 0; e < xfes->GetNE(); e++)
338 auto factor = compute_element_coefficient(mesh, e);
340 error_estimates(e) = sqrt(factor * error_estimates(e));
343 total_error = error_estimates.
Norml2();
353 ParGridFunction *pflux =
dynamic_cast<ParGridFunction*
>(flux);
354 MFEM_VERIFY(pflux,
"flux is not a ParGridFunction pointer");
356 ParMesh *pmesh =
dynamic_cast<ParMesh*
>(mesh);
357 MFEM_VERIFY(pmesh,
"mesh is not a ParMesh pointer");
359 pflux->ExchangeFaceNbrData();
361 for (
int sf = 0; sf < pmesh->GetNSharedFaces(); sf++)
363 auto FT = pmesh->GetSharedFaceTransformations(sf,
true);
364 if (attributes.
Size() &&
365 (attributes.
FindSorted(FT->Elem1->Attribute) == -1
366 || attributes.
FindSorted(FT->Elem2->Attribute) == -1))
371 auto &int_rule =
IntRules.
Get(FT->FaceGeom, 2 * xfes->GetFaceOrder(0));
379 for (
int i = 0; i < nip; i++)
382 auto &fip = int_rule.IntPoint(i);
384 FT->Loc1.Transform(fip, ip);
386 Vector val(flux_space->
GetVDim());
389 Vector normal(mesh->SpaceDimension());
390 FT->Face->SetIntPoint(&fip);
391 if (mesh->Dimension() == mesh->SpaceDimension())
397 Vector ref_normal(mesh->Dimension());
398 FT->Loc1.Transf.SetIntPoint(&fip);
399 CalcOrtho(FT->Loc1.Transf.Jacobian(), ref_normal);
400 auto &e1 = FT->GetElement1Transformation();
401 e1.AdjugateJacobian().MultTranspose(ref_normal, normal);
402 normal /= e1.Weight();
405 jumps(i) = val * normal * fip.weight * FT->Face->Weight();
410 for (
int i = 0; i < nip; i++)
413 auto &fip = int_rule.IntPoint(i);
415 FT->Loc2.Transform(fip, ip);
417 Vector val(flux_space->
GetVDim());
421 Vector normal(mesh->SpaceDimension());
422 FT->Face->SetIntPoint(&fip);
423 if (mesh->Dimension() == mesh->SpaceDimension())
429 Vector ref_normal(mesh->Dimension());
430 CalcOrtho(FT->Loc1.Transf.Jacobian(), ref_normal);
431 auto &e1 = FT->GetElement1Transformation();
432 e1.AdjugateJacobian().MultTranspose(ref_normal, normal);
433 normal /= e1.Weight();
436 jumps(i) -= val * normal * fip.weight * FT->Face->Weight();
440 for (
int i = 0; i < nip; i++)
442 jumps(i) *= jumps(i);
444 auto h_k_face = compute_face_coefficient(mesh, sf,
true);
445 double jump_integral = h_k_face*jumps.Sum();
447 error_estimates(FT->Elem1No) += jump_integral;
455 for (
int e = 0; e < xfes->GetNE(); e++)
457 auto factor = compute_element_coefficient(mesh, e);
459 error_estimates(e) = sqrt(factor * error_estimates(e));
463 auto pfes =
dynamic_cast<ParFiniteElementSpace*
>(xfes);
464 MFEM_VERIFY(pfes,
"xfes is not a ParFiniteElementSpace pointer");
466 double process_local_error = pow(error_estimates.
Norml2(),2.0);
467 MPI_Allreduce(&process_local_error, &total_error, 1, MPI_DOUBLE,
468 MPI_SUM, pfes->GetComm());
469 total_error = sqrt(total_error);
470 #endif // MFEM_USE_MPI
475 MFEM_VERIFY(
coef != NULL ||
vcoef != NULL,
476 "LpErrorEstimator has no coefficient! Call SetCoef first.");
493 MPI_Allreduce(&process_local_error, &
total_error, 1, MPI_DOUBLE,
494 MPI_SUM, pfes->GetComm());
496 #endif // MFEM_USE_MPI
int GetNPoints() const
Returns the number of the points in the integration rule.
int Size() const
Return the logical size of the array.
Class for grid function - Vector with associated FE space.
virtual void Update(bool want_transform=true)
Reflect changes in the mesh: update number of DOFs, etc. Also, calculate GridFunction transformation ...
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void ComputeEstimates()
Compute the element error estimates.
BilinearFormIntegrator & integ
virtual void ComputeElementLpErrors(const double p, Coefficient &exsol, Vector &error, Coefficient *weight=NULL, const IntegrationRule *irs[]=NULL) const
virtual void GetVectorValue(int i, const IntegrationPoint &ip, Vector &val) const
ParFiniteElementSpace * flux_space
Ownership based on the flag own_flux_fes. Its Update() method is called automatically by this class w...
void SetSize(int s)
Resize the vector to size s.
VectorCoefficient * vcoef
double Norml2() const
Returns the l2 norm of the vector.
virtual void Update(bool want_transform=true)
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...
int GetNE() const
Returns number of elements.
Abstract parallel finite element space.
void CalcOrtho(const DenseMatrix &J, Vector &n)
virtual FaceElementTransformations * GetFaceElementTransformations(int FaceNo, int mask=31)
const IntegrationRule * GetVertices(int GeomType)
Return an IntegrationRule consisting of all vertices of the given Geometry::Type, GeomType...
void ComputeEstimates()
Compute the element error estimates.
BilinearFormIntegrator & integ
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
double f(const Vector &xvec)
Mesh * GetMesh() const
Returns the mesh.
void ResetCoefficientFunctions()
Change the coefficients back to default as described above.
int local_norm_p
Local L_p norm to use, default is 1.
void Sort()
Sorts the array in ascending order. This requires operator< to be defined for T.
int FindSorted(const T &el) const
Do bisection search for 'el' in a sorted array; return -1 if not found.
int GetVDim() const
Returns vector dimension.
FiniteElementSpace * FESpace()
KellyErrorEstimator(BilinearFormIntegrator &di_, GridFunction &sol_, FiniteElementSpace &flux_fes_, const Array< int > &attributes_=Array< int >())
Construct a new KellyErrorEstimator object for a scalar field.
BilinearFormIntegrator & integ
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
int SpaceDimension() const
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
double L2ZZErrorEstimator(BilinearFormIntegrator &flux_integrator, const ParGridFunction &x, ParFiniteElementSpace &smooth_flux_fes, ParFiniteElementSpace &flux_fes, Vector &errors, int norm_p, double solver_tol, int solver_max_it)
ParFiniteElementSpace * smooth_flux_space
Ownership based on the flag own_flux_fes. Its Update() method is called automatically by this class w...
ParGridFunction & solution
double ZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, GridFunction &flux, Vector &error_estimates, Array< int > *aniso_flags, int with_subdomains, bool with_coeff)
FiniteElementSpace * flux_space
Ownership based on own_flux_fes. Its Update() method is called automatically by this class when neede...
double LSZZErrorEstimator(BilinearFormIntegrator &blfi, GridFunction &u, Vector &error_estimates, bool subdomain_reconstruction, bool with_coeff, double tichonov_coeff)
A ``true'' ZZ error estimator that uses face-based patches for flux reconstruction.
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
void ComputeEstimates()
Compute the element error estimates.
void ComputeEstimates()
Compute the element error estimates.
Class for parallel grid function.
bool subdomain_reconstruction
Class for parallel meshes.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
double Sum() const
Return the sum of the vector entries.