Reference documentation for deal.II version 9.4.1
|
#include <deal.II/grid/grid_tools_cache.h>
Public Member Functions | |
Cache (const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >())) | |
~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 RTree< std::pair< BoundingBox< spacedim >, typename Triangulation< dim, spacedim >::active_cell_iterator > > & | get_locally_owned_cell_bounding_boxes_rtree () const |
const std::vector< std::set< unsigned int > > & | get_vertex_to_neighbor_subdomain () 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 unsigned int level=0) const |
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 |
std::map< unsigned int, RTree< std::pair< BoundingBox< spacedim >, unsigned int > > > | covering_rtree |
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 |
RTree< std::pair< BoundingBox< spacedim >, typename Triangulation< dim, spacedim >::active_cell_iterator > > | locally_owned_cell_bounding_boxes_rtree |
std::vector< std::set< unsigned int > > | vertex_to_neighbor_subdomain |
boost::signals2::connection | tria_signal |
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) |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
static std::mutex | mutex |
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 |
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.
Definition at line 65 of file grid_tools_cache.h.
GridTools::Cache< dim, spacedim >::Cache | ( | const Triangulation< dim, spacedim > & | tria, |
const Mapping< dim, spacedim > & | mapping = (ReferenceCells::get_hypercube<dim>() .template get_default_linear_mapping<dim, spacedim>() ) |
||
) |
Constructor.
If you provide the optional mapping
argument, then this is used whenever a mapping is required.
tria | The triangulation for which to store information |
mapping | The mapping to use when computing cached objects |
Definition at line 28 of file grid_tools_cache.cc.
|
override |
Destructor.
Definition at line 39 of file grid_tools_cache.cc.
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.
flags | What to mark for update |
Definition at line 51 of file grid_tools_cache.cc.
const std::vector< std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > > & GridTools::Cache< dim, spacedim >::get_vertex_to_cell_map |
Return the cached vertex_to_cell_map as computed by GridTools::vertex_to_cell_map().
Definition at line 61 of file grid_tools_cache.cc.
const std::vector< std::vector< Tensor< 1, spacedim > > > & GridTools::Cache< dim, spacedim >::get_vertex_to_cell_centers_directions |
Return the cached vertex_to_cell_centers_directions as computed by GridTools::vertex_to_cell_centers_directions().
Definition at line 75 of file grid_tools_cache.cc.
const std::map< unsigned int, Point< spacedim > > & GridTools::Cache< dim, spacedim >::get_used_vertices |
Return the cached map of used vertices, as computed by GridTools::extract_used_vertices().
Definition at line 90 of file grid_tools_cache.cc.
const RTree< std::pair< Point< spacedim >, unsigned int > > & GridTools::Cache< dim, spacedim >::get_used_vertices_rtree |
Return the cached RTree object for the vertices, constructed using the used vertices of the triangulation.
Definition at line 104 of file grid_tools_cache.cc.
const RTree< std::pair< BoundingBox< spacedim >, typename Triangulation< dim, spacedim >::active_cell_iterator > > & GridTools::Cache< dim, spacedim >::get_cell_bounding_boxes_rtree |
Return the cached RTree object of the cell bounding boxes, constructed using the active cell iterators of the stored triangulation. For parallel::distributed::Triangulation objects, this function will return also the bounding boxes of ghost and artificial cells.
Definition at line 126 of file grid_tools_cache.cc.
const RTree< std::pair< BoundingBox< spacedim >, typename Triangulation< dim, spacedim >::active_cell_iterator > > & GridTools::Cache< dim, spacedim >::get_locally_owned_cell_bounding_boxes_rtree |
Return the cached RTree object of bounding boxes containing locally owned active cells, constructed using the active cell iterators of the stored triangulation.
In contrast to the previous function, this function builds the RTree using only the locally owned cells, i.e., not including ghost or artificial cells. The two functions return the same result in serial computations.
Definition at line 150 of file grid_tools_cache.cc.
const std::vector< std::set< unsigned int > > & GridTools::Cache< dim, spacedim >::get_vertex_to_neighbor_subdomain |
Returns the vector of set of integer containing the subdomain id to which each vertex is connected to. This feature is used extensively in the particle_handler to detect on which processors ghost particles must be built.
Definition at line 203 of file grid_tools_cache.cc.
|
inline |
Return a reference to the stored triangulation.
Definition at line 306 of file grid_tools_cache.h.
|
inline |
Return a reference to the stored mapping.
Definition at line 315 of file grid_tools_cache.h.
const RTree< std::pair< BoundingBox< spacedim >, unsigned int > > & GridTools::Cache< dim, spacedim >::get_covering_rtree | ( | const unsigned int | level = 0 | ) | 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 the cells falling within the given bounding box.
Given a point on a parallel::TriangulationBase, 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.
The local bounding boxes are constructed by extracting the specified level
from the rtree object returned by the get_locally_owned_cell_bounding_boxes_rtree(). Notice that level
in this context does not refer to the level of the triangulation, but refer to the level of the RTree object (see, e.g., https://en.wikipedia.org/wiki/R-tree).
Definition at line 174 of file grid_tools_cache.cc.
|
mutableprivate |
Keep track of what needs to be updated next.
Definition at line 219 of file grid_tools_cache.h.
|
private |
A pointer to the Triangulation.
Definition at line 224 of file grid_tools_cache.h.
|
private |
Mapping to use when computing on the Triangulation.
Definition at line 229 of file grid_tools_cache.h.
|
mutableprivate |
Store vertex to cell map information, as generated by GridTools::vertex_to_cell_map()
Definition at line 238 of file grid_tools_cache.h.
|
mutableprivate |
Store vertex to cell center directions, as generated by GridTools::vertex_to_cell_centers_directions().
Definition at line 245 of file grid_tools_cache.h.
|
mutableprivate |
A collection of rtree objects covering the whole mesh.
Each entry of the map is constructed from the function extract_rtree_level() applied to the rtree returned by the function get_locally_owned_cell_bounding_boxes_rtree(), with the specified level in input.
Definition at line 257 of file grid_tools_cache.h.
|
mutableprivate |
Store the used vertices of the Triangulation, as generated by GridTools::extract_used_vertices().
Definition at line 263 of file grid_tools_cache.h.
|
mutableprivate |
Store an RTree object, containing the used vertices of the triangulation.
Definition at line 268 of file grid_tools_cache.h.
|
mutableprivate |
Store an RTree object, containing the bounding boxes of the cells of the triangulation.
Definition at line 277 of file grid_tools_cache.h.
|
mutableprivate |
Store an RTree object, containing the bounding boxes of the locally owned cells of the triangulation.
Definition at line 286 of file grid_tools_cache.h.
|
mutableprivate |
Store an std::vector of std::set of integer containing the id of all subdomain to which a vertex is connected to.
Definition at line 293 of file grid_tools_cache.h.
|
private |
Storage for the status of the triangulation signal.
Definition at line 298 of file grid_tools_cache.h.