Reference documentation for deal.II version 9.5.0
|
#include <deal.II/grid/grid_tools.h>
Public Types | |
using | IntersectionType = std::array<::Point< spacedim >, structdim+1 > |
Public Member Functions | |
template<int dim> | |
GridTools::internal::DistributedComputePointLocationsInternal< dim, spacedim > | convert_to_distributed_compute_point_locations_internal (const unsigned int n_points_1D, const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping, const bool consistent_numbering_of_sender_and_receiver=false) const |
template<int dim> | |
DistributedComputePointLocationsInternal< dim, spacedim > | convert_to_distributed_compute_point_locations_internal (const unsigned int n_points_1D, const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping, const bool consistent_numbering_of_sender_and_receiver) const |
Public Attributes | |
std::vector< std::tuple< std::pair< int, int >, unsigned int, unsigned int, IntersectionType > > | send_components |
std::vector< std::tuple< unsigned int, unsigned int, IntersectionType > > | recv_components |
std::vector< unsigned int > | recv_ptrs |
Private Member Functions | |
std::map< unsigned int, std::vector< unsigned int > > | communicate_indices (const std::vector< std::tuple< unsigned int, unsigned int, unsigned int > > &point_recv_components, const MPI_Comm comm) const |
Data structure returned by GridTools::internal::distributed_compute_intersection_locations(). It can be converted to GridTools::internal::DistributedComputePointLocationsInternal, which can be used to reinit Utilities::MPI::RemotePointEvaluation.
Definition at line 1325 of file grid_tools.h.
using GridTools::internal::DistributedComputeIntersectionLocationsInternal< structdim, spacedim >::IntersectionType = std::array<::Point<spacedim>, structdim + 1> |
Intersections are assumed to be simplices (as, e.g., provided by CGAL)
Definition at line 1330 of file grid_tools.h.
GridTools::internal::DistributedComputePointLocationsInternal< dim, spacedim > GridTools::internal::DistributedComputeIntersectionLocationsInternal< structdim, spacedim >::convert_to_distributed_compute_point_locations_internal | ( | const unsigned int | n_points_1D, |
const Triangulation< dim, spacedim > & | tria, | ||
const Mapping< dim, spacedim > & | mapping, | ||
const bool | consistent_numbering_of_sender_and_receiver = false |
||
) | const |
Distribute quadrature points according to QGaussSimplex<structdim>(n_points_1D) on found intersections and construct GridTools::internal::DistributedComputePointLocationsInternal from class members. This can be done without searching for points again since all information is locally known.
The parameter consistent_numbering_of_sender_and_receiver
can be used to ensure points on sender and receiver side are numbered consistently. This parameter is optional if DistributedComputePointLocationsInternal is used to set up RemotePointEvaluation, but might be helpful for debugging or other usage of DistributedComputePointLocationsInternal. Note that setting this parameter true requires an additional communication step during the setup phase.
|
private |
Helper function for convert_to_distributed_compute_point_locations_internal(). It sends the indices associated to quadrature points at the receiver side to the sender side, where the information is needed to build GridTools::internal::DistributedComputePointLocationsInternal::send_components
Definition at line 6450 of file grid_tools.cc.
DistributedComputePointLocationsInternal< dim, spacedim > GridTools::internal::DistributedComputeIntersectionLocationsInternal< structdim, spacedim >::convert_to_distributed_compute_point_locations_internal | ( | const unsigned int | n_points_1D, |
const Triangulation< dim, spacedim > & | tria, | ||
const Mapping< dim, spacedim > & | mapping, | ||
const bool | consistent_numbering_of_sender_and_receiver | ||
) | const |
Definition at line 6343 of file grid_tools.cc.
std::vector<std::tuple<std::pair<int, int>, unsigned int, unsigned int, IntersectionType> > GridTools::internal::DistributedComputeIntersectionLocationsInternal< structdim, spacedim >::send_components |
Information of each intersection on sending/evaluation side. The elements of the tuple are as follows: 0) cell level and index, 1) rank of the owning process, 2) local index of the owning process, 3) found intersection.
Definition at line 1345 of file grid_tools.h.
std::vector<std::tuple<unsigned int, unsigned int, IntersectionType> > GridTools::internal::DistributedComputeIntersectionLocationsInternal< structdim, spacedim >::recv_components |
Information of each received data value. The elements of the tuple are as follows: 0) rank of sender, 1) local index, 2) found intersections.
Definition at line 1356 of file grid_tools.h.
std::vector<unsigned int> GridTools::internal::DistributedComputeIntersectionLocationsInternal< structdim, spacedim >::recv_ptrs |
Pointers of ranges to found intersections for requested intersection.
Definition at line 1361 of file grid_tools.h.