Reference documentation for deal.II version 9.0.0
partitioner.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2011 - 2018 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_partitioner_h
17 #define dealii_partitioner_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/mpi.h>
21 #include <deal.II/base/types.h>
22 #include <deal.II/base/utilities.h>
23 #include <deal.II/base/memory_consumption.h>
24 #include <deal.II/base/index_set.h>
25 #include <deal.II/base/array_view.h>
26 
27 #include <deal.II/lac/communication_pattern_base.h>
28 #include <deal.II/lac/vector.h>
29 #include <deal.II/lac/vector_operation.h>
30 
31 #include <limits>
32 
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 namespace Utilities
37 {
38  namespace MPI
39  {
125  {
126  public:
130  Partitioner ();
131 
136  Partitioner (const unsigned int size);
137 
145  Partitioner (const IndexSet &locally_owned_indices,
146  const IndexSet &ghost_indices_in,
147  const MPI_Comm communicator_in);
148 
156  Partitioner (const IndexSet &locally_owned_indices,
157  const MPI_Comm communicator_in);
158 
166  virtual void reinit(const IndexSet &vector_space_vector_index_set,
167  const IndexSet &read_write_vector_index_set,
168  const MPI_Comm &communicator);
169 
173  void set_owned_indices (const IndexSet &locally_owned_indices);
174 
186  const IndexSet &larger_ghost_index_set = IndexSet());
187 
192 
197  unsigned int local_size() const;
198 
205  const IndexSet &locally_owned_range() const;
206 
212  std::pair<types::global_dof_index,types::global_dof_index>
213  local_range() const;
214 
219  bool in_local_range (const types::global_dof_index global_index) const;
220 
230  unsigned int
231  global_to_local (const types::global_dof_index global_index) const;
232 
241  local_to_global (const unsigned int local_index) const;
242 
248  bool is_ghost_entry (const types::global_dof_index global_index) const;
249 
253  const IndexSet &ghost_indices() const;
254 
259  unsigned int n_ghost_indices() const;
260 
271  const std::vector<std::pair<unsigned int, unsigned int> > &
273 
279  const std::vector<std::pair<unsigned int, unsigned int> > &
280  ghost_targets() const;
281 
291  const std::vector<std::pair<unsigned int, unsigned int> > &
292  import_indices() const;
293 
298  unsigned int n_import_indices() const;
299 
309  const std::vector<std::pair<unsigned int, unsigned int> > &
310  import_targets() const;
311 
321  bool is_compatible (const Partitioner &part) const;
322 
336  bool is_globally_compatible (const Partitioner &part) const;
337 
342  unsigned int this_mpi_process () const;
343 
348  unsigned int n_mpi_processes () const;
349 
353  DEAL_II_DEPRECATED
354  const MPI_Comm &get_communicator() const;
355 
359  virtual const MPI_Comm &get_mpi_communicator() const;
360 
366  bool ghost_indices_initialized() const;
367 
368 #ifdef DEAL_II_WITH_MPI
369 
403  template <typename Number>
404  void
405  export_to_ghosted_array_start(const unsigned int communication_channel,
406  const ArrayView<const Number> &locally_owned_array,
407  const ArrayView<Number> &temporary_storage,
408  const ArrayView<Number> &ghost_array,
409  std::vector<MPI_Request> &requests) const;
410 
428  template <typename Number>
429  void
431  std::vector<MPI_Request> &requests) const;
432 
469  template <typename Number>
470  void
472  const unsigned int communication_channel,
473  const ArrayView<Number> &ghost_array,
474  const ArrayView<Number> &temporary_storage,
475  std::vector<MPI_Request> &requests) const;
476 
511  template <typename Number>
512  void
514  const ArrayView<const Number> &temporary_storage,
515  const ArrayView<Number> &locally_owned_storage,
516  const ArrayView<Number> &ghost_array,
517  std::vector<MPI_Request> &requests) const;
518 #endif
519 
523  std::size_t memory_consumption() const;
524 
530  unsigned int,
531  << "Global index " << arg1
532  << " neither owned nor ghost on proc " << arg2 << ".");
533 
538  unsigned int,
539  unsigned int,
540  unsigned int,
541  << "The size of the ghost index array (" << arg1
542  << ") must either equal the number of ghost in the "
543  << "partitioner (" << arg2
544  << ") or be equal in size to a more comprehensive index"
545  << "set which contains " << arg3
546  << " elements for this partitioner.");
547 
548  private:
553 
558 
563  std::pair<types::global_dof_index,types::global_dof_index> local_range_data;
564 
570 
575  unsigned int n_ghost_indices_data;
576 
581  std::vector<std::pair<unsigned int, unsigned int> > ghost_targets_data;
582 
589  std::vector<std::pair<unsigned int, unsigned int> > import_indices_data;
590 
595  unsigned int n_import_indices_data;
596 
601  std::vector<std::pair<unsigned int, unsigned int> > import_targets_data;
602 
607  std::vector<unsigned int> import_indices_chunks_by_rank_data;
608 
614 
619  std::vector<unsigned int> ghost_indices_subset_chunks_by_rank_data;
620 
626  std::vector<std::pair<unsigned int, unsigned int> > ghost_indices_subset_data;
627 
631  unsigned int my_pid;
632 
636  unsigned int n_procs;
637 
641  MPI_Comm communicator;
642 
647  };
648 
649 
650 
651  /*----------------------- Inline functions ----------------------------------*/
652 
653 #ifndef DOXYGEN
654 
655  inline
657  {
658  return global_size;
659  }
660 
661 
662 
663  inline
665  {
667  }
668 
669 
670 
671  inline
672  std::pair<types::global_dof_index,types::global_dof_index>
674  {
675  return local_range_data;
676  }
677 
678 
679 
680  inline
681  unsigned int
682  Partitioner::local_size () const
683  {
685  Assert(size<=std::numeric_limits<unsigned int>::max(),
687  return static_cast<unsigned int>(size);
688  }
689 
690 
691 
692  inline
693  bool
694  Partitioner::in_local_range (const types::global_dof_index global_index) const
695  {
696  return (local_range_data.first <= global_index &&
697  global_index < local_range_data.second);
698  }
699 
700 
701 
702  inline
703  bool
704  Partitioner::is_ghost_entry (const types::global_dof_index global_index) const
705  {
706  // if the index is in the global range, it is trivially not a ghost
707  if (in_local_range(global_index) == true)
708  return false;
709  else
710  return ghost_indices().is_element(global_index);
711  }
712 
713 
714 
715  inline
716  unsigned int
717  Partitioner::global_to_local (const types::global_dof_index global_index) const
718  {
719  Assert(in_local_range(global_index) || is_ghost_entry (global_index),
720  ExcIndexNotPresent(global_index, my_pid));
721  if (in_local_range(global_index))
722  return static_cast<unsigned int>(global_index - local_range_data.first);
723  else if (is_ghost_entry (global_index))
724  return (local_size() +
725  static_cast<unsigned int>(ghost_indices_data.index_within_set (global_index)));
726  else
727  // should only end up here in optimized mode, when we use this large
728  // number to trigger a segfault when using this method for array
729  // access
731  }
732 
733 
734 
735  inline
737  Partitioner::local_to_global (const unsigned int local_index) const
738  {
740  if (local_index < local_size())
741  return local_range_data.first + types::global_dof_index(local_index);
742  else
743  return ghost_indices_data.nth_index_in_set (local_index-local_size());
744  }
745 
746 
747 
748  inline
749  const IndexSet &Partitioner::ghost_indices() const
750  {
751  return ghost_indices_data;
752  }
753 
754 
755 
756  inline
757  unsigned int
759  {
760  return n_ghost_indices_data;
761  }
762 
763 
764 
765  inline
766  const std::vector<std::pair<unsigned int, unsigned int> > &
768  {
770  }
771 
772 
773 
774  inline
775  const std::vector<std::pair<unsigned int, unsigned int> > &
777  {
778  return ghost_targets_data;
779  }
780 
781 
782  inline
783  const std::vector<std::pair<unsigned int, unsigned int> > &
785  {
786  return import_indices_data;
787  }
788 
789 
790 
791  inline
792  unsigned int
794  {
795  return n_import_indices_data;
796  }
797 
798 
799 
800  inline
801  const std::vector<std::pair<unsigned int, unsigned int> > &
803  {
804  return import_targets_data;
805  }
806 
807 
808 
809  inline
810  unsigned int
812  {
813  // return the id from the variable stored in this class instead of
814  // Utilities::MPI::this_mpi_process() in order to make this query also
815  // work when MPI is not initialized.
816  return my_pid;
817  }
818 
819 
820 
821  inline
822  unsigned int
824  {
825  // return the number of MPI processes from the variable stored in this
826  // class instead of Utilities::MPI::n_mpi_processes() in order to make
827  // this query also work when MPI is not initialized.
828  return n_procs;
829  }
830 
831 
832 
833  inline
834  const MPI_Comm &
836  {
837  return communicator;
838  }
839 
840 
841 
842  inline
843  const MPI_Comm &
845  {
846  return communicator;
847  }
848 
849 
850 
851  inline
852  bool
854  {
855  return have_ghost_indices;
856  }
857 
858 #endif // ifndef DOXYGEN
859 
860  } // end of namespace MPI
861 
862 } // end of namespace Utilities
863 
864 
865 DEAL_II_NAMESPACE_CLOSE
866 
867 #endif
std::vector< std::pair< unsigned int, unsigned int > > import_indices_data
Definition: partitioner.h:589
const IndexSet & locally_owned_range() const
bool is_globally_compatible(const Partitioner &part) const
Definition: partitioner.cc:454
static const unsigned int invalid_unsigned_int
Definition: types.h:173
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:358
virtual void reinit(const IndexSet &vector_space_vector_index_set, const IndexSet &read_write_vector_index_set, const MPI_Comm &communicator)
Definition: partitioner.cc:98
bool is_compatible(const Partitioner &part) const
Definition: partitioner.cc:428
static ::ExceptionBase & ExcIndexNotPresent(types::global_dof_index arg1, unsigned int arg2)
unsigned int local_size() const
const std::vector< std::pair< unsigned int, unsigned int > > & import_targets() const
size_type nth_index_in_set(const unsigned int local_index) const
Definition: index_set.h:1769
const std::vector< std::pair< unsigned int, unsigned int > > & import_indices() const
virtual const MPI_Comm & get_mpi_communicator() const
types::global_dof_index global_size
Definition: partitioner.h:552
types::global_dof_index local_to_global(const unsigned int local_index) const
#define AssertIndexRange(index, range)
Definition: exceptions.h:1284
unsigned int n_ghost_indices_data
Definition: partitioner.h:575
const std::vector< std::pair< unsigned int, unsigned int > > & ghost_indices_within_larger_ghost_set() 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
void import_from_ghosted_array_start(const VectorOperation::values vector_operation, const unsigned int communication_channel, const ArrayView< Number > &ghost_array, const ArrayView< Number > &temporary_storage, std::vector< MPI_Request > &requests) const
void set_owned_indices(const IndexSet &locally_owned_indices)
Definition: partitioner.cc:111
std::vector< unsigned int > ghost_indices_subset_chunks_by_rank_data
Definition: partitioner.h:619
std::vector< unsigned int > import_indices_chunks_by_rank_data
Definition: partitioner.h:607
void export_to_ghosted_array_finish(const ArrayView< Number > &ghost_array, std::vector< MPI_Request > &requests) const
bool ghost_indices_initialized() const
unsigned int this_mpi_process() const
const IndexSet & ghost_indices() const
const MPI_Comm & get_communicator() const
unsigned int n_ghost_indices() const
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:1142
size_type index_within_set(const size_type global_index) const
Definition: index_set.h:1811
unsigned int n_mpi_processes() const
std::vector< std::pair< unsigned int, unsigned int > > ghost_indices_subset_data
Definition: partitioner.h:626
void export_to_ghosted_array_start(const unsigned int communication_channel, const ArrayView< const Number > &locally_owned_array, const ArrayView< Number > &temporary_storage, const ArrayView< Number > &ghost_array, std::vector< MPI_Request > &requests) const
const std::vector< std::pair< unsigned int, unsigned int > > & ghost_targets() const
Definition: cuda.h:31
types::global_dof_index size() const
std::size_t memory_consumption() const
Definition: partitioner.cc:463
void import_from_ghosted_array_finish(const VectorOperation::values vector_operation, const ArrayView< const Number > &temporary_storage, const ArrayView< Number > &locally_owned_storage, const ArrayView< Number > &ghost_array, std::vector< MPI_Request > &requests) const
std::vector< std::pair< unsigned int, unsigned int > > import_targets_data
Definition: partitioner.h:601
std::vector< std::pair< unsigned int, unsigned int > > ghost_targets_data
Definition: partitioner.h:581
static ::ExceptionBase & ExcNotImplemented()
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
bool is_element(const size_type index) const
Definition: index_set.h:1645
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:370
std::pair< types::global_dof_index, types::global_dof_index > local_range_data
Definition: partitioner.h:563
std::pair< types::global_dof_index, types::global_dof_index > local_range() const
unsigned int n_import_indices() const
unsigned int n_import_indices_data
Definition: partitioner.h:595
bool is_ghost_entry(const types::global_dof_index global_index) const
static ::ExceptionBase & ExcGhostIndexArrayHasWrongSize(unsigned int arg1, unsigned int arg2, unsigned int arg3)
unsigned int n_ghost_indices_in_larger_set
Definition: partitioner.h:613
void set_ghost_indices(const IndexSet &ghost_indices, const IndexSet &larger_ghost_index_set=IndexSet())
Definition: partitioner.cc:147