ErrorEstimator Class Reference
#include <error_estimator.h>

Public Types | |
| typedef std::map< std::pair < const System *, unsigned int > , ErrorVector * > | ErrorMap |
Public Member Functions | |
| ErrorEstimator () | |
| virtual | ~ErrorEstimator () |
| virtual void | estimate_error (const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=NULL, 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=NULL, 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=NULL, bool estimate_parent_error=false) |
Public Attributes | |
| SystemNorm | error_norm |
Protected Member Functions | |
| void | reduce_error (std::vector< float > &error_per_cell) 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 that in general the computed errors are stored as floats rather than doubles since the required precision is low.
Definition at line 50 of file error_estimator.h.
Member Typedef Documentation
| typedef std::map<std::pair<const System*, unsigned int>, ErrorVector*> ErrorEstimator::ErrorMap |
When calculating many error vectors at once, we need a data structure to hold them all
Definition at line 105 of file error_estimator.h.
Constructor & Destructor Documentation
| ErrorEstimator::ErrorEstimator | ( | ) | [inline] |
Constructor. Empty. Derived classes should reset error_norm as appropriate.
Definition at line 58 of file error_estimator.h.
00058 : error_norm() {}
| virtual ErrorEstimator::~ErrorEstimator | ( | ) | [inline, virtual] |
Member Function Documentation
| virtual void ErrorEstimator::estimate_error | ( | const System & | system, | |
| ErrorVector & | error_per_cell, | |||
| const NumericVector< Number > * | solution_vector = NULL, |
|||
| 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 NULL, 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 AdjointResidualErrorEstimator, ExactErrorEstimator, JumpErrorEstimator, PatchRecoveryErrorEstimator, and UniformRefinementEstimator.
Referenced by estimate_errors().
| void ErrorEstimator::estimate_errors | ( | const EquationSystems & | equation_systems, | |
| ErrorMap & | errors_per_cell, | |||
| const std::map< const System *, const NumericVector< Number > * > * | solution_vectors = NULL, |
|||
| 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.
Definition at line 92 of file error_estimator.C.
References error_norm, estimate_error(), EquationSystems::get_system(), EquationSystems::n_systems(), System::n_vars(), n_vars, and SystemNorm::type().
00096 { 00097 SystemNorm old_error_norm = this->error_norm; 00098 00099 // Find the requested error values from each system 00100 for (unsigned int s = 0; s != equation_systems.n_systems(); ++s) 00101 { 00102 const System &sys = equation_systems.get_system(s); 00103 00104 unsigned int n_vars = sys.n_vars(); 00105 00106 for (unsigned int v = 0; v != n_vars; ++v) 00107 { 00108 // Only fill in ErrorVectors the user asks for 00109 if (errors_per_cell.find(std::make_pair(&sys, v)) == 00110 errors_per_cell.end()) 00111 continue; 00112 00113 // Calculate error in only one variable 00114 std::vector<Real> weights(n_vars, 0.0); 00115 weights[v] = 1.0; 00116 this->error_norm = 00117 SystemNorm(std::vector<FEMNormType>(n_vars, old_error_norm.type(0)), 00118 weights); 00119 00120 const NumericVector<Number>* solution_vector = NULL; 00121 if (solution_vectors && 00122 solution_vectors->find(&sys) != solution_vectors->end()) 00123 solution_vector = solution_vectors->find(&sys)->second; 00124 00125 this->estimate_error 00126 (sys, *errors_per_cell[std::make_pair(&sys, v)], 00127 solution_vector, estimate_parent_error); 00128 } 00129 } 00130 00131 // Restore our old state before returning 00132 this->error_norm = old_error_norm; 00133 }
| void 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 = NULL, |
|||
| 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 UniformRefinementEstimator.
Definition at line 46 of file error_estimator.C.
References error_norm, estimate_error(), EquationSystems::get_system(), and EquationSystems::n_systems().
00051 { 00052 SystemNorm old_error_norm = this->error_norm; 00053 00054 // Sum the error values from each system 00055 for (unsigned int s = 0; s != equation_systems.n_systems(); ++s) 00056 { 00057 ErrorVector system_error_per_cell; 00058 const System &sys = equation_systems.get_system(s); 00059 if (error_norms.find(&sys) == error_norms.end()) 00060 this->error_norm = old_error_norm; 00061 else 00062 this->error_norm = error_norms.find(&sys)->second; 00063 00064 const NumericVector<Number>* solution_vector = NULL; 00065 if (solution_vectors && 00066 solution_vectors->find(&sys) != solution_vectors->end()) 00067 solution_vector = solution_vectors->find(&sys)->second; 00068 00069 this->estimate_error(sys, system_error_per_cell, 00070 solution_vector, estimate_parent_error); 00071 00072 if (s) 00073 { 00074 libmesh_assert(error_per_cell.size() == system_error_per_cell.size()); 00075 for (unsigned int i=0; i != error_per_cell.size(); ++i) 00076 error_per_cell[i] += system_error_per_cell[i]; 00077 } 00078 else 00079 error_per_cell = system_error_per_cell; 00080 } 00081 00082 // Restore our old state before returning 00083 this->error_norm = old_error_norm; 00084 }
| void ErrorEstimator::reduce_error | ( | std::vector< float > & | error_per_cell | ) | 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.
Referenced by PatchRecoveryErrorEstimator::estimate_error(), and JumpErrorEstimator::estimate_error().
Member Data Documentation
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 137 of file error_estimator.h.
Referenced by KellyErrorEstimator::boundary_side_integration(), DiscontinuityMeasure::boundary_side_integration(), DiscontinuityMeasure::DiscontinuityMeasure(), JumpErrorEstimator::estimate_error(), estimate_errors(), ExactErrorEstimator::ExactErrorEstimator(), ExactErrorEstimator::find_squared_element_error(), KellyErrorEstimator::internal_side_integration(), LaplacianErrorEstimator::internal_side_integration(), DiscontinuityMeasure::internal_side_integration(), KellyErrorEstimator::KellyErrorEstimator(), LaplacianErrorEstimator::LaplacianErrorEstimator(), PatchRecoveryErrorEstimator::EstimateError::operator()(), PatchRecoveryErrorEstimator::PatchRecoveryErrorEstimator(), and UniformRefinementEstimator::UniformRefinementEstimator().
The documentation for this class was generated from the following files: