Reference documentation for deal.II version 9.4.1
|
#include <deal.II/grid/intergrid_map.h>
Public Types | |
using | cell_iterator = typename MeshType::cell_iterator |
Public Member Functions | |
InterGridMap () | |
void | make_mapping (const MeshType &source_grid, const MeshType &destination_grid) |
cell_iterator | operator[] (const cell_iterator &source_cell) const |
void | clear () |
const MeshType & | get_source_grid () const |
const MeshType & | get_destination_grid () const |
std::size_t | memory_consumption () const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcInvalidKey (cell_iterator arg1) |
static ::ExceptionBase & | ExcIncompatibleGrids () |
Private Member Functions | |
void | set_mapping (const cell_iterator &src_cell, const cell_iterator &dst_cell) |
void | set_entries_to_cell (const cell_iterator &src_cell, const cell_iterator &dst_cell) |
Private Attributes | |
std::vector< std::vector< cell_iterator > > | mapping |
SmartPointer< const MeshType, InterGridMap< MeshType > > | source_grid |
SmartPointer< const MeshType, InterGridMap< MeshType > > | destination_grid |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
void | check_no_subscribers () const noexcept |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
static std::mutex | mutex |
This class provides a map between two grids which are derived from the same coarse grid. For each cell iterator of the source map, it provides the respective cell iterator on the destination map, through its operator []
.
Usually, the two grids will be refined differently. Then, the value returned for an iterator on the source grid will be either:
Keys for this map are all cells on the source grid, whether active or not.
For example, consider these two one-dimensional grids:
* Grid 1: * x--x--x-----x-----------x * 1 2 3 4 * * Grid 2: * x-----x-----x-----x-----x * 1 2 3 4 *
(Cell numbers are only given as an example and will not correspond to real cell iterator's indices.) The mapping from grid 1 to grid 2 will then be as follows:
* Cell on grid 1 Cell on grid 2 * 1 ------------------> 1 * 2 ------------------> 1 * 3 ------------------> 2 * 4 ------------------> mother cell of cells 3 and 4 * (a non-active cell, not shown here) *
Besides the mappings shown here, the non-active cells on grid 1 are also valid keys. For example, the mapping for the mother cell of cells 1 and 2 on the first grid will point to cell 1 on the second grid.
MeshType | This class may be used with any class that satisfies the MeshType concept. The extension to other classes offering iterator functions and some minor additional requirements is simple. |
Note that this class could in principle be based on the C++ std::map<Key,Value>
data type. Instead, it uses another data format which is more effective both in terms of computing time for access as well as with regard to memory consumption.
In practice, use of this class is as follows:
Note that the template parameters to this class have to be given as InterGridMap<DoFHandler<2> >
, which here is DoFHandler (and could equally well be Triangulation or PersistentTriangulation).
Definition at line 114 of file intergrid_map.h.
using InterGridMap< MeshType >::cell_iterator = typename MeshType::cell_iterator |
Typedef to the iterator type of the grid class under consideration.
Definition at line 120 of file intergrid_map.h.
InterGridMap< MeshType >::InterGridMap |
Constructor setting the class name arguments in the SmartPointer members.
Definition at line 37 of file intergrid_map.cc.
void InterGridMap< MeshType >::make_mapping | ( | const MeshType & | source_grid, |
const MeshType & | destination_grid | ||
) |
Create the mapping between the two grids.
Definition at line 46 of file intergrid_map.cc.
InterGridMap< MeshType >::cell_iterator InterGridMap< MeshType >::operator[] | ( | const cell_iterator & | source_cell | ) | const |
Access operator: give a cell on the source grid and receive the respective cell on the other grid, or if that does not exist, the most refined cell of which the source cell would be created if it were further refined.
Definition at line 157 of file intergrid_map.cc.
void InterGridMap< MeshType >::clear |
Delete all data of this class.
Definition at line 174 of file intergrid_map.cc.
const MeshType & InterGridMap< MeshType >::get_source_grid |
Return a reference to the source grid.
Definition at line 185 of file intergrid_map.cc.
const MeshType & InterGridMap< MeshType >::get_destination_grid |
Return a reference to the destination grid.
Definition at line 194 of file intergrid_map.cc.
std::size_t InterGridMap< MeshType >::memory_consumption |
Determine an estimate for the memory consumption (in bytes) of this object.
Definition at line 203 of file intergrid_map.cc.
|
private |
Set the mapping for the pair of cells given. These shall match in level of refinement and all other properties.
Definition at line 102 of file intergrid_map.cc.
|
private |
Set the value of the key src_cell
to dst_cell
. Do so as well for all the children and their children of src_cell
. This function is used for cells which are more refined on src_grid
than on dst_grid
; then all values of the hierarchy of cells and their children point to one cell on the dst_grid
.
Definition at line 141 of file intergrid_map.cc.
|
private |
The actual data. Hold one iterator for each cell on each level.
Definition at line 183 of file intergrid_map.h.
|
private |
Store a pointer to the source grid.
Definition at line 188 of file intergrid_map.h.
|
private |
Likewise for the destination grid.
Definition at line 193 of file intergrid_map.h.