Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
ug::SideAndElemErrEstData< TDomain > Class Template Reference

Error estimator data class storing a number vector per side and per element. More...

#include <err_est_data.h>

+ Inheritance diagram for ug::SideAndElemErrEstData< TDomain >:

Classes

struct  GetQuadRules
 helper struct for getting quadrature rules by use of mpl::lists More...
 

Public Types

typedef Attachment< std::vector< number > > attachment_type
 attachment type
 
typedef TDomain domain_type
 domain type
 
typedef domain_traits< dim >::element_type elem_type
 
typedef domain_traits< dim >::side_type side_type
 type of the sides (face, edge) and the elems (volume, face)
 
typedef SideAndElemErrEstData< TDomain > this_type
 this class
 

Public Member Functions

MathVector< TDomain::dim > * all_side_global_ips (GridObject *elem, const MathVector< dim > vCornerCoords[])
 get all global side integration points
 
virtual void alloc_err_est_data (ConstSmartPtr< SurfaceView > spSV, const GridLevel &gl)
 virtual function to allocate data structures for the error estimator
 
MathVector< TDomain::dim > * elem_global_ips (GridObject *elem, const MathVector< dim > vCornerCoords[])
 get the global elem integration points for a specific roid
 
template<int refDim>
const MathVector< refDim > * elem_local_ips (const ReferenceObjectID roid)
 get the local elem integration points for a specific roid
 
size_t elem_order () const
 getting the elem integration order
 
size_t first_side_ips (const ReferenceObjectID roid, const size_t side)
 get number of first IP belonging to a specific side
 
virtual number get_elem_error_indicator (GridObject *elem, const MathVector< dim > vCornerCoords[])
 calculate L2 integrals
 
size_t num_all_side_ips (const ReferenceObjectID roid)
 get number of side IPs
 
size_t num_elem_ips (const ReferenceObjectID roid)
 get number of elem IPs
 
size_t num_side_ips (const ReferenceObjectID roid)
 get number of side IPs of a specific side type
 
size_t num_side_ips (const side_type *pSide)
 get number of side IPs of a specific side
 
numberoperator() (elem_type *pElem, size_t ip)
 get the data reference for a given elem and ip
 
numberoperator() (side_type *pSide, size_t ip)
 get the data reference for a given side and ip
 
virtual void release_err_est_data ()
 virtual function to release data structures of the error estimator
 
void set_type (int type)
 select L2/H1 Estimator
 
MathVector< TDomain::dim > * side_global_ips (GridObject *elem, const MathVector< dim > vCornerCoords[])
 get the global side integration points for a specific side roid
 
size_t side_ip_index (const ReferenceObjectID roid, const size_t side, const size_t ip)
 get index of specific side IP in sideIP array returned by side_local_ips
 
template<int refDim>
const MathVector< refDim > * side_local_ips (const ReferenceObjectID roid)
 get the local side integration points for a specific roid
 
size_t side_order () const
 getting the side integration order
 
 SideAndElemErrEstData (size_t _sideOrder, size_t _elemOrder, const char *subsets)
 constructors
 
 SideAndElemErrEstData (size_t _sideOrder, size_t _elemOrder, std::vector< std::string > subsets=std::vector< std::string >(0))
 
virtual void summarize_err_est_data (SmartPtr< TDomain > spDomain)
 virtual function called after the computation of the error estimator data in all the elements
 
ConstSmartPtr< SurfaceView > & surface_view ()
 get the surface view
 
virtual ~SideAndElemErrEstData ()
 virtual class destructor
 
- Public Member Functions inherited from ug::IErrEstData< TDomain >
bool consider_me () const
 virtual function granting get access to the m_consider member
 
 IErrEstData ()
 class constructor
 
number scaling_factor ()
 get scaling factor
 
void set_consider_me (bool b)
 whether or not this instance is to be considered by domainDisc
 
void set_scaling_factor (number scale)
 set scaling factor for final error calculation
 
virtual ~IErrEstData ()
 virtual class destructor
 

Static Public Attributes

static const int dim = TDomain::dim
 world dimension
 
static const int MAX_NUM_SIDES = 8
 maximal number of sides of any element
 
- Static Public Attributes inherited from ug::IErrEstData< TDomain >
static const int dim = TDomain::dim
 world dimension
 

Protected Member Functions

void init_quadrature ()
 initialization of quadrature (to be called during construction)
 

Private Types

enum  type { L2_ERROR_TYPE =0 , H1_ERROR_TYPE }
 set type More...
 

Private Attributes

size_t elemOrder
 
MultiGrid::AttachmentAccessor< elem_type, attachment_typem_aaElem
 vector of elem attachment accessors
 
MultiGrid::AttachmentAccessor< side_type, attachment_typem_aaSide
 vector of side attachment accessors
 
attachment_type m_aElem
 vector of attachments for elems
 
attachment_type m_aSide
 vector of attachments for sides
 
std::vector< MathVector< TDomain::dim > > m_elemGlobalIPcoords
 
GridLevel m_errEstGL
 Finest grid level.
 
std::vector< MathVector< TDomain::dim > > m_sideGlobalIPcoords
 storage for global elem and side IPs
 
std::vector< MathVector< TDomain::dim > > m_SideIPcoords [NUM_REFERENCE_OBJECTS]
 extra storage for local side IPs (elem IPs are contained in elem quad rules)
 
size_t m_sideIPsStartIndex [NUM_REFERENCE_OBJECTS][MAX_NUM_SIDES]
 the first index for IPs of a specific side in the sideIP series for a roid
 
std::vector< MathVector< TDomain::dim > > m_singleSideGlobalIPcoords
 
ConstSmartPtr< SurfaceViewm_spSV
 Grid for the attachment.
 
SubsetGroup m_ssg
 
int m_type
 
std::vector< std::string > m_vSs
 the subsets this error estimator will produce values for
 
QuadratureRule< dim > * quadRuleElem [NUM_REFERENCE_OBJECTS]
 
QuadratureRule< dim-1 > * quadRuleSide [NUM_REFERENCE_OBJECTS]
 storage for integration rules
 
size_t sideOrder
 order of side and elem function approximations for integrating
 

Detailed Description

template<typename TDomain>
class ug::SideAndElemErrEstData< TDomain >

Error estimator data class storing a number vector per side and per element.

This class represents an H1 error estimator. It can integrate expressions on elements and their sides with arbitrary order. A vector of values at defined integration points is attached to any element and side to that end.

RECOMMENDED (INTENDED) USAGE The data will typically consist of the values of certain functions at integration points (IP) on the sides and the element. A pointer to an object of this class can be handed to any element discretization involved in a discretization. They will access the attachments in their method compute_err_est_elem and add their respective parts of the function to be integrated for the error estimator. Exactly one of them (or maybe some other object) then has to do the actual integration using the given values at the IPs and add up side and element terms according to the error estimator formula used.

Template Parameters
TDomaindomain type

Member Typedef Documentation

◆ attachment_type

template<typename TDomain >
typedef Attachment<std::vector<number> > ug::SideAndElemErrEstData< TDomain >::attachment_type

attachment type

◆ domain_type

template<typename TDomain >
typedef TDomain ug::SideAndElemErrEstData< TDomain >::domain_type

domain type

◆ elem_type

template<typename TDomain >
typedef domain_traits<dim>::element_type ug::SideAndElemErrEstData< TDomain >::elem_type

◆ side_type

template<typename TDomain >
typedef domain_traits<dim>::side_type ug::SideAndElemErrEstData< TDomain >::side_type

type of the sides (face, edge) and the elems (volume, face)

◆ this_type

template<typename TDomain >
typedef SideAndElemErrEstData<TDomain> ug::SideAndElemErrEstData< TDomain >::this_type

this class

Member Enumeration Documentation

◆ type

template<typename TDomain >
enum ug::SideAndElemErrEstData::type
private

set type

Enumerator
L2_ERROR_TYPE 
H1_ERROR_TYPE 

Constructor & Destructor Documentation

◆ SideAndElemErrEstData() [1/2]

template<typename TDomain >
ug::SideAndElemErrEstData< TDomain >::SideAndElemErrEstData ( size_t  _sideOrder,
size_t  _elemOrder,
const char *  subsets 
)

◆ SideAndElemErrEstData() [2/2]

template<typename TDomain >
ug::SideAndElemErrEstData< TDomain >::SideAndElemErrEstData ( size_t  _sideOrder,
size_t  _elemOrder,
std::vector< std::string >  subsets = std::vector<std::string>(0) 
)

◆ ~SideAndElemErrEstData()

template<typename TDomain >
virtual ug::SideAndElemErrEstData< TDomain >::~SideAndElemErrEstData ( )
inlinevirtual

virtual class destructor

Member Function Documentation

◆ all_side_global_ips()

template<typename TDomain >
MathVector< TDomain::dim > * ug::SideAndElemErrEstData< TDomain >::all_side_global_ips ( GridObject elem,
const MathVector< dim vCornerCoords[] 
)

◆ alloc_err_est_data()

template<typename TDomain >
void ug::SideAndElemErrEstData< TDomain >::alloc_err_est_data ( ConstSmartPtr< SurfaceView spSV,
const GridLevel gl 
)
virtual

virtual function to allocate data structures for the error estimator

Allocates data structures for the error estimator.

Implements ug::IErrEstData< TDomain >.

References ug::GridLevel::SURFACE, ug::GridLevel::type(), and UG_COND_THROW.

◆ elem_global_ips()

template<typename TDomain >
MathVector< TDomain::dim > * ug::SideAndElemErrEstData< TDomain >::elem_global_ips ( GridObject elem,
const MathVector< dim vCornerCoords[] 
)

get the global elem integration points for a specific roid

References ug::DimReferenceMapping< TDim, TWorldDim >::local_to_global(), ug::GridObject::reference_object_id(), and UG_THROW.

◆ elem_local_ips()

template<typename TDomain >
template<int refDim>
const MathVector< refDim > * ug::SideAndElemErrEstData< TDomain >::elem_local_ips ( const ReferenceObjectID  roid)

get the local elem integration points for a specific roid

References UG_COND_THROW.

◆ elem_order()

template<typename TDomain >
size_t ug::SideAndElemErrEstData< TDomain >::elem_order ( ) const
inline

getting the elem integration order

References ug::SideAndElemErrEstData< TDomain >::elemOrder.

◆ first_side_ips()

template<typename TDomain >
size_t ug::SideAndElemErrEstData< TDomain >::first_side_ips ( const ReferenceObjectID  roid,
const size_t  side 
)

get number of first IP belonging to a specific side

◆ get_elem_error_indicator()

◆ init_quadrature()

◆ num_all_side_ips()

template<typename TDomain >
size_t ug::SideAndElemErrEstData< TDomain >::num_all_side_ips ( const ReferenceObjectID  roid)

get number of side IPs

◆ num_elem_ips()

template<typename TDomain >
size_t ug::SideAndElemErrEstData< TDomain >::num_elem_ips ( const ReferenceObjectID  roid)

get number of elem IPs

References UG_COND_THROW.

◆ num_side_ips() [1/2]

template<typename TDomain >
size_t ug::SideAndElemErrEstData< TDomain >::num_side_ips ( const ReferenceObjectID  roid)

get number of side IPs of a specific side type

References UG_COND_THROW.

◆ num_side_ips() [2/2]

template<typename TDomain >
size_t ug::SideAndElemErrEstData< TDomain >::num_side_ips ( const side_type pSide)

get number of side IPs of a specific side

Referenced by ug::FV1InnerBoundaryElemDisc< TImpl, TDomain >::compute_err_est_A_elem().

◆ operator()() [1/2]

template<typename TDomain >
number & ug::SideAndElemErrEstData< TDomain >::operator() ( elem_type pElem,
size_t  ip 
)

get the data reference for a given elem and ip

Parameters
pElempointer to the elem
ipintegration point id on the elem

References UG_CATCH_THROW, and UG_THROW.

◆ operator()() [2/2]

template<typename TDomain >
number & ug::SideAndElemErrEstData< TDomain >::operator() ( side_type pSide,
size_t  ip 
)

get the data reference for a given side and ip

Parameters
pSidepointer to the side
ipintegration point id on the side

References UG_CATCH_THROW, and UG_THROW.

◆ release_err_est_data()

template<typename TDomain >
void ug::SideAndElemErrEstData< TDomain >::release_err_est_data ( )
virtual

virtual function to release data structures of the error estimator

Releases data structures of the error estimator.

Implements ug::IErrEstData< TDomain >.

◆ set_type()

template<typename TDomain >
void ug::SideAndElemErrEstData< TDomain >::set_type ( int  type)
inline

◆ side_global_ips()

template<typename TDomain >
MathVector< TDomain::dim > * ug::SideAndElemErrEstData< TDomain >::side_global_ips ( GridObject elem,
const MathVector< dim vCornerCoords[] 
)

◆ side_ip_index()

template<typename TDomain >
size_t ug::SideAndElemErrEstData< TDomain >::side_ip_index ( const ReferenceObjectID  roid,
const size_t  side,
const size_t  ip 
)

get index of specific side IP in sideIP array returned by side_local_ips

◆ side_local_ips()

template<typename TDomain >
template<int refDim>
const MathVector< refDim > * ug::SideAndElemErrEstData< TDomain >::side_local_ips ( const ReferenceObjectID  roid)

get the local side integration points for a specific roid

References UG_THROW.

◆ side_order()

template<typename TDomain >
size_t ug::SideAndElemErrEstData< TDomain >::side_order ( ) const
inline

getting the side integration order

References ug::SideAndElemErrEstData< TDomain >::sideOrder.

◆ summarize_err_est_data()

template<typename TDomain >
void ug::SideAndElemErrEstData< TDomain >::summarize_err_est_data ( SmartPtr< TDomain >  spDomain)
virtual

virtual function called after the computation of the error estimator data in all the elements

Called after the computation of the error estimator data in all the elements.

Because of the multigrid hierarchy, the sides at the rim of a multigrid level only have one of the two flux terms (cis or trans). In order to calculate a jump, we need to add up the resp. attachments on parent and child for that side.

Implements ug::IErrEstData< TDomain >.

References ug::CollectCornerCoordinates(), pcl::InterfaceCommunicator< TLayout >::communicate(), ug::Grid::distributed_grid_manager(), pcl::InterfaceCommunicator< TLayout >::exchange_data(), ug::ComPol_CopyAttachment< TLayout, TAttachment >::extract_on_constrained_elems_only(), ug::DistributedGridManager::grid_layout_map(), ug::INT_H_MASTER, ug::INT_H_SLAVE, ug::INT_V_MASTER, ug::INT_V_SLAVE, PCL_RO_SUM, ug::SurfaceView::SHADOW_RIM, UG_THROW, ug::VecLengthSq(), and ug::VecSubtract().

◆ surface_view()

template<typename TDomain >
ConstSmartPtr< SurfaceView > & ug::SideAndElemErrEstData< TDomain >::surface_view ( )
inline

get the surface view

References ug::SideAndElemErrEstData< TDomain >::m_spSV.

Member Data Documentation

◆ dim

template<typename TDomain >
const int ug::SideAndElemErrEstData< TDomain >::dim = TDomain::dim
static

world dimension

