Domi
Multi-dimensional, distributed data structures
 All Classes Files Functions Variables Typedefs Friends
Domi_MDMap.hpp
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Domi: Multi-dimensional Distributed Linear Algebra Services
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia
8 // Corporation, the U.S. Government retains certain rights in this
9 // software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact William F. Spotz (wfspotz@sandia.gov)
39 //
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef DOMI_MDMAP_HPP
44 #define DOMI_MDMAP_HPP
45 
46 // System includes
47 #include <algorithm>
48 #include <limits>
49 #include <sstream>
50 
51 // Domi includes
52 #include "Domi_ConfigDefs.hpp"
53 #include "Domi_Utils.hpp"
54 #include "Domi_getValidParameters.hpp"
55 #include "Domi_Slice.hpp"
56 #include "Domi_MDComm.hpp"
57 #include "Domi_MDArray.hpp"
58 
59 // Teuchos includes
60 #include "Teuchos_Comm.hpp"
61 #include "Teuchos_DefaultComm.hpp"
62 #include "Teuchos_CommHelpers.hpp"
63 #include "Teuchos_Tuple.hpp"
64 
65 // Kokkos includes
66 #include "Kokkos_Core.hpp"
67 
68 #ifdef HAVE_EPETRA
69 #include "Epetra_Map.h"
70 #endif
71 
72 #ifdef HAVE_TPETRA
73 #include "Tpetra_Map.hpp"
74 #include "Tpetra_Util.hpp"
75 typedef typename Tpetra::Map<>::local_ordinal_type TpetraLOType;
76 typedef typename Tpetra::Map<>::global_ordinal_type TpetraGOType;
77 typedef typename Tpetra::Map<>::node_type TpetraNodeType;
78 #endif
79 
80 namespace Domi
81 {
82 
145 class MDMap
146 {
147 public:
148 
151 
185  MDMap(const Teuchos::RCP< const MDComm > mdComm,
186  const Teuchos::ArrayView< const dim_type > & dimensions,
187  const Teuchos::ArrayView< const int > & commPad =
188  Teuchos::ArrayView< const int >(),
189  const Teuchos::ArrayView< const int > & bndryPad =
190  Teuchos::ArrayView< const int >(),
191  const Teuchos::ArrayView< const int > & replicatedBoundary =
192  Teuchos::ArrayView< const int >(),
193  const Layout layout = DEFAULT_ORDER);
194 
206  MDMap(Teuchos::ParameterList & plist);
207 
220  MDMap(const Teuchos::RCP< const Teuchos::Comm< int > > teuchosComm,
221  Teuchos::ParameterList & plist);
222 
235  MDMap(const Teuchos::RCP< const MDComm > mdComm,
236  Teuchos::ParameterList & plist);
237 
267  MDMap(const Teuchos::RCP< const MDComm > mdComm,
268  const Teuchos::ArrayView< Slice > & myGlobalBounds,
269  const Teuchos::ArrayView< padding_type > & padding =
270  Teuchos::ArrayView< padding_type >(),
271  const Teuchos::ArrayView< const int > & replicatedBoundary =
272  Teuchos::ArrayView< const int >(),
273  const Layout layout = DEFAULT_ORDER);
274 
279  MDMap(const MDMap & source);
280 
294  MDMap(const MDMap & parent,
295  int axis,
296  dim_type index);
297 
319  MDMap(const MDMap & parent,
320  int axis,
321  const Slice & slice,
322  int bndryPad = 0);
323 
337  MDMap(const MDMap & parent,
338  const Teuchos::ArrayView< Slice > & slices,
339  const Teuchos::ArrayView< int > & bndryPad =
340  Teuchos::ArrayView< int >());
341 
344  ~MDMap();
345 
347 
350 
355  MDMap & operator=(const MDMap & source);
356 
358 
361 
368  inline Teuchos::RCP< const Teuchos::Comm< int > > getTeuchosComm() const;
369 
375  inline Teuchos::RCP< const MDComm > getMDComm() const;
376 
383  inline bool onSubcommunicator() const;
384 
391  inline int numDims() const;
392 
399  inline Teuchos::Array< int > getCommDims() const;
400 
410  inline int getCommDim(int axis) const;
411 
421  inline bool isPeriodic(int axis) const;
422 
432  inline int getCommIndex(int axis) const;
433 
450  inline int getLowerNeighbor(int axis) const;
451 
469  inline int getUpperNeighbor(int axis) const;
470 
472 
475 
479  inline Teuchos::Array< dim_type > getGlobalDims() const;
480 
489  dim_type getGlobalDim(int axis,
490  bool withBndryPad=false) const;
491 
500  Slice getGlobalBounds(int axis,
501  bool withBndryPad=false) const;
502 
516  Slice getGlobalRankBounds(int axis,
517  bool withBndryPad=false) const;
518 
521  Teuchos::Array< dim_type > getLocalDims() const;
522 
531  dim_type getLocalDim(int axis,
532  bool withPad=false) const;
533 
541  Teuchos::ArrayView< const Slice > getLocalBounds() const;
542 
556  Slice getLocalBounds(int axis,
557  bool withPad=false) const;
558 
576  Slice getLocalInteriorBounds(int axis) const;
577 
579 
582 
593  bool hasPadding() const;
594 
603  int getLowerPadSize(int axis) const;
604 
613  int getUpperPadSize(int axis) const;
614 
624  int getCommPadSize(int axis) const;
625 
632  int getLowerBndryPad(int axis) const;
633 
640  int getUpperBndryPad(int axis) const;
641 
648  Teuchos::Array< int > getBndryPadSizes() const;
649 
659  int getBndryPadSize(int axis) const;
660 
661  /* \brief Return whether given local index is in the padding
662  *
663  * \param index [in] An array of indexes ((i) for 1D, (i,j) for 2D,
664  * (i,j,k) for 3D, etc)
665  */
666  bool isPad(const Teuchos::ArrayView< dim_type > & index) const;
667 
668  /* \brief Return whether given local index is in the communication
669  * padding
670  *
671  * \param index [in] An array of indexes ((i) for 1D, (i,j) for 2D,
672  * (i,j,k) for 3D, etc)
673  */
674  bool isCommPad(const Teuchos::ArrayView< dim_type > & index) const;
675 
676  /* \brief Return whether given local index is in the boundary
677  * padding
678  *
679  * \param index [in] An array of indexes ((i) for 1D, (i,j) for 2D,
680  * (i,j,k) for 3D, etc)
681  */
682  bool isBndryPad(const Teuchos::ArrayView< dim_type > & index) const;
683 
686  bool isReplicatedBoundary(int axis) const;
687 
690  Layout getLayout() const;
691 
693 
696 
701  Teuchos::ArrayView< Teuchos::RCP< const Domi::MDMap > >
702  getAxisMaps() const;
703 
708  Teuchos::RCP< const Domi::MDMap > getAxisMap(int axis) const;
709 
726  Teuchos::RCP< const MDMap >
727  getAugmentedMDMap(const dim_type leadingDim,
728  const dim_type trailingDim=0) const;
729 
730 #ifdef HAVE_EPETRA
731 
740  Teuchos::RCP< const Epetra_Map >
741  getEpetraMap(bool withCommPad=true) const;
742 
753  Teuchos::RCP< const Epetra_Map >
754  getEpetraAxisMap(int axis,
755  bool withCommPad=true) const;
756 
757 #endif
758 
759 #ifdef HAVE_TPETRA
760 
770  template< class LocalOrdinal = TpetraLOType,
771  class GlobalOrdinal = TpetraGOType,
772  class Node = TpetraNodeType>
773  Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > >
774  getTpetraMap(bool withCommPad=true) const;
775 
788  template< class LocalOrdinal = TpetraLOType,
789  class GlobalOrdinal = TpetraGOType,
790  class Node = TpetraNodeType >
791  Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > >
792  getTpetraAxisMap(int axis,
793  bool withCommPad=true) const;
794 
795 #endif
796 
798 
801 
806  Teuchos::Array< dim_type >
807  getGlobalIndex(size_type globalID) const;
808 
813  Teuchos::Array< dim_type >
814  getLocalIndex(size_type localID) const;
815 
820  size_type getGlobalID(size_type localID) const;
821 
826  size_type
827  getGlobalID(const Teuchos::ArrayView< dim_type > & globalIndex) const;
828 
836  size_type getLocalID(size_type globalID) const;
837 
842  size_type
843  getLocalID(const Teuchos::ArrayView< dim_type > & localIndex) const;
844 
846 
849 
864  bool isCompatible(const MDMap & mdMap) const;
865 
887  bool isSameAs(const MDMap & mdMap,
888  const int verbose = 0) const;
889 
899  bool isContiguous() const;
900 
902 
903 private:
904 
905  // A private method for computing the bounds and local dimensions,
906  // after the global dimensions, communication and boundary padding
907  // have been properly assigned.
908  void computeBounds();
909 
910  // The underlying multi-dimensional communicator.
911  Teuchos::RCP< const MDComm > _mdComm;
912 
913  // The size of the global dimensions along each axis. This includes
914  // the values of the boundary padding.
915  Teuchos::Array< dim_type > _globalDims;
916 
917  // Store the start and stop indexes of the global problem along each
918  // axis. This includes the values of the boundary padding.
919  Teuchos::Array< Slice > _globalBounds;
920 
921  // A double array of Slices that stores the starting and stopping
922  // indexes for each axis processor along each axis. This data
923  // structure allows any processor to know the map bounds on any
924  // other processor. These bounds do NOT include the boundary or
925  // communication padding. The outer array will have a size equal to
926  // the number of dimensions. Each inner array will have a size
927  // equal to the number of axis processors along that axis.
928  Teuchos::Array< Teuchos::Array< Slice > > _globalRankBounds;
929 
930  // The global stride between adjacent IDs. This quantity is
931  // "virtual" as it does not describe actual memory layout, but it is
932  // useful for index conversion algorithms.
933  Teuchos::Array< size_type > _globalStrides;
934 
935  // The minumum global ID of the global data structure, including
936  // boundary padding. This will only be non-zero on a sub-map.
937  size_type _globalMin;
938 
939  // The maximum global ID of the global data structure, including
940  // boundary padding.
941  size_type _globalMax;
942 
943  // The size of the local dimensions along each axis. This includes
944  // the values of the padding.
945  Teuchos::Array< dim_type > _localDims;
946 
947  // The local loop bounds along each axis, stored as an array of
948  // Slices. These bounds DO include the padding. By definition, the
949  // start() attribute of these Slices will all be zero by definition.
950  // This may seem wasteful (we could store an array of dim_type
951  // rather than an array of Slice), but this storage is convenient
952  // when it comes to the getLocalBounds() method.
953  Teuchos::Array< Slice > _localBounds;
954 
955  // The local stride between adjacent elements in memory.
956  Teuchos::Array< size_type > _localStrides;
957 
958  // The minimum local ID of the local data structure, including
959  // padding.
960  size_type _localMin;
961 
962  // The maximum local ID of the local data structure, including
963  // padding.
964  size_type _localMax;
965 
966  // The communication padding that was specified at construction, one
967  // value along each axis.
968  Teuchos::Array< int > _commPadSizes;
969 
970  // The actual padding stored on this processor along each axis. The
971  // padding can be either communication padding or boundary padding
972  // based on the processor position on the boundary of a domain.
973  Teuchos::Array< padding_type > _pad;
974 
975  // The size of the boundary padding along each axis.
976  Teuchos::Array< int > _bndryPadSizes;
977 
978  // The actual boundary padding stored along each axis. For a full
979  // communicator, the lower and upper boundary padding are both the
980  // same as the corresponding value in _bndryPadSizes. However, the
981  // introduction of sub-maps creates the possibility of different
982  // upper and lower boundary padding values. If an axis is periodic,
983  // then these values will be set to the communication padding.
984  Teuchos::Array< padding_type > _bndryPad;
985 
986  // An array of flags, one for each axis, that specifies whether the
987  // endpoints of a periodic boundary are replicated (1 or true) or
988  // unique (0 or false). These flags only have meaning for axes that
989  // are periodic. The arrray type is int, beacause Array< bool > is
990  // specialized and sometimes difficult to work with.
991  Teuchos::Array< int > _replicatedBoundary;
992 
993  // The storage order
994  Layout _layout;
995 
996  // An array of axis maps that correspond to the full MDMap. An axis
997  // map describes the map along a single axis. This member is mutable
998  // because it is logically const but does not get constructed until
999  // requested by the user.
1000  mutable Teuchos::Array< Teuchos::RCP< const MDMap > > _axisMaps;
1001 
1002 #ifdef HAVE_EPETRA
1003  // An RCP pointing to an Epetra_Map that is equivalent to this
1004  // MDMap, including communication padding. It is mutable because we
1005  // do not construct it until it asked for by a get method that is
1006  // const.
1007  mutable Teuchos::RCP< const Epetra_Map > _epetraMap;
1008 
1009  // An RCP pointing to an Epetra_Map that is equivalent to this
1010  // MDMap, excluding communication padding. The returned Epetra_Map
1011  // thus indicates what IDs are owned by this processor. It is
1012  // mutable because we do not construct it until it asked for by a
1013  // get method that is const.
1014  mutable Teuchos::RCP< const Epetra_Map > _epetraOwnMap;
1015 
1016  // An ArrayRCP that stores Epetra_Maps that represent the
1017  // distribution of indexes along each axis, including communication
1018  // padding. It is mutable because we do not construct it until it
1019  // asked for by a get method that is const.
1020  mutable Teuchos::Array< Teuchos::RCP< const Epetra_Map > > _epetraAxisMaps;
1021 
1022  // An ArrayRCP that stores Epetra_Maps that represent the
1023  // distribution of indexes along each axis, excluding communication
1024  // padding. It is mutable because we do not construct it until it
1025  // asked for by a get method that is const.
1026  mutable Teuchos::Array< Teuchos::RCP< const Epetra_Map > > _epetraAxisOwnMaps;
1027 #endif
1028 
1029 };
1030 
1032 // Inline implementations //
1034 
1036 
1037 Teuchos::RCP< const Teuchos::Comm< int > >
1039 {
1040  return _mdComm->getTeuchosComm();
1041 }
1042 
1044 
1045 Teuchos::RCP< const MDComm >
1047 {
1048  return _mdComm;
1049 }
1050 
1052 
1053 bool
1055 {
1056  return _mdComm->onSubcommunicator();
1057 }
1058 
1060 
1061 int
1063 {
1064  return _mdComm->numDims();
1065 }
1066 
1068 
1069 Teuchos::Array< int >
1071 {
1072  return _mdComm->getCommDims();
1073 }
1074 
1076 
1077 int
1078 MDMap::getCommDim(int axis) const
1079 {
1080  return _mdComm->getCommDim(axis);
1081 }
1082 
1084 
1085 bool
1086 MDMap::isPeriodic(int axis) const
1087 {
1088  return _mdComm->isPeriodic(axis);
1089 }
1090 
1092 
1093 int
1094 MDMap::getCommIndex(int axis) const
1095 {
1096  return _mdComm->getCommIndex(axis);
1097 }
1098 
1100 
1101 int
1103 {
1104  return _mdComm->getLowerNeighbor(axis);
1105 }
1106 
1108 
1109 int
1111 {
1112  return _mdComm->getUpperNeighbor(axis);
1113 }
1114 
1116 
1117 Teuchos::Array< dim_type >
1118 MDMap::
1120 {
1121  return _globalDims;
1122 }
1123 
1125 
1127 // Templated method implementations //
1129 
1131 
1132 #ifdef HAVE_TPETRA
1133 
1134 template< class LocalOrdinal,
1135  class GlobalOrdinal,
1136  class Node >
1137 Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > >
1138 MDMap::getTpetraMap(bool withCommPad) const
1139 {
1140  if (withCommPad)
1141  {
1142  // Allocate the elementsMDArray and the index array
1143  int num_dims = numDims();
1144  Teuchos::Array<dim_type> localDims(num_dims);
1145  for (int axis = 0; axis < num_dims; ++axis)
1146  localDims[axis] = _localDims[axis];
1147  MDArray< GlobalOrdinal > elementMDArray(localDims);
1148  Teuchos::Array< LocalOrdinal > index(num_dims);
1149 
1150  // Iterate over the local MDArray and assign global IDs
1151  for (typename MDArray< GlobalOrdinal >::iterator it = elementMDArray.begin();
1152  it != elementMDArray.end(); ++it)
1153  {
1154  GlobalOrdinal globalID = 0;
1155  for (int axis = 0; axis < num_dims; ++axis)
1156  {
1157  int axisRank = getCommIndex(axis);
1158  GlobalOrdinal start = _globalRankBounds[axis][axisRank].start() -
1159  _pad[axis][0];
1160  globalID += (start + it.index(axis)) * _globalStrides[axis];
1161  }
1162  *it = globalID;
1163  }
1164 
1165  // Return the Tpetra::Map
1166  const Teuchos::Array< GlobalOrdinal > & myElements =
1167  elementMDArray.array();
1168  Teuchos::RCP< const Teuchos::Comm< int > > teuchosComm =
1169  _mdComm->getTeuchosComm();
1170  return
1171  Teuchos::rcp(new Tpetra::Map< LocalOrdinal,
1172  GlobalOrdinal,
1173  Node >(Teuchos::OrdinalTraits< Tpetra::global_size_t >::invalid(),
1174  myElements(),
1175  0,
1176  teuchosComm));
1177  }
1178  else
1179  {
1180  // Allocate the elementMDArray MDArray and the index array
1181  int num_dims = numDims();
1182  Teuchos::Array< LocalOrdinal > index(num_dims);
1183  Teuchos::Array< dim_type > myDims(num_dims);
1184  for (int axis = 0; axis < num_dims; ++axis)
1185  {
1186  myDims[axis] =
1187  _localDims[axis] - _pad[axis][0] - _pad[axis][1];
1188  int axisRank = getCommIndex(axis);
1189  if (axisRank == 0)
1190  myDims[axis] += _bndryPad[axis][0];
1191  if (axisRank == getCommDim(axis)-1)
1192  myDims[axis] += _bndryPad[axis][1];
1193  }
1194  MDArray< GlobalOrdinal > elementMDArray(myDims());
1195 
1196  // Iterate over the local MDArray and assign global IDs
1197  for (typename MDArray< GlobalOrdinal >::iterator it = elementMDArray.begin();
1198  it != elementMDArray.end(); ++it)
1199  {
1200  GlobalOrdinal globalID = 0;
1201  for (int axis = 0; axis < num_dims; ++axis)
1202  {
1203  int axisRank = getCommIndex(axis);
1204  GlobalOrdinal start = _globalRankBounds[axis][axisRank].start();
1205  if (axisRank == 0)
1206  start -= _bndryPad[axis][0];
1207  if (axisRank == getCommDim(axis)-1)
1208  start += _bndryPad[axis][1];
1209  globalID += (start + it.index(axis)) * _globalStrides[axis];
1210  }
1211  }
1212 
1213  // Return the Tpetra::Map
1214  const Teuchos::Array< GlobalOrdinal> & myElements =
1215  elementMDArray.array();
1216  Teuchos::RCP< const Teuchos::Comm< int > > teuchosComm =
1217  _mdComm->getTeuchosComm();
1218  return
1219  Teuchos::rcp(new Tpetra::Map< LocalOrdinal,
1220  GlobalOrdinal,
1221  Node >(Teuchos::OrdinalTraits< Tpetra::global_size_t>::invalid(),
1222  myElements(),
1223  0,
1224  teuchosComm));
1225  }
1226 }
1227 #endif
1228 
1230 
1231 #ifdef HAVE_TPETRA
1232 
1233 template< class LocalOrdinal,
1234  class GlobalOrdinal,
1235  class Node >
1236 Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > >
1237 MDMap::
1238 getTpetraAxisMap(int axis,
1239  bool withCommPad) const
1240 {
1241 #ifdef HAVE_DOMI_ARRAY_BOUNDSCHECK
1242  TEUCHOS_TEST_FOR_EXCEPTION(
1243  ((axis < 0) || (axis >= numDims())),
1244  RangeError,
1245  "invalid axis index = " << axis << " (number of dimensions = " <<
1246  numDims() << ")");
1247 #endif
1248  Teuchos::RCP< const Teuchos::Comm< int > > teuchosComm =
1249  _mdComm->getTeuchosComm();
1250  Teuchos::Array< GlobalOrdinal > elements(getLocalDim(axis,withCommPad));
1251  GlobalOrdinal start = getGlobalRankBounds(axis,true).start();
1252  if (withCommPad && (getCommIndex(axis) != 0)) start -= _pad[axis][0];
1253  for (LocalOrdinal i = 0; i < elements.size(); ++i)
1254  elements[i] = i + start;
1255  return Teuchos::rcp(new Tpetra::Map< LocalOrdinal,
1256  GlobalOrdinal,
1257  Node >(Teuchos::OrdinalTraits< Tpetra::global_size_t>::invalid(),
1258  elements,
1259  0,
1260  teuchosComm));
1261 }
1262 #endif
1263 
1265 
1266 }
1267 
1268 #endif
size_type getGlobalID(size_type localID) const
Convert a local ID to a global ID.
Teuchos::Array< int > getCommDims() const
Get the communicator sizes along each axis.
Definition: Domi_MDMap.hpp:1070
Slice getGlobalBounds(int axis, bool withBndryPad=false) const
Get the bounds of the global problem.
Slice getLocalInteriorBounds(int axis) const
Get the local interior loop bounds along the specified axis.
Iterator class suitable for multi-dimensional arrays.
Definition: Domi_MDIterator.hpp:100
int getCommPadSize(int axis) const
Get the communication padding size along the given axis.
Slice getGlobalRankBounds(int axis, bool withBndryPad=false) const
Get the global loop bounds along the specified axis.
bool isCompatible(const MDMap &mdMap) const
True if two MDMaps are &quot;compatible&quot;.
int getBndryPadSize(int axis) const
Get the boundary padding size along the given axis.
int getUpperPadSize(int axis) const
Get the size of the upper padding along the given axis.
Teuchos::Array< dim_type > getLocalDims() const
Get an array of the local dimensions, including padding.
Teuchos::ArrayView< Teuchos::RCP< const Domi::MDMap > > getAxisMaps() const
Return an array of axis maps.
int getCommIndex(int axis) const
Get the axis rank of this processor.
Definition: Domi_MDMap.hpp:1094
~MDMap()
Destructor.
int getLowerPadSize(int axis) const
Get the size of the lower padding along the given axis.
A Slice defines a subset of a container.
int getUpperNeighbor(int axis) const
Get the rank of the upper neighbor.
Definition: Domi_MDMap.hpp:1110
size_type getLocalID(size_type globalID) const
Convert a global ID to a local ID.
A Slice contains a start, stop, and step index, describing a subset of an ordered container...
Definition: Domi_Slice.hpp:137
bool onSubcommunicator() const
Query whether this processor is on the sub-communicator.
Definition: Domi_MDMap.hpp:1054
int numDims() const
Get the number of dimensions.
Definition: Domi_MDMap.hpp:1062
int getLowerNeighbor(int axis) const
Get the rank of the lower neighbor.
Definition: Domi_MDMap.hpp:1102
Teuchos::Array< dim_type > getGlobalIndex(size_type globalID) const
Convert a global ID to a global index.
Teuchos::RCP< const Domi::MDMap > getAxisMap(int axis) const
Return an axis map for the given axis.
dim_type getGlobalDim(int axis, bool withBndryPad=false) const
Get the global dimension along the specified axis.
MDMap & operator=(const MDMap &source)
Assignment operator.
bool isSameAs(const MDMap &mdMap, const int verbose=0) const
True if two MDMaps are &quot;identical&quot;.
Teuchos::RCP< const MDComm > getMDComm() const
Access the underlying MDComm object.
Definition: Domi_MDMap.hpp:1046
Teuchos::RCP< const MDMap > getAugmentedMDMap(const dim_type leadingDim, const dim_type trailingDim=0) const
Return an RCP to a new MDMap that is a simple augmentation of this MDMap.
const dim_type start() const
Start index.
Definition: Domi_Slice.hpp:431
bool isReplicatedBoundary(int axis) const
Return whether the given axis supports a replicated boundary.
bool isContiguous() const
True if there are no stride gaps on any processor.
Teuchos::RCP< const Teuchos::Comm< int > > getTeuchosComm() const
Get the Teuchos communicator.
Definition: Domi_MDMap.hpp:1038
int getLowerBndryPad(int axis) const
Get the size of the lower boundary padding along the given axis.
Teuchos::Array< int > getBndryPadSizes() const
Get the boundary padding sizes along each axis.
Teuchos::Array< dim_type > getLocalIndex(size_type localID) const
Convert a local ID to a local index.
Teuchos::Array< dim_type > getGlobalDims() const
Get an array of the the global dimensions, including boundary padding.
Definition: Domi_MDMap.hpp:1119
bool hasPadding() const
Return true if there is any padding stored locally.
int getCommDim(int axis) const
Get the communicator size along the given axis.
Definition: Domi_MDMap.hpp:1078
Multi-dimensional map.
Definition: Domi_MDMap.hpp:145
Memory-safe templated multi-dimensional array class.
Definition: Domi_MDArray.hpp:65
dim_type getLocalDim(int axis, bool withPad=false) const
Get the local dimension along the specified axis.
int getUpperBndryPad(int axis) const
Get the size of the upper boundary padding along the given axis.
bool isPeriodic(int axis) const
Return the periodic flag for the given axis.
Definition: Domi_MDMap.hpp:1086
Layout getLayout() const
Get the storage order.
MDMap(const Teuchos::RCP< const MDComm > mdComm, const Teuchos::ArrayView< const dim_type > &dimensions, const Teuchos::ArrayView< const int > &commPad=Teuchos::ArrayView< const int >(), const Teuchos::ArrayView< const int > &bndryPad=Teuchos::ArrayView< const int >(), const Teuchos::ArrayView< const int > &replicatedBoundary=Teuchos::ArrayView< const int >(), const Layout layout=DEFAULT_ORDER)
Main constructor.
Teuchos::ArrayView< const Slice > getLocalBounds() const
Get the local loop bounds along every axis.

Generated on Fri Apr 19 2024 09:28:43 for Domi by doxygen 1.8.5