Reference documentation for deal.II version GIT relicensing-487-ge9eb5ab491 2024-04-25 07:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | Private Member Functions | List of all members
GridTools::internal::DistributedComputeIntersectionLocationsInternal< structdim, spacedim > Struct Template Reference

#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, std::vector< Quadrature< spacedim > > *mapped_quadratures_recv_comp=nullptr, const bool consistent_numbering_of_sender_and_receiver=false) 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 intrecv_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
 
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, std::vector< Quadrature< spacedim > > *mapped_quadratures_recv_comp, const bool consistent_numbering_of_sender_and_receiver) const
 

Detailed Description

template<int structdim, int spacedim>
struct GridTools::internal::DistributedComputeIntersectionLocationsInternal< structdim, spacedim >

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 843 of file grid_tools.h.

Member Typedef Documentation

◆ IntersectionType

template<int structdim, int spacedim>
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 848 of file grid_tools.h.

Member Function Documentation

◆ convert_to_distributed_compute_point_locations_internal() [1/2]

template<int structdim, int spacedim>
template<int dim>
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,
std::vector< Quadrature< spacedim > > *  mapped_quadratures_recv_comp = nullptr,
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.

mapped_quadratures_recv_comp is a pointer to an empty vector of mapped quadratures. By default it is a nullptr and the parameter is ignored. Otherwise, the vector is filled with the mapped quadrature rules (in real coordinates) corresponding to recv_components.

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.

◆ communicate_indices()

template<int structdim, int spacedim>
std::map< unsigned int, std::vector< unsigned int > > GridTools::internal::DistributedComputeIntersectionLocationsInternal< structdim, spacedim >::communicate_indices ( const std::vector< std::tuple< unsigned int, unsigned int, unsigned int > > &  point_recv_components,
const MPI_Comm  comm 
) const
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 4221 of file grid_tools.cc.

◆ convert_to_distributed_compute_point_locations_internal() [2/2]

template<int structdim, int spacedim>
template<int dim>
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,
std::vector< Quadrature< spacedim > > *  mapped_quadratures_recv_comp,
const bool  consistent_numbering_of_sender_and_receiver 
) const
private

Definition at line 4103 of file grid_tools.cc.

Member Data Documentation

◆ send_components

template<int structdim, int spacedim>
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.

Note
The vector is sorted according to 1), 2).

Definition at line 863 of file grid_tools.h.

◆ recv_components

template<int structdim, int spacedim>
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.

Note
The vector is sorted according to 1), 0), 2).
Multiple intersections between cells can be found

Definition at line 874 of file grid_tools.h.

◆ recv_ptrs

template<int structdim, int spacedim>
std::vector<unsigned int> GridTools::internal::DistributedComputeIntersectionLocationsInternal< structdim, spacedim >::recv_ptrs

Pointers of ranges to found intersections for requested intersection.

Definition at line 879 of file grid_tools.h.


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