◆ elemOrder

template<typename TDomain >
size_t ug::SideAndElemErrEstData< TDomain >::elemOrder
private

◆ m_aaElem

template<typename TDomain >
MultiGrid::AttachmentAccessor<elem_type, attachment_type > ug::SideAndElemErrEstData< TDomain >::m_aaElem
private

vector of elem attachment accessors

◆ m_aaSide

template<typename TDomain >
MultiGrid::AttachmentAccessor<side_type, attachment_type > ug::SideAndElemErrEstData< TDomain >::m_aaSide
private

vector of side attachment accessors

◆ m_aElem

template<typename TDomain >
attachment_type ug::SideAndElemErrEstData< TDomain >::m_aElem
private

vector of attachments for elems

◆ m_aSide

template<typename TDomain >
attachment_type ug::SideAndElemErrEstData< TDomain >::m_aSide
private

vector of attachments for sides

◆ m_elemGlobalIPcoords

template<typename TDomain >
std::vector<MathVector<TDomain::dim> > ug::SideAndElemErrEstData< TDomain >::m_elemGlobalIPcoords
private

◆ m_errEstGL

template<typename TDomain >
GridLevel ug::SideAndElemErrEstData< TDomain >::m_errEstGL
private

Finest grid level.

◆ m_sideGlobalIPcoords

template<typename TDomain >
std::vector<MathVector<TDomain::dim> > ug::SideAndElemErrEstData< TDomain >::m_sideGlobalIPcoords
private

storage for global elem and side IPs

◆ m_SideIPcoords

template<typename TDomain >
std::vector<MathVector<TDomain::dim> > ug::SideAndElemErrEstData< TDomain >::m_SideIPcoords[NUM_REFERENCE_OBJECTS]
private

extra storage for local side IPs (elem IPs are contained in elem quad rules)

◆ m_sideIPsStartIndex

template<typename TDomain >
size_t ug::SideAndElemErrEstData< TDomain >::m_sideIPsStartIndex[NUM_REFERENCE_OBJECTS][MAX_NUM_SIDES]
private

the first index for IPs of a specific side in the sideIP series for a roid

◆ m_singleSideGlobalIPcoords

template<typename TDomain >
std::vector<MathVector<TDomain::dim> > ug::SideAndElemErrEstData< TDomain >::m_singleSideGlobalIPcoords
private

◆ m_spSV

template<typename TDomain >
ConstSmartPtr<SurfaceView> ug::SideAndElemErrEstData< TDomain >::m_spSV
private

Grid for the attachment.

Referenced by ug::SideAndElemErrEstData< TDomain >::surface_view().

◆ m_ssg

template<typename TDomain >
SubsetGroup ug::SideAndElemErrEstData< TDomain >::m_ssg
private

◆ m_type

template<typename TDomain >
int ug::SideAndElemErrEstData< TDomain >::m_type
private

◆ m_vSs

template<typename TDomain >
std::vector<std::string> ug::SideAndElemErrEstData< TDomain >::m_vSs
private

◆ MAX_NUM_SIDES

template<typename TDomain >
const int ug::SideAndElemErrEstData< TDomain >::MAX_NUM_SIDES = 8
static

maximal number of sides of any element

◆ quadRuleElem

template<typename TDomain >
QuadratureRule<dim>* ug::SideAndElemErrEstData< TDomain >::quadRuleElem[NUM_REFERENCE_OBJECTS]
private

◆ quadRuleSide

template<typename TDomain >
QuadratureRule<dim-1>* ug::SideAndElemErrEstData< TDomain >::quadRuleSide[NUM_REFERENCE_OBJECTS]
private

storage for integration rules

◆ sideOrder

template<typename TDomain >
size_t ug::SideAndElemErrEstData< TDomain >::sideOrder
private

order of side and elem function approximations for integrating

Referenced by ug::SideAndElemErrEstData< TDomain >::side_order().


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