libMesh
Public Types | Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
libMesh::JumpErrorEstimator Class Referenceabstract

This abstract base class implements utility functions for error estimators which are based on integrated jumps between elements. More...

#include <jump_error_estimator.h>

Inheritance diagram for libMesh::JumpErrorEstimator:
[legend]

Public Types

typedef std::map< std::pair< const System *, unsigned int >, ErrorVector * > ErrorMap
 When calculating many error vectors at once, we need a data structure to hold them all. More...
 

Public Member Functions

 JumpErrorEstimator ()
 Constructor. More...
 
 JumpErrorEstimator (const JumpErrorEstimator &)=delete
 This class cannot be (default) copy constructed/assigned because it has unique_ptr members. More...
 
JumpErrorEstimatoroperator= (const JumpErrorEstimator &)=delete
 
 JumpErrorEstimator (JumpErrorEstimator &&)=default
 Defaulted move ctor, move assignment operator, and destructor. More...
 
JumpErrorEstimatoroperator= (JumpErrorEstimator &&)=default
 
virtual ~JumpErrorEstimator ()=default
 
virtual void estimate_error (const System &system, ErrorVector &error_per_cell, const NumericVector< Number > *solution_vector=nullptr, bool estimate_parent_error=false) override
 This function uses the derived class's jump error estimate formula to estimate the error on each cell. More...
 
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)
 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. More...
 
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)
 This virtual function can be redefined in derived classes, but by default it calls estimate_error repeatedly to calculate the requested error vectors. More...
 
virtual ErrorEstimatorType type () const =0
 

Public Attributes

bool scale_by_n_flux_faces
 This boolean flag allows you to scale the error indicator result for each element by the number of "flux faces" the element actually has. More...
 
SystemNorm 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. More...
 

Protected Member Functions

void reinit_sides ()
 A utility function to reinit the finite element data on elements sharing a side. More...
 
float coarse_n_flux_faces_increment ()
 A utility function to correctly increase n_flux_faces for the coarse element. More...
 
virtual void init_context (FEMContext &c)
 An initialization function, to give derived classes a chance to request specific data from the FE objects. More...
 
virtual void internal_side_integration ()=0
 The function, to be implemented by derived classes, which calculates an error term on an internal side. More...
 
virtual bool boundary_side_integration ()
 The function, to be implemented by derived classes, which calculates an error term on a boundary side. More...
 
void reduce_error (std::vector< ErrorVectorReal > &error_per_cell, const Parallel::Communicator &comm) const
 This method takes the local error contributions in error_per_cell from each processor and combines them to get the global error vector. More...
 

Protected Attributes

bool integrate_boundary_sides
 A boolean flag, by default false, to be set to true if integrations with boundary_side_integration() should be performed. More...
 
std::unique_ptr< FEMContextfine_context
 Context objects for integrating on the fine and coarse elements sharing a face. More...
 
std::unique_ptr< FEMContextcoarse_context
 
Real fine_error
 The fine and coarse error values to be set by each side_integration();. More...
 
Real coarse_error
 
unsigned int var
 The variable number currently being evaluated. More...
 

Detailed Description

This abstract base class implements utility functions for error estimators which are based on integrated jumps between elements.

Author
Roy H. Stogner
Date
2006

Definition at line 49 of file jump_error_estimator.h.

Member Typedef Documentation

◆ ErrorMap

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 124 of file error_estimator.h.

Constructor & Destructor Documentation

◆ JumpErrorEstimator() [1/3]

libMesh::JumpErrorEstimator::JumpErrorEstimator ( )

Constructor.

Definition at line 56 of file jump_error_estimator.h.

57  : ErrorEstimator(),
58  scale_by_n_flux_faces(false),
60  fine_context(),
62  fine_error(0),
63  coarse_error(0) {}
std::unique_ptr< FEMContext > fine_context
Context objects for integrating on the fine and coarse elements sharing a face.
std::unique_ptr< FEMContext > coarse_context
ErrorEstimator()=default
Constructor.
Real fine_error
The fine and coarse error values to be set by each side_integration();.
bool scale_by_n_flux_faces
This boolean flag allows you to scale the error indicator result for each element by the number of "f...
bool integrate_boundary_sides
A boolean flag, by default false, to be set to true if integrations with boundary_side_integration() ...

◆ JumpErrorEstimator() [2/3]

libMesh::JumpErrorEstimator::JumpErrorEstimator ( const JumpErrorEstimator )
delete

