libMesh::ErrorEstimator Class Referenceabstract

#include <error_estimator.h>

Inheritance diagram for libMesh::ErrorEstimator:

Public Types

typedef std::map< std::pair< const System *, unsigned int >, ErrorVector * > ErrorMap
 

Public Member Functions

 ErrorEstimator ()=default
 
 ErrorEstimator (const ErrorEstimator &)=default
 
 ErrorEstimator (ErrorEstimator &&)=default
 
ErrorEstimatoroperator= (const ErrorEstimator &)=default
 
ErrorEstimatoroperator= (ErrorEstimator &&)=default
 
virtual ~ErrorEstimator ()=default
 
virtual void estimate_error (const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false)=0
 
virtual void estimate_errors (const EquationSystems &equation_systems, ErrorVector &error_per_cell, const std::map< const System *, SystemNorm > &error_norms, const std::map< const System *, const NumericVector< Number > *> *solution_vectors=nullptr, bool estimate_parent_error=false)
 
virtual void estimate_errors (const EquationSystems &equation_systems, ErrorMap &errors_per_cell, const std::map< const System *, const NumericVector< Number > *> *solution_vectors=nullptr, bool estimate_parent_error=false)
 
virtual ErrorEstimatorType type () const =0
 

Public Attributes

SystemNorm error_norm
 

Protected Member Functions

void reduce_error (std::vector< ErrorVectorReal > &error_per_cell, const Parallel::Communicator &comm) const
 

Detailed Description

This class holds functions that will estimate the error in a finite element solution on a given mesh. These error estimates can be useful in their own right, or may be used to guide adaptive mesh refinement.

Note
The computed errors are stored as floats rather than doubles since the required precision is low.
Author
Benjamin S. Kirk
Date
2003

Definition at line 64 of file error_estimator.h.

Member Typedef Documentation

◆ ErrorMap

typedef std::map<std::pair<const System *, unsigned int>, ErrorVector *> libMesh::ErrorEstimator::ErrorMap

When calculating many error vectors at once, we need a data structure to hold them all

Definition at line 124 of file error_estimator.h.

Constructor & Destructor Documentation

◆ ErrorEstimator() [1/3]

libMesh::ErrorEstimator::ErrorEstimator ( )
default

Constructor. Empty. Derived classes should reset error_norm as appropriate.

◆ ErrorEstimator() [2/3]

libMesh::ErrorEstimator::ErrorEstimator ( const ErrorEstimator )
default

Copy/move ctor, copy/move assignment operator, and destructor are all explicitly defaulted for this simple class.

◆ ErrorEstimator() [3/3]

libMesh::ErrorEstimator::ErrorEstimator ( ErrorEstimator &&  )
default

◆ ~ErrorEstimator()

virtual libMesh::ErrorEstimator::~ErrorEstimator ( )
virtualdefault

Member Function Documentation

◆ estimate_error()

virtual void libMesh::ErrorEstimator::estimate_error ( const System system,
ErrorVector error_per_cell,
const NumericVector< Number > *  solution_vector = nullptr,
bool  estimate_parent_error = false 
)
pure virtual

This pure virtual function must be redefined in derived classes to compute the error for each active element and place it in the "error_per_cell" vector.

If solution_vector is not nullptr, the estimator will (if able) attempt to estimate an error in that field instead of in system.solution.

If estimate_parent_error is not false, the estimator will (if able) attempt to give a consistent estimate of errors in parent elements that would be generated by coarsening.

Implemented in libMesh::ExactErrorEstimator, libMesh::AdjointResidualErrorEstimator, libMesh::AdjointRefinementEstimator, libMesh::JumpErrorEstimator, libMesh::PatchRecoveryErrorEstimator, libMesh::UniformRefinementEstimator, and libMesh::WeightedPatchRecoveryErrorEstimator.

Referenced by estimate_errors().

◆ estimate_errors() [1/2]

void libMesh::ErrorEstimator::estimate_errors ( const EquationSystems equation_systems,
ErrorVector error_per_cell,
const std::map< const System *, SystemNorm > &  error_norms,
const std::map< const System *, const NumericVector< Number > *> *  solution_vectors = nullptr,
bool  estimate_parent_error = false 
)
virtual

This virtual function can be redefined in derived classes, but by default computes the sum of the error_per_cell for each system in the equation_systems.

Currently this function ignores the error_norm member variable, and uses the function argument error_norms instead.

This function is named estimate_errors instead of estimate_error because otherwise C++ can get confused.

Reimplemented in libMesh::UniformRefinementEstimator.

Definition at line 47 of file error_estimator.C.

References error_norm, estimate_error(), libMesh::EquationSystems::get_system(), libMesh::index_range(), and libMesh::EquationSystems::n_systems().

52 {
53  SystemNorm old_error_norm = this->error_norm;
54 
55  // Sum the error values from each system
56  for (unsigned int s = 0; s != equation_systems.n_systems(); ++s)
57  {
58  ErrorVector system_error_per_cell;
59  const System & sys = equation_systems.get_system(s);
60  if (error_norms.find(&sys) == error_norms.end())
61  this->error_norm = old_error_norm;
62  else
63  this->error_norm = error_norms.find(&sys)->second;
64 
65  const NumericVector<Number> * solution_vector = nullptr;
66  if (solution_vectors &&
67  solution_vectors->find(&sys) != solution_vectors->end())
68  solution_vector = solution_vectors->find(&sys)->second;
69 
70  this->estimate_error(sys, system_error_per_cell,
71  solution_vector, estimate_parent_error);
72 
73  if (s)
74  {
75  libmesh_assert_equal_to (error_per_cell.size(), system_error_per_cell.size());
76  for (auto i : index_range(error_per_cell))
77  error_per_cell[i] += system_error_per_cell[i];
78  }
79  else
80  error_per_cell = system_error_per_cell;
81  }
82 
83  // Restore our old state before returning
84  this->error_norm = old_error_norm;
85 }
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Definition: int_range.h:104
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false)=0

