Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
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"
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
30namespace KokkosClassic {
36 ReadWrite = 0,
37 OverwriteAll = 1
38};
39} // namespace KokkosClassic
40
41namespace Tpetra {
42
125template <class DistObjectType>
126void 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
167template <class DistObjectType>
168void removeEmptyProcessesInPlace(Teuchos::RCP<DistObjectType>& input);
169
281template <class Packet,
282 class LocalOrdinal,
283 class GlobalOrdinal,
284 class Node>
285class 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
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
381 const CombineMode CM,
382 const bool restrictedMode = false);
383
411 void
414 const CombineMode CM,
415 const bool restrictedMode = false);
416
445 void
448 const CombineMode CM,
449 const bool restrictedMode = false);
450
479 void
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
640 protected:
650 DoForward, //*!< Perform the transfer in forward mode.
651 DoReverse //*!< Perform the transfer in reverse mode.
652 };
653
670 virtual size_t constantNumberOfPackets() const;
671
691 virtual void
692 doTransfer(const SrcDistObject& src,
693 const ::Tpetra::Details::Transfer<local_ordinal_type, global_ordinal_type, node_type>& transfer,
694 const char modeString[],
695 const ReverseOption revOp,
696 const CombineMode CM,
697 const bool restrictedMode);
698
713 virtual bool
715 const size_t numImportLIDs);
716
720 ::Tpetra::Details::DefaultTypes::comm_buffer_memory_space<device_type>;
721
722 public:
734 Kokkos::Device<typename device_type::execution_space,
736
737 protected:
743 void beginTransfer(const SrcDistObject& src,
744 const ::Tpetra::Details::Transfer<local_ordinal_type, global_ordinal_type, node_type>& transfer,
745 const char modeString[],
746 const ReverseOption revOp,
747 const CombineMode CM,
748 const bool restrictedMode);
749
750 void endTransfer(const SrcDistObject& src,
751 const ::Tpetra::Details::Transfer<local_ordinal_type, global_ordinal_type, node_type>& transfer,
752 const char modeString[],
753 const ReverseOption revOp,
754 const CombineMode CM,
755 const bool restrictedMode);
756
757 void doPostRecvs(const Details::DistributorPlan& distributorPlan,
758 size_t constantNumPackets,
759 bool commOnHost,
760 std::shared_ptr<std::string> prefix,
761 const bool canTryAliasing,
762 const CombineMode CM);
763
764 void doPostSends(const Details::DistributorPlan& distributorPlan,
765 size_t constantNumPackets,
766 bool commOnHost,
767 std::shared_ptr<std::string> prefix);
768
769 void doPackAndPrepare(const SrcDistObject& src,
770 const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>& exportLIDs,
771 size_t& constantNumPackets,
772 const execution_space& space);
773
774 void doUnpackAndCombine(const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>& remoteLIDs,
775 size_t constantNumPackets,
777 const execution_space& space);
778
789
790
794 virtual bool
796
826 virtual void
828 const size_t numSameIDs,
829 const Kokkos::DualView<const local_ordinal_type*,
831 const Kokkos::DualView<const local_ordinal_type*,
833 const CombineMode CM);
834
837 virtual void copyAndPermute(
838 const SrcDistObject& source, const size_t numSameIDs,
839 const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>& permuteToLIDs,
840 const Kokkos::DualView<const local_ordinal_type*, buffer_device_type>& permuteFromLIDs,
841 const CombineMode CM, const execution_space& space);
842
880 virtual void
882 const Kokkos::DualView<const local_ordinal_type*,
884 Kokkos::DualView<packet_type*,
885 buffer_device_type>& exports,
886 Kokkos::DualView<size_t*,
889 size_t& constantNumPackets);
890
893 virtual void
895 const Kokkos::DualView<const local_ordinal_type*,
897 Kokkos::DualView<packet_type*,
898 buffer_device_type>& exports,
899 Kokkos::DualView<size_t*,
902 size_t& constantNumPackets,
903 const execution_space& space);
904
944 virtual void
945 unpackAndCombine(const Kokkos::DualView<const local_ordinal_type*,
947 Kokkos::DualView<packet_type*,
949 imports,
950 Kokkos::DualView<size_t*,
953 const size_t constantNumPackets,
955
956 virtual void
957 unpackAndCombine(const Kokkos::DualView<const local_ordinal_type*,
959 Kokkos::DualView<packet_type*,
961 imports,
962 Kokkos::DualView<size_t*,
965 const size_t constantNumPackets,
967 const execution_space& space);
968
970 Teuchos::RCP<const map_type> map_;
971
972 protected:
973 std::unique_ptr<std::string>
974 createPrefix(const char className[],
975 const char methodName[]) const;
976
983 Kokkos::DualView<packet_type*, buffer_device_type> imports_;
984
1003 virtual bool
1004 reallocImportsIfNeeded(const size_t newSize,
1005 const bool verbose,
1006 const std::string* prefix,
1007 const bool remoteLIDsContiguous = false,
1008 const CombineMode CM = INSERT);
1009
1023 Kokkos::DualView<size_t*, buffer_device_type> numImportPacketsPerLID_;
1024
1030 Kokkos::DualView<packet_type*, buffer_device_type> exports_;
1031
1045 Kokkos::DualView<size_t*, buffer_device_type> numExportPacketsPerLID_;
1046
1047 private:
1049
1050 Details::DistributorActor distributorActor_;
1051
1052}; // class DistObject
1053} // namespace Tpetra
1054
1055#endif // TPETRA_DISTOBJECT_DECL_HPP
ReadWriteOption
Read/write options for non-const views.
Forward declaration of Tpetra::DistObject.
Abstract base class for sources of an Import or Export.
Struct that holds views of the contents of a CrsMatrix.
Base class for distributed Tpetra objects that support data redistribution.
DistObject(const DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &)=default
Copy constructor (default).
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.
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.
::Tpetra::Details::DefaultTypes::comm_buffer_memory_space< device_type > buffer_memory_space
Kokkos memory space for communication buffers.
DistObject(DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &&)=default
Move constructor (default).
Kokkos::DualView< size_t *, buffer_device_type > numImportPacketsPerLID_
Number of packets to receive for each receive operation.
Kokkos::DualView< packet_type *, buffer_device_type > exports_
Buffer from which packed data are exported (sent to other processes).
Kokkos::DualView< packet_type *, buffer_device_type > imports_
Buffer into which packed data are imported (received from other processes).
virtual bool reallocArraysForNumPacketsPerLid(const size_t numExportLIDs, const size_t numImportLIDs)
Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary.
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 ("forward mode").
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.
bool transferArrived() const
Whether the data from an import/export operation has arrived, and is ready for the unpack and combine...
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).
DistObject & operator=(const DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &)=default
Assignment operator (default).
DistObject & operator=(DistObject< Packet, LocalOrdinal, GlobalOrdinal, Node > &&)=default
Move assignment (default).
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type
Kokkos::Device specialization for communication buffers.
LocalOrdinal local_ordinal_type
The type of local indices.
typename Node::device_type device_type
The Kokkos Device type.
typename ::Kokkos::ArithTraits< Packet >::val_type packet_type
The type of each datum being sent or received in an Import or Export.
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.
virtual bool checkSizes(const SrcDistObject &source)=0
Compare the source and target (this) objects for compatibility.
typename device_type::execution_space execution_space
The Kokkos execution space.
void print(std::ostream &os) const
Print this object to the given output stream.
GlobalOrdinal global_ordinal_type
The type of global indices.
Kokkos::DualView< size_t *, buffer_device_type > numExportPacketsPerLID_
Number of packets to send for each send operation.
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 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.
Teuchos::RCP< const map_type > map_
The Map over which this object is distributed.
ReverseOption
Whether the data transfer should be performed in forward or reverse mode.
Node node_type
The Node type. If you don't know what this is, don't use it.
virtual size_t constantNumberOfPackets() const
Whether the implementation's instance promises always to have a constant number of packets per LID (l...
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 ("forward mode").
virtual std::string description() const
One-line descriptiion of this object.
virtual Teuchos::RCP< const map_type > getMap() const
The Map describing the parallel distribution of this object.
virtual ~DistObject()=default
Destructor (virtual for memory safety of derived classes).
virtual void removeEmptyProcessesInPlace(const Teuchos::RCP< const map_type > &newMap)
Remove processes which contain no entries in this object's Map.
bool isDistributed() const
Whether this is a globally distributed object.
A parallel distribution of indices over processes.
Abstract base class for objects that can be the source of an Import or Export operation.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
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's Map.
CombineMode
Rule for combining data in an Import or Export.
@ INSERT
Insert new values that don't currently exist.