MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
estimators.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#ifndef MFEM_ERROR_ESTIMATORS
13#define MFEM_ERROR_ESTIMATORS
14
15#include <functional>
16
17#include "../config/config.hpp"
18#include "../linalg/vector.hpp"
19#include "bilinearform.hpp"
20#ifdef MFEM_USE_MPI
21#include "pgridfunc.hpp"
22#endif
23
24namespace mfem
25{
26
27/** @brief Base class for all error estimators.
28 */
30{
31public:
33};
34
35
36/** @brief Base class for all element based error estimators.
37
38 At a minimum, an ErrorEstimator must be able compute one non-negative real
39 (double) number for each element in the Mesh.
40 */
42{
43public:
44 /// Return the total error from the last error estimate.
45 /** @note This method is optional for derived classes to override and the
46 base class implementation simply returns 0. */
47 virtual real_t GetTotalError() const { return 0.0; }
48
49 /// Get a Vector with all element errors.
50 virtual const Vector &GetLocalErrors() = 0;
51
52 /// Force recomputation of the estimates on the next call to GetLocalErrors.
53 virtual void Reset() = 0;
54
55 /// Destruct the error estimator
56 virtual ~ErrorEstimator() { }
57};
58
59
60/** @brief The AnisotropicErrorEstimator class is the base class for all error
61 estimators that compute one non-negative real (double) number and an
62 anisotropic flag for every element in the Mesh.
63 */
65{
66public:
67 /** @brief Get an Array<int> with anisotropic flags for all mesh elements.
68 @return An empty array when anisotropic estimates are not available or
69 enabled. */
70 virtual const Array<int> &GetAnisotropicFlags() = 0;
71};
72
73
74/** @brief The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu
75 error estimation procedure.
76
77 [1] Zienkiewicz, O.C. and Zhu, J.Z., The superconvergent patch recovery
78 and a posteriori error estimates. Part 1: The recovery technique.
79 Int. J. Num. Meth. Engng. 33, 1331-1364 (1992).
80
81 [2] Zienkiewicz, O.C. and Zhu, J.Z., The superconvergent patch recovery
82 and a posteriori error estimates. Part 2: Error estimates and adaptivity.
83 Int. J. Num. Meth. Engng. 33, 1365-1382 (1992).
84
85 The required BilinearFormIntegrator must implement the methods
86 ComputeElementFlux() and ComputeFluxEnergy().
87 */
89{
90protected:
96 int flux_averaging; // see SetFluxAveraging()
97
100
101 FiniteElementSpace *flux_space; /**< @brief Ownership based on own_flux_fes.
102 Its Update() method is called automatically by this class when needed. */
104 bool own_flux_fes; ///< Ownership flag for flux_space.
105
106 /// Check if the mesh of the solution was modified.
108 {
109 long mesh_sequence = solution.FESpace()->GetMesh()->GetSequence();
110 MFEM_ASSERT(mesh_sequence >= current_sequence, "");
111 return (mesh_sequence > current_sequence);
112 }
113
114 /// Compute the element error estimates.
115 void ComputeEstimates();
116
117public:
118 /** @brief Construct a new ZienkiewiczZhuEstimator object.
119 @param integ This BilinearFormIntegrator must implement the methods
120 ComputeElementFlux() and ComputeFluxEnergy().
121 @param sol The solution field whose error is to be estimated.
122 @param flux_fes The ZienkiewiczZhuEstimator assumes ownership of this
123 FiniteElementSpace and will call its Update() method when
124 needed.*/
126 FiniteElementSpace *flux_fes)
127 : current_sequence(-1),
128 total_error(),
129 anisotropic(false),
131 integ(integ),
132 solution(sol),
133 flux_space(flux_fes),
134 with_coeff(false),
135 own_flux_fes(true)
136 { }
137
138 /** @brief Construct a new ZienkiewiczZhuEstimator object.
139 @param integ This BilinearFormIntegrator must implement the methods
140 ComputeElementFlux() and ComputeFluxEnergy().
141 @param sol The solution field whose error is to be estimated.
142 @param flux_fes The ZienkiewiczZhuEstimator does NOT assume ownership of
143 this FiniteElementSpace; will call its Update() method
144 when needed. */
146 FiniteElementSpace &flux_fes)
147 : current_sequence(-1),
148 total_error(),
149 anisotropic(false),
151 integ(integ),
152 solution(sol),
153 flux_space(&flux_fes),
154 with_coeff(false),
155 own_flux_fes(false)
156 { }
157
158 /** @brief Consider the coefficient in BilinearFormIntegrator to calculate
159 the fluxes for the error estimator.*/
160 void SetWithCoeff(bool w_coeff = true) { with_coeff = w_coeff; }
161
162 /** @brief Enable/disable anisotropic estimates. To enable this option, the
163 BilinearFormIntegrator must support the 'd_energy' parameter in its
164 ComputeFluxEnergy() method. */
165 void SetAnisotropic(bool aniso = true) { anisotropic = aniso; }
166
167 /** @brief Set the way the flux is averaged (smoothed) across elements.
168
169 When @a fa is zero (default), averaging is performed across interfaces
170 between different mesh attributes. When @a fa is non-zero, the flux is
171 not averaged across interfaces between different mesh attributes. */
172 void SetFluxAveraging(int fa) { flux_averaging = fa; }
173
174 /// Return the total error from the last error estimate.
175 virtual real_t GetTotalError() const override { return total_error; }
176
177 /// Get a Vector with all element errors.
178 virtual const Vector &GetLocalErrors() override
179 {
180 if (MeshIsModified()) { ComputeEstimates(); }
181 return error_estimates;
182 }
183
184 /** @brief Get an Array<int> with anisotropic flags for all mesh elements.
185 Return an empty array when anisotropic estimates are not available or
186 enabled. */
187 virtual const Array<int> &GetAnisotropicFlags() override
188 {
189 if (MeshIsModified()) { ComputeEstimates(); }
190 return aniso_flags;
191 }
192
193 /// Reset the error estimator.
194 virtual void Reset() override { current_sequence = -1; }
195
196 /** @brief Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the
197 FiniteElementSpace, flux_space. */
199 {
200 if (own_flux_fes) { delete flux_space; }
201 }
202};
203
204
205/** @brief The LSZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu
206 error estimation procedure [1,2] using face-based patches [3].
207
208 [1] Zienkiewicz, O.C. and Zhu, J.Z., The superconvergent patch recovery
209 and a posteriori error estimates. Part 1: The recovery technique.
210 Int. J. Num. Meth. Engng. 33, 1331-1364 (1992).
211
212 [2] Zienkiewicz, O.C. and Zhu, J.Z., The superconvergent patch recovery
213 and a posteriori error estimates. Part 2: Error estimates and adaptivity.
214 Int. J. Num. Meth. Engng. 33, 1365-1382 (1992).
215
216 [3] Bartels, S. and Carstensen, C., Each averaging technique yields reliable
217 a posteriori error control in FEM on unstructured grids. Part II: Higher
218 order FEM. Math. Comp. 71(239), 971-994 (2002)
219
220 The required BilinearFormIntegrator must implement the method
221 ComputeElementFlux().
222
223 @note
224 - The present implementation ignores all single-element patches corresponding
225 to boundary faces. This is appropriate for Dirichlet boundaries, but
226 suboptimal for Neumann boundaries. Reference 3 shows that a constrained
227 least-squares problem, where the reconstructed flux is constrained by the
228 Neumann boundary data, is appropriate to handle this case.
229 NOTE THAT THIS CONSTRAINED LS PROBLEM IS NOT YET IMPLEMENTED, so it is
230 possible that the local error estimates for elements on a Neumann boundary
231 are suboptimal.
232 - The global polynomial basis used for the flux reconstruction is, by default,
233 aligned with the physical Cartesian axis. For patches with 2D elements, this
234 has been improved on so that the basis is aligned with the physical patch
235 orientation. Reorientation of the flux reconstruction basis is helpful to
236 maintain symmetry in the refinement pattern and could be extended to 3D.
237 - This estimator is ONLY implemented IN SERIAL.
238 - Anisotropic refinement is NOT YET SUPPORTED.
239
240 */
242{
243protected:
249
253
254 /// Check if the mesh of the solution was modified.
256 {
257 long mesh_sequence = solution.FESpace()->GetMesh()->GetSequence();
258 MFEM_ASSERT(mesh_sequence >= current_sequence, "");
259 return (mesh_sequence > current_sequence);
260 }
261
262 /// Compute the element error estimates.
263 void ComputeEstimates();
264
265public:
266 /** @brief Construct a new LSZienkiewiczZhuEstimator object.
267 @param integ This BilinearFormIntegrator must implement only the
268 method ComputeElementFlux().
269 @param sol The solution field whose error is to be estimated.
270 */
280
281 /** @brief Consider the coefficient in BilinearFormIntegrator to calculate
282 the fluxes for the error estimator.*/
283 void SetWithCoeff(bool w_coeff = true) { with_coeff = w_coeff; }
284
285 /** @brief Disable reconstructing the flux in patches spanning different
286 * subdomains. */
288
289 /** @brief Solve a Tichonov-regularized least-squares problem for the
290 * reconstructed fluxes. This is especially helpful for when not
291 * using tensor product elements, which typically require fewer
292 * integration points and, therefore, may lead to an
293 * ill-conditioned linear system. */
294 void SetTichonovRegularization(real_t tcoeff = 1.0e-8)
295 {
296 MFEM_VERIFY(tcoeff >= 0.0, "Tichonov coefficient cannot be negative");
297 tichonov_coeff = tcoeff;
298 }
299
300 /// Return the total error from the last error estimate.
301 virtual real_t GetTotalError() const override { return total_error; }
302
303 /// Get a Vector with all element errors.
304 virtual const Vector &GetLocalErrors() override
305 {
306 if (MeshIsModified()) { ComputeEstimates(); }
307 return error_estimates;
308 }
309
310 /// Reset the error estimator.
311 virtual void Reset() override { current_sequence = -1; }
312
314};
315
316
317#ifdef MFEM_USE_MPI
318
319/** @brief The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu
320 error estimation procedure where the flux averaging is replaced by a global
321 L2 projection (requiring a mass matrix solve).
322
323 The required BilinearFormIntegrator must implement the methods
324 ComputeElementFlux() and ComputeFluxEnergy().
325
326 Implemented for the parallel case only.
327 */
329{
330protected:
332 int local_norm_p; ///< Local L_p norm to use, default is 1.
335
338
339 ParFiniteElementSpace *flux_space; /**< @brief Ownership based on the flag
340 own_flux_fes. Its Update() method is called automatically by this class
341 when needed. */
342 ParFiniteElementSpace *smooth_flux_space; /**< @brief Ownership based on the
343 flag own_flux_fes. Its Update() method is called automatically by this
344 class when needed.*/
345 bool own_flux_fes; ///< Ownership flag for flux_space and smooth_flux_space.
346
347 /// Check if the mesh of the solution was modified.
349 {
350 long mesh_sequence = solution.FESpace()->GetMesh()->GetSequence();
351 MFEM_ASSERT(mesh_sequence >= current_sequence, "");
352 return (mesh_sequence > current_sequence);
353 }
354
355 /// Compute the element error estimates.
356 void ComputeEstimates();
357
358public:
359 /** @brief Construct a new L2ZienkiewiczZhuEstimator object.
360 @param integ This BilinearFormIntegrator must implement the methods
361 ComputeElementFlux() and ComputeFluxEnergy().
362 @param sol The solution field whose error is to be estimated.
363 @param flux_fes The L2ZienkiewiczZhuEstimator assumes ownership of this
364 FiniteElementSpace and will call its Update() method when
365 needed.
366 @param smooth_flux_fes
367 The L2ZienkiewiczZhuEstimator assumes ownership of this
368 FiniteElementSpace and will call its Update() method when
369 needed. */
371 ParGridFunction &sol,
372 ParFiniteElementSpace *flux_fes,
373 ParFiniteElementSpace *smooth_flux_fes)
374 : current_sequence(-1),
375 local_norm_p(1),
376 total_error(0.0),
377 integ(integ),
378 solution(sol),
379 flux_space(flux_fes),
380 smooth_flux_space(smooth_flux_fes),
381 own_flux_fes(true)
382 { }
383
384 /** @brief Construct a new L2ZienkiewiczZhuEstimator object.
385 @param integ This BilinearFormIntegrator must implement the methods
386 ComputeElementFlux() and ComputeFluxEnergy().
387 @param sol The solution field whose error is to be estimated.
388 @param flux_fes The L2ZienkiewiczZhuEstimator does NOT assume ownership
389 of this FiniteElementSpace; will call its Update() method
390 when needed.
391 @param smooth_flux_fes
392 The L2ZienkiewiczZhuEstimator does NOT assume ownership
393 of this FiniteElementSpace; will call its Update() method
394 when needed. */
396 ParGridFunction &sol,
397 ParFiniteElementSpace &flux_fes,
398 ParFiniteElementSpace &smooth_flux_fes)
399 : current_sequence(-1),
400 local_norm_p(1),
401 total_error(0.0),
402 integ(integ),
403 solution(sol),
404 flux_space(&flux_fes),
405 smooth_flux_space(&smooth_flux_fes),
406 own_flux_fes(false)
407 { }
408
409 /** @brief Set the exponent, p, of the Lp norm used for computing the local
410 element errors. Default value is 1. */
412
413 /// Return the total error from the last error estimate.
414 virtual real_t GetTotalError() const override { return total_error; }
415
416 /// Get a Vector with all element errors.
417 virtual const Vector &GetLocalErrors() override
418 {
419 if (MeshIsModified()) { ComputeEstimates(); }
420 return error_estimates;
421 }
422
423 /// Reset the error estimator.
424 virtual void Reset() override { current_sequence = -1; }
425
426 /** @brief Destroy a L2ZienkiewiczZhuEstimator object. Destroys, if owned,
427 the FiniteElementSpace, flux_space. */
429 {
430 if (own_flux_fes) { delete flux_space; delete smooth_flux_space; }
431 }
432};
433
434#endif // MFEM_USE_MPI
435
436
437/** @brief The LpErrorEstimator class compares the solution to a known
438 coefficient.
439
440 This class can be used, for example, to adapt a mesh to a non-trivial
441 initial condition in a time-dependent simulation. It can also be used to
442 force refinement in the neighborhood of small features before switching to a
443 more traditional error estimator.
444
445 The LpErrorEstimator supports either scalar or vector coefficients and works
446 both in serial and in parallel.
447*/
449{
450protected:
454
456
460
461 /// Check if the mesh of the solution was modified.
463 {
464 long mesh_sequence = sol->FESpace()->GetMesh()->GetSequence();
465 MFEM_ASSERT(mesh_sequence >= current_sequence, "");
466 return (mesh_sequence > current_sequence);
467 }
468
469 /// Compute the element error estimates.
470 void ComputeEstimates();
471
472public:
473 /** @brief Construct a new LpErrorEstimator object for a scalar field.
474 @param p Integer which selects which Lp norm to use.
475 @param sol The GridFunction representation of the scalar field.
476 Note: the coefficient must be set before use with the SetCoef method.
477 */
481
482 /** @brief Construct a new LpErrorEstimator object for a scalar field.
483 @param p Integer which selects which Lp norm to use.
484 @param coef The scalar Coefficient to compare to the solution.
485 @param sol The GridFunction representation of the scalar field.
486 */
490
491 /** @brief Construct a new LpErrorEstimator object for a vector field.
492 @param p Integer which selects which Lp norm to use.
493 @param coef The vector VectorCoefficient to compare to the solution.
494 @param sol The GridFunction representation of the vector field.
495 */
499
500 /** @brief Set the exponent, p, of the Lp norm used for computing the local
501 element errors. */
503
504 void SetCoef(Coefficient &A) { coef = &A; }
505 void SetCoef(VectorCoefficient &A) { vcoef = &A; }
506
507 /// Reset the error estimator.
508 virtual void Reset() override { current_sequence = -1; }
509
510 /// Get a Vector with all element errors.
511 virtual const Vector &GetLocalErrors() override
512 {
513 if (MeshIsModified()) { ComputeEstimates(); }
514 return error_estimates;
515 }
516
517 /// Destructor
518 virtual ~LpErrorEstimator() {}
519};
520
521
522/** @brief The KellyErrorEstimator class provides a fast error indication
523 strategy for smooth scalar parallel problems.
524
525 The Kelly error indicator is based on the following papers:
526
527 [1] Kelly, D. W., et al. "A posteriori error analysis and adaptive processes
528 in the finite element method: Part I—Error analysis." International journal
529 for numerical methods in engineering 19.11 (1983): 1593-1619.
530
531 [2] De SR Gago, J. P., et al. "A posteriori error analysis and adaptive
532 processes in the finite element method: Part II—Adaptive mesh refinement."
533 International journal for numerical methods in engineering 19.11 (1983):
534 1621-1656.
535
536 It can be roughly described by:
537 ||∇(u-uₕ)||ₑ ≅ √( C hₑ ∑ₖ (hₖ ∫ |J[∇uₕ]|²) dS )
538 where "e" denotes an element, ||⋅||ₑ the corresponding local norm and k the
539 corresponding faces. u is the analytic solution and uₕ the discretized
540 solution. hₖ and hₑ are factors dependent on the face and element geometry.
541 J is the jump function, i.e. the difference between the limits at each point
542 for each side of the face. A custom method to compute hₖ can be provided. It
543 is also possible to estimate the error only on a subspace by feeding this
544 class an attribute array describing the subspace.
545
546 @note This algorithm is appropriate only for problems with scalar diffusion
547 coefficients (e.g. Poisson problems), because it measures only the flux of
548 the gradient of the discrete solution. The current implementation also does
549 not include the volume term present in Equation 75 of Kelly et al [1].
550 Instead, it includes only the flux term recommended in Equation 82. The
551 current implementation also does not include the constant factor "C". It
552 furthermore assumes that the approximation error at the boundary is small
553 enough, as the implementation ignores boundary faces.
554*/
556{
557public:
558 /// Function type to compute the local coefficient hₑ of an element.
560 std::function<real_t(Mesh*, const int)>;
561 /** @brief Function type to compute the local coefficient hₖ of a face. The
562 third argument is true for shared faces and false for local faces. */
564 std::function<real_t(Mesh*, const int, const bool)>;
565
566private:
567 int current_sequence = -1;
568
569 Vector error_estimates;
570
571 real_t total_error = 0.0;
572
573 Array<int> attributes;
574
575 /** @brief A method to compute hₑ on per-element basis.
576
577 This method weights the error approximation on the element level.
578
579 Defaults to hₑ=1.0.
580 */
581 ElementCoefficientFunction compute_element_coefficient;
582
583 /** @brief A method to compute hₖ on per-face basis.
584
585 This method weights the error approximation on the face level. The
586 background here is that classical Kelly error estimator implementations
587 approximate the geometrical characteristic hₖ with the face diameter,
588 which should be also be a possibility in this implementation.
589
590 Defaults to hₖ=diameter/2p.
591 */
592 FaceCoefficientFunction compute_face_coefficient;
593
594 BilinearFormIntegrator* flux_integrator; ///< Not owned.
595 GridFunction* solution; ///< Not owned.
596
598 flux_space; /**< @brief Ownership based on own_flux_fes. */
599 bool own_flux_fespace; ///< Ownership flag for flux_space.
600
601#ifdef MFEM_USE_MPI
602 const bool isParallel;
603#endif
604
605 /// Check if the mesh of the solution was modified.
606 bool MeshIsModified()
607 {
608 long mesh_sequence = solution->FESpace()->GetMesh()->GetSequence();
609 MFEM_ASSERT(mesh_sequence >= current_sequence,
610 "improper mesh update sequence");
611 return (mesh_sequence > current_sequence);
612 }
613
614 /** @brief Compute the element error estimates.
615
616 Algorithm outline:
617 1. Compute flux field for each element
618 2. Add error contribution from local interior faces
619 3. Add error contribution from shared interior faces
620 4. Finalize by computing hₖ and scale errors.
621 */
622 void ComputeEstimates();
623
624public:
625 /** @brief Construct a new KellyErrorEstimator object for a scalar field.
626 @param di_ The bilinearform to compute the interface flux.
627 @param sol_ The solution field whose error is to be estimated.
628 @param flux_fes_ The finite element space for the interface flux.
629 @param attributes_ The attributes of the subdomain(s) for which the
630 error should be estimated. An empty array results in
631 estimating the error over the complete domain.
632 */
634 FiniteElementSpace& flux_fes_,
635 const Array<int> &attributes_ = Array<int>());
636
637 /** @brief Construct a new KellyErrorEstimator object for a scalar field.
638 @param di_ The bilinearform to compute the interface flux.
639 @param sol_ The solution field whose error is to be estimated.
640 @param flux_fes_ The finite element space for the interface flux.
641 @param attributes_ The attributes of the subdomain(s) for which the
642 error should be estimated. An empty array results in
643 estimating the error over the complete domain.
644 */
646 FiniteElementSpace* flux_fes_,
647 const Array<int> &attributes_ = Array<int>());
648
650
651 /// Get a Vector with all element errors.
652 const Vector& GetLocalErrors() override
653 {
654 if (MeshIsModified())
655 {
656 ComputeEstimates();
657 }
658 return error_estimates;
659 }
660
661 /// Reset the error estimator.
662 void Reset() override { current_sequence = -1; };
663
664 virtual real_t GetTotalError() const override { return total_error; }
665
666 /** @brief Change the method to compute hₑ on a per-element basis.
667 @param compute_element_coefficient_
668 A function taking a mesh and an element index to
669 compute the local hₑ for the element.
670 */
672 compute_element_coefficient_)
673 {
674 compute_element_coefficient = compute_element_coefficient_;
675 }
676
677 /** @brief Change the method to compute hₖ on a per-element basis.
678 @param compute_face_coefficient_
679 A function taking a mesh and a face index to
680 compute the local hₖ for the face.
681 */
684 compute_face_coefficient_)
685 {
686 compute_face_coefficient = compute_face_coefficient_;
687 }
688
689 /// Change the coefficients back to default as described above.
691};
692
693} // namespace mfem
694
695#endif // MFEM_ERROR_ESTIMATORS
Base class for all error estimators.
The AnisotropicErrorEstimator class is the base class for all error estimators that compute one non-n...
virtual const Array< int > & GetAnisotropicFlags()=0
Get an Array<int> with anisotropic flags for all mesh elements.
Abstract base class BilinearFormIntegrator.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Base class for all element based error estimators.
virtual real_t GetTotalError() const
Return the total error from the last error estimate.
virtual void Reset()=0
Force recomputation of the estimates on the next call to GetLocalErrors.
virtual ~ErrorEstimator()
Destruct the error estimator.
virtual const Vector & GetLocalErrors()=0
Get a Vector with all element errors.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
Mesh * GetMesh() const
Returns the mesh.
Definition fespace.hpp:559
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
FiniteElementSpace * FESpace()
Definition gridfunc.hpp:696
The KellyErrorEstimator class provides a fast error indication strategy for smooth scalar parallel pr...
void ResetCoefficientFunctions()
Change the coefficients back to default as described above.
void Reset() override
Reset the error estimator.
virtual real_t GetTotalError() const override
Return the total error from the last error estimate.
const Vector & GetLocalErrors() override
Get a Vector with all element errors.
std::function< real_t(Mesh *, const int)> ElementCoefficientFunction
Function type to compute the local coefficient hₑ of an element.
void SetElementCoefficientFunction(ElementCoefficientFunction compute_element_coefficient_)
Change the method to compute hₑ on a per-element basis.
KellyErrorEstimator(BilinearFormIntegrator &di_, GridFunction &sol_, FiniteElementSpace &flux_fes_, const Array< int > &attributes_=Array< int >())
Construct a new KellyErrorEstimator object for a scalar field.
std::function< real_t(Mesh *, const int, const bool)> FaceCoefficientFunction
Function type to compute the local coefficient hₖ of a face. The third argument is true for shared fa...
void SetFaceCoefficientFunction(FaceCoefficientFunction compute_face_coefficient_)
Change the method to compute hₖ on a per-element basis.
The L2ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure where t...
L2ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace *flux_fes, ParFiniteElementSpace *smooth_flux_fes)
Construct a new L2ZienkiewiczZhuEstimator object.
void ComputeEstimates()
Compute the element error estimates.
void SetLocalErrorNormP(int p)
Set the exponent, p, of the Lp norm used for computing the local element errors. Default value is 1.
virtual ~L2ZienkiewiczZhuEstimator()
Destroy a L2ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace,...
virtual real_t GetTotalError() const override
Return the total error from the last error estimate.
ParFiniteElementSpace * flux_space
Ownership based on the flag own_flux_fes. Its Update() method is called automatically by this class w...
ParFiniteElementSpace * smooth_flux_space
Ownership based on the flag own_flux_fes. Its Update() method is called automatically by this class w...
virtual void Reset() override
Reset the error estimator.
L2ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, ParGridFunction &sol, ParFiniteElementSpace &flux_fes, ParFiniteElementSpace &smooth_flux_fes)
Construct a new L2ZienkiewiczZhuEstimator object.
bool own_flux_fes
Ownership flag for flux_space and smooth_flux_space.
bool MeshIsModified()
Check if the mesh of the solution was modified.
int local_norm_p
Local L_p norm to use, default is 1.
virtual const Vector & GetLocalErrors() override
Get a Vector with all element errors.
BilinearFormIntegrator & integ
The LSZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure [1,...
virtual real_t GetTotalError() const override
Return the total error from the last error estimate.
virtual const Vector & GetLocalErrors() override
Get a Vector with all element errors.
virtual void Reset() override
Reset the error estimator.
LSZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, GridFunction &sol)
Construct a new LSZienkiewiczZhuEstimator object.
BilinearFormIntegrator & integ
void DisableReconstructionAcrossSubdomains()
Disable reconstructing the flux in patches spanning different subdomains.
bool MeshIsModified()
Check if the mesh of the solution was modified.
void ComputeEstimates()
Compute the element error estimates.
void SetTichonovRegularization(real_t tcoeff=1.0e-8)
Solve a Tichonov-regularized least-squares problem for the reconstructed fluxes. This is especially h...
void SetWithCoeff(bool w_coeff=true)
Consider the coefficient in BilinearFormIntegrator to calculate the fluxes for the error estimator.
The LpErrorEstimator class compares the solution to a known coefficient.
LpErrorEstimator(int p, Coefficient &coef, GridFunction &sol)
Construct a new LpErrorEstimator object for a scalar field.
LpErrorEstimator(int p, GridFunction &sol)
Construct a new LpErrorEstimator object for a scalar field.
void SetCoef(VectorCoefficient &A)
LpErrorEstimator(int p, VectorCoefficient &coef, GridFunction &sol)
Construct a new LpErrorEstimator object for a vector field.
void ComputeEstimates()
Compute the element error estimates.
void SetLocalErrorNormP(int p)
Set the exponent, p, of the Lp norm used for computing the local element errors.
bool MeshIsModified()
Check if the mesh of the solution was modified.
virtual const Vector & GetLocalErrors() override
Get a Vector with all element errors.
void SetCoef(Coefficient &A)
virtual ~LpErrorEstimator()
Destructor.
virtual void Reset() override
Reset the error estimator.
VectorCoefficient * vcoef
Mesh data type.
Definition mesh.hpp:56
long GetSequence() const
Definition mesh.hpp:2242
Abstract parallel finite element space.
Definition pfespace.hpp:29
Class for parallel grid function.
Definition pgridfunc.hpp:33
Base class for vector Coefficients that optionally depend on time and space.
Vector data type.
Definition vector.hpp:80
The ZienkiewiczZhuEstimator class implements the Zienkiewicz-Zhu error estimation procedure.
virtual const Vector & GetLocalErrors() override
Get a Vector with all element errors.
virtual ~ZienkiewiczZhuEstimator()
Destroy a ZienkiewiczZhuEstimator object. Destroys, if owned, the FiniteElementSpace,...
BilinearFormIntegrator & integ
bool MeshIsModified()
Check if the mesh of the solution was modified.
ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace *flux_fes)
Construct a new ZienkiewiczZhuEstimator object.
virtual real_t GetTotalError() const override
Return the total error from the last error estimate.
bool own_flux_fes
Ownership flag for flux_space.
void ComputeEstimates()
Compute the element error estimates.
void SetWithCoeff(bool w_coeff=true)
Consider the coefficient in BilinearFormIntegrator to calculate the fluxes for the error estimator.
ZienkiewiczZhuEstimator(BilinearFormIntegrator &integ, GridFunction &sol, FiniteElementSpace &flux_fes)
Construct a new ZienkiewiczZhuEstimator object.
FiniteElementSpace * flux_space
Ownership based on own_flux_fes. Its Update() method is called automatically by this class when neede...
void SetAnisotropic(bool aniso=true)
Enable/disable anisotropic estimates. To enable this option, the BilinearFormIntegrator must support ...
void SetFluxAveraging(int fa)
Set the way the flux is averaged (smoothed) across elements.
virtual const Array< int > & GetAnisotropicFlags() override
Get an Array<int> with anisotropic flags for all mesh elements. Return an empty array when anisotropi...
virtual void Reset() override
Reset the error estimator.
float real_t
Definition config.hpp:43
real_t p(const Vector &x, real_t t)