◆ estimate_errors() [2/2]

void libMesh::ErrorEstimator::estimate_errors ( const EquationSystems equation_systems,
ErrorMap errors_per_cell,
const std::map< const System *, const NumericVector< Number > *> *  solution_vectors = nullptr,
bool  estimate_parent_error = false 
)
virtual

This virtual function can be redefined in derived classes, but by default it calls estimate_error repeatedly to calculate the requested error vectors.

Currently this function ignores the error_norm.weight() values because it calculates each variable's error individually, unscaled.

The user selects which errors get computed by filling a map with error vectors: If errors_per_cell[&system][v] exists, it will be filled with the error values in variable v of system

FIXME: This is a default implementation - derived classes should reimplement it for efficiency.

Reimplemented in libMesh::UniformRefinementEstimator.

Definition at line 93 of file error_estimator.C.

References error_norm, estimate_error(), libMesh::EquationSystems::get_system(), libMesh::EquationSystems::n_systems(), n_vars, libMesh::System::n_vars(), and libMesh::SystemNorm::type().

97 {
98  SystemNorm old_error_norm = this->error_norm;
99 
100  // Find the requested error values from each system
101  for (unsigned int s = 0; s != equation_systems.n_systems(); ++s)
102  {
103  const System & sys = equation_systems.get_system(s);
104 
105  unsigned int n_vars = sys.n_vars();
106 
107  for (unsigned int v = 0; v != n_vars; ++v)
108  {
109  // Only fill in ErrorVectors the user asks for
110  if (errors_per_cell.find(std::make_pair(&sys, v)) ==
111  errors_per_cell.end())
112  continue;
113 
114  // Calculate error in only one variable
115  std::vector<Real> weights(n_vars, 0.0);
116  weights[v] = 1.0;
117  this->error_norm =
118  SystemNorm(std::vector<FEMNormType>(n_vars, old_error_norm.type(v)),
119  weights);
120 
121  const NumericVector<Number> * solution_vector = nullptr;
122  if (solution_vectors &&
123  solution_vectors->find(&sys) != solution_vectors->end())
124  solution_vector = solution_vectors->find(&sys)->second;
125 
126  this->estimate_error
127  (sys, *errors_per_cell[std::make_pair(&sys, v)],
128  solution_vector, estimate_parent_error);
129  }
130  }
131 
132  // Restore our old state before returning
133  this->error_norm = old_error_norm;
134 }
virtual void estimate_error(const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false)=0
const unsigned int n_vars
Definition: tecplot_io.C:69

◆ operator=() [1/2]

ErrorEstimator& libMesh::ErrorEstimator::operator= ( const ErrorEstimator )
default

◆ operator=() [2/2]

ErrorEstimator& libMesh::ErrorEstimator::operator= ( ErrorEstimator &&  )
default

◆ reduce_error()

void libMesh::ErrorEstimator::reduce_error ( std::vector< ErrorVectorReal > &  error_per_cell,
const Parallel::Communicator comm 
) const
protected

This method takes the local error contributions in error_per_cell from each processor and combines them to get the global error vector.

Definition at line 32 of file error_estimator.C.

References libMesh::Parallel::Communicator::sum().

Referenced by libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), and libMesh::ExactErrorEstimator::estimate_error().

34 {
35  // This function must be run on all processors at once
36  // parallel_object_only();
37 
38  // Each processor has now computed the error contributions
39  // for its local elements. We may need to sum the vector to
40  // recover the error for each element.
41 
42  comm.sum(error_per_cell);
43 }

◆ type()

Member Data Documentation

◆ error_norm

SystemNorm libMesh::ErrorEstimator::error_norm

When estimating the error in a single system, the error_norm is used to control the scaling and norm choice for each variable. Not all estimators will support all norm choices. The default scaling is for all variables to be weighted equally. The default norm choice depends on the error estimator.

Part of this functionality was supported via component_scale and sobolev_order in older libMesh versions, and a small part was supported via component_mask in even older versions. Hopefully the encapsulation here will allow us to avoid changing this API again.

Definition at line 161 of file error_estimator.h.

Referenced by libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::AdjointRefinementEstimator::AdjointRefinementEstimator(), libMesh::DiscontinuityMeasure::boundary_side_integration(), libMesh::KellyErrorEstimator::boundary_side_integration(), libMesh::DiscontinuityMeasure::DiscontinuityMeasure(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointResidualErrorEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), estimate_errors(), libMesh::ExactErrorEstimator::ExactErrorEstimator(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::LaplacianErrorEstimator::init_context(), libMesh::DiscontinuityMeasure::init_context(), libMesh::KellyErrorEstimator::init_context(), libMesh::LaplacianErrorEstimator::internal_side_integration(), libMesh::DiscontinuityMeasure::internal_side_integration(), libMesh::KellyErrorEstimator::internal_side_integration(), libMesh::KellyErrorEstimator::KellyErrorEstimator(), libMesh::LaplacianErrorEstimator::LaplacianErrorEstimator(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::PatchRecoveryErrorEstimator(), libMesh::JumpErrorEstimator::reinit_sides(), and libMesh::UniformRefinementEstimator::UniformRefinementEstimator().


The documentation for this class was generated from the following files: