Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Member Functions | Private Attributes | List of all members
GridTools::Cache< dim, spacedim > Class Template Reference

#include <deal.II/grid/grid_tools_cache.h>

Inheritance diagram for GridTools::Cache< dim, spacedim >:
[legend]

Public Member Functions

 Cache (const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping=StaticMappingQ1< dim, spacedim >::mapping)
 
 ~Cache () override
 
void mark_for_update (const CacheUpdateFlags &flags=update_all)
 
const std::vector< std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > > & get_vertex_to_cell_map () const
 
const std::vector< std::vector< Tensor< 1, spacedim > > > & get_vertex_to_cell_centers_directions () const
 
const std::map< unsigned int, Point< spacedim > > & get_used_vertices () const
 
const RTree< std::pair< Point< spacedim >, unsigned int > > & get_used_vertices_rtree () const
 
const RTree< std::pair< BoundingBox< spacedim >, typename Triangulation< dim, spacedim >::active_cell_iterator > > & get_cell_bounding_boxes_rtree () const
 
const Triangulation< dim, spacedim > & get_triangulation () const
 
const Mapping< dim, spacedim > & get_mapping () const
 
const RTree< std::pair< BoundingBox< spacedim >, unsigned int > > & get_covering_rtree () const
 
const KDTree< spacedim > & get_vertex_kdtree () const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
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)
 

Private Attributes

CacheUpdateFlags update_flags
 
SmartPointer< const Triangulation< dim, spacedim >, Cache< dim, spacedim > > tria
 
SmartPointer< const Mapping< dim, spacedim >, Cache< dim, spacedim > > mapping
 
std::vector< std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > > vertex_to_cells
 
std::vector< std::vector< Tensor< 1, spacedim > > > vertex_to_cell_centers
 
RTree< std::pair< BoundingBox< spacedim >, unsigned int > > covering_rtree
 
KDTree< spacedim > vertex_kdtree
 
std::map< unsigned int, Point< spacedim > > used_vertices
 
RTree< std::pair< Point< spacedim >, unsigned int > > used_vertices_rtree
 
RTree< std::pair< BoundingBox< spacedim >, typename Triangulation< dim, spacedim >::active_cell_iterator > > cell_bounding_boxes_rtree
 
boost::signals2::connection tria_signal
 

Additional Inherited Members

- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Detailed Description

template<int dim, int spacedim = dim>
class GridTools::Cache< dim, spacedim >

A class that caches computationally intensive information about a Triangulation.

This class attaches a signal to the Triangulation passed at construction time to keep track about refinement changes, and allows the user to query some of the data structures constructed using functions in the GridTools namespace which are computed only once, and then cached inside this class for faster access whenever the triangulation has not changed.

Notice that this class only notices if the underlying Triangulation has changed due to a Triangulation::Signals::any_change() signal being triggered.

If the triangulation changes for other reasons, for example because you use it in conjunction with a MappingQEulerian object that sees the vertices through its own transformation, or because you manually change some vertex locations, then some of the structures in this class become obsolete, and you will have to mark them as outdated, by calling the method mark_for_update() manually.

Author
Luca Heltai, 2017.

Definition at line 125 of file grid_tools.h.

Constructor & Destructor Documentation

◆ Cache()

template<int dim, int spacedim>
GridTools::Cache< dim, spacedim >::Cache ( const Triangulation< dim, spacedim > &  tria,
const Mapping< dim, spacedim > &  mapping = StaticMappingQ1<dim, spacedim>::mapping 
)

Constructor.

If you provide the optional mapping argument, then this is used whenever a mapping is required.

Parameters
triaThe triangulation for which to store information
mappingThe mapping to use when computing cached objects

Definition at line 30 of file grid_tools_cache.cc.

◆ ~Cache()

template<int dim, int spacedim>
GridTools::Cache< dim, spacedim >::~Cache ( )
override

Destructor.

Definition at line 41 of file grid_tools_cache.cc.

Member Function Documentation

◆ mark_for_update()

template<int dim, int spacedim>
void GridTools::Cache< dim, spacedim >::mark_for_update ( const CacheUpdateFlags flags = update_all)

Make sure that the objects marked for update are recomputed during subsequent calls to the get_* functions defined in this class.

Notice that no work is performed when you call this function. The actual data structures are computed the next time you call the corresponding get_* method.

Parameters
flagsWhat to mark for update

Definition at line 53 of file grid_tools_cache.cc.

◆ get_vertex_to_cell_map()

template<int dim, int spacedim>
const std::vector< std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > > & GridTools::Cache< dim, spacedim >::get_vertex_to_cell_map ( ) const

Return the cached vertex_to_cell_map as computed by GridTools::vertex_to_cell_map().

Definition at line 63 of file grid_tools_cache.cc.

◆ get_vertex_to_cell_centers_directions()

template<int dim, int spacedim>
const std::vector< std::vector< Tensor< 1, spacedim > > > & GridTools::Cache< dim, spacedim >::get_vertex_to_cell_centers_directions ( ) const

Return the cached vertex_to_cell_centers_directions as computed by GridTools::vertex_to_cell_centers_directions().

Definition at line 77 of file grid_tools_cache.cc.

◆ get_used_vertices()

template<int dim, int spacedim>
const std::map< unsigned int, Point< spacedim > > & GridTools::Cache< dim, spacedim >::get_used_vertices ( ) const

Return the cached map of used vertices, as computed by GridTools::extract_used_vertices().

Definition at line 92 of file grid_tools_cache.cc.

◆ get_used_vertices_rtree()

template<int dim, int spacedim>
const RTree< std::pair< Point< spacedim >, unsigned int > > & GridTools::Cache< dim, spacedim >::get_used_vertices_rtree ( ) const

Return the cached RTree object for the vertices, constructed using the used vertices of the triangulation.

Definition at line 106 of file grid_tools_cache.cc.

◆ get_cell_bounding_boxes_rtree()

template<int dim, int spacedim>
const RTree< std::pair< BoundingBox< spacedim >, typename Triangulation< dim, spacedim >::active_cell_iterator > > & GridTools::Cache< dim, spacedim >::get_cell_bounding_boxes_rtree ( ) const

Return the cached RTree object of the cell bounging boxes, constructed using the active cell iterators of the stored triangulation.

Definition at line 127 of file grid_tools_cache.cc.

◆ get_triangulation()

template<int dim, int spacedim>
const Triangulation< dim, spacedim > & GridTools::Cache< dim, spacedim >::get_triangulation ( ) const
inline

Return a reference to the stored triangulation.

Definition at line 265 of file grid_tools_cache.h.

◆ get_mapping()

template<int dim, int spacedim>
const Mapping< dim, spacedim > & GridTools::Cache< dim, spacedim >::get_mapping ( ) const
inline

Return a reference to the stored mapping.

Definition at line 274 of file grid_tools_cache.h.

◆ get_covering_rtree()

template<int dim, int spacedim>
const RTree< std::pair< BoundingBox< spacedim >, unsigned int > > & GridTools::Cache< dim, spacedim >::get_covering_rtree ( ) const

This function returns an object that allows identifying which process(es) in a parallel computation may own the cell that surrounds a given point. The elements of this object – an Rtree – are pairs of bounding boxes denoting areas that cover all or parts of the local portion of a parallel triangulation, and an unsigned int representing the process or subdomain that owns these cells. Given a point on a parallel::Triangulation, this tree allows to identify one, or few candidate processes, for which the point lies on a locally owned cell.

Constructing or updating the rtree requires a call to GridTools::build_global_description_tree(), which exchanges bounding boxes between all processes using Utilities::MPI::all_gather(), a collective operation. Therefore this function must be called by all processes at the same time.

While each box may only cover part of a process's locally owned part of the triangulation, the boxes associated with each process jointly cover the entire local portion.

Definition at line 163 of file grid_tools_cache.cc.

◆ get_vertex_kdtree()

template<int dim, int spacedim>
const KDTree< spacedim > & GridTools::Cache< dim, spacedim >::get_vertex_kdtree ( ) const

Return the cached vertex_kdtree object, constructed with the vertices of the stored triangulation.

Definition at line 148 of file grid_tools_cache.cc.

Member Data Documentation

◆ update_flags

template<int dim, int spacedim = dim>
CacheUpdateFlags GridTools::Cache< dim, spacedim >::update_flags
mutableprivate

Keep track of what needs to be updated next.

Definition at line 193 of file grid_tools_cache.h.

◆ tria

template<int dim, int spacedim = dim>
SmartPointer<const Triangulation<dim, spacedim>, Cache<dim, spacedim> > GridTools::Cache< dim, spacedim >::tria
private

A pointer to the Triangulation.

Definition at line 198 of file grid_tools_cache.h.

◆ mapping

template<int dim, int spacedim = dim>
SmartPointer<const Mapping<dim, spacedim>, Cache<dim, spacedim> > GridTools::Cache< dim, spacedim >::mapping
private

Mapping to use when computing on the Triangulation.

Definition at line 203 of file grid_tools_cache.h.

◆ vertex_to_cells

template<int dim, int spacedim = dim>
std::vector< std::set<typename Triangulation<dim, spacedim>::active_cell_iterator> > GridTools::Cache< dim, spacedim >::vertex_to_cells
mutableprivate

Store vertex to cell map information, as generated by GridTools::vertex_to_cell_map()

Definition at line 212 of file grid_tools_cache.h.

◆ vertex_to_cell_centers

template<int dim, int spacedim = dim>
std::vector<std::vector<Tensor<1, spacedim> > > GridTools::Cache< dim, spacedim >::vertex_to_cell_centers
mutableprivate

Store vertex to cell center directions, as generated by GridTools::vertex_to_cell_centers_directions().

Definition at line 219 of file grid_tools_cache.h.

◆ covering_rtree

template<int dim, int spacedim = dim>
RTree<std::pair<BoundingBox<spacedim>, unsigned int> > GridTools::Cache< dim, spacedim >::covering_rtree
mutableprivate

An rtree object covering the whole mesh.

Definition at line 225 of file grid_tools_cache.h.

◆ vertex_kdtree

template<int dim, int spacedim = dim>
KDTree<spacedim> GridTools::Cache< dim, spacedim >::vertex_kdtree
mutableprivate

A KDTree object constructed with the vertices of the triangulation.

Definition at line 231 of file grid_tools_cache.h.

◆ used_vertices

template<int dim, int spacedim = dim>
std::map<unsigned int, Point<spacedim> > GridTools::Cache< dim, spacedim >::used_vertices
mutableprivate

Store the used vertices of the Triangulation, as generated by GridTools::extract_used_vertices().

Definition at line 238 of file grid_tools_cache.h.

◆ used_vertices_rtree

template<int dim, int spacedim = dim>
RTree<std::pair<Point<spacedim>, unsigned int> > GridTools::Cache< dim, spacedim >::used_vertices_rtree
mutableprivate

Store an RTree object, containing the used vertices of the triangulation.

Definition at line 243 of file grid_tools_cache.h.

◆ cell_bounding_boxes_rtree

template<int dim, int spacedim = dim>
RTree< std::pair<BoundingBox<spacedim>, typename Triangulation<dim, spacedim>::active_cell_iterator> > GridTools::Cache< dim, spacedim >::cell_bounding_boxes_rtree
mutableprivate

Store an RTree object, containing the bounding boxes of the cells of the triangulation.

Definition at line 252 of file grid_tools_cache.h.

◆ tria_signal

template<int dim, int spacedim = dim>
boost::signals2::connection GridTools::Cache< dim, spacedim >::tria_signal
private

Storage for the status of the triangulation signal.

Definition at line 257 of file grid_tools_cache.h.


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