Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_DistObject_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Tpetra: Templated Linear Algebra Services Package
4 //
5 // Copyright 2008 NTESS and the Tpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef TPETRA_DISTOBJECT_DECL_HPP
11 #define TPETRA_DISTOBJECT_DECL_HPP
12 
15 
16 #include "Tpetra_Details_DistributorActor.hpp"
17 #include "Tpetra_Map.hpp"
18 #include "Tpetra_Import.hpp"
19 #include "Tpetra_Export.hpp"
20 #include "Tpetra_SrcDistObject.hpp"
22 #if KOKKOS_VERSION >= 40799
23 #include "KokkosKernels_ArithTraits.hpp"
24 #else
25 #include "Kokkos_ArithTraits.hpp"
26 #endif
27 #include <memory>
28 #include <type_traits>
29 
30 namespace KokkosClassic {
35 enum ReadWriteOption {
36  ReadWrite = 0,
37  OverwriteAll = 1
38 };
39 } // namespace KokkosClassic
40 
41 namespace Tpetra {
42 
125 template <class DistObjectType>
126 void removeEmptyProcessesInPlace(Teuchos::RCP<DistObjectType>& input,
127  const Teuchos::RCP<const Map<typename DistObjectType::local_ordinal_type,
128  typename DistObjectType::global_ordinal_type,
129  typename DistObjectType::node_type> >& newMap);
130 
167 template <class DistObjectType>
168 void removeEmptyProcessesInPlace(Teuchos::RCP<DistObjectType>& input);
169 
281 template <class Packet,
282  class LocalOrdinal,
283  class GlobalOrdinal,
284  class Node>
285 class DistObject : virtual public SrcDistObject,
286  virtual public Teuchos::Describable {
287  public:
289 
290 
295 #if KOKKOS_VERSION >= 40799
296  using packet_type = typename ::KokkosKernels::ArithTraits<Packet>::val_type;
297 #else
298  using packet_type = typename ::Kokkos::ArithTraits<Packet>::val_type;
299 #endif
300  using local_ordinal_type = LocalOrdinal;
303  using global_ordinal_type = GlobalOrdinal;
305  using node_type = Node;
306 
308  using device_type = typename Node::device_type;
310  using execution_space = typename device_type::execution_space;
311 
314 
316 
318 
322  explicit DistObject(const Teuchos::RCP<const map_type>& map);
323 
326 
329 
332 
335 
345  virtual ~DistObject() = default;
346 
348 
350 
378  void
379  doImport(const SrcDistObject& source,
381  const CombineMode CM,
382  const bool restrictedMode = false);
383 
411  void
412  doExport(const SrcDistObject& source,
414  const CombineMode CM,
415  const bool restrictedMode = false);
416 
445  void
446  doImport(const SrcDistObject& source,
448  const CombineMode CM,
449  const bool restrictedMode = false);
450 
479  void
480  doExport(const SrcDistObject& source,
482  const CombineMode CM,
483  const bool restrictedMode = false);
484 
485  void
486  beginImport(const SrcDistObject& source,
488  const CombineMode CM,
489  const bool restrictedMode = false);
490 
491  void
492  beginExport(const SrcDistObject& source,
494  const CombineMode CM,
495  const bool restrictedMode = false);
496 
497  void
498  beginImport(const SrcDistObject& source,
500  const CombineMode CM,
501  const bool restrictedMode = false);
502 
503  void
504  beginExport(const SrcDistObject& source,
506  const CombineMode CM,
507  const bool restrictedMode = false);
508 
509  void
510  endImport(const SrcDistObject& source,
512  const CombineMode CM,
513  const bool restrictedMode = false);
514 
515  void
516  endExport(const SrcDistObject& source,
518  const CombineMode CM,
519  const bool restrictedMode = false);
520 
521  void
522  endImport(const SrcDistObject& source,
524  const CombineMode CM,
525  const bool restrictedMode = false);
526 
527  void
528  endExport(const SrcDistObject& source,
530  const CombineMode CM,
531  const bool restrictedMode = false);
532 
535  bool transferArrived() const;
536 
538 
540 
546  bool isDistributed() const;
547 
554  virtual Teuchos::RCP<const map_type> getMap() const { return map_; }
555 
557 
559 
564  void print(std::ostream& os) const;
565 
567 
569 
574  virtual std::string description() const;
575 
580  virtual void
581  describe(Teuchos::FancyOStream& out,
582  const Teuchos::EVerbosityLevel verbLevel =
583  Teuchos::Describable::verbLevel_default) const;
584 
586 
588 
635  virtual void
636  removeEmptyProcessesInPlace(const Teuchos::RCP<const map_type>& newMap);
637 
639  const Details::DistributorActor& getActor() { return distributorActor_; }
640 
642 
643  protected:
653  DoForward, //*!< Perform the transfer in forward mode.
654  DoReverse //*!< Perform the transfer in reverse mode.
655  };
656 
673  virtual size_t constantNumberOfPackets() const;
674 
694  virtual void
695  doTransfer(const SrcDistObject& src,
696  const ::Tpetra::Details::Transfer<local_ordinal_type, global_ordinal_type, node_type>& transfer,
697  const char modeString[],
698  const ReverseOption revOp,
699  const CombineMode CM,
700  const bool restrictedMode);
701 
716  virtual bool
717  reallocArraysForNumPacketsPerLid(const size_t numExportLIDs,
718  const size_t numImportLIDs);
719 
722  using buffer_memory_space =
723  ::Tpetra::Details::DefaultTypes::comm_buffer_memory_space<device_type>;
724 
725  public:
736  using buffer_device_type =
737  Kokkos::Device<typename device_type::execution_space,
739 
740  protected:
746  void beginTransfer(const SrcDistObject& src,
747  const ::Tpetra::Details::Transfer<local_ordinal_type, global_ordinal_type, node_type>& transfer,
748  const char modeString[],
749  const ReverseOption revOp,
750  const CombineMode CM,
751  const bool restrictedMode);
752 
753  void endTransfer(const SrcDistObject& src,
754  const ::Tpetra::Details::Transfer<local_ordinal_type, global_ordinal_type, node_type>& transfer,
755  const char modeString[],
756  const ReverseOption revOp,
757  const CombineMode CM,
758  const bool restrictedMode);
759 
760  void doPostRecvs(const Details::DistributorPlan& distributorPlan,
761  size_t constantNumPackets,
762  bool commOnHost,
763  std::shared_ptr<std::string> prefix,
764  const bool canTryAliasing,
765  const CombineMode CM);
766 
767  void doPostSends(const Details::DistributorPlan& distributorPlan,
768  size_t constantNumPackets,
769  bool commOnHost,
770  std::shared_ptr<std::string> prefix);
771 
772  void doPackAndPrepare(const SrcDistObject& src,
773  const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>& exportLIDs,
774  size_t& constantNumPackets,
775  const execution_space& space);
776 
777  void doUnpackAndCombine(const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>& remoteLIDs,
778  size_t constantNumPackets,
779  CombineMode CM,
780  const execution_space& space);
781 
792 
793 
797  virtual bool
798  checkSizes(const SrcDistObject& source) = 0;
799 
829  virtual void
830  copyAndPermute(const SrcDistObject& source,
831  const size_t numSameIDs,
832  const Kokkos::DualView<const local_ordinal_type*,
833  buffer_device_type>& permuteToLIDs,
834  const Kokkos::DualView<const local_ordinal_type*,
835  buffer_device_type>& permuteFromLIDs,
836  const CombineMode CM);
837 
840  virtual void copyAndPermute(
841  const SrcDistObject& source, const size_t numSameIDs,
842  const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>& permuteToLIDs,
843  const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>& permuteFromLIDs,
844  const CombineMode CM, const execution_space& space);
845 
883  virtual void
884  packAndPrepare(const SrcDistObject& source,
885  const Kokkos::DualView<const local_ordinal_type*,
886  buffer_device_type>& exportLIDs,
887  Kokkos::DualView<packet_type*,
888  buffer_device_type>& exports,
889  Kokkos::DualView<size_t*,
891  numPacketsPerLID,
892  size_t& constantNumPackets);
893 
896  virtual void
897  packAndPrepare(const SrcDistObject& source,
898  const Kokkos::DualView<const local_ordinal_type*,
899  buffer_device_type>& exportLIDs,
900  Kokkos::DualView<packet_type*,
901  buffer_device_type>& exports,
902  Kokkos::DualView<size_t*,
904  numPacketsPerLID,
905  size_t& constantNumPackets,
906  const execution_space& space);
907 
947  virtual void
948  unpackAndCombine(const Kokkos::DualView<const local_ordinal_type*,
949  buffer_device_type>& importLIDs,
950  Kokkos::DualView<packet_type*,
952  imports,
953  Kokkos::DualView<size_t*,
955  numPacketsPerLID,
956  const size_t constantNumPackets,
957  const CombineMode combineMode);
958 
959  virtual void
960  unpackAndCombine(const Kokkos::DualView<const local_ordinal_type*,
961  buffer_device_type>& importLIDs,
962  Kokkos::DualView<packet_type*,
964  imports,
965  Kokkos::DualView<size_t*,
967  numPacketsPerLID,
968  const size_t constantNumPackets,
969  const CombineMode combineMode,
970  const execution_space& space);
971 
973  Teuchos::RCP<const map_type> map_;
974 
975  protected:
976  std::unique_ptr<std::string>
977  createPrefix(const char className[],
978  const char methodName[]) const;
979 
986  Kokkos::DualView<packet_type*, buffer_device_type> imports_;
987 
1006  virtual bool
1007  reallocImportsIfNeeded(const size_t newSize,
1008  const bool verbose,
1009  const std::string* prefix,
1010  const bool remoteLIDsContiguous = false,
1011  const CombineMode CM = INSERT);
1012 
1026  Kokkos::DualView<size_t*, buffer_device_type> numImportPacketsPerLID_;
1027 
1033  Kokkos::DualView<packet_type*, buffer_device_type> exports_;
1034 
1048  Kokkos::DualView<size_t*, buffer_device_type> numExportPacketsPerLID_;
1049 
1050  private:
1052 
1053  Details::DistributorActor distributorActor_;
1054 
1055 }; // class DistObject
1056 } // namespace Tpetra
1057 
1058 #endif // TPETRA_DISTOBJECT_DECL_HPP
Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distr...
Node node_type
The Node type. If you don&#39;t know what this is, don&#39;t use it.
Kokkos::DualView< size_t *, buffer_device_type > numExportPacketsPerLID_
Number of packets to send for each send operation.
virtual void copyAndPermute(const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM)
Perform copies and permutations that are local to the calling (MPI) process.
void doImport(const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
Import data into this object using an Import object (&quot;forward mode&quot;).
typename::Kokkos::ArithTraits< Packet >::val_type packet_type
The type of each datum being sent or received in an Import or Export.
void print(std::ostream &os) const
Print this object to the given output stream.
virtual bool reallocArraysForNumPacketsPerLid(const size_t numExportLIDs, const size_t numImportLIDs)
Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary.
Abstract base class for sources of an Import or Export.
bool isDistributed() const
Whether this is a globally distributed object.
Forward declaration of Tpetra::DistObject.
void removeEmptyProcessesInPlace(Teuchos::RCP< DistObjectType > &input, const Teuchos::RCP< const Map< typename DistObjectType::local_ordinal_type, typename DistObjectType::global_ordinal_type, typename DistObjectType::node_type > > &newMap)
Remove processes which contain no elements in this object&#39;s Map.
virtual void unpackAndCombine(const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode)
Perform any unpacking and combining after communication.
virtual void doTransfer(const SrcDistObject &src, const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &transfer, const char modeString[], const ReverseOption revOp, const CombineMode CM, const bool restrictedMode)
Redistribute data across (MPI) processes.
void beginTransfer(const SrcDistObject &src, const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &transfer, const char modeString[], const ReverseOption revOp, const CombineMode CM, const bool restrictedMode)
Implementation detail of doTransfer.
virtual ~DistObject()=default
Destructor (virtual for memory safety of derived classes).
Kokkos::DualView< size_t *, buffer_device_type > numImportPacketsPerLID_
Number of packets to receive for each receive operation.
Insert new values that don&#39;t currently exist.
typename device_type::execution_space execution_space
The Kokkos execution space.
Kokkos::DualView< packet_type *, buffer_device_type > imports_
Buffer into which packed data are imported (received from other processes).
Kokkos::DualView< packet_type *, buffer_device_type > exports_
Buffer from which packed data are exported (sent to other processes).
::Tpetra::Details::DefaultTypes::comm_buffer_memory_space< device_type > buffer_memory_space
Kokkos memory space for communication buffers.
const Details::DistributorActor & getActor()
Getter for the DistributorActor.
Communication plan for data redistribution from a (possibly) multiply-owned to a uniquely-owned distr...
virtual void packAndPrepare(const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets)
Pack data and metadata for communication (sends).
CombineMode
Rule for combining data in an Import or Export.
Teuchos::RCP< const map_type > map_
The Map over which this object is distributed.
Abstract base class for objects that can be the source of an Import or Export operation.
DistObject & operator=(const DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &)=default
Assignment operator (default).
DistObject(const Teuchos::RCP< const map_type > &map)
Constructor.
virtual std::string description() const
One-line descriptiion of this object.
bool transferArrived() const
Whether the data from an import/export operation has arrived, and is ready for the unpack and combine...
virtual size_t constantNumberOfPackets() const
Whether the implementation&#39;s instance promises always to have a constant number of packets per LID (l...
virtual bool reallocImportsIfNeeded(const size_t newSize, const bool verbose, const std::string *prefix, const bool remoteLIDsContiguous=false, const CombineMode CM=INSERT)
Reallocate imports_ if needed.
A parallel distribution of indices over processes.
void doExport(const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
Export data into this object using an Export object (&quot;forward mode&quot;).
virtual Teuchos::RCP< const map_type > getMap() const
The Map describing the parallel distribution of this object.
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print a descriptiion of this object to the given output stream.
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type
Base class for distributed Tpetra objects that support data redistribution.
virtual void removeEmptyProcessesInPlace(const Teuchos::RCP< const map_type > &newMap)
Remove processes which contain no entries in this object&#39;s Map.
virtual bool checkSizes(const SrcDistObject &source)=0
Compare the source and target (this) objects for compatibility.