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_TPETRABLOCKMAP_HPP
00047 #define XPETRA_TPETRABLOCKMAP_HPP
00048
00049 #include "Xpetra_ConfigDefs.hpp"
00050
00051 #ifndef HAVE_XPETRA_TPETRA
00052 #error This file should be included only if HAVE_XPETRA_TPETRA is defined.
00053 #endif
00054
00055 #include "Xpetra_Map.hpp"
00056 #include "Xpetra_BlockMap.hpp"
00057
00058 #include "Tpetra_BlockMap.hpp"
00059 #include "Tpetra_Map.hpp"
00060
00061 #include "Xpetra_Exceptions.hpp"
00062
00063 namespace Xpetra {
00064
00065 template <class LocalOrdinal, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType>
00066 class TpetraBlockMap
00067 : public BlockMap<LocalOrdinal,GlobalOrdinal,Node>
00068 {
00069 public:
00070
00071
00072 typedef TpetraMap<LocalOrdinal, GlobalOrdinal, Node> TpetraMapClass;
00073
00075
00076
00079 TpetraBlockMap(global_size_t numGlobalBlocks,
00080 LocalOrdinal blockSize,
00081 GlobalOrdinal indexBase,
00082 const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
00083 const Teuchos::RCP<Node> &node = Kokkos::DefaultNode::getDefaultNode()) : map_(rcp(new Tpetra::BlockMap<LocalOrdinal, GlobalOrdinal, Node>(numGlobalBlocks, blockSize, indexBase, comm, node))) { }
00084
00087 TpetraBlockMap(global_size_t numGlobalBlocks,
00088 size_t numLocalBlocks,
00089 LocalOrdinal blockSize,
00090 GlobalOrdinal indexBase,
00091 const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
00092 const Teuchos::RCP<Node> &node = Kokkos::DefaultNode::getDefaultNode()) : map_(rcp(new Tpetra::BlockMap<LocalOrdinal, GlobalOrdinal, Node>(numGlobalBlocks, numLocalBlocks, blockSize, indexBase, comm, node))) { }
00093
00096 TpetraBlockMap(global_size_t numGlobalBlocks,
00097 const Teuchos::ArrayView<const GlobalOrdinal>& myGlobalBlockIDs,
00098 const Teuchos::ArrayView<const GlobalOrdinal>& myFirstGlobalPointInBlocks,
00099 const Teuchos::ArrayView<const LocalOrdinal>& myBlockSizes,
00100 GlobalOrdinal indexBase,
00101 const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
00102 const Teuchos::RCP<Node> &node = Kokkos::DefaultNode::getDefaultNode()) : map_(rcp(new Tpetra::BlockMap<LocalOrdinal, GlobalOrdinal, Node>(numGlobalBlocks, myGlobalBlockIDs, myFirstGlobalPointInBlocks, myBlockSizes, indexBase, comm, node))) { }
00103
00110 TpetraBlockMap(const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& pointMap,
00111 const Teuchos::ArrayView<const GlobalOrdinal>& myGlobalBlockIDs,
00112 const Teuchos::ArrayView<const LocalOrdinal>& myBlockSizes,
00113 const Teuchos::RCP<Node> &node = Kokkos::DefaultNode::getDefaultNode()) {
00114
00115 XPETRA_RCP_DYNAMIC_CAST(const TpetraMapClass, pointMap, tPointMap, "Xpetra::TpetraBlockMap constructors only accept Xpetra::TpetraMap as input arguments.");
00116 map_ = rcp(new Tpetra::BlockMap<LocalOrdinal, GlobalOrdinal, Node>(tPointMap->getTpetra_Map(), myGlobalBlockIDs, myBlockSizes, node));
00117 }
00118
00119 TpetraBlockMap(const Teuchos::RCP<const Tpetra::BlockMap<LocalOrdinal, GlobalOrdinal, Node> > &map) : map_(map) { }
00120
00122 virtual ~TpetraBlockMap(){ }
00123
00125
00127
00128
00129 global_size_t getGlobalNumBlocks() const { return map_->getGlobalNumBlocks(); }
00130
00132 size_t getNodeNumBlocks() const { return map_->getNodeNumBlocks(); }
00133
00134 Teuchos::ArrayView<const GlobalOrdinal> getNodeBlockIDs() const { return map_->getNodeBlockIDs(); }
00135
00136 bool isBlockSizeConstant() const { return map_->isBlockSizeConstant(); }
00137
00139 Teuchos::ArrayRCP<const LocalOrdinal> getNodeFirstPointInBlocks() const { return map_->getNodeFirstPointInBlocks(); }
00140
00142
00145 Teuchos::ArrayRCP<const LocalOrdinal> getNodeFirstPointInBlocks_Device() const { return map_->getNodeFirstPointInBlocks_Device(); }
00146
00148
00150 GlobalOrdinal getGlobalBlockID(LocalOrdinal localBlockID) const { return map_->getGlobalBlockID(localBlockID); }
00151
00153
00155 LocalOrdinal getLocalBlockID(GlobalOrdinal globalBlockID) const { return map_->getLocalBlockID(globalBlockID); }
00156
00158
00161 LocalOrdinal getLocalBlockSize(LocalOrdinal localBlockID) const { return map_->getLocalBlockSize(localBlockID); }
00162
00164
00167 LocalOrdinal getFirstLocalPointInLocalBlock(LocalOrdinal localBlockID) const { return map_->getFirstLocalPointInLocalBlock(localBlockID); }
00168
00170
00173 GlobalOrdinal getFirstGlobalPointInLocalBlock(LocalOrdinal localBlockID) const { return map_->getFirstGlobalPointInLocalBlock(localBlockID); }
00174
00176
00177 RCP< const Tpetra::BlockMap<LocalOrdinal, GlobalOrdinal, Node> > getTpetra_BlockMap() const { return map_; }
00178
00179 private:
00180 RCP< const Tpetra::BlockMap<LocalOrdinal, GlobalOrdinal, Node> > map_;
00181
00182
00183 };
00184
00185 }
00186
00187 #define XPETRA_TPETRABLOCKMAP_SHORT
00188 #endif