Reference documentation for deal.II version 9.4.1
\(\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
partitioner.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2011 - 2022 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_partitioner_h
17#define dealii_partitioner_h
18
19#include <deal.II/base/config.h>
20
25#include <deal.II/base/mpi.h>
26#include <deal.II/base/types.h>
27
29
30#include <limits>
31
32
34
35namespace Utilities
36{
37 namespace MPI
38 {
192 {
193 public:
198
203 Partitioner(const unsigned int size);
204
220 const types::global_dof_index ghost_size,
221 const MPI_Comm & communicator);
222
230 Partitioner(const IndexSet &locally_owned_indices,
231 const IndexSet &ghost_indices_in,
232 const MPI_Comm &communicator_in);
233
241 Partitioner(const IndexSet &locally_owned_indices,
242 const MPI_Comm &communicator_in);
243
251 virtual void
252 reinit(const IndexSet &vector_space_vector_index_set,
253 const IndexSet &read_write_vector_index_set,
254 const MPI_Comm &communicator) override;
255
259 void
260 set_owned_indices(const IndexSet &locally_owned_indices);
261
272 void
274 const IndexSet &larger_ghost_index_set = IndexSet());
275
280 size() const;
281
292 unsigned int
293 local_size() const;
294
301 unsigned int
303
310 const IndexSet &
312
318 std::pair<types::global_dof_index, types::global_dof_index>
319 local_range() const;
320
325 bool
326 in_local_range(const types::global_dof_index global_index) const;
327
338 unsigned int
339 global_to_local(const types::global_dof_index global_index) const;
340
350 local_to_global(const unsigned int local_index) const;
351
357 bool
358 is_ghost_entry(const types::global_dof_index global_index) const;
359
363 const IndexSet &
365
370 unsigned int
372
383 const std::vector<std::pair<unsigned int, unsigned int>> &
385
391 const std::vector<std::pair<unsigned int, unsigned int>> &
393
402 const std::vector<std::pair<unsigned int, unsigned int>> &
404
409 unsigned int
411
421 const std::vector<std::pair<unsigned int, unsigned int>> &
423
434 bool
435 is_compatible(const Partitioner &part) const;
436
453 bool
455
460 unsigned int
462
467 unsigned int
469
473 virtual const MPI_Comm &
474 get_mpi_communicator() const override;
475
481 bool
483
484#ifdef DEAL_II_WITH_MPI
520 template <typename Number, typename MemorySpaceType = MemorySpace::Host>
521 void
523 const unsigned int communication_channel,
524 const ArrayView<const Number, MemorySpaceType> &locally_owned_array,
525 const ArrayView<Number, MemorySpaceType> & temporary_storage,
526 const ArrayView<Number, MemorySpaceType> & ghost_array,
527 std::vector<MPI_Request> & requests) const;
528
546 template <typename Number, typename MemorySpaceType = MemorySpace::Host>
547 void
549 const ArrayView<Number, MemorySpaceType> &ghost_array,
550 std::vector<MPI_Request> & requests) const;
551
589 template <typename Number, typename MemorySpaceType = MemorySpace::Host>
590 void
592 const VectorOperation::values vector_operation,
593 const unsigned int communication_channel,
594 const ArrayView<Number, MemorySpaceType> &ghost_array,
595 const ArrayView<Number, MemorySpaceType> &temporary_storage,
596 std::vector<MPI_Request> & requests) const;
597
632 template <typename Number, typename MemorySpaceType = MemorySpace::Host>
633 void
635 const VectorOperation::values vector_operation,
636 const ArrayView<const Number, MemorySpaceType> &temporary_storage,
637 const ArrayView<Number, MemorySpaceType> & locally_owned_storage,
638 const ArrayView<Number, MemorySpaceType> & ghost_array,
639 std::vector<MPI_Request> & requests) const;
640#endif
641
645 std::size_t
647
653 unsigned int,
654 << "Global index " << arg1
655 << " neither owned nor ghost on proc " << arg2 << '.');
656
661 unsigned int,
662 unsigned int,
663 unsigned int,
664 << "The size of the ghost index array (" << arg1
665 << ") must either equal the number of ghost in the "
666 << "partitioner (" << arg2
667 << ") or be equal in size to a more comprehensive index"
668 << "set which contains " << arg3
669 << " elements for this partitioner.");
670
671 private:
676 void
678
683
688
693 std::pair<types::global_dof_index, types::global_dof_index>
695
701
707
712 std::vector<std::pair<unsigned int, unsigned int>> ghost_targets_data;
713
720 std::vector<std::pair<unsigned int, unsigned int>> import_indices_data;
721
728 // The variable is mutable to enable lazy initialization in
729 // export_to_ghosted_array_start(). This way partitioner does not have to
730 // be templated on the MemorySpaceType.
731 mutable std::vector<
732 std::pair<std::unique_ptr<unsigned int[], void (*)(unsigned int *)>,
733 unsigned int>>
735
741
746 std::vector<std::pair<unsigned int, unsigned int>> import_targets_data;
747
752 std::vector<unsigned int> import_indices_chunks_by_rank_data;
753
759
765
771 std::vector<std::pair<unsigned int, unsigned int>>
773
777 unsigned int my_pid;
778
782 unsigned int n_procs;
783
788
793 };
794
795
796
797 /*--------------------- Inline functions --------------------------------*/
798
799#ifndef DOXYGEN
800
802 Partitioner::size() const
803 {
804 return global_size;
805 }
806
807
808
809 inline const IndexSet &
811 {
813 }
814
815
816
817 inline std::pair<types::global_dof_index, types::global_dof_index>
819 {
820 return local_range_data;
821 }
822
823
824
825 inline unsigned int
827 {
828 return locally_owned_size();
829 }
830
831
832
833 inline unsigned int
835 {
837 local_range_data.second - local_range_data.first;
838 Assert(size <= std::numeric_limits<unsigned int>::max(),
840 return static_cast<unsigned int>(size);
841 }
842
843
844
845 inline bool
847 const types::global_dof_index global_index) const
848 {
849 return (local_range_data.first <= global_index &&
850 global_index < local_range_data.second);
851 }
852
853
854
855 inline bool
857 const types::global_dof_index global_index) const
858 {
859 // if the index is in the global range, it is trivially not a ghost
860 if (in_local_range(global_index) == true)
861 return false;
862 else
863 return ghost_indices().is_element(global_index);
864 }
865
866
867
868 inline unsigned int
870 const types::global_dof_index global_index) const
871 {
872 Assert(in_local_range(global_index) || is_ghost_entry(global_index),
873 ExcIndexNotPresent(global_index, my_pid));
874 if (in_local_range(global_index))
875 return static_cast<unsigned int>(global_index - local_range_data.first);
876 else if (is_ghost_entry(global_index))
877 return (locally_owned_size() +
878 static_cast<unsigned int>(
879 ghost_indices_data.index_within_set(global_index)));
880 else
881 // should only end up here in optimized mode, when we use this large
882 // number to trigger a segfault when using this method for array
883 // access
885 }
886
887
888
890 Partitioner::local_to_global(const unsigned int local_index) const
891 {
892 AssertIndexRange(local_index,
894 if (local_index < locally_owned_size())
895 return local_range_data.first + types::global_dof_index(local_index);
896 else
897 return ghost_indices_data.nth_index_in_set(local_index -
899 }
900
901
902
903 inline const IndexSet &
905 {
906 return ghost_indices_data;
907 }
908
909
910
911 inline unsigned int
913 {
915 }
916
917
918
919 inline const std::vector<std::pair<unsigned int, unsigned int>> &
921 {
923 }
924
925
926
927 inline const std::vector<std::pair<unsigned int, unsigned int>> &
929 {
930 return ghost_targets_data;
931 }
932
933
934 inline const std::vector<std::pair<unsigned int, unsigned int>> &
936 {
937 return import_indices_data;
938 }
939
940
941
942 inline unsigned int
944 {
946 }
947
948
949
950 inline const std::vector<std::pair<unsigned int, unsigned int>> &
952 {
953 return import_targets_data;
954 }
955
956
957
958 inline unsigned int
960 {
961 // return the id from the variable stored in this class instead of
962 // Utilities::MPI::this_mpi_process() in order to make this query also
963 // work when MPI is not initialized.
964 return my_pid;
965 }
966
967
968
969 inline unsigned int
971 {
972 // return the number of MPI processes from the variable stored in this
973 // class instead of Utilities::MPI::n_mpi_processes() in order to make
974 // this query also work when MPI is not initialized.
975 return n_procs;
976 }
977
978
979
980 inline const MPI_Comm &
982 {
983 return communicator;
984 }
985
986
987
988 inline bool
990 {
991 return have_ghost_indices;
992 }
993
994#endif // ifndef DOXYGEN
995
996 } // end of namespace MPI
997
998} // end of namespace Utilities
999
1000
1002
1003#endif
size_type index_within_set(const size_type global_index) const
Definition: index_set.h:1923
bool is_element(const size_type index) const
Definition: index_set.h:1767
size_type nth_index_in_set(const size_type local_index) const
Definition: index_set.h:1882
unsigned int n_ghost_indices_data
Definition: partitioner.h:706
std::size_t memory_consumption() const
unsigned int global_to_local(const types::global_dof_index global_index) const
bool in_local_range(const types::global_dof_index global_index) const
std::vector< std::pair< unsigned int, unsigned int > > import_indices_data
Definition: partitioner.h:720
bool is_globally_compatible(const Partitioner &part) const
Partitioner(const IndexSet &locally_owned_indices, const IndexSet &ghost_indices_in, const MPI_Comm &communicator_in)
void set_owned_indices(const IndexSet &locally_owned_indices)
Partitioner(const types::global_dof_index local_size, const types::global_dof_index ghost_size, const MPI_Comm &communicator)
const IndexSet & locally_owned_range() const
unsigned int n_ghost_indices_in_larger_set
Definition: partitioner.h:758
unsigned int n_ghost_indices() const
std::vector< std::pair< unsigned int, unsigned int > > ghost_targets_data
Definition: partitioner.h:712
unsigned int local_size() const
virtual void reinit(const IndexSet &vector_space_vector_index_set, const IndexSet &read_write_vector_index_set, const MPI_Comm &communicator) override
const std::vector< std::pair< unsigned int, unsigned int > > & ghost_indices_within_larger_ghost_set() const
const IndexSet & ghost_indices() const
const std::vector< std::pair< unsigned int, unsigned int > > & ghost_targets() const
void export_to_ghosted_array_start(const unsigned int communication_channel, const ArrayView< const Number, MemorySpaceType > &locally_owned_array, const ArrayView< Number, MemorySpaceType > &temporary_storage, const ArrayView< Number, MemorySpaceType > &ghost_array, std::vector< MPI_Request > &requests) const
unsigned int this_mpi_process() const
Partitioner(const IndexSet &locally_owned_indices, const MPI_Comm &communicator_in)
const std::vector< std::pair< unsigned int, unsigned int > > & import_targets() const
std::vector< std::pair< std::unique_ptr< unsigned int[], void(*)(unsigned int *)>, unsigned int > > import_indices_plain_dev
Definition: partitioner.h:734
void import_from_ghosted_array_finish(const VectorOperation::values vector_operation, const ArrayView< const Number, MemorySpaceType > &temporary_storage, const ArrayView< Number, MemorySpaceType > &locally_owned_storage, const ArrayView< Number, MemorySpaceType > &ghost_array, std::vector< MPI_Request > &requests) const
unsigned int n_mpi_processes() const
void initialize_import_indices_plain_dev() const
Partitioner(const unsigned int size)
types::global_dof_index global_size
Definition: partitioner.h:682
unsigned int locally_owned_size() const
types::global_dof_index local_to_global(const unsigned int local_index) const
std::pair< types::global_dof_index, types::global_dof_index > local_range_data
Definition: partitioner.h:694
std::vector< unsigned int > import_indices_chunks_by_rank_data
Definition: partitioner.h:752
unsigned int n_import_indices_data
Definition: partitioner.h:740
void set_ghost_indices(const IndexSet &ghost_indices, const IndexSet &larger_ghost_index_set=IndexSet())
void export_to_ghosted_array_finish(const ArrayView< Number, MemorySpaceType > &ghost_array, std::vector< MPI_Request > &requests) const
std::pair< types::global_dof_index, types::global_dof_index > local_range() const
unsigned int n_import_indices() const
const std::vector< std::pair< unsigned int, unsigned int > > & import_indices() const
std::vector< unsigned int > ghost_indices_subset_chunks_by_rank_data
Definition: partitioner.h:764
void import_from_ghosted_array_start(const VectorOperation::values vector_operation, const unsigned int communication_channel, const ArrayView< Number, MemorySpaceType > &ghost_array, const ArrayView< Number, MemorySpaceType > &temporary_storage, std::vector< MPI_Request > &requests) const
virtual const MPI_Comm & get_mpi_communicator() const override
std::vector< std::pair< unsigned int, unsigned int > > import_targets_data
Definition: partitioner.h:746
types::global_dof_index size() const
bool is_ghost_entry(const types::global_dof_index global_index) const
bool is_compatible(const Partitioner &part) const
std::vector< std::pair< unsigned int, unsigned int > > ghost_indices_subset_data
Definition: partitioner.h:772
bool ghost_indices_initialized() const
#define DEAL_II_DEPRECATED
Definition: config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1473
static ::ExceptionBase & ExcIndexNotPresent(types::global_dof_index arg1, unsigned int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:532
static ::ExceptionBase & ExcGhostIndexArrayHasWrongSize(unsigned int arg1, unsigned int arg2, unsigned int arg3)
#define AssertIndexRange(index, range)
Definition: exceptions.h:1732
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:555
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
static const unsigned int invalid_unsigned_int
Definition: types.h:201
unsigned int global_dof_index
Definition: types.h:76