libMesh::UniformRefinementEstimator Class Reference
#include <uniform_refinement_estimator.h>

Public Types | |
| typedef std::map< std::pair < const System *, unsigned int > , ErrorVector * > | ErrorMap |
Public Member Functions | |
| UniformRefinementEstimator () | |
| ~UniformRefinementEstimator () | |
| virtual void | estimate_error (const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=NULL, bool estimate_parent_error=false) |
| 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 | |
| unsigned char | number_h_refinements |
| unsigned char | number_p_refinements |
| SystemNorm | error_norm |
Protected Member Functions | |
| virtual void | _estimate_error (const EquationSystems *equation_systems, const System *system, ErrorVector *error_per_cell, std::map< std::pair< const System *, unsigned int >, ErrorVector * > *errors_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) |
| void | reduce_error (std::vector< float > &error_per_cell) const |
Detailed Description
This class implements a ``brute force'' error estimator which integrates differences between the current solution and the solution on a uniformly refined (in h and/or p, for an arbitrary number of levels) grid.
Definition at line 44 of file uniform_refinement_estimator.h.
Member Typedef Documentation
typedef std::map<std::pair<const System*, unsigned int>, ErrorVector*> libMesh::ErrorEstimator::ErrorMap [inherited] |
When calculating many error vectors at once, we need a data structure to hold them all
Definition at line 107 of file error_estimator.h.
Constructor & Destructor Documentation
| libMesh::UniformRefinementEstimator::UniformRefinementEstimator | ( | ) | [inline] |
Constructor. Sets the most common default parameter values.
Definition at line 51 of file uniform_refinement_estimator.h.
References libMesh::ErrorEstimator::error_norm, and libMeshEnums::H1.
00051 : number_h_refinements(1), 00052 number_p_refinements(0) 00053 { error_norm = H1; }
| libMesh::UniformRefinementEstimator::~UniformRefinementEstimator | ( | ) | [inline] |
Member Function Documentation
| virtual void libMesh::UniformRefinementEstimator::_estimate_error | ( | const EquationSystems * | equation_systems, | |
| const System * | system, | |||
| ErrorVector * | error_per_cell, | |||
| std::map< std::pair< const System *, unsigned int >, ErrorVector * > * | errors_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 | |||
| ) | [protected, virtual] |
The code for estimate_error and both estimate_errors versions is very similar, so we use the same function for all three
Referenced by estimate_error(), and estimate_errors().
| void libMesh::UniformRefinementEstimator::estimate_error | ( | const System & | system, | |
| ErrorVector & | error_per_cell, | |||
| const NumericVector< Number > * | solution_vector = NULL, |
|||
| bool | estimate_parent_error = false | |||
| ) | [virtual] |
This function does uniform refinements and a solve to get an improved solution on each cell, then estimates the error by integrating differences between the coarse and fine solutions.
system.solve() must be called, and so should have no side effects.
Only the provided system is solved on the refined mesh; for problems decoupled into multiple systems, use of estimate_errors() should be more reliable.
The estimated error is output in the vector error_per_cell
Implements libMesh::ErrorEstimator.
Definition at line 49 of file uniform_refinement_estimator.C.
References _estimate_error().
00053 { 00054 START_LOG("estimate_error()", "UniformRefinementEstimator"); 00055 std::map<const System*, const NumericVector<Number>* > solution_vectors; 00056 solution_vectors[&_system] = solution_vector; 00057 this->_estimate_error (NULL, &_system, &error_per_cell, NULL, NULL, 00058 &solution_vectors, estimate_parent_error); 00059 STOP_LOG("estimate_error()", "UniformRefinementEstimator"); 00060 }
| void libMesh::UniformRefinementEstimator::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] |
Currently this function ignores the component_scale member variable, because it calculates each 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
Reimplemented from libMesh::ErrorEstimator.
Definition at line 75 of file uniform_refinement_estimator.C.
References _estimate_error().
00079 { 00080 START_LOG("estimate_errors()", "UniformRefinementEstimator"); 00081 this->_estimate_error (&_es, NULL, NULL, &errors_per_cell, NULL, 00082 solution_vectors, estimate_parent_error); 00083 STOP_LOG("estimate_errors()", "UniformRefinementEstimator"); 00084 }
| void libMesh::UniformRefinementEstimator::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] |
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 from libMesh::ErrorEstimator.
Definition at line 62 of file uniform_refinement_estimator.C.
References _estimate_error().
00067 { 00068 START_LOG("estimate_errors()", "UniformRefinementEstimator"); 00069 this->_estimate_error (&_es, NULL, &error_per_cell, NULL, 00070 &error_norms, solution_vectors, 00071 estimate_parent_error); 00072 STOP_LOG("estimate_errors()", "UniformRefinementEstimator"); 00073 }
| void libMesh::ErrorEstimator::reduce_error | ( | std::vector< float > & | error_per_cell | ) | const [protected, inherited] |
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 libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), and libMesh::JumpErrorEstimator::estimate_error().
Member Data Documentation
SystemNorm libMesh::ErrorEstimator::error_norm [inherited] |
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 139 of file error_estimator.h.
Referenced by libMesh::AdjointRefinementEstimator::AdjointRefinementEstimator(), libMesh::KellyErrorEstimator::boundary_side_integration(), libMesh::DiscontinuityMeasure::boundary_side_integration(), libMesh::DiscontinuityMeasure::DiscontinuityMeasure(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointResidualErrorEstimator::estimate_error(), libMesh::ErrorEstimator::estimate_errors(), libMesh::ExactErrorEstimator::ExactErrorEstimator(), libMesh::ExactErrorEstimator::find_squared_element_error(), libMesh::KellyErrorEstimator::internal_side_integration(), libMesh::LaplacianErrorEstimator::internal_side_integration(), libMesh::DiscontinuityMeasure::internal_side_integration(), libMesh::KellyErrorEstimator::KellyErrorEstimator(), libMesh::LaplacianErrorEstimator::LaplacianErrorEstimator(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::PatchRecoveryErrorEstimator(), and UniformRefinementEstimator().
How many h refinements to perform to get the fine grid
Definition at line 111 of file uniform_refinement_estimator.h.
How many p refinements to perform to get the fine grid
Definition at line 116 of file uniform_refinement_estimator.h.
The documentation for this class was generated from the following files:
Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:42 UTC
Hosted By: