00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046 #ifndef XPETRA_TPETRAIMPORT_HPP
00047 #define XPETRA_TPETRAIMPORT_HPP
00048
00049
00050
00051 #include "Xpetra_TpetraConfigDefs.hpp"
00052
00053 #include "Xpetra_Import.hpp"
00054 #include "Xpetra_Exceptions.hpp"
00055
00056 #include "Xpetra_TpetraMap.hpp"
00057 #include "Tpetra_Import.hpp"
00058
00059 namespace Xpetra {
00060
00061
00062 template <class LocalOrdinal, class GlobalOrdinal, class Node>
00063 const Tpetra::Import<LocalOrdinal,GlobalOrdinal,Node> & toTpetra(const Import<LocalOrdinal,GlobalOrdinal,Node> &);
00064
00065 template <class LocalOrdinal, class GlobalOrdinal, class Node>
00066 RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > toXpetra(const RCP< const Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node > > &);
00067
00068
00069 template <class LocalOrdinal, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType>
00070 class TpetraImport
00071 : public Import<LocalOrdinal, GlobalOrdinal, Node>
00072 {
00073
00074 public:
00075
00077 typedef Map<LocalOrdinal,GlobalOrdinal,Node> map_type;
00078
00080
00081
00083 TpetraImport(const Teuchos::RCP< const map_type > &source, const Teuchos::RCP< const map_type > &target)
00084 : import_(Teuchos::rcp(new Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >(toTpetra(source), toTpetra(target)))) { }
00085
00087 TpetraImport(const Teuchos::RCP< const map_type > &source, const Teuchos::RCP< const map_type > &target, const Teuchos::RCP< Teuchos::ParameterList > &plist)
00088 : import_(Teuchos::rcp(new Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >(toTpetra(source), toTpetra(target), plist))) { }
00089
00091 TpetraImport(const Import< LocalOrdinal, GlobalOrdinal, Node > &import)
00092 : import_(Teuchos::rcp(new Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >(toTpetra(import)))) { }
00093
00095 ~TpetraImport() { }
00096
00098
00100
00101
00103 size_t getNumSameIDs() const { XPETRA_MONITOR("TpetraImport::getNumSameIDs"); return import_->getNumSameIDs(); }
00104
00106 size_t getNumPermuteIDs() const { XPETRA_MONITOR("TpetraImport::getNumPermuteIDs"); return import_->getNumPermuteIDs(); }
00107
00109 ArrayView< const LocalOrdinal > getPermuteFromLIDs() const { XPETRA_MONITOR("TpetraImport::getPermuteFromLIDs"); return import_->getPermuteFromLIDs(); }
00110
00112 ArrayView< const LocalOrdinal > getPermuteToLIDs() const { XPETRA_MONITOR("TpetraImport::getPermuteToLIDs"); return import_->getPermuteToLIDs(); }
00113
00115 size_t getNumRemoteIDs() const { XPETRA_MONITOR("TpetraImport::getNumRemoteIDs"); return import_->getNumRemoteIDs(); }
00116
00118 ArrayView< const LocalOrdinal > getRemoteLIDs() const { XPETRA_MONITOR("TpetraImport::getRemoteLIDs"); return import_->getRemoteLIDs(); }
00119
00121 size_t getNumExportIDs() const { XPETRA_MONITOR("TpetraImport::getNumExportIDs"); return import_->getNumExportIDs(); }
00122
00124 ArrayView< const LocalOrdinal > getExportLIDs() const { XPETRA_MONITOR("TpetraImport::getExportLIDs"); return import_->getExportLIDs(); }
00125
00127 ArrayView< const int > getExportImageIDs() const { XPETRA_MONITOR("TpetraImport::getExportImageIDs"); return import_->getExportImageIDs(); }
00128
00130 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getSourceMap() const { XPETRA_MONITOR("TpetraImport::getSourceMap"); return toXpetra(import_->getSourceMap()); }
00131
00133 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getTargetMap() const { XPETRA_MONITOR("TpetraImport::getTargetMap"); return toXpetra(import_->getTargetMap()); }
00134
00136
00138
00139
00141 void print(std::ostream &os) const { XPETRA_MONITOR("TpetraImport::print"); import_->print(os); }
00142
00144
00146
00147
00149 TpetraImport(const RCP<const Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node > > &import) : import_(import) { }
00150
00151 RCP< const Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node > > getTpetra_Import() const { return import_; }
00152
00154
00155 private:
00156
00157 RCP<const Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node > > import_;
00158
00159 };
00160
00161
00162 template <class LocalOrdinal, class GlobalOrdinal, class Node>
00163 const Tpetra::Import<LocalOrdinal,GlobalOrdinal,Node> & toTpetra(const Import<LocalOrdinal,GlobalOrdinal,Node> &import) {
00164
00165 const TpetraImport<LocalOrdinal,GlobalOrdinal,Node> & tpetraImport = dynamic_cast<const TpetraImport<LocalOrdinal,GlobalOrdinal,Node> &>(import);
00166 return *tpetraImport.getTpetra_Import();
00167 }
00168
00169 template <class LocalOrdinal, class GlobalOrdinal, class Node>
00170 RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > toXpetra(const RCP< const Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node > > &import) {
00171 return rcp( new TpetraImport<LocalOrdinal, GlobalOrdinal, Node>(import));
00172 }
00173
00174
00175 }
00176
00177 #define XPETRA_TPETRAIMPORT_SHORT
00178 #endif // XPETRA_TPETRAIMPORT_HPP