This class cannot be (default) copy constructed/assigned because it has unique_ptr members.

Explicitly deleting these functions is the best way to document this fact.

◆ JumpErrorEstimator() [3/3]

libMesh::JumpErrorEstimator::JumpErrorEstimator ( JumpErrorEstimator &&  )
default

Defaulted move ctor, move assignment operator, and destructor.

◆ ~JumpErrorEstimator()

virtual libMesh::JumpErrorEstimator::~JumpErrorEstimator ( )
virtualdefault

Member Function Documentation

◆ boundary_side_integration()

virtual bool libMesh::JumpErrorEstimator::boundary_side_integration ( )
protectedvirtual

The function, to be implemented by derived classes, which calculates an error term on a boundary side.

Returns
true if the flux bc function is in fact defined on the current side.

Reimplemented in libMesh::KellyErrorEstimator, and libMesh::DiscontinuityMeasure.

Definition at line 132 of file jump_error_estimator.h.

132 { return false; }

◆ coarse_n_flux_faces_increment()

float libMesh::JumpErrorEstimator::coarse_n_flux_faces_increment ( )
protected

A utility function to correctly increase n_flux_faces for the coarse element.

◆ estimate_error()

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

This function uses the derived class's jump error estimate formula to estimate the error on each cell.

The estimated error is output in the vector error_per_cell

Implements libMesh::ErrorEstimator.

◆ estimate_errors() [1/2]

virtual 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 
)
virtualinherited

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.

◆ estimate_errors() [2/2]

virtual 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 
)
virtualinherited

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

Reimplemented in libMesh::UniformRefinementEstimator.

◆ init_context()

virtual void libMesh::JumpErrorEstimator::init_context ( FEMContext c)
protectedvirtual

An initialization function, to give derived classes a chance to request specific data from the FE objects.

Reimplemented in libMesh::KellyErrorEstimator, libMesh::DiscontinuityMeasure, and libMesh::LaplacianErrorEstimator.

◆ internal_side_integration()

virtual void libMesh::JumpErrorEstimator::internal_side_integration ( )
protectedpure virtual

The function, to be implemented by derived classes, which calculates an error term on an internal side.

Implemented in libMesh::KellyErrorEstimator, libMesh::DiscontinuityMeasure, and libMesh::LaplacianErrorEstimator.

◆ operator=() [1/2]

JumpErrorEstimator& libMesh::JumpErrorEstimator::operator= ( const JumpErrorEstimator )
delete

◆ operator=() [2/2]

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

◆ reduce_error()

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

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

◆ reinit_sides()

void libMesh::JumpErrorEstimator::reinit_sides ( )
protected

A utility function to reinit the finite element data on elements sharing a side.

◆ type()

virtual ErrorEstimatorType libMesh::ErrorEstimator::type ( ) const
pure virtualinherited

Member Data Documentation

◆ coarse_context

std::unique_ptr<FEMContext> libMesh::JumpErrorEstimator::coarse_context
protected

Definition at line 144 of file jump_error_estimator.h.

◆ coarse_error

Real libMesh::JumpErrorEstimator::coarse_error
protected

Definition at line 149 of file jump_error_estimator.h.

◆ error_norm

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 161 of file error_estimator.h.

◆ fine_context

std::unique_ptr<FEMContext> libMesh::JumpErrorEstimator::fine_context
protected

Context objects for integrating on the fine and coarse elements sharing a face.

Definition at line 144 of file jump_error_estimator.h.

◆ fine_error

Real libMesh::JumpErrorEstimator::fine_error
protected

The fine and coarse error values to be set by each side_integration();.

Definition at line 149 of file jump_error_estimator.h.

◆ integrate_boundary_sides

bool libMesh::JumpErrorEstimator::integrate_boundary_sides
protected

A boolean flag, by default false, to be set to true if integrations with boundary_side_integration() should be performed.

Definition at line 138 of file jump_error_estimator.h.

◆ scale_by_n_flux_faces

bool libMesh::JumpErrorEstimator::scale_by_n_flux_faces

This boolean flag allows you to scale the error indicator result for each element by the number of "flux faces" the element actually has.

This tends to weight more evenly cells which are on the boundaries and thus have fewer contributions to their flux. The value is initialized to false, simply set it to true if you want to use the feature.

Definition at line 99 of file jump_error_estimator.h.

◆ var

unsigned int libMesh::JumpErrorEstimator::var
protected

The variable number currently being evaluated.

Definition at line 154 of file jump_error_estimator.h.


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