ug4
|
describes the ansatz spaces on a domain More...
#include <approximation_space.h>
Public Types | |
typedef MultiGrid | grid_type |
Type of Subset Handler. | |
typedef MGSubsetHandler | subset_handler_type |
Type of Subset Handler. | |
Public Member Functions | |
void | add (const char *name, const char *type) |
adds function using string to indicate finite element type | |
void | add (const char *name, const char *type, const char *subsets) |
adds function using string to indicate finite element type | |
void | add (const char *name, const char *type, int order) |
adds function using string to indicate finite element type | |
void | add (const char *name, const char *type, int order, const char *subsets) |
adds function using string to indicate finite element type | |
void | add (const std::vector< std::string > &vName, const char *type) |
adds function using string to indicate finite element type | |
void | add (const std::vector< std::string > &vName, const char *type, const std::vector< std::string > &vSubsets) |
adds function using string to indicate finite element type | |
void | add (const std::vector< std::string > &vName, const char *type, int order) |
adds function using string to indicate finite element type | |
void | add (const std::vector< std::string > &vName, const char *type, int order, const std::vector< std::string > &vSubsets) |
adds function using string to indicate finite element type | |
void | add (const std::vector< std::string > &vName, LFEID id) |
add single solutions of LocalShapeFunctionSetID to the entire domain | |
void | add (const std::vector< std::string > &vName, LFEID id, const std::vector< std::string > &vSubset) |
add single solutions of LocalShapeFunctionSetID to selected subsets | |
void | clear () |
clears functions | |
std::vector< SmartPtr< DoFDistribution > > | dof_distributions () const |
returns all currently created dof distributions | |
bool | grouped () const |
returns if dofs are grouped | |
IApproximationSpace (SmartPtr< subset_handler_type > spMGSH, SmartPtr< grid_type > spMG) | |
Constructor. | |
IApproximationSpace (SmartPtr< subset_handler_type > spMGSH, SmartPtr< grid_type >, const AlgebraType &algebraType) | |
Constructor setting the grouping flag. | |
void | init_levels () |
initializes all level dof distributions | |
void | init_surfaces () |
initializes all surface dof distributions | |
void | init_top_surface () |
initializes all top surface dof distributions | |
bool | might_contain_ghosts () const |
returns if ghosts might be present on any level | |
bool | might_contain_ghosts (int lvl) const |
returns if ghosts might be present on a level | |
size_t | num_levels () const |
returns the number of level | |
void | print_layout_statistic () const |
prints statistic on layouts | |
void | print_statistic () const |
prints statistic about DoF Distribution | |
void | print_statistic (std::string flags) const |
prints statistic about DoF Distribution | |
const RevisionCounter & | revision () const |
returns the current revision | |
ConstSmartPtr< MGSubsetHandler > | subset_handler () const |
get underlying subset handler | |
ConstSmartPtr< SurfaceView > | surface_view () const |
returns the approximation space | |
~IApproximationSpace () | |
Destructor. | |
SmartPtr< DoFDistribution > | dof_distribution (const GridLevel &gl, bool bCreate=true) |
SmartPtr< DoFDistribution > | dd (const GridLevel &gl, bool bCreate=true) |
ConstSmartPtr< DoFDistribution > | dof_distribution (const GridLevel &gl, bool bCreate=true) const |
ConstSmartPtr< DoFDistribution > | dd (const GridLevel &gl, bool bCreate=true) const |
ConstSmartPtr< DoFDistributionInfo > | ddinfo () const |
ConstSmartPtr< DoFDistributionInfo > | dof_distribution_info () const |
![]() | |
SubsetGroup | all_subsets_grp_except_for (const char *names) const |
returns subset group of all the subsets except for (by name) | |
int | dim (size_t fct) const |
returns the dimension in which solution lives | |
int | dim_subset (int si) const |
returns dimension of subset | |
ConstSmartPtr< DoFDistributionInfo > | dof_distribution_info () const |
returns underlying info | |
DoFDistributionInfoProvider () | |
constructor | |
DoFDistributionInfoProvider (ConstSmartPtr< DoFDistributionInfo > spDDI) | |
constructor | |
FunctionGroup | fct_grp_by_name (const char *names) const |
returns a function group to a string of functions | |
size_t | fct_id_by_name (const char *name) const |
returns fct id by name | |
ConstSmartPtr< FunctionPattern > | function_pattern () const |
returns the function pattern | |
bool | is_def_everywhere (size_t fct) const |
returns true if the discrete function nr_fct is defined everywhere | |
bool | is_def_in_subset (size_t fct, int si) const |
returns if a function is defined on a subset | |
size_t | max_dofs (const GridBaseObjectId gbo) const |
returns the maximum number of dofs on a grid base object type | |
size_t | max_dofs (const GridBaseObjectId gbo, const int si) const |
returns the maximum number of dofs on a grid base object on a subset | |
size_t | max_dofs (const int dim) const |
returns the maximum number of dofs on grid objects in a dimension | |
size_t | max_dofs (const int dim, const int si) const |
returns the maximum number of dofs in a dimension on a subset | |
size_t | max_dofs (const ReferenceObjectID roid) const |
returns the maximum number of dofs on reference object type | |
size_t | max_fct_dofs (const size_t fct, const GridBaseObjectId gbo) const |
returns the maximal number of dofs on a base object type for a function component | |
size_t | max_fct_dofs (const size_t fct, const GridBaseObjectId gbo, const int si) const |
returns the maximum number of dofs on a grid base object on a subset for a function component | |
size_t | max_fct_dofs (const size_t fct, const int dim, const int si) const |
returns the maximum number of dofs in a dimension on a subset for a function component | |
size_t | max_fct_dofs (const size_t fct, const ReferenceObjectID roid) const |
returns the number of dofs on a reference object for a function component | |
size_t | max_fct_dofs (const size_t fct, int dim) const |
returns the maximal number of dofs on a dimension for a function component | |
std::string | name (size_t fct) const |
returns the name of the discrete function nr_fct | |
std::vector< std::string > | names () const |
returns the names of the discrete functions | |
size_t | num_dofs (const ReferenceObjectID roid, const int si) const |
returns the number of dofs on a Reference Object on a subset | |
size_t | num_fct () const |
number of discrete functions on subset si | |
size_t | num_fct (int si) const |
number of discrete functions on subset si | |
size_t | num_fct_dofs (const size_t fct, const ReferenceObjectID roid, const int si) const |
returns the number of dofs on a Reference Object on a subset for a function component | |
int | num_subsets () const |
returns number of subsets | |
size_t | offset (const ReferenceObjectID roid, const int si, const size_t fct) const |
returns the offset for reference element, subset and function | |
void | print_local_dof_statistic () const |
prints statistic on local dof distribution | |
void | print_local_dof_statistic (int verboseLev) const |
prints informations | |
void | set_dof_distribution_info (ConstSmartPtr< DoFDistributionInfo > spDDI) |
sets the dd info | |
SubsetGroup | subset_grp_by_name (const char *names) const |
returns subset group by name | |
ConstSmartPtr< ISubsetHandler > | subset_handler () const |
returns the subset handler | |
int | subset_id_by_name (const char *name) const |
returns the subset id | |
std::string | subset_name (int si) const |
returns subset name | |
const LFEID & | local_finite_element_id (size_t fct) const |
const LFEID & | lfeid (size_t fct) const |
Protected Member Functions | |
void | create_dof_distribution (const GridLevel &gl) |
creates a dof distribution | |
void | dof_distribution_info_required () |
create dof distribution info | |
void | grid_changed_callback (const GridMessage_Adaption &msg) |
void | grid_distribution_callback (const GridMessage_Distribution &msg) |
called during parallel redistribution | |
void | init (SmartPtr< subset_handler_type > spMGSH, SmartPtr< grid_type > spMG, const AlgebraType &algebraType) |
initializing | |
void | register_at_adaption_msg_hub () |
registers at message hub for grid adaption | |
void | reinit () |
reinits all data after grid adaption | |
void | surface_view_required () |
creates surface SurfaceView if needed | |
Protected Attributes | |
AlgebraType | m_algebraType |
suitable algebra type for the index distribution pattern | |
bool | m_bAdaptionIsActive |
bool | m_bGrouped |
flag if DoFs should be grouped | |
RevisionCounter | m_RevCnt |
revision counter | |
SmartPtr< DoFDistributionInfo > | m_spDoFDistributionInfo |
DofDistributionInfo. | |
MessageHub::SPCallbackId | m_spGridAdaptionCallbackID |
message hub id | |
MessageHub::SPCallbackId | m_spGridDistributionCallbackID |
SmartPtr< MultiGrid > | m_spMG |
multigrid, where elements are stored | |
SmartPtr< MGSubsetHandler > | m_spMGSH |
subsethandler, where elements are stored | |
SmartPtr< SurfaceView > | m_spSurfaceView |
Surface View. | |
std::vector< SmartPtr< DoFDistribution > > | m_vDD |
MG Level DoF Distribution. | |
SmartPtr< DoFIndexStorage > | m_spDoFIndexStrgForLevelNoGhost |
SmartPtr< DoFIndexStorage > | m_spDoFIndexStrgForLevelWithGhost |
![]() | |
ConstSmartPtr< DoFDistributionInfo > | m_spDDI |
Function Pattern. | |
describes the ansatz spaces on a domain
This class provides grid function spaces on a domain.
The Domain defines a partition of the Grid/Multigrid in terms of subsets. The user can add discrete functions on this subsets or unions of them.
Once finalized, this function pattern is fixed. Internally DoF indices are created. Using this Approximation Space the user can create GridFunctions of the following types:
Type of Subset Handler.
Type of Subset Handler.
ug::IApproximationSpace::IApproximationSpace | ( | SmartPtr< subset_handler_type > | spMGSH, |
SmartPtr< grid_type > | spMG | ||
) |
Constructor.
References ug::DefaultAlgebra::get(), and init().
ug::IApproximationSpace::IApproximationSpace | ( | SmartPtr< subset_handler_type > | spMGSH, |
SmartPtr< grid_type > | spMG, | ||
const AlgebraType & | algebraType | ||
) |
Constructor setting the grouping flag.
References init().
ug::IApproximationSpace::~IApproximationSpace | ( | ) |
Destructor.
References m_spSurfaceView.
void ug::IApproximationSpace::add | ( | const char * | name, |
const char * | type | ||
) |
adds function using string to indicate finite element type
References add(), name, and ug::TokenizeTrimString().
void ug::IApproximationSpace::add | ( | const char * | name, |
const char * | type, | ||
const char * | subsets | ||
) |
adds function using string to indicate finite element type
[in] | name | name(s) of single solution (comma separated) |
[in] | type | type of local finite element space |
[in] | subsets | Subsets separated by ',' |
References add(), name, and ug::TokenizeTrimString().
void ug::IApproximationSpace::add | ( | const char * | name, |
const char * | type, | ||
int | order | ||
) |
adds function using string to indicate finite element type
References add(), name, and ug::TokenizeTrimString().
void ug::IApproximationSpace::add | ( | const char * | name, |
const char * | type, | ||
int | order, | ||
const char * | subsets | ||
) |
adds function using string to indicate finite element type
[in] | name | name(s) of single solution (comma separated) |
[in] | type | type of local finite element space |
[in] | order | order of local finite element space |
[in] | subsets | Subsets separated by ',' |
References add(), name, and ug::TokenizeTrimString().
void ug::IApproximationSpace::add | ( | const std::vector< std::string > & | vName, |
const char * | type | ||
) |
adds function using string to indicate finite element type
References add(), ug::ConvertStringToLFEID(), ug::DoFDistributionInfoProvider::dim(), ug::DIM_SUBSET_EMPTY_GRID, ug::DimensionOfSubsets(), m_spMGSH, and UG_THROW.
void ug::IApproximationSpace::add | ( | const std::vector< std::string > & | vName, |
const char * | type, | ||
const std::vector< std::string > & | vSubsets | ||
) |
adds function using string to indicate finite element type
References add(), ug::ConvertStringToLFEID(), ug::DoFDistributionInfoProvider::dim(), ug::DIM_SUBSET_EMPTY_GRID, ug::SubsetGroup::get_highest_subset_dimension(), m_spMGSH, and UG_THROW.
void ug::IApproximationSpace::add | ( | const std::vector< std::string > & | vName, |
const char * | type, | ||
int | order | ||
) |
adds function using string to indicate finite element type
References add(), ug::ConvertStringToLFEID(), ug::DoFDistributionInfoProvider::dim(), ug::DIM_SUBSET_EMPTY_GRID, ug::DimensionOfSubsets(), m_spMGSH, and UG_THROW.
void ug::IApproximationSpace::add | ( | const std::vector< std::string > & | vName, |
const char * | type, | ||
int | order, | ||
const std::vector< std::string > & | vSubsets | ||
) |
adds function using string to indicate finite element type
References add(), ug::ConvertStringToLFEID(), ug::DoFDistributionInfoProvider::dim(), ug::DIM_SUBSET_EMPTY_GRID, ug::SubsetGroup::get_highest_subset_dimension(), m_spMGSH, and UG_THROW.
|
inline |
|
inline |
add single solutions of LocalShapeFunctionSetID to selected subsets
[in] | name | name(s) of single solution (comma separated) |
[in] | id | Shape Function set id |
[in] | subsets | Subsets separated by ',' |
References m_spDoFDistributionInfo.
|
inline |
clears functions
References m_spDoFDistributionInfo.
|
protected |
creates a dof distribution
References dof_distribution_info_required(), ug::GridLevel::ghosts(), ug::GridLevel::is_level(), m_bGrouped, m_spDoFDistributionInfo, m_spDoFIndexStrgForLevelNoGhost, m_spDoFIndexStrgForLevelWithGhost, m_spMG, m_spMGSH, m_spSurfaceView, m_vDD, ug::SortDD(), and surface_view_required().
Referenced by dof_distribution().
SmartPtr< DoFDistribution > ug::IApproximationSpace::dd | ( | const GridLevel & | gl, |
bool | bCreate = true |
||
) |
returns dof distribution for a grid level
References dof_distribution().
ConstSmartPtr< DoFDistribution > ug::IApproximationSpace::dd | ( | const GridLevel & | gl, |
bool | bCreate = true |
||
) | const |
returns dof distribution for a grid level
References dof_distribution().
|
inline |
returns dof distribution info
References m_spDoFDistributionInfo.
Referenced by dof_distribution_info().
SmartPtr< DoFDistribution > ug::IApproximationSpace::dof_distribution | ( | const GridLevel & | gl, |
bool | bCreate = true |
||
) |
returns dof distribution for a grid level
References create_dof_distribution(), dof_distribution(), m_vDD, and UG_THROW.
Referenced by ug::ApproximationSpace< TGridFunction::domain_type >::ApproximationSpace(), dd(), dd(), dof_distribution(), dof_distribution(), init_levels(), and init_top_surface().
ConstSmartPtr< DoFDistribution > ug::IApproximationSpace::dof_distribution | ( | const GridLevel & | gl, |
bool | bCreate = true |
||
) | const |
returns dof distribution for a grid level
References dof_distribution().
|
inline |
returns dof distribution info
References ddinfo().
|
protected |
create dof distribution info
References ug::AlgebraType::blocksize(), m_algebraType, m_bGrouped, ug::DoFDistributionInfoProvider::m_spDDI, m_spDoFDistributionInfo, ug::NUM_REFERENCE_OBJECTS, UG_THROW, and ug::AlgebraType::VariableBlockSize.
Referenced by create_dof_distribution().
std::vector< SmartPtr< DoFDistribution > > ug::IApproximationSpace::dof_distributions | ( | ) | const |
returns all currently created dof distributions
References m_vDD.
Referenced by ug::OrderCuthillMcKee(), ug::OrderDownwind(), and ug::OrderLex().
|
protected |
this callback is called by the message hub, when a grid change has been performed. It will call all necessary actions in order to keep the grid correct for computations.
References ug::GridMessage_Adaption::adaption_begins(), ug::GridMessage_Adaption::adaption_ends(), m_bAdaptionIsActive, m_spMG, m_spSurfaceView, pcl::ProcRank(), reinit(), ug::SaveParallelGridLayout(), ug::SaveSurfaceViewTransformed(), UG_LOG, and UG_THROW.
Referenced by register_at_adaption_msg_hub().
|
protected |
called during parallel redistribution
References ug::GMDT_DISTRIBUTION_STARTS, ug::GMDT_DISTRIBUTION_STOPS, m_spMG, m_spSurfaceView, ug::GridMessage_Distribution::msg(), pcl::ProcRank(), PROFILE_FUNC, reinit(), ug::SaveParallelGridLayout(), ug::SaveSurfaceViewTransformed(), and UG_LOG.
Referenced by register_at_adaption_msg_hub().
|
inline |
returns if dofs are grouped
References m_bGrouped.
|
protected |
void ug::IApproximationSpace::init_levels | ( | ) |
initializes all level dof distributions
References dof_distribution(), ug::GridLevel::LEVEL, num_levels(), and PROFILE_FUNC.
void ug::IApproximationSpace::init_surfaces | ( | ) |
initializes all surface dof distributions
void ug::IApproximationSpace::init_top_surface | ( | ) |
initializes all top surface dof distributions
References dof_distribution(), PROFILE_FUNC, ug::GridLevel::SURFACE, and ug::GridLevel::TOP.
Referenced by ug::ApproximationSpace< TGridFunction::domain_type >::ApproximationSpace().
bool ug::IApproximationSpace::might_contain_ghosts | ( | ) | const |
returns if ghosts might be present on any level
References might_contain_ghosts(), and num_levels().
Referenced by might_contain_ghosts().
bool ug::IApproximationSpace::might_contain_ghosts | ( | int | lvl | ) | const |
returns if ghosts might be present on a level
References ug::EDGE, ug::FACE, ug::DistributedGridManager::grid_layout_map(), m_spMG, ug::DoFDistributionInfoProvider::max_dofs(), num_levels(), UG_THROW, ug::VERTEX, and ug::VOLUME.
|
inline |
returns the number of level
References m_spMGSH.
Referenced by ug::ApproximationSpace< TGridFunction::domain_type >::ApproximationSpace(), init_levels(), might_contain_ghosts(), and might_contain_ghosts().
void ug::IApproximationSpace::print_layout_statistic | ( | ) | const |
prints statistic on layouts
References ug::LogAssistant::get_output_process(), ug::GetLogAssistant(), m_vDD, ug::NumIndices(), pcl::NumProcs(), ug::repeat(), and UG_LOG.
void ug::IApproximationSpace::print_statistic | ( | ) | const |
prints statistic about DoF Distribution
References print_statistic().
Referenced by print_statistic().
void ug::IApproximationSpace::print_statistic | ( | std::string | flags | ) | const |
prints statistic about DoF Distribution
References ug::AlgebraType::blocksize(), ug::DefaultAlgebra::get(), ug::GetLogAssistant(), m_vDD, pcl::NumProcs(), ug::PrintDoFCount(), pcl::ProcRank(), PROFILE_BEGIN, PROFILE_END, PROFILE_FUNC, ug::repeat(), ug::ToLower(), UG_LOG, and ug::AlgebraType::VariableBlockSize.
|
protected |
registers at message hub for grid adaption
References grid_changed_callback(), grid_distribution_callback(), m_spGridAdaptionCallbackID, m_spGridDistributionCallbackID, and m_spMGSH.
Referenced by init().
|
protected |
reinits all data after grid adaption
References m_RevCnt, m_spSurfaceView, m_vDD, and PROFILE_FUNC.
Referenced by grid_changed_callback(), and grid_distribution_callback().
|
inline |
returns the current revision
References m_RevCnt.
|
inline |
get underlying subset handler
References m_spMGSH.
|
inline |
returns the approximation space
References m_spSurfaceView.
|
protected |
creates surface SurfaceView if needed
References m_spMGSH, and m_spSurfaceView.
Referenced by create_dof_distribution().
|
protected |
suitable algebra type for the index distribution pattern
Referenced by dof_distribution_info_required(), and init().
|
protected |
Referenced by grid_changed_callback(), and init().
|
protected |
flag if DoFs should be grouped
Referenced by create_dof_distribution(), dof_distribution_info_required(), grouped(), and init().
|
protected |
revision counter
Referenced by init(), reinit(), and revision().
|
protected |
DofDistributionInfo.
Referenced by add(), add(), clear(), create_dof_distribution(), ddinfo(), dof_distribution_info_required(), and init().
|
protected |
Index Storage for Level (ghost / noghost)
Referenced by create_dof_distribution().
|
protected |
Index Storage for Level (ghost / noghost)
Referenced by create_dof_distribution().
|
protected |
message hub id
Referenced by register_at_adaption_msg_hub().
|
protected |
Referenced by register_at_adaption_msg_hub().
multigrid, where elements are stored
Referenced by create_dof_distribution(), grid_changed_callback(), grid_distribution_callback(), init(), and might_contain_ghosts().
|
protected |
subsethandler, where elements are stored
Referenced by add(), add(), add(), add(), ug::ApproximationSpace< TDomain >::ApproximationSpace(), ug::ApproximationSpace< TDomain >::ApproximationSpace(), create_dof_distribution(), init(), num_levels(), register_at_adaption_msg_hub(), subset_handler(), and surface_view_required().
|
protected |
Surface View.
Referenced by create_dof_distribution(), grid_changed_callback(), grid_distribution_callback(), reinit(), surface_view(), surface_view_required(), and ~IApproximationSpace().
|
protected |
MG Level DoF Distribution.
Referenced by create_dof_distribution(), dof_distribution(), dof_distributions(), print_layout_statistic(), print_statistic(), and reinit().