Reference documentation for deal.II version GIT relicensing-478-g3275795167 2024-04-24 07:10: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
Classes | Public Member Functions | Private Attributes | List of all members
Utilities::MPI::RemotePointEvaluation< dim, spacedim > Class Template Reference

#include <deal.II/base/mpi_remote_point_evaluation.h>

Inheritance diagram for Utilities::MPI::RemotePointEvaluation< dim, spacedim >:
Inheritance graph
[legend]

Classes

struct  AdditionalData
 
class  CellData
 

Public Member Functions

 RemotePointEvaluation (const AdditionalData &additional_data=AdditionalData())
 
 RemotePointEvaluation (const double tolerance, const bool enforce_unique_mapping=false, const unsigned int rtree_level=0, const std::function< std::vector< bool >()> &marked_vertices={})
 
 ~RemotePointEvaluation ()
 
void reinit (const std::vector< Point< spacedim > > &points, const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping)
 
void reinit (const GridTools::internal::DistributedComputePointLocationsInternal< dim, spacedim > &data, const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping)
 
const CellDataget_cell_data () const
 
template<typename T >
void evaluate_and_process (std::vector< T > &output, std::vector< T > &buffer, const std::function< void(const ArrayView< T > &, const CellData &)> &evaluation_function) const
 
template<typename T >
std::vector< T > evaluate_and_process (const std::function< void(const ArrayView< T > &, const CellData &)> &evaluation_function) const
 
template<typename T >
void process_and_evaluate (const std::vector< T > &input, std::vector< T > &buffer, const std::function< void(const ArrayView< const T > &, const CellData &)> &evaluation_function) const
 
template<typename T >
void process_and_evaluate (const std::vector< T > &input, const std::function< void(const ArrayView< const T > &, const CellData &)> &evaluation_function) const
 
const std::vector< unsigned int > & get_point_ptrs () const
 
bool is_map_unique () const
 
bool all_points_found () const
 
bool point_found (const unsigned int i) const
 
const Triangulation< dim, spacedim > & get_triangulation () const
 
const Mapping< dim, spacedim > & get_mapping () const
 
bool is_ready () const
 

Private Attributes

const AdditionalData additional_data
 
boost::signals2::connection tria_signal
 
bool ready_flag
 
SmartPointer< const Triangulation< dim, spacedim > > tria
 
SmartPointer< const Mapping< dim, spacedim > > mapping
 
bool unique_mapping
 
bool all_points_found_flag
 
std::vector< unsigned intpoint_ptrs
 
std::vector< unsigned intrecv_permutation
 
std::vector< unsigned intrecv_ptrs
 
std::vector< unsigned intrecv_ranks
 
std::unique_ptr< CellDatacell_data
 
std::vector< unsigned intsend_permutation
 
std::vector< unsigned intsend_ranks
 
std::vector< unsigned intsend_ptrs
 

Detailed Description

template<int dim, int spacedim = dim>
class Utilities::MPI::RemotePointEvaluation< dim, spacedim >

Helper class to access values on non-matching grids.

Note
The name of the fields are chosen with the method evaluate_and_process() in mind. Here, quantities are computed at specified arbitrary positioned points (and even on remote processes in the MPI universe) cell by cell and these values are sent to requesting processes, which receive the result and resort the result according to the points.

Definition at line 52 of file mpi_remote_point_evaluation.h.

Constructor & Destructor Documentation

◆ RemotePointEvaluation() [1/2]

template<int dim, int spacedim>
Utilities::MPI::RemotePointEvaluation< dim, spacedim >::RemotePointEvaluation ( const AdditionalData additional_data = AdditionalData())

Constructor.

Parameters
additional_dataConfigure options for RemotePointEvaluation.

Definition at line 52 of file mpi_remote_point_evaluation.cc.

◆ RemotePointEvaluation() [2/2]

template<int dim, int spacedim>
Utilities::MPI::RemotePointEvaluation< dim, spacedim >::RemotePointEvaluation ( const double  tolerance,
const bool  enforce_unique_mapping = false,
const unsigned int  rtree_level = 0,
const std::function< std::vector< bool >()> &  marked_vertices = {} 
)

Constructor. This constructor is deprecated. Use the other constructor taking AdditionalData instead.

Parameters
toleranceTolerance in terms of unit cell coordinates for determining all cells around a point passed to the class during reinit(). Depending on the problem, it might be necessary to adjust the tolerance in order to be able to identify a cell. Floating point arithmetic implies that a point will, in general, not lie exactly on a vertex, edge, or face.
enforce_unique_mappingEnforce unique mapping, i.e., (one-to-one) relation of points and cells.
rtree_levelRTree level to be used during the construction of the bounding boxes.
marked_verticesFunction that marks relevant vertices to make search of active cells around point more efficient.
Deprecated:

Definition at line 61 of file mpi_remote_point_evaluation.cc.

◆ ~RemotePointEvaluation()

template<int dim, int spacedim>
Utilities::MPI::RemotePointEvaluation< dim, spacedim >::~RemotePointEvaluation ( )

Destructor.

Definition at line 76 of file mpi_remote_point_evaluation.cc.

Member Function Documentation

◆ reinit() [1/2]

template<int dim, int spacedim>
void Utilities::MPI::RemotePointEvaluation< dim, spacedim >::reinit ( const std::vector< Point< spacedim > > &  points,
const Triangulation< dim, spacedim > &  tria,
const Mapping< dim, spacedim > &  mapping 
)

Set up internal data structures and communication pattern based on a list of points points and mesh description (tria and mapping).

Warning
This is a collective call that needs to be executed by all processors in the communicator.
Note
If you want to be sure that all points have been found, call all_points_found() after calling this function.

Definition at line 86 of file mpi_remote_point_evaluation.cc.

◆ reinit() [2/2]

template<int dim, int spacedim>
void Utilities::MPI::RemotePointEvaluation< dim, spacedim >::reinit ( const GridTools::internal::DistributedComputePointLocationsInternal< dim, spacedim > &  data,
const Triangulation< dim, spacedim > &  tria,
const Mapping< dim, spacedim > &  mapping 
)

Set up internal data structures and communication pattern based on GridTools::internal::DistributedComputePointLocationsInternal.

This function is called internally by the reinit() function above. Having it as a separate function makes it possible to set up the class if it is known in which cells corresponding reference points are located (e.g. if intersections of cells are known).

Definition at line 137 of file mpi_remote_point_evaluation.cc.

◆ get_cell_data()

template<int dim, int spacedim>
const RemotePointEvaluation< dim, spacedim >::CellData & Utilities::MPI::RemotePointEvaluation< dim, spacedim >::get_cell_data ( ) const

Return internal data structure storing the data of points positioned in cells.

Definition at line 283 of file mpi_remote_point_evaluation.cc.

◆ evaluate_and_process() [1/2]

template<int dim, int spacedim>
template<typename T >
void Utilities::MPI::RemotePointEvaluation< dim, spacedim >::evaluate_and_process ( std::vector< T > &  output,
std::vector< T > &  buffer,
const std::function< void(const ArrayView< T > &, const CellData &)> &  evaluation_function 
) const

Evaluate function evaluation_function in the given points and triangulation. The result is stored in output.

Note
If the map of points to cells is not a one-to-one relation (is_map_unique()==false), the result needs to be processed with the help of get_point_ptrs(). This might be the case if a point coincides with a geometric entity (e.g., vertex) that is shared by multiple cells or a point is outside of the computational domain.
Warning
This is a collective call that needs to be executed by all processors in the communicator.

Definition at line 472 of file mpi_remote_point_evaluation.h.

◆ evaluate_and_process() [2/2]

template<int dim, int spacedim>
template<typename T >
std::vector< T > Utilities::MPI::RemotePointEvaluation< dim, spacedim >::evaluate_and_process ( const std::function< void(const ArrayView< T > &, const CellData &)> &  evaluation_function) const

Same as above but with the result provided as return value and without external allocation of a user-provided buffer.

Definition at line 631 of file mpi_remote_point_evaluation.h.

◆ process_and_evaluate() [1/2]

template<int dim, int spacedim>
template<typename T >
void Utilities::MPI::RemotePointEvaluation< dim, spacedim >::process_and_evaluate ( const std::vector< T > &  input,
std::vector< T > &  buffer,
const std::function< void(const ArrayView< const T > &, const CellData &)> &  evaluation_function 
) const

This method is the inverse of the method evaluate_and_process(). It makes the data at the points, provided by input, available in the function evaluation_function.

Warning
This is a collective call that needs to be executed by all processors in the communicator.

Definition at line 648 of file mpi_remote_point_evaluation.h.

◆ process_and_evaluate() [2/2]

template<int dim, int spacedim>
template<typename T >
void Utilities::MPI::RemotePointEvaluation< dim, spacedim >::process_and_evaluate ( const std::vector< T > &  input,
const std::function< void(const ArrayView< const T > &, const CellData &)> &  evaluation_function 
) const

Same as above but without external allocation of a user-provided buffer.

Definition at line 827 of file mpi_remote_point_evaluation.h.

◆ get_point_ptrs()

template<int dim, int spacedim>
const std::vector< unsigned int > & Utilities::MPI::RemotePointEvaluation< dim, spacedim >::get_point_ptrs ( ) const

Return a CRS-like data structure to determine the position of the result corresponding a point and the amount.

Definition at line 292 of file mpi_remote_point_evaluation.cc.

◆ is_map_unique()

template<int dim, int spacedim>
bool Utilities::MPI::RemotePointEvaluation< dim, spacedim >::is_map_unique ( ) const

Return if points and cells have a one-to-one relation. This is not the case if a point is not owned by any cell (the point is outside of the domain) or if multiple cells own the point (the point is positioned on a geometric entity shared by neighboring cells).

Definition at line 301 of file mpi_remote_point_evaluation.cc.

◆ all_points_found()

template<int dim, int spacedim>
bool Utilities::MPI::RemotePointEvaluation< dim, spacedim >::all_points_found ( ) const

Return if all points could be found in the domain.

Definition at line 310 of file mpi_remote_point_evaluation.cc.

◆ point_found()

template<int dim, int spacedim>
bool Utilities::MPI::RemotePointEvaluation< dim, spacedim >::point_found ( const unsigned int  i) const

Return if point i could be found in the domain.

Definition at line 319 of file mpi_remote_point_evaluation.cc.

◆ get_triangulation()

template<int dim, int spacedim>
const Triangulation< dim, spacedim > & Utilities::MPI::RemotePointEvaluation< dim, spacedim >::get_triangulation ( ) const

Return the Triangulation object used during reinit().

Definition at line 334 of file mpi_remote_point_evaluation.cc.

◆ get_mapping()

template<int dim, int spacedim>
const Mapping< dim, spacedim > & Utilities::MPI::RemotePointEvaluation< dim, spacedim >::get_mapping ( ) const

Return the Mapping object used during reinit().

Definition at line 343 of file mpi_remote_point_evaluation.cc.

◆ is_ready()

template<int dim, int spacedim>
bool Utilities::MPI::RemotePointEvaluation< dim, spacedim >::is_ready ( ) const

Return if the internal data structures have been set up and if yes whether they are still valid (and not invalidated due to changes of the Triangulation).

Definition at line 352 of file mpi_remote_point_evaluation.cc.

Member Data Documentation

◆ additional_data

template<int dim, int spacedim = dim>
const AdditionalData Utilities::MPI::RemotePointEvaluation< dim, spacedim >::additional_data
private

Additional data with basic settings.

Definition at line 373 of file mpi_remote_point_evaluation.h.

◆ tria_signal

template<int dim, int spacedim = dim>
boost::signals2::connection Utilities::MPI::RemotePointEvaluation< dim, spacedim >::tria_signal
private

Storage for the status of the triangulation signal.

Definition at line 378 of file mpi_remote_point_evaluation.h.

◆ ready_flag

template<int dim, int spacedim = dim>
bool Utilities::MPI::RemotePointEvaluation< dim, spacedim >::ready_flag
private

Flag indicating if the reinit() function has been called and if yes the triangulation has not been modified since then (potentially invalidating the communication pattern).

Definition at line 385 of file mpi_remote_point_evaluation.h.

◆ tria

template<int dim, int spacedim = dim>
SmartPointer<const Triangulation<dim, spacedim> > Utilities::MPI::RemotePointEvaluation< dim, spacedim >::tria
private

Reference to the Triangulation object used during reinit().

Definition at line 390 of file mpi_remote_point_evaluation.h.

◆ mapping

template<int dim, int spacedim = dim>
SmartPointer<const Mapping<dim, spacedim> > Utilities::MPI::RemotePointEvaluation< dim, spacedim >::mapping
private

Reference to the Mapping object used during reinit().

Definition at line 395 of file mpi_remote_point_evaluation.h.

◆ unique_mapping

template<int dim, int spacedim = dim>
bool Utilities::MPI::RemotePointEvaluation< dim, spacedim >::unique_mapping
private

(One-to-one) relation of points and cells.

Definition at line 400 of file mpi_remote_point_evaluation.h.

◆ all_points_found_flag

template<int dim, int spacedim = dim>
bool Utilities::MPI::RemotePointEvaluation< dim, spacedim >::all_points_found_flag
private

Cache if all points passed in during reinit() have been found.

Definition at line 405 of file mpi_remote_point_evaluation.h.

◆ point_ptrs

template<int dim, int spacedim = dim>
std::vector<unsigned int> Utilities::MPI::RemotePointEvaluation< dim, spacedim >::point_ptrs
private

Since for each point multiple or no results can be available, the pointers in this vector indicate the first and last entry associated with a point in a CRS-like fashion.

Definition at line 412 of file mpi_remote_point_evaluation.h.

◆ recv_permutation

template<int dim, int spacedim = dim>
std::vector<unsigned int> Utilities::MPI::RemotePointEvaluation< dim, spacedim >::recv_permutation
private

Permutation index within a recv buffer.

Definition at line 417 of file mpi_remote_point_evaluation.h.

◆ recv_ptrs

template<int dim, int spacedim = dim>
std::vector<unsigned int> Utilities::MPI::RemotePointEvaluation< dim, spacedim >::recv_ptrs
private

Pointers of ranges within a receive buffer that are filled by ranks specified by recv_ranks.

Definition at line 423 of file mpi_remote_point_evaluation.h.

◆ recv_ranks

template<int dim, int spacedim = dim>
std::vector<unsigned int> Utilities::MPI::RemotePointEvaluation< dim, spacedim >::recv_ranks
private

Ranks from where data is received.

Definition at line 428 of file mpi_remote_point_evaluation.h.

◆ cell_data

template<int dim, int spacedim = dim>
std::unique_ptr<CellData> Utilities::MPI::RemotePointEvaluation< dim, spacedim >::cell_data
private

Point data sorted according to cells so that evaluation (incl. reading of degrees of freedoms) needs to be performed only once per cell.

Definition at line 434 of file mpi_remote_point_evaluation.h.

◆ send_permutation

template<int dim, int spacedim = dim>
std::vector<unsigned int> Utilities::MPI::RemotePointEvaluation< dim, spacedim >::send_permutation
private

Permutation index within a send buffer.

Definition at line 439 of file mpi_remote_point_evaluation.h.

◆ send_ranks

template<int dim, int spacedim = dim>
std::vector<unsigned int> Utilities::MPI::RemotePointEvaluation< dim, spacedim >::send_ranks
private

Ranks to send to.

Definition at line 444 of file mpi_remote_point_evaluation.h.

◆ send_ptrs

template<int dim, int spacedim = dim>
std::vector<unsigned int> Utilities::MPI::RemotePointEvaluation< dim, spacedim >::send_ptrs
private

Pointers of ranges within a send buffer to be sent to the ranks specified by send_ranks.

Definition at line 450 of file mpi_remote_point_evaluation.h.


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