InterGridMap< GridClass > Class Template Reference
[Grid classes]

List of all members.

Classes

class  ExcIncompatibleGrids
class  ExcInvalidKey

Public Types

typedef GridClass::cell_iterator cell_iterator

Public Member Functions

 InterGridMap ()
void make_mapping (const GridClass &source_grid, const GridClass &destination_grid)
cell_iterator operator[] (const cell_iterator &source_cell) const
void clear ()
const GridClass & get_source_grid () const
const GridClass & get_destination_grid () const
unsigned int memory_consumption () const

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 GridClass > source_grid
SmartPointer< const GridClass > destination_grid


Detailed Description

template<class GridClass>
class InterGridMap< GridClass >

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.

The implementation of this class is such that not only cell iterators into triangulations can be mapped, but also iterators into objects of type DoFHandler, hp::DoFHandler and MGDoFHandler. 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++ 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 consumpion.

Usage

In practice, use of this class is as follows:

 *                   // have two grids, which are derived from the
 *                   // same coarse grid
 *   Triangulation<dim> tria1, tria2;
 *   DoFHandler<dim> dof_handler_1(tria1), dof_handler_2(tria2);
 *   ...
 *                   // do something with these objects, e.g.
 *                   // refine the triangulations differently,
 *                   // distribute degrees of freedom, etc
 *   ...
 *                   // create the mapping
 *   InterGridMap<DoFHandler<dim> > grid_1_to_2_map;
 *   grid_1_to_2_map.make_mapping (dof_handler_1,
 *                                 dof_handler_2);
 *   ...
 *   typename DoFHandler<dim>::cell_iterator cell = dof_handler_1.begin(),
 *                                           endc = dof_handler_1.end();
 *   for (; cell!=endc; ++cell)
 *                    // now do something with the cell of dof_handler_2
 *                    // corresponding to @p cell (which is one of
 *                    // dof_handler_1
 *     f( grid_1_to_2_map[cell]);
 * 

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, PersistentTriangulation, hp::DoFHandler or MGDoFHandler).

Author:
Wolfgang Bangerth, 1999

Member Typedef Documentation

template<class GridClass >
typedef GridClass::cell_iterator InterGridMap< GridClass >::cell_iterator

Typedef to the iterator type of the grid class under consideration.


Constructor & Destructor Documentation

template<class GridClass >
InterGridMap< GridClass >::InterGridMap (  ) 

Constructor setting the class name arguments in the SmartPointer members.


Member Function Documentation

template<class GridClass >
void InterGridMap< GridClass >::make_mapping ( const GridClass &  source_grid,
const GridClass &  destination_grid 
)

Create the mapping between the two grids.

template<class GridClass >
cell_iterator InterGridMap< GridClass >::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.

template<class GridClass >
void InterGridMap< GridClass >::clear (  ) 

Delete all data of this class.

template<class GridClass >
const GridClass& InterGridMap< GridClass >::get_source_grid (  )  const

Return a pointer to the source grid.

template<class GridClass >
const GridClass& InterGridMap< GridClass >::get_destination_grid (  )  const

Return a pointer to the destination grid.

template<class GridClass >
unsigned int InterGridMap< GridClass >::memory_consumption (  )  const

Determine an estimate for the memory consumption (in bytes) of this object.

template<class GridClass >
void InterGridMap< GridClass >::set_mapping ( const cell_iterator src_cell,
const cell_iterator dst_cell 
) [private]

Set the mapping for the pair of cells given. These shall match in level of refinement and all other properties.

template<class GridClass >
void InterGridMap< GridClass >::set_entries_to_cell ( const cell_iterator src_cell,
const cell_iterator dst_cell 
) [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.


Member Data Documentation

template<class GridClass >
std::vector<std::vector<cell_iterator> > InterGridMap< GridClass >::mapping [private]

The actual data. Hold one iterator for each cell on each level.

template<class GridClass >
SmartPointer<const GridClass > InterGridMap< GridClass >::source_grid [private]

Store a pointer to the source grid.

template<class GridClass >
SmartPointer<const GridClass > InterGridMap< GridClass >::destination_grid [private]

Likewise for the destination grid.


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

deal.II documentation generated on Sat Aug 15 16:52:03 2009 by doxygen 1.5.9