deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
|
#include <deal.II/base/mpi_remote_point_evaluation.h>
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::Cache< dim, spacedim > &cache, const std::vector< Point< spacedim > > &points) |
void | reinit (const GridTools::internal::DistributedComputePointLocationsInternal< dim, spacedim > &data, const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping) |
const CellData & | get_cell_data () const |
template<typename T , unsigned int n_components = 1> | |
void | evaluate_and_process (std::vector< T > &output, std::vector< T > &buffer, const std::function< void(const ArrayView< T > &, const CellData &)> &evaluation_function, const bool sort_data=true) const |
template<typename T , unsigned int n_components = 1> | |
std::vector< T > | evaluate_and_process (const std::function< void(const ArrayView< T > &, const CellData &)> &evaluation_function, const bool sort_data=true) const |
template<typename T , unsigned int n_components = 1> | |
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 bool sort_data=true) const |
template<typename T , unsigned int n_components = 1> | |
void | process_and_evaluate (const std::vector< T > &input, const std::function< void(const ArrayView< const T > &, const CellData &)> &evaluation_function, const bool sort_data=true) 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 |
const std::vector< unsigned int > & | get_send_permutation () const |
const std::vector< unsigned int > & | get_inverse_recv_permutation () const |
Private Attributes | |
const AdditionalData | additional_data |
boost::signals2::connection | tria_signal |
bool | ready_flag |
ObserverPointer< const Triangulation< dim, spacedim > > | tria |
ObserverPointer< const Mapping< dim, spacedim > > | mapping |
bool | unique_mapping |
bool | all_points_found_flag |
std::vector< unsigned int > | point_ptrs |
std::vector< unsigned int > | recv_permutation |
std::vector< unsigned int > | recv_permutation_inv |
std::vector< unsigned int > | recv_ptrs |
std::vector< unsigned int > | recv_ranks |
std::unique_ptr< CellData > | cell_data |
std::vector< unsigned int > | send_permutation |
std::vector< unsigned int > | send_permutation_inv |
std::vector< unsigned int > | send_ranks |
std::vector< unsigned int > | send_ptrs |
unsigned int | buffer_size_with_sorting |
unsigned int | buffer_size_without_sorting |
Helper class to access values on non-matching grids.
Definition at line 55 of file mpi_remote_point_evaluation.h.
Utilities::MPI::RemotePointEvaluation< dim, spacedim >::RemotePointEvaluation | ( | const AdditionalData & | additional_data = AdditionalData() | ) |
Constructor.
additional_data | Configure options for RemotePointEvaluation. |
Definition at line 52 of file mpi_remote_point_evaluation.cc.
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.
tolerance | Tolerance 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_mapping | Enforce unique mapping, i.e., (one-to-one) relation of points and cells. |
rtree_level | RTree level to be used during the construction of the bounding boxes. |
marked_vertices | Function that marks relevant vertices to make search of active cells around point more efficient. |
Definition at line 61 of file mpi_remote_point_evaluation.cc.
Utilities::MPI::RemotePointEvaluation< dim, spacedim >::~RemotePointEvaluation | ( | ) |
Destructor.
Definition at line 76 of file mpi_remote_point_evaluation.cc.
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
).
Definition at line 86 of file mpi_remote_point_evaluation.cc.
void Utilities::MPI::RemotePointEvaluation< dim, spacedim >::reinit | ( | const GridTools::Cache< dim, spacedim > & | cache, |
const std::vector< Point< spacedim > > & | points | ||
) |
Set up internal data structures and communication pattern based on an existing Cache and a list of points points
.
Definition at line 100 of file mpi_remote_point_evaluation.cc.
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 140 of file mpi_remote_point_evaluation.cc.
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 311 of file mpi_remote_point_evaluation.cc.
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 bool | sort_data = true |
||
) | const |
Evaluate function evaluation_function
in the given points and triangulation. The result is stored in output
.
n_components
. This allows to provide unrolled tensors, which is useful, e.g., if its dimension and its rank is not known at compile time. Definition at line 709 of file mpi_remote_point_evaluation.h.
std::vector< T > Utilities::MPI::RemotePointEvaluation< dim, spacedim >::evaluate_and_process | ( | const std::function< void(const ArrayView< T > &, const CellData &)> & | evaluation_function, |
const bool | sort_data = true |
||
) | const |
Same as above but with the result provided as return value and without external allocation of a user-provided buffer.
Definition at line 900 of file mpi_remote_point_evaluation.h.
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 bool | sort_data = true |
||
) | 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
.
n_components
. This allows to provide unrolled tensors, which is useful, e.g., if its dimension and its rank is not known at compile time. Definition at line 921 of file mpi_remote_point_evaluation.h.
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 bool | sort_data = true |
||
) | const |
Same as above but without external allocation of a user-provided buffer.
Definition at line 1124 of file mpi_remote_point_evaluation.h.
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 320 of file mpi_remote_point_evaluation.cc.
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 329 of file mpi_remote_point_evaluation.cc.
bool Utilities::MPI::RemotePointEvaluation< dim, spacedim >::all_points_found | ( | ) | const |
Return if all points could be found in the domain.
Definition at line 338 of file mpi_remote_point_evaluation.cc.
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 347 of file mpi_remote_point_evaluation.cc.
const Triangulation< dim, spacedim > & Utilities::MPI::RemotePointEvaluation< dim, spacedim >::get_triangulation | ( | ) | const |
Return the Triangulation object used during reinit().
Definition at line 362 of file mpi_remote_point_evaluation.cc.
const Mapping< dim, spacedim > & Utilities::MPI::RemotePointEvaluation< dim, spacedim >::get_mapping | ( | ) | const |
Return the Mapping object used during reinit().
Definition at line 371 of file mpi_remote_point_evaluation.cc.
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 380 of file mpi_remote_point_evaluation.cc.
const std::vector< unsigned int > & Utilities::MPI::RemotePointEvaluation< dim, spacedim >::get_send_permutation | ( | ) | const |
Return permutation needed for sending. If this is applied, evaluation results do not have to be sorted according to the receiving ranks during evaluate_and_process().
Definition at line 389 of file mpi_remote_point_evaluation.cc.
const std::vector< unsigned int > & Utilities::MPI::RemotePointEvaluation< dim, spacedim >::get_inverse_recv_permutation | ( | ) | const |
Data is contiguous rank by rank directly after MPI communication. If no sorting is requested during evaluate_and_process(), this permutation allows to access all values corresponding to an evaluation point (sorted into different cells).
Definition at line 398 of file mpi_remote_point_evaluation.cc.
|
private |
Additional data with basic settings.
Definition at line 424 of file mpi_remote_point_evaluation.h.
|
private |
Storage for the status of the triangulation signal.
Definition at line 429 of file mpi_remote_point_evaluation.h.
|
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 436 of file mpi_remote_point_evaluation.h.
|
private |
Reference to the Triangulation object used during reinit().
Definition at line 441 of file mpi_remote_point_evaluation.h.
|
private |
Reference to the Mapping object used during reinit().
Definition at line 446 of file mpi_remote_point_evaluation.h.
|
private |
(One-to-one) relation of points and cells.
Definition at line 451 of file mpi_remote_point_evaluation.h.
|
private |
Cache if all points passed in during reinit() have been found.
Definition at line 456 of file mpi_remote_point_evaluation.h.
|
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 463 of file mpi_remote_point_evaluation.h.
|
private |
Permutation index within a recv buffer.
Definition at line 468 of file mpi_remote_point_evaluation.h.
|
private |
Inverse of permutation index within a recv buffer.
Definition at line 473 of file mpi_remote_point_evaluation.h.
|
private |
Pointers of ranges within a receive buffer that are filled by ranks specified by recv_ranks.
Definition at line 479 of file mpi_remote_point_evaluation.h.
|
private |
Ranks from where data is received.
Definition at line 484 of file mpi_remote_point_evaluation.h.
|
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 490 of file mpi_remote_point_evaluation.h.
|
private |
Permutation index within a send buffer.
Definition at line 495 of file mpi_remote_point_evaluation.h.
|
private |
Inverse of permutation index within a send buffer.
Definition at line 500 of file mpi_remote_point_evaluation.h.
|
private |
Ranks to send to.
Definition at line 505 of file mpi_remote_point_evaluation.h.
|
private |
Pointers of ranges within a send buffer to be sent to the ranks specified by send_ranks.
Definition at line 511 of file mpi_remote_point_evaluation.h.
|
private |
Buffer size (if internal sorting is requested) determined as max of:
needed during evaluate_and_process()
needed during process_and_evaluate()
Definition at line 526 of file mpi_remote_point_evaluation.h.
|
private |
Buffer size (if sorting is not requested); corresponds to the number of evaluation points.
Definition at line 532 of file mpi_remote_point_evaluation.h.