43 #ifndef DOMI_MDVECTOR_HPP
44 #define DOMI_MDVECTOR_HPP
52 #include "Domi_ConfigDefs.hpp"
53 #include "Domi_MDMap.hpp"
54 #include "Domi_MDArrayRCP.hpp"
57 #include "Teuchos_DataAccess.hpp"
58 #include "Teuchos_Describable.hpp"
59 #include "Teuchos_ScalarTraitsDecl.hpp"
60 #include "Teuchos_Comm.hpp"
61 #include "Teuchos_CommHelpers.hpp"
64 #include "Epetra_IntVector.h"
65 #include "Epetra_Vector.h"
69 #include "Tpetra_Vector.hpp"
80 struct ompi_datatype_t {};
174 template<
class Scalar >
189 MDVector(
const Teuchos::RCP< const MDMap > & mdMap,
190 bool zeroOut =
true);
213 MDVector(
const Teuchos::RCP< const MDMap > & mdMap,
214 const dim_type leadingDim,
215 const dim_type trailingDim = 0,
216 bool zeroOut =
true);
225 MDVector(
const Teuchos::RCP< const MDMap > & mdMap,
235 MDVector(
const Teuchos::RCP< const MDMap > & mdMap,
251 Teuchos::DataAccess access = Teuchos::View);
266 MDVector(
const Teuchos::RCP<
const Teuchos::Comm< int > > teuchosComm,
267 Teuchos::ParameterList & plist);
281 MDVector(
const Teuchos::RCP< const MDComm > mdComm,
282 Teuchos::ParameterList & plist);
331 const Teuchos::ArrayView< Slice > & slices,
332 const Teuchos::ArrayView< int > & bndryPad =
333 Teuchos::ArrayView< int >());
353 inline const Teuchos::RCP< const MDMap >
370 inline Teuchos::RCP< const Teuchos::Comm< int > >
getTeuchosComm()
const;
458 inline dim_type
getGlobalDim(
int axis,
bool withBndryPad=
false)
const;
494 inline dim_type
getLocalDim(
int axis,
bool withCommPad=
false)
const;
668 Teuchos::RCP< Epetra_IntVector > getEpetraIntVectorView()
const;
685 Teuchos::RCP< Epetra_Vector > getEpetraVectorView()
const;
707 Teuchos::RCP< Epetra_MultiVector > getEpetraMultiVectorView()
const;
718 Teuchos::RCP< Epetra_IntVector > getEpetraIntVectorCopy()
const;
729 Teuchos::RCP< Epetra_Vector > getEpetraVectorCopy()
const;
747 Teuchos::RCP< Epetra_MultiVector > getEpetraMultiVectorCopy()
const;
763 template<
class LocalOrdinal,
764 class GlobalOrdinal = LocalOrdinal >
765 Teuchos::RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal > >
766 getTpetraVectorView()
const;
785 template <
class LocalOrdinal,
786 class GlobalOrdinal = LocalOrdinal >
787 Teuchos::RCP< Tpetra::MultiVector< Scalar,
790 getTpetraMultiVectorView()
const;
797 template<
class LocalOrdinal,
798 class GlobalOrdinal = LocalOrdinal >
799 Teuchos::RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal > >
800 getTpetraVectorCopy()
const;
814 template <
class LocalOrdinal,
815 class GlobalOrdinal = LocalOrdinal >
816 Teuchos::RCP< Tpetra::MultiVector< Scalar,
819 getTpetraMultiVectorCopy()
const;
884 typename Teuchos::ScalarTraits< Scalar >::magnitudeType
norm1()
const;
888 typename Teuchos::ScalarTraits< Scalar >::magnitudeType
norm2()
const;
892 typename Teuchos::ScalarTraits< Scalar >::magnitudeType
normInf()
const;
898 typename Teuchos::ScalarTraits< Scalar >::magnitudeType
922 describe(Teuchos::FancyOStream &out,
923 const Teuchos::EVerbosityLevel verbLevel =
924 Teuchos::Describable::verbLevel_default)
const;
939 bool includePadding =
true);
1039 bool includeBndryPad =
false)
const;
1048 void readBinary(
const std::string & filename,
1049 bool includeBndryPad =
false);
1058 Teuchos::RCP< const Teuchos::Comm< int > > _teuchosComm;
1062 Teuchos::RCP< const MDMap > _mdMap;
1083 Teuchos::Array< MPI_Request > _requests;
1096 Teuchos::RCP< MPI_Datatype > datatype;
1110 Teuchos::Array< Teuchos::Tuple< MessageInfo, 2 > > _sendMessages;
1115 Teuchos::Array< Teuchos::Tuple< MessageInfo, 2 > > _recvMessages;
1119 void initializeMessages();
1131 Teuchos::Array< dim_type > fileShape;
1132 Teuchos::Array< dim_type > bufferShape;
1133 Teuchos::Array< dim_type > dataShape;
1134 Teuchos::Array< dim_type > fileStart;
1135 Teuchos::Array< dim_type > dataStart;
1137 Teuchos::RCP< MPI_Datatype > filetype;
1138 Teuchos::RCP< MPI_Datatype > datatype;
1146 mutable Teuchos::RCP< FileInfo > _fileInfo;
1152 mutable Teuchos::RCP< FileInfo > _fileInfoWithBndry;
1160 Teuchos::RCP< FileInfo > & computeFileInfo(
bool includeBndryPad)
const;
1170 template<
class Scalar >
1174 _teuchosComm(mdMap->getTeuchosComm()),
1185 setObjectLabel(
"Domi::MDVector");
1188 int numDims = _mdMap->numDims();
1189 Teuchos::Array< dim_type > dims(numDims);
1190 for (
int axis = 0; axis <
numDims; ++axis)
1191 dims[axis] = _mdMap->getLocalDim(axis,
true);
1194 _mdArrayRcp.
resize(dims);
1195 _mdArrayView = _mdArrayRcp();
1200 template<
class Scalar >
1203 const dim_type leadingDim,
1204 const dim_type trailingDim,
1206 _teuchosComm(mdMap->getTeuchosComm()),
1207 _mdMap(mdMap->getAugmentedMDMap(leadingDim, trailingDim)),
1217 setObjectLabel(
"Domi::MDVector");
1220 int numDims = _mdMap->numDims();
1221 Teuchos::Array< dim_type > dims(numDims);
1222 for (
int axis = 0; axis <
numDims; ++axis)
1223 dims[axis] = _mdMap->getLocalDim(axis,
true);
1226 _mdArrayRcp.
resize(dims);
1227 _mdArrayView = _mdArrayRcp();
1232 template<
class Scalar >
1237 _mdArrayRcp(source),
1238 _mdArrayView(_mdArrayRcp()),
1246 setObjectLabel(
"Domi::MDVector");
1247 int numDims = _mdMap->numDims();
1248 TEUCHOS_TEST_FOR_EXCEPTION(
1249 numDims != _mdArrayRcp.
numDims(),
1251 "MDMap and source array do not have the same number of dimensions");
1253 for (
int axis = 0; axis <
numDims; ++axis)
1255 TEUCHOS_TEST_FOR_EXCEPTION(
1256 _mdMap->getLocalDim(axis) != _mdArrayRcp.
dimension(axis),
1258 "Axis " << axis <<
": MDMap dimension = " << _mdMap->getLocalDim(axis)
1259 <<
", MDArray dimension = " << _mdArrayRcp.
dimension(axis));
1265 template<
class Scalar >
1270 _mdArrayRcp(mdArrayRcp),
1271 _mdArrayView(_mdArrayRcp()),
1279 #ifdef DOMI_MDVECTOR_VERBOSE
1280 cout <<
"_mdArrayRcp = " << _mdArrayRcp << endl;
1281 cout <<
"_mdArrayView.getRawPtr() = " << _mdArrayView.
getRawPtr()
1282 <<
" (constructor)" << endl;
1283 cout <<
"_mdArrayView = " << _mdArrayView << endl;
1285 setObjectLabel(
"Domi::MDVector");
1286 int numDims = _mdMap->numDims();
1287 TEUCHOS_TEST_FOR_EXCEPTION(
1288 numDims != _mdArrayRcp.
numDims(),
1290 "MDMap and source array do not have the same number of dimensions");
1292 for (
int axis = 0; axis <
numDims; ++axis)
1294 TEUCHOS_TEST_FOR_EXCEPTION(
1295 _mdMap->getLocalDim(axis) != _mdArrayRcp.
dimension(axis),
1297 "Axis " << axis <<
": MDMap dimension = " << _mdMap->getLocalDim(axis)
1298 <<
", MDArray dimension = " << _mdArrayRcp.
dimension(axis));
1304 template<
class Scalar >
1307 Teuchos::DataAccess access) :
1308 _teuchosComm(source.getMDMap()->getTeuchosComm()),
1309 _mdMap(source.getMDMap()),
1310 _mdArrayRcp(source._mdArrayRcp),
1311 _mdArrayView(source._mdArrayView),
1319 setObjectLabel(
"Domi::MDVector");
1321 if (access == Teuchos::Copy)
1323 #ifdef DOMI_MDVECTOR_VERBOSE
1324 cout <<
"Inside MDVector copy constructor with copy access" << endl;
1327 int numDims = _mdMap->numDims();
1328 Teuchos::Array< dim_type > dims(numDims);
1329 for (
int axis = 0; axis <
numDims; ++axis)
1330 dims[axis] = _mdMap->getLocalDim(axis,
true);
1334 _mdArrayView = _mdArrayRcp();
1339 const_iterator src = source.
getData().begin();
1340 for (iterator trg = _mdArrayView.
begin(); trg != _mdArrayView.
end(); ++trg)
1346 #ifdef DOMI_MDVECTOR_VERBOSE
1349 cout <<
"Inside MDVector copy constructor with view access"
1357 template<
class Scalar >
1359 MDVector(
const Teuchos::RCP<
const Teuchos::Comm< int > > teuchosComm,
1360 Teuchos::ParameterList & plist) :
1361 _teuchosComm(teuchosComm),
1372 setObjectLabel(
"Domi::MDVector");
1375 MDMap * myMdMap =
new MDMap(teuchosComm, plist);
1376 dim_type leadingDim = plist.get(
"leading dimension" , 0);
1377 dim_type trailingDim = plist.get(
"trailing dimension", 0);
1378 if (leadingDim + trailingDim > 0)
1384 _mdMap = Teuchos::rcp(myMdMap);
1387 int numDims = _mdMap->numDims();
1388 Teuchos::Array< dim_type > dims(numDims);
1389 for (
int axis = 0; axis <
numDims; ++axis)
1390 dims[axis] = _mdMap->getLocalDim(axis,
true);
1393 _mdArrayRcp.
resize(dims);
1394 _mdArrayView = _mdArrayRcp();
1399 template<
class Scalar >
1402 Teuchos::ParameterList & plist) :
1403 _teuchosComm(mdComm->getTeuchosComm()),
1404 _mdMap(Teuchos::rcp(new
MDMap(mdComm, plist))),
1414 setObjectLabel(
"Domi::MDVector");
1418 dim_type leadingDim = plist.get(
"leading dimension" , 0);
1419 dim_type trailingDim = plist.get(
"trailing dimension", 0);
1420 if (leadingDim + trailingDim > 0)
1426 _mdMap = Teuchos::rcp(myMdMap);
1429 int numDims = _mdMap->numDims();
1430 Teuchos::Array< dim_type > dims(numDims);
1431 for (
int axis = 0; axis <
numDims; ++axis)
1432 dims[axis] = _mdMap->getLocalDim(axis,
true);
1435 _mdArrayRcp.
resize(dims);
1436 _mdArrayView = _mdArrayRcp();
1441 template<
class Scalar >
1445 dim_type globalIndex) :
1446 _teuchosComm(parent._teuchosComm),
1448 _mdArrayRcp(parent._mdArrayRcp),
1449 _mdArrayView(parent._mdArrayView),
1457 setObjectLabel(
"Domi::MDVector");
1460 Teuchos::RCP< const MDMap > parentMdMap = parent.
getMDMap();
1463 _mdMap = Teuchos::rcp(
new MDMap(*parentMdMap,
1468 if (_mdMap->onSubcommunicator())
1475 dim_type origin = parentMdMap->getGlobalRankBounds(axis,
false).start() -
1476 parentMdMap->getLowerPadSize(axis);
1481 dim_type localIndex = globalIndex - origin;
1485 _mdArrayView = newView;
1491 _mdArrayRcp.
clear();
1498 template<
class Scalar >
1502 const Slice & slice,
1506 _mdArrayRcp(parent._mdArrayRcp),
1507 _mdArrayView(parent._mdArrayView),
1515 #ifdef DOMI_MDVECTOR_VERBOSE
1516 cout <<
"slice axis " << axis << endl;
1517 cout <<
" _mdArrayRcp @ " << _mdArrayRcp.
getRawPtr() << endl;
1518 cout <<
" _mdArrayView @ " << _mdArrayView.
getRawPtr() << endl;
1521 setObjectLabel(
"Domi::MDVector");
1524 Teuchos::RCP< const MDMap > parentMdMap = parent.
getMDMap();
1527 _mdMap = Teuchos::rcp(
new MDMap(*parentMdMap,
1531 _teuchosComm = _mdMap->getTeuchosComm();
1534 if (_mdMap->onSubcommunicator())
1537 Slice bounds = slice.
bounds(parentMdMap->getGlobalDim(axis,
true));
1544 dim_type origin = parentMdMap->getGlobalRankBounds(axis,
false).
start() -
1545 parentMdMap->getLowerPadSize(axis);
1552 dim_type start = std::max(0, bounds.
start() - origin - bndryPad);
1559 dim_type stop = std::min(bounds.
stop() - origin + bndryPad,
1560 parentMdMap->getLocalDim(axis,
true));
1564 _mdArrayView = newView;
1570 _mdArrayRcp.
clear();
1573 #ifdef DOMI_MDVECTOR_VERBOSE
1574 cout <<
" _mdArrayView @ " << _mdArrayView.
getRawPtr() << endl;
1575 cout <<
" offset = " << int(_mdArrayView.
getRawPtr() -
1582 template<
class Scalar >
1585 const Teuchos::ArrayView< Slice > & slices,
1586 const Teuchos::ArrayView< int > & bndryPad)
1588 setObjectLabel(
"Domi::MDVector");
1591 int numDims = parent.
numDims();
1594 TEUCHOS_TEST_FOR_EXCEPTION(
1595 (slices.size() != numDims),
1597 "number of slices = " << slices.size() <<
" != parent MDVector number of "
1598 "dimensions = " << numDims);
1602 for (
int axis = 0; axis < numDims; ++axis)
1604 int bndryPadding = (axis < bndryPad.size()) ? bndryPad[axis] : 0;
1609 tempMDVector1 = tempMDVector2;
1611 *
this = tempMDVector1;
1616 template<
class Scalar >
1621 _teuchosComm = source._teuchosComm;
1622 _mdMap = source._mdMap;
1623 _mdArrayRcp = source._mdArrayRcp;
1624 _mdArrayView = source._mdArrayView;
1625 _nextAxis = source._nextAxis;
1627 _requests = source._requests;
1629 _sendMessages = source._sendMessages;
1630 _recvMessages = source._recvMessages;
1636 template<
class Scalar >
1643 template<
class Scalar >
1644 const Teuchos::RCP< const MDMap >
1653 template<
class Scalar >
1658 return _mdMap->onSubcommunicator();
1663 template<
class Scalar >
1664 Teuchos::RCP< const Teuchos::Comm< int > >
1668 return _mdMap->getTeuchosComm();
1673 template<
class Scalar >
1678 return _mdMap->numDims();
1683 template<
class Scalar >
1688 return _mdMap->getCommDim(axis);
1693 template<
class Scalar >
1698 return _mdMap->isPeriodic(axis);
1703 template<
class Scalar >
1708 return _mdMap->getCommIndex(axis);
1713 template<
class Scalar >
1718 return _mdMap->getLowerNeighbor(axis);
1723 template<
class Scalar >
1728 return _mdMap->getUpperNeighbor(axis);
1733 template<
class Scalar >
1738 return _mdMap->getGlobalDim(axis, withBndryPad);
1743 template<
class Scalar >
1748 return _mdMap->getGlobalBounds(axis, withBndryPad);
1753 template<
class Scalar >
1758 return _mdMap->getGlobalRankBounds(axis, withBndryPad);
1763 template<
class Scalar >
1768 return _mdMap->getLocalDim(axis, withCommPad);
1773 template<
class Scalar >
1774 Teuchos::ArrayView< const Slice >
1778 return _mdMap->getLocalBounds();
1783 template<
class Scalar >
1788 return _mdMap->getLocalBounds(axis, withCommPad);
1793 template<
class Scalar >
1798 return _mdMap->getLocalInteriorBounds(axis);
1803 template<
class Scalar >
1808 return _mdMap->hasPadding();
1813 template<
class Scalar >
1818 return _mdMap->getLowerPadSize(axis);
1823 template<
class Scalar >
1828 return _mdMap->getUpperPadSize(axis);
1833 template<
class Scalar >
1838 return _mdMap->getCommPadSize(axis);
1843 template<
class Scalar >
1848 return _mdMap->getLowerBndryPad(axis);
1853 template<
class Scalar >
1858 return _mdMap->getUpperBndryPad(axis);
1863 template<
class Scalar >
1868 return _mdMap->getBndryPadSize(axis);
1873 template<
class Scalar >
1885 template<
class Scalar >
1897 template<
class Scalar >
1902 return _mdMap->isReplicatedBoundary(axis);
1907 template<
class Scalar >
1912 return _mdMap->getLayout();
1917 template<
class Scalar >
1922 return _mdMap->isContiguous();
1936 template<
class Scalar >
1937 Teuchos::RCP< Epetra_IntVector >
1942 TEUCHOS_TEST_FOR_EXCEPTION(
1943 typeid(Scalar) !=
typeid(
int),
1945 "MDVector is of scalar type '" <<
typeid(Scalar).name() <<
"', but "
1946 "Epetra_IntVector requires scalar type 'int'");
1949 TEUCHOS_TEST_FOR_EXCEPTION(
1952 "This MDVector's MDMap is non-contiguous. This can happen when you take "
1953 "a slice of a parent MDVector.");
1956 Teuchos::RCP< const Epetra_Map > epetraMap = _mdMap->getEpetraMap(
true);
1961 void * buffer = (
void*) _mdArrayView.getRawPtr();
1962 return Teuchos::rcp(
new Epetra_IntVector(View,
1976 template<
class Scalar >
1977 Teuchos::RCP< Epetra_Vector >
1978 MDVector< Scalar >::
1979 getEpetraVectorView()
const
1982 TEUCHOS_TEST_FOR_EXCEPTION(
1983 typeid(Scalar) !=
typeid(
double),
1985 "MDVector is of scalar type '" <<
typeid(Scalar).name() <<
"', but "
1986 "Epetra_Vector requires scalar type 'double'");
1989 TEUCHOS_TEST_FOR_EXCEPTION(
1991 MDMapNoncontiguousError,
1992 "This MDVector's MDMap is non-contiguous. This can happen when you take "
1993 "a slice of a parent MDVector.");
1996 Teuchos::RCP< const Epetra_Map > epetraMap = _mdMap->getEpetraMap(
true);
2001 void * buffer = (
void*) _mdArrayView.getRawPtr();
2002 return Teuchos::rcp(
new Epetra_Vector(View,
2016 template<
class Scalar >
2017 Teuchos::RCP< Epetra_MultiVector >
2018 MDVector< Scalar >::
2019 getEpetraMultiVectorView()
const
2022 TEUCHOS_TEST_FOR_EXCEPTION(
2023 typeid(Scalar) !=
typeid(
double),
2025 "MDVector is of scalar type '" <<
typeid(Scalar).name() <<
"', but "
2026 "Epetra_Vector requires scalar type 'double'");
2029 int vectorAxis = (getLayout() == C_ORDER) ? 0 : numDims()-1;
2030 int padding = getLowerPadSize(vectorAxis) + getUpperPadSize(vectorAxis);
2031 int commDim = getCommDim(vectorAxis);
2032 int numVectors = getGlobalDim(vectorAxis);
2035 Teuchos::RCP< const MDMap > newMdMap;
2036 if (padding == 0 && commDim == 1)
2037 newMdMap = Teuchos::rcp(
new MDMap(*_mdMap, vectorAxis, 0));
2043 TEUCHOS_TEST_FOR_EXCEPTION(
2044 ! newMdMap->isContiguous(),
2045 MDMapNoncontiguousError,
2046 "This MDVector's MDMap is non-contiguous. This can happen when you take "
2047 "a slice of a parent MDVector.");
2052 size_type stride = newMdMap->getLocalDim(0,
true);
2053 for (
int axis = 1; axis < newMdMap->numDims(); ++axis)
2054 stride *= newMdMap->getLocalDim(axis,
true);
2055 TEUCHOS_TEST_FOR_EXCEPTION(
2056 stride*numVectors > Teuchos::OrdinalTraits<int>::max(),
2058 "Buffer size " << stride*numVectors <<
" is too large for Epetra int "
2060 int lda = (int)stride;
2063 Teuchos::RCP< const Epetra_Map > epetraMap = newMdMap->getEpetraMap(
true);
2068 void * buffer = (
void*) _mdArrayView.getRawPtr();
2069 return Teuchos::rcp(
new Epetra_MultiVector(View,
2078 template<
class Scalar >
2079 Teuchos::RCP< Epetra_IntVector >
2080 MDVector< Scalar >::
2081 getEpetraIntVectorCopy()
const
2083 typedef typename MDArrayView< const Scalar >::iterator iterator;
2086 Teuchos::RCP< const Epetra_Map > epetraMap = _mdMap->getEpetraMap(
true);
2089 Teuchos::RCP< Epetra_IntVector > result =
2090 Teuchos::rcp(
new Epetra_IntVector(*epetraMap));
2095 for (iterator it = _mdArrayView.begin(); it != _mdArrayView.end(); ++it)
2096 result->operator[](ii++) = (int) *it;
2104 template<
class Scalar >
2105 Teuchos::RCP< Epetra_Vector >
2106 MDVector< Scalar >::
2107 getEpetraVectorCopy()
const
2109 typedef typename MDArrayView< const Scalar >::iterator iterator;
2112 Teuchos::RCP< const Epetra_Map > epetraMap = _mdMap->getEpetraMap(
true);
2115 Teuchos::RCP< Epetra_Vector > result =
2116 Teuchos::rcp(
new Epetra_Vector(*epetraMap));
2121 for (iterator it = _mdArrayView.begin(); it != _mdArrayView.end(); ++it)
2122 result->operator[](ii++) = (double) *it;
2130 template<
class Scalar >
2131 Teuchos::RCP< Epetra_MultiVector >
2132 MDVector< Scalar >::
2133 getEpetraMultiVectorCopy()
const
2135 typedef typename MDArrayView< Scalar >::iterator iterator;
2136 typedef typename MDArrayView< const Scalar >::iterator citerator;
2139 int vectorAxis = (getLayout() == C_ORDER) ? 0 : numDims()-1;
2140 int padding = getLowerPadSize(vectorAxis) + getUpperPadSize(vectorAxis);
2141 int commDim = getCommDim(vectorAxis);
2142 int numVectors = getGlobalDim(vectorAxis);
2145 Teuchos::RCP< const MDMap > newMdMap;
2146 if (padding == 0 && commDim == 1)
2147 newMdMap = Teuchos::rcp(
new MDMap(*_mdMap, vectorAxis, 0));
2155 Teuchos::RCP< const Epetra_Map > epetraMap = newMdMap->getEpetraMap(
true);
2158 Teuchos::RCP< Epetra_MultiVector > result =
2159 Teuchos::rcp(
new Epetra_MultiVector(*epetraMap, numVectors));
2164 if (numVectors == 1)
2166 for (citerator it = _mdArrayView.begin(); it != _mdArrayView.end(); ++it)
2167 result->operator[](0)[ii++] = (double) *it;
2171 for (
int iv = 0; iv < numVectors; ++iv)
2175 for (iterator it = subArray.begin(); it != subArray.end(); ++it)
2176 result->operator[](iv)[ii++] = (double) *it;
2190 template<
class Scalar >
2191 template<
class LocalOrdinal,
2192 class GlobalOrdinal >
2193 Teuchos::RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal > >
2194 MDVector< Scalar >::
2195 getTpetraVectorView()
const
2198 TEUCHOS_TEST_FOR_EXCEPTION(
2200 MDMapNoncontiguousError,
2201 "This MDVector's MDMap is non-contiguous. This can happen when you take "
2202 "a slice of a parent MDVector.");
2205 Teuchos::RCP<
const Tpetra::Map< LocalOrdinal,
2206 GlobalOrdinal > > tpetraMap =
2207 _mdMap->template getTpetraMap< LocalOrdinal, GlobalOrdinal >(
true);
2210 return Teuchos::rcp(
new Tpetra::Vector< Scalar,
2212 GlobalOrdinal >(tpetraMap,
2213 _mdArrayView.arrayView()));
2218 template<
class Scalar >
2219 template<
class LocalOrdinal,
2220 class GlobalOrdinal >
2221 Teuchos::RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal > >
2222 MDVector< Scalar >::
2223 getTpetraVectorCopy()
const
2225 typedef typename MDArrayView< const Scalar >::iterator iterator;
2228 Teuchos::RCP<
const Tpetra::Map< LocalOrdinal,
2229 GlobalOrdinal > > tpetraMap =
2230 _mdMap->template getTpetraMap< LocalOrdinal, GlobalOrdinal >(
true);
2233 Teuchos::RCP< Tpetra::Vector< Scalar,
2235 GlobalOrdinal > > result =
2236 Teuchos::rcp(
new Tpetra::Vector< Scalar,
2238 GlobalOrdinal >(tpetraMap) );
2242 Teuchos::ArrayRCP< Scalar > tpetraVectorArray = result->getDataNonConst();
2244 for (iterator it = _mdArrayView.begin(); it != _mdArrayView.end(); ++it)
2245 tpetraVectorArray[ii++] = *it;
2253 template<
class Scalar >
2254 template <
class LocalOrdinal,
2255 class GlobalOrdinal >
2256 Teuchos::RCP< Tpetra::MultiVector< Scalar,
2259 MDVector< Scalar >::
2260 getTpetraMultiVectorView()
const
2263 int vectorAxis = (getLayout() == C_ORDER) ? 0 : numDims()-1;
2264 int padding = getLowerPadSize(vectorAxis) + getUpperPadSize(vectorAxis);
2265 int commDim = getCommDim(vectorAxis);
2266 size_t numVectors = getGlobalDim(vectorAxis);
2269 Teuchos::RCP< const MDMap > newMdMap;
2270 if (padding == 0 && commDim == 1)
2271 newMdMap = Teuchos::rcp(
new MDMap(*_mdMap, vectorAxis, 0));
2277 TEUCHOS_TEST_FOR_EXCEPTION(
2278 ! newMdMap->isContiguous(),
2279 MDMapNoncontiguousError,
2280 "This MDVector's MDMap is non-contiguous. This can happen when you take "
2281 "a slice of a parent MDVector.");
2286 size_type stride = newMdMap->getLocalDim(0,
true);
2287 for (
int axis = 1; axis < newMdMap->numDims(); ++axis)
2288 stride *= newMdMap->getLocalDim(axis,
true);
2289 TEUCHOS_TEST_FOR_EXCEPTION(
2290 (
long long int)(stride*numVectors) > Teuchos::OrdinalTraits<GlobalOrdinal>::max(),
2292 "Buffer size " << stride*numVectors <<
" is too large for Tpetra "
2293 "GlobalOrdinal = " <<
typeid(GlobalOrdinal).name() );
2294 size_t lda = (size_t)stride;
2297 Teuchos::RCP<
const Tpetra::Map< LocalOrdinal,
2298 GlobalOrdinal > > tpetraMap =
2299 newMdMap->template getTpetraMap< LocalOrdinal, GlobalOrdinal >(
true);
2302 return Teuchos::rcp(
new Tpetra::MultiVector< Scalar,
2304 GlobalOrdinal >(tpetraMap,
2305 _mdArrayView.arrayView(),
2312 template<
class Scalar >
2313 template <
class LocalOrdinal,
2314 class GlobalOrdinal >
2315 Teuchos::RCP< Tpetra::MultiVector< Scalar,
2318 MDVector< Scalar >::
2319 getTpetraMultiVectorCopy()
const
2321 typedef typename MDArrayView< Scalar >::iterator iterator;
2322 typedef typename MDArrayView< const Scalar >::iterator citerator;
2325 int vectorAxis = (getLayout() == C_ORDER) ? 0 : numDims()-1;
2326 int padding = getLowerPadSize(vectorAxis) + getUpperPadSize(vectorAxis);
2327 int commDim = getCommDim(vectorAxis);
2328 int numVectors = getGlobalDim(vectorAxis);
2331 Teuchos::RCP< const MDMap > newMdMap;
2332 if (padding == 0 && commDim == 1)
2333 newMdMap = Teuchos::rcp(
new MDMap(*_mdMap, vectorAxis, 0));
2341 Teuchos::RCP<
const Tpetra::Map< LocalOrdinal,
2342 GlobalOrdinal > > tpetraMap =
2343 newMdMap->template getTpetraMap< LocalOrdinal, GlobalOrdinal >(
true);
2346 Teuchos::RCP< Tpetra::MultiVector< Scalar,
2348 GlobalOrdinal > > result =
2349 Teuchos::rcp(
new Tpetra::MultiVector< Scalar,
2351 GlobalOrdinal >(tpetraMap, numVectors));
2356 if (numVectors == 1)
2358 Teuchos::ArrayRCP< Scalar > tpetraVectorArray = result->getDataNonConst(0);
2359 for (citerator it = _mdArrayView.begin(); it != _mdArrayView.end(); ++it)
2360 tpetraVectorArray[ii++] = (
double) *it;
2364 for (
int iv = 0; iv < numVectors; ++iv)
2367 Teuchos::ArrayRCP< Scalar > tpetraVectorArray =
2368 result->getDataNonConst(iv);
2370 for (iterator it = subArray.begin(); it != subArray.end(); ++it)
2371 tpetraVectorArray[ii++] = (
double) *it;
2383 template<
class Scalar >
2388 #ifdef DOMI_MDVECTOR_VERBOSE
2389 cout <<
"_mdArrayView.getRawPtr() = " << _mdArrayView.
getRawPtr()
2391 cout <<
"_mdArrayView = " << _mdArrayView << endl;
2394 return _mdArrayView;
2396 for (
int axis = 0; axis < numDims(); ++axis)
2398 int lo = getLowerPadSize(axis);
2399 int hi = getLocalDim(axis,
true) - getUpperPadSize(axis);
2407 template<
class Scalar >
2413 return _mdArrayView.getConst();
2415 for (
int axis = 0; axis < numDims(); ++axis)
2417 int lo = getLowerPadSize(axis);
2418 int hi = getLocalDim(axis,
true) - getUpperPadSize(axis);
2426 template<
class Scalar >
2433 Slice(getLowerPadSize(axis)));
2434 return newArrayView;
2439 template<
class Scalar >
2446 Slice(getLowerPadSize(axis)));
2447 return newArrayView;
2452 template<
class Scalar >
2457 dim_type n = getLocalDim(axis,
true);
2458 int pad = getUpperPadSize(axis);
2460 if (pad) slice =
Slice(n-pad,n);
2461 else slice =
Slice(n-1,n-1);
2465 return newArrayView;
2470 template<
class Scalar >
2475 dim_type n = getLocalDim(axis,
true);
2476 int pad = getUpperPadSize(axis);
2478 if (pad) slice =
Slice(n-pad,n);
2479 else slice =
Slice(n-1,n-1);
2483 return newArrayView;
2488 template<
class Scalar >
2495 TEUCHOS_TEST_FOR_EXCEPTION(
2496 ! _mdMap->isCompatible(*(a._mdMap)),
2498 "MDMap of calling MDVector and argument 'a' are incompatible");
2501 Scalar local_dot = 0;
2502 iterator a_it = aView.begin();
2503 for (iterator it = _mdArrayView.begin(); it != _mdArrayView.end();
2505 local_dot += *it * *a_it;
2506 Scalar global_dot = 0;
2507 Teuchos::reduceAll(*_teuchosComm,
2508 Teuchos::REDUCE_SUM,
2517 template<
class Scalar >
2518 typename Teuchos::ScalarTraits< Scalar >::magnitudeType
2522 typedef typename Teuchos::ScalarTraits< Scalar >::magnitudeType mag;
2525 mag local_norm1 = 0;
2526 for (iterator it = _mdArrayView.begin(); it != _mdArrayView.end(); ++it)
2527 local_norm1 += std::abs(*it);
2528 mag global_norm1 = 0;
2529 Teuchos::reduceAll(*_teuchosComm,
2530 Teuchos::REDUCE_SUM,
2534 return global_norm1;
2539 template<
class Scalar >
2540 typename Teuchos::ScalarTraits< Scalar >::magnitudeType
2544 typedef typename Teuchos::ScalarTraits< Scalar >::magnitudeType mag;
2545 mag norm2 = dot(*
this);
2546 return Teuchos::ScalarTraits<mag>::squareroot(norm2);
2551 template<
class Scalar >
2552 typename Teuchos::ScalarTraits< Scalar >::magnitudeType
2556 typedef typename Teuchos::ScalarTraits< Scalar >::magnitudeType mag;
2559 mag local_normInf = 0;
2560 for (iterator it = _mdArrayView.begin(); it != _mdArrayView.end(); ++it)
2561 local_normInf = std::max(local_normInf, std::abs(*it));
2562 mag global_normInf = 0;
2563 Teuchos::reduceAll(*_teuchosComm,
2564 Teuchos::REDUCE_MAX,
2568 return global_normInf;
2573 template<
class Scalar >
2574 typename Teuchos::ScalarTraits< Scalar >::magnitudeType
2578 typedef typename Teuchos::ScalarTraits< Scalar >::magnitudeType mag;
2581 TEUCHOS_TEST_FOR_EXCEPTION(
2582 ! _mdMap->isCompatible(*(weights._mdMap)),
2584 "MDMap of calling MDVector and argument 'weights' are incompatible");
2587 mag local_wNorm = 0;
2588 iterator w_it = wView.begin();
2589 for (iterator it = _mdArrayView.begin(); it != _mdArrayView.end();
2591 local_wNorm += *it * *it * *w_it;
2592 mag global_wNorm = 0;
2593 Teuchos::reduceAll(*_teuchosComm,
2594 Teuchos::REDUCE_SUM,
2598 Teuchos::Array< dim_type > dimensions(numDims());
2599 for (
int i = 0; i < numDims(); ++i)
2600 dimensions[i] = _mdMap->getGlobalDim(i);
2601 size_type n = computeSize(dimensions);
2602 if (n == 0)
return 0;
2603 return Teuchos::ScalarTraits<mag>::squareroot(global_wNorm / n);
2608 template<
class Scalar >
2613 typedef typename Teuchos::ScalarTraits< Scalar >::magnitudeType mag;
2617 for (iterator it = _mdArrayView.begin(); it != _mdArrayView.end(); ++it)
2620 Teuchos::reduceAll(*_teuchosComm,
2621 Teuchos::REDUCE_SUM,
2625 Teuchos::Array< dim_type > dimensions(numDims());
2626 for (
int i = 0; i < numDims(); ++i)
2627 dimensions[i] = _mdMap->getGlobalDim(i);
2628 size_type n = computeSize(dimensions);
2629 if (n == 0)
return 0;
2630 return global_sum / n;
2635 template<
class Scalar >
2640 using Teuchos::TypeNameTraits;
2642 Teuchos::Array< dim_type > dims(numDims());
2643 for (
int axis = 0; axis < numDims(); ++axis)
2644 dims[axis] = getGlobalDim(axis,
true);
2646 std::ostringstream oss;
2647 oss <<
"\"Domi::MDVector\": {"
2648 <<
"Template parameters: {"
2649 <<
"Scalar: " << TypeNameTraits<Scalar>::name()
2651 if (this->getObjectLabel() !=
"")
2652 oss <<
", Label: \"" << this->getObjectLabel () <<
"\", ";
2653 oss <<
"Global dimensions: " << dims <<
" }";
2659 template<
class Scalar >
2663 const Teuchos::EVerbosityLevel verbLevel)
const
2666 using Teuchos::Comm;
2668 using Teuchos::TypeNameTraits;
2669 using Teuchos::VERB_DEFAULT;
2670 using Teuchos::VERB_NONE;
2671 using Teuchos::VERB_LOW;
2672 using Teuchos::VERB_MEDIUM;
2673 using Teuchos::VERB_HIGH;
2674 using Teuchos::VERB_EXTREME;
2676 const Teuchos::EVerbosityLevel vl =
2677 (verbLevel == VERB_DEFAULT) ? VERB_LOW : verbLevel;
2679 const MDMap & mdMap = *(getMDMap());
2680 Teuchos::RCP< const Teuchos::Comm< int > > comm = mdMap.
getTeuchosComm();
2681 const int myImageID = comm->getRank();
2682 const int numImages = comm->getSize();
2683 Teuchos::OSTab tab0(out);
2685 if (vl != VERB_NONE)
2689 out <<
"\"Domi::MDVector\":" << endl;
2691 Teuchos::OSTab tab1(out);
2694 out <<
"Template parameters:";
2696 Teuchos::OSTab tab2(out);
2697 out <<
"Scalar: " << TypeNameTraits<Scalar>::name() << endl;
2700 if (this->getObjectLabel() !=
"")
2702 out <<
"Label: \"" << getObjectLabel() <<
"\"" << endl;
2704 Teuchos::Array< dim_type > globalDims(numDims());
2705 for (
int axis = 0; axis < numDims(); ++axis)
2706 globalDims[axis] = getGlobalDim(axis,
true);
2707 out <<
"Global dimensions: " << globalDims << endl;
2709 for (
int imageCtr = 0; imageCtr < numImages; ++imageCtr)
2711 if (myImageID == imageCtr)
2716 out <<
"Process: " << myImageID << endl;
2717 Teuchos::OSTab tab2(out);
2718 Teuchos::Array< dim_type > localDims(numDims());
2719 for (
int axis = 0; axis < numDims(); ++axis)
2720 localDims[axis] = getLocalDim(axis,
true);
2721 out <<
"Local dimensions: " << localDims << endl;
2734 template<
class Scalar >
2738 bool includePadding)
2743 for (iterator it = newArray.
begin(); it != newArray.
end(); ++it)
2749 template<
class Scalar >
2756 Teuchos::ScalarTraits< Scalar >::seedrandom(time(NULL));
2757 for (iterator it = _mdArrayView.begin(); it != _mdArrayView.end(); ++it)
2758 *it = Teuchos::ScalarTraits< Scalar >::random();
2763 template<
class Scalar >
2768 for (
int axis = 0; axis < numDims(); ++axis)
2770 updateCommPad(axis);
2776 template<
class Scalar >
2781 startUpdateCommPad(axis);
2782 endUpdateCommPad(axis);
2787 template<
class Scalar >
2796 if (_sendMessages.empty()) initializeMessages();
2799 int rank = _teuchosComm->getRank();
2800 int numProc = _teuchosComm->getSize();
2806 Teuchos::RCP< const Teuchos::MpiComm< int > > mpiComm =
2807 Teuchos::rcp_dynamic_cast<
const Teuchos::MpiComm< int > >(_teuchosComm);
2808 const Teuchos::OpaqueWrapper< MPI_Comm > & communicator =
2809 *(mpiComm->getRawMpiComm());
2812 MPI_Request request;
2813 for (
int boundary = 0; boundary < 2; ++boundary)
2815 MessageInfo message = _sendMessages[axis][boundary];
2816 if (message.proc >= 0)
2818 tag = 2 * (rank * numProc + message.proc) + boundary;
2820 #ifdef DOMI_MDVECTOR_OUTPUT_UPDATECOMMPAD
2821 cout << rank <<
": post send for axis " << axis <<
", boundary "
2822 << boundary <<
", buffer = " << message.buffer <<
", proc = "
2823 << message.proc <<
", tag = " << tag << endl;
2826 if (MPI_Isend(message.buffer,
2828 *(message.datatype),
2833 throw std::runtime_error(
"Domi::MDVector: Error in MPI_Isend");
2834 _requests.push_back(request);
2839 for (
int boundary = 0; boundary < 2; ++boundary)
2841 MessageInfo message = _recvMessages[axis][boundary];
2842 if (message.proc >= 0)
2844 tag = 2 * (message.proc * numProc + rank) + (1-boundary);
2846 #ifdef DOMI_MDVECTOR_OUTPUT_UPDATECOMMPAD
2847 cout << rank <<
": post recv for axis " << axis <<
", boundary "
2848 << boundary <<
", buffer = " << message.buffer <<
", proc = "
2849 << message.proc <<
", tag = " << tag << endl;
2852 if (MPI_Irecv(message.buffer,
2854 *(message.datatype),
2859 throw std::runtime_error(
"Domi::MDVector: Error in MPI_Irecv");
2860 _requests.push_back(request);
2867 if (isPeriodic(axis))
2869 for (
int sendBndry = 0; sendBndry < 2; ++sendBndry)
2871 int recvBndry = 1 - sendBndry;
2881 for ( ; it_recv != recvView.
end(); ++it_recv, ++it_send)
2882 *it_recv = *it_send;
2890 template<
class Scalar >
2896 if (_requests.size() > 0)
2898 Teuchos::Array< MPI_Status > status(_requests.size());
2899 if (MPI_Waitall(_requests.size(),
2902 throw std::runtime_error(
"Domi::MDVector: Error in MPI_Waitall");
2910 template<
class Scalar >
2918 int newAxis = _nextAxis + 1;
2919 if (newAxis >= result.
numDims()) newAxis = 0;
2920 result._nextAxis = newAxis;
2926 template<
class Scalar >
2934 int newAxis = _nextAxis + 1;
2935 if (newAxis >= result.
numDims()) newAxis = 0;
2936 result._nextAxis = newAxis;
2942 template<
class Scalar >
2949 int ndims = numDims();
2950 Teuchos::Array<int> sizes(ndims);
2951 Teuchos::Array<int> subsizes(ndims);
2952 Teuchos::Array<int> starts(ndims);
2953 MessageInfo messageInfo;
2955 _sendMessages.resize(ndims);
2956 _recvMessages.resize(ndims);
2958 #ifdef DOMI_MDVECTOR_MESSAGE_INITIALIZE
2959 std::stringstream msg;
2960 int rank = _teuchosComm->getRank();
2964 int order = mpiOrder(getLayout());
2965 MPI_Datatype datatype = mpiType< Scalar >();
2969 for (
int msgAxis = 0; msgAxis < ndims; ++msgAxis)
2972 for (
int axis = 0; axis < ndims; ++axis)
2974 sizes[axis] = _mdArrayRcp.dimension(axis);
2975 subsizes[axis] = _mdArrayView.dimension(axis);
2983 int proc = getLowerNeighbor(msgAxis);
2985 #ifdef DOMI_MDVECTOR_MESSAGE_INITIALIZE
2986 msg << endl <<
"P" << rank <<
": axis " << msgAxis <<
", lower neighbor = "
2991 subsizes[msgAxis] = getLowerPadSize(msgAxis);
2993 if (subsizes[msgAxis] == 0) proc = -1;
2995 messageInfo.buffer = (
void*) getData().getRawPtr();
2996 messageInfo.proc = proc;
2997 messageInfo.axis = msgAxis;
3007 #ifdef DOMI_MDVECTOR_MESSAGE_INITIALIZE
3008 msg << endl <<
"P" << rank <<
": axis " << msgAxis
3009 <<
", Lower receive message:" << endl <<
" ndims = " << ndims
3010 << endl <<
" sizes = " << sizes << endl <<
" subsizes = "
3011 << subsizes << endl <<
" starts = " << starts << endl
3012 <<
" order = " << order << endl;
3014 Teuchos::RCP< MPI_Datatype > commPad = Teuchos::rcp(
new MPI_Datatype);
3015 MPI_Type_create_subarray(ndims,
3022 MPI_Type_commit(commPad.get());
3023 messageInfo.datatype = commPad;
3025 messageInfo.dataview = _mdArrayView;
3026 for (
int axis = 0; axis < numDims(); ++axis)
3028 Slice slice(starts[axis], starts[axis] + subsizes[axis]);
3036 _recvMessages[msgAxis][0] = messageInfo;
3042 starts[msgAxis] = getLowerPadSize(msgAxis);
3043 if (isReplicatedBoundary(msgAxis) && getCommIndex(msgAxis) == 0)
3044 starts[msgAxis] += 1;
3050 #ifdef DOMI_MDVECTOR_MESSAGE_INITIALIZE
3051 msg << endl <<
"P" << rank <<
": axis " << msgAxis
3052 <<
", Lower send message:" << endl <<
" ndims = " << ndims
3053 << endl <<
" sizes = " << sizes << endl <<
" subsizes = "
3054 << subsizes << endl <<
" starts = " << starts << endl
3055 <<
" order = " << order << endl;
3057 Teuchos::RCP< MPI_Datatype > commPad = Teuchos::rcp(
new MPI_Datatype);
3058 MPI_Type_create_subarray(ndims,
3065 MPI_Type_commit(commPad.get());
3066 messageInfo.datatype = commPad;
3068 messageInfo.dataview = _mdArrayView;
3069 for (
int axis = 0; axis < numDims(); ++axis)
3071 Slice slice(starts[axis], starts[axis] + subsizes[axis]);
3079 _sendMessages[msgAxis][0] = messageInfo;
3085 proc = getUpperNeighbor(msgAxis);
3087 #ifdef DOMI_MDVECTOR_MESSAGE_INITIALIZE
3088 msg << endl <<
"P" << rank <<
": axis " << msgAxis <<
", upper neighbor = "
3092 subsizes[msgAxis] = getUpperPadSize(msgAxis);
3093 starts[msgAxis] = _mdArrayView.dimension(msgAxis) -
3094 getUpperPadSize(msgAxis);
3095 if (subsizes[msgAxis] == 0) proc = -1;
3096 messageInfo.proc = proc;
3105 #ifdef DOMI_MDVECTOR_MESSAGE_INITIALIZE
3106 msg << endl <<
"P" << rank <<
": axis " << msgAxis
3107 <<
", Upper receive message:" << endl <<
" ndims = " << ndims
3108 << endl <<
" sizes = " << sizes << endl <<
" subsizes = "
3109 << subsizes << endl <<
" starts = " << starts << endl
3110 <<
" order = " << order << endl;
3112 Teuchos::RCP< MPI_Datatype > commPad = Teuchos::rcp(
new MPI_Datatype);
3113 MPI_Type_create_subarray(ndims,
3120 MPI_Type_commit(commPad.get());
3121 messageInfo.datatype = commPad;
3123 messageInfo.dataview = _mdArrayView;
3124 for (
int axis = 0; axis < numDims(); ++axis)
3126 Slice slice(starts[axis], starts[axis] + subsizes[axis]);
3133 _recvMessages[msgAxis][1] = messageInfo;
3139 starts[msgAxis] -= getUpperPadSize(msgAxis);
3140 if (isReplicatedBoundary(msgAxis) &&
3141 getCommIndex(msgAxis) == getCommDim(msgAxis)-1)
3142 starts[msgAxis] -= 1;
3148 #ifdef DOMI_MDVECTOR_MESSAGE_INITIALIZE
3149 msg << endl <<
"P" << rank <<
": axis " << msgAxis
3150 <<
", Upper send message:" << endl <<
" ndims = " << ndims
3151 << endl <<
" sizes = " << sizes << endl <<
" subsizes = "
3152 << subsizes << endl <<
" starts = " << starts << endl
3153 <<
" order = " << order << endl;
3155 Teuchos::RCP< MPI_Datatype > commPad = Teuchos::rcp(
new MPI_Datatype);
3156 MPI_Type_create_subarray(ndims,
3163 MPI_Type_commit(commPad.get());
3164 messageInfo.datatype = commPad;
3166 messageInfo.dataview = _mdArrayView;
3167 for (
int axis = 0; axis < numDims(); ++axis)
3169 Slice slice(starts[axis], starts[axis] + subsizes[axis]);
3177 _sendMessages[msgAxis][1] = messageInfo;
3180 #ifdef DOMI_MDVECTOR_MESSAGE_INITIALIZE
3181 for (
int proc = 0; proc < _teuchosComm->getSize(); ++proc)
3188 _teuchosComm->barrier();
3196 template<
class Scalar >
3200 bool includeBndryPad)
const
3207 int pid = _teuchosComm->getRank();
3210 datafile = fopen(filename.c_str(),
"w");
3213 _teuchosComm->barrier();
3216 const Scalar * buffer = getData(
true).getRawPtr();
3220 Teuchos::RCP< FileInfo > & fileInfo = computeFileInfo(includeBndryPad);
3229 Teuchos::RCP< const Teuchos::MpiComm< int > > mpiComm =
3230 Teuchos::rcp_dynamic_cast<
const Teuchos::MpiComm< int > >(_teuchosComm);
3231 const Teuchos::OpaqueWrapper< MPI_Comm > & communicator =
3232 *(mpiComm->getRawMpiComm());
3235 int access = MPI_MODE_WRONLY | MPI_MODE_CREATE;
3240 char * cstr =
new char[filename.size()+1];
3241 std::strcpy(cstr, filename.c_str());
3246 char datarep[7] =
"native";
3247 MPI_File_open(communicator(), cstr, access, MPI_INFO_NULL, &mpiFile);
3248 MPI_File_set_view(mpiFile, 0, mpiType< Scalar >(),
3249 *(fileInfo->filetype), datarep, MPI_INFO_NULL);
3250 MPI_File_write_all(mpiFile, (
void*)buffer, 1, *(fileInfo->datatype),
3252 MPI_File_close(&mpiFile);
3264 datafile = fopen(filename.c_str(),
"w");
3271 for (const_iterator it = mdArrayView.begin(); it != mdArrayView.end(); ++it)
3273 fwrite((
const void *) &(*it),
sizeof(Scalar), 1, datafile);
3285 template<
class Scalar >
3289 bool includeBndryPad)
3292 const Scalar * buffer = getDataNonConst(
true).getRawPtr();
3296 Teuchos::RCP< FileInfo > & fileInfo = computeFileInfo(includeBndryPad);
3305 Teuchos::RCP< const Teuchos::MpiComm< int > > mpiComm =
3306 Teuchos::rcp_dynamic_cast<
const Teuchos::MpiComm< int > >(_teuchosComm);
3307 const Teuchos::OpaqueWrapper< MPI_Comm > & communicator =
3308 *(mpiComm->getRawMpiComm());
3311 int access = MPI_MODE_RDONLY;
3316 char * cstr =
new char[filename.size()+1];
3317 std::strcpy(cstr, filename.c_str());
3322 char datarep[7] =
"native";
3323 MPI_File_open(communicator(), cstr, access, MPI_INFO_NULL, &mpiFile);
3324 MPI_File_set_view(mpiFile, 0, mpiType< Scalar >(),
3325 *(fileInfo->filetype), datarep, MPI_INFO_NULL);
3326 MPI_File_read_all(mpiFile, (
void*)buffer, 1, *(fileInfo->datatype),
3328 MPI_File_close(&mpiFile);
3337 int ndims = _mdMap->numDims();
3341 datafile = fopen(filename.c_str(),
"r");
3348 Teuchos::Array< Ordinal > index(3);
3349 for (
int axis = 0; axis < ndims; ++axis)
3350 index[axis] = fileInfo->dataStart[axis];
3354 for (iterator it = mdArrayView.begin(); it != mdArrayView.end(); ++it)
3356 fread(&(*it),
sizeof(Scalar), 1, datafile);
3368 template<
class Scalar >
3369 Teuchos::RCP< typename MDVector< Scalar >::FileInfo > &
3375 Teuchos::RCP< MDVector< Scalar >::FileInfo > & fileInfo =
3376 includeBndryPad ? _fileInfoWithBndry : _fileInfo;
3379 if (!fileInfo.is_null())
return fileInfo;
3382 int ndims = _mdMap->
numDims();
3384 fileInfo->fileShape.resize(ndims);
3385 fileInfo->bufferShape.resize(ndims);
3386 fileInfo->dataShape.resize(ndims);
3387 fileInfo->fileStart.resize(ndims);
3388 fileInfo->dataStart.resize(ndims);
3391 for (
int axis = 0; axis < ndims; ++axis)
3395 fileInfo->fileShape[axis] = getGlobalDim(axis,includeBndryPad);
3396 fileInfo->bufferShape[axis] = getLocalDim(axis,
true );
3397 fileInfo->dataShape[axis] = getLocalDim(axis,
false);
3398 fileInfo->fileStart[axis] = getGlobalRankBounds(axis,includeBndryPad).start();
3399 fileInfo->dataStart[axis] = getLocalBounds(axis).start();
3401 if (includeBndryPad)
3403 int commIndex = _mdMap->getCommIndex(axis);
3406 int pad = getLowerBndryPad(axis);
3407 fileInfo->dataShape[axis] += pad;
3408 fileInfo->dataStart[axis] -= pad;
3410 if (commIndex == _mdMap->getCommDim(axis)-1)
3412 fileInfo->dataShape[axis] += getUpperBndryPad(axis);
3417 #ifdef DOMI_MDVECTOR_DEBUG_IO
3418 cout << pid <<
": fileShape = " << fileInfo->fileShape() << endl;
3419 cout << pid <<
": bufferShape = " << fileInfo->bufferShape() << endl;
3420 cout << pid <<
": dataShape = " << fileInfo->dataShape() << endl;
3421 cout << pid <<
": fileStart = " << fileInfo->fileStart() << endl;
3422 cout << pid <<
": dataStart = " << fileInfo->dataStart() << endl;
3426 int mpi_order = getLayout() == C_ORDER ? MPI_ORDER_C : MPI_ORDER_FORTRAN;
3428 fileInfo->filetype = Teuchos::rcp(
new MPI_Datatype);
3429 MPI_Type_create_subarray(ndims,
3430 fileInfo->fileShape.getRawPtr(),
3431 fileInfo->dataShape.getRawPtr(),
3432 fileInfo->fileStart.getRawPtr(),
3434 mpiType< Scalar >(),
3435 fileInfo->filetype.get());
3436 MPI_Type_commit(fileInfo->filetype.get());
3439 fileInfo->datatype = Teuchos::rcp(
new MPI_Datatype);
3440 MPI_Type_create_subarray(ndims,
3441 fileInfo->bufferShape.getRawPtr(),
3442 fileInfo->dataShape.getRawPtr(),
3443 fileInfo->dataStart.getRawPtr(),
3445 mpiType< Scalar >(),
3446 fileInfo->datatype.get());
3447 MPI_Type_commit(fileInfo->datatype.get());
3448 #endif // DGM_PARALLEL
const Teuchos::RCP< const MDMap > getMDMap() const
MDMap accessor method.
Definition: Domi_MDVector.hpp:1646
const_pointer getRawPtr() const
Return a const raw pointer to the beginning of the MDArrayRCP or NULL if unsized. ...
Definition: Domi_MDArrayRCP.hpp:1718
int numDims() const
Return the number of dimensions.
Definition: Domi_MDArrayRCP.hpp:999
Teuchos::ArrayView< const Slice > getLocalBounds() const
Get the local loop bounds along every axis.
Definition: Domi_MDVector.hpp:1776
Slice getGlobalBounds(int axis, bool withBndryPadding=false) const
Get the bounds of the global problem.
Definition: Domi_MDVector.hpp:1746
int getLowerPadSize(int axis) const
Get the size of the lower padding along the given axis.
Definition: Domi_MDVector.hpp:1816
bool hasPadding() const
Return true if there is any padding stored locally.
Definition: Domi_MDVector.hpp:1806
Layout getLayout() const
Get the storage order.
Definition: Domi_MDVector.hpp:1910
void resize(const Teuchos::ArrayView< dim_type > &dims)
Resize the MDArrayRCP based on the given dimensions.
Definition: Domi_MDArrayRCP.hpp:1684
MDArrayView< const Scalar > getData(bool includePadding=true) const
Get a const view of the data as an MDArrayView.
Definition: Domi_MDVector.hpp:2410
int getCommIndex(int axis) const
Get the axis rank of this processor.
Definition: Domi_MDVector.hpp:1706
void endUpdateCommPad(int axis)
Complete an asyncronous update of the communication padding.
Definition: Domi_MDVector.hpp:2893
void assign(const T &value)
Assign a value to all elements of the MDArrayView
Definition: Domi_MDArrayView.hpp:1413
void updateCommPad()
Sum values of a locally replicated multivector across all processes.
Definition: Domi_MDVector.hpp:2766
Teuchos::ScalarTraits< Scalar >::magnitudeType normWeighted(const MDVector< Scalar > &weights) const
Compute the weighted norm of this.
Definition: Domi_MDVector.hpp:2576
void putScalar(const Scalar &value, bool includePadding=true)
Set all values in the multivector with the given value.
Definition: Domi_MDVector.hpp:2737
MDArrayView< Scalar > getLowerPadDataNonConst(int axis)
Get a non-const view of the lower padding data along the given axis as an MDArrayView.
Definition: Domi_MDVector.hpp:2429
dim_type dimension(int axis) const
Return the dimension of the given axis.
Definition: Domi_MDArrayRCP.hpp:1017
iterator begin()
Return the beginning iterator.
Definition: Domi_MDArrayView.hpp:960
int getBndryPadSize(int axis) const
Get the boundary padding size along the given axis.
Definition: Domi_MDVector.hpp:1866
Iterator class suitable for multi-dimensional arrays.
Definition: Domi_MDIterator.hpp:100
iterator end()
Return the ending iterator.
Definition: Domi_MDArrayView.hpp:969
void readBinary(const std::string &filename, bool includeBndryPad=false)
Read the MDVector from a binary file.
Definition: Domi_MDVector.hpp:3288
Scalar meanValue() const
Compute the mean (average) value of this MDVector.
Definition: Domi_MDVector.hpp:2611
int getUpperPadSize(int axis) const
Get the size of the upper padding along the given axis.
Definition: Domi_MDVector.hpp:1826
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to a FancyOStream.
Definition: Domi_MDVector.hpp:2662
int getLowerNeighbor(int axis) const
Get the rank of the lower neighbor.
Definition: Domi_MDVector.hpp:1716
const dim_type stop() const
Stop index.
Definition: Domi_Slice.hpp:438
Teuchos::ScalarTraits< Scalar >::magnitudeType norm1() const
Compute the 1-norm of this MDVector.
Definition: Domi_MDVector.hpp:2520
int getCommPadSize(int axis) const
Get the communication padding size along the given axis.
Definition: Domi_MDVector.hpp:1836
void setUpperPad(int axis, const Scalar value)
Assign all elements of the upper pad a constant value.
Definition: Domi_MDVector.hpp:1888
MDVector< Scalar > operator[](dim_type index) const
Sub-vector access operator.
Definition: Domi_MDVector.hpp:2913
int getUpperNeighbor(int axis) const
Get the rank of the upper neighbor.
Definition: Domi_MDVector.hpp:1726
Teuchos::RCP< const Teuchos::Comm< int > > getTeuchosComm() const
Get the Teuchos communicator.
Definition: Domi_MDVector.hpp:1666
virtual Slice bounds(dim_type size) const
Return a Slice with concrete start and stop values.
dim_type getGlobalDim(int axis, bool withBndryPad=false) const
Get the global dimension along the specified axis.
Definition: Domi_MDVector.hpp:1736
MDArrayView< const Scalar > getUpperPadData(int axis) const
Get a const view of the upper padding data along the given axis as an MDArrayView.
Definition: Domi_MDVector.hpp:2473
MDArrayView< const T > getConst() const
Return an MDArrayView< const T > of an MDArrayView< T > object.
Definition: Domi_MDArrayView.hpp:1055
void startUpdateCommPad(int axis)
Start an asyncronous update of the communication padding.
Definition: Domi_MDVector.hpp:2790
Teuchos::ScalarTraits< Scalar >::magnitudeType norm2() const
Compute the 2-norm of this MDVector.
Definition: Domi_MDVector.hpp:2542
iterator begin()
Return the beginning iterator.
Definition: Domi_MDArrayRCP.hpp:1065
A Slice contains a start, stop, and step index, describing a subset of an ordered container...
Definition: Domi_Slice.hpp:137
void setLowerPad(int axis, const Scalar value)
Assign all elements of the lower pad a constant value.
Definition: Domi_MDVector.hpp:1876
int getUpperBndryPad(int axis) const
Get the size of the upper boundary padding along the given axis.
Definition: Domi_MDVector.hpp:1856
virtual std::string description() const
A simple one-line description of this MDVector.
Definition: Domi_MDVector.hpp:2638
MDArrayView< Scalar > getUpperPadDataNonConst(int axis)
Get a non-const view of the upper padding data along the given axis as an MDArrayView.
Definition: Domi_MDVector.hpp:2455
void clear()
Clear the MDArrayRCP
Definition: Domi_MDArrayRCP.hpp:1652
void randomize()
Set all values in the multivector to pseudorandom numbers.
Definition: Domi_MDVector.hpp:2752
bool isContiguous() const
True if there are no stride gaps on any processor.
Definition: Domi_MDVector.hpp:1920
bool isReplicatedBoundary(int axis) const
Return whether the given axis supports a replicated boundary.
Definition: Domi_MDVector.hpp:1900
MDVector(const Teuchos::RCP< const MDMap > &mdMap, bool zeroOut=true)
Main constructor.
Definition: Domi_MDVector.hpp:1172
MDMap Error exception type.
Definition: Domi_Exceptions.hpp:114
Type Error exception type.
Definition: Domi_Exceptions.hpp:126
bool onSubcommunicator() const
Query whether this processor is on the sub-communicator.
Definition: Domi_MDVector.hpp:1656
Teuchos::ScalarTraits< Scalar >::magnitudeType normInf() const
Compute the infinity-norm of this MDVector.
Definition: Domi_MDVector.hpp:2554
int numDims() const
Get the number of dimensions.
Definition: Domi_MDVector.hpp:1676
bool isPeriodic(int axis) const
Return the periodic flag for the given axis.
Definition: Domi_MDVector.hpp:1696
Multi-dimensional distributed vector.
Definition: Domi_MDVector.hpp:175
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.
MDArrayView< const Scalar > getLowerPadData(int axis) const
Get a const view of the lower padding data along the given axis as an MDArrayView.
Definition: Domi_MDVector.hpp:2442
void writeBinary(const std::string &filename, bool includeBndryPad=false) const
Write the MDVector to a binary file.
Definition: Domi_MDVector.hpp:3199
const dim_type start() const
Start index.
Definition: Domi_Slice.hpp:431
virtual ~MDVector()
Destructor.
Definition: Domi_MDVector.hpp:1637
Teuchos::RCP< const Teuchos::Comm< int > > getTeuchosComm() const
Get the Teuchos communicator.
Definition: Domi_MDMap.hpp:1038
iterator end()
Return the ending iterator.
Definition: Domi_MDArrayRCP.hpp:1074
Scalar dot(const MDVector< Scalar > &a) const
Compute the dot product of this MDVector and MDVector a.
Definition: Domi_MDVector.hpp:2491
Multi-dimensional map.
Definition: Domi_MDMap.hpp:145
Slice getLocalInteriorBounds(int axis) const
Get the local interior looping bounds along the specified axis.
Definition: Domi_MDVector.hpp:1796
dim_type getLocalDim(int axis, bool withCommPad=false) const
Get the local dimension along the specified axis.
Definition: Domi_MDVector.hpp:1766
const T * getRawPtr() const
Return a const raw pointer to the beginning of the MDArrayView or NULL if unsized.
Definition: Domi_MDArrayView.hpp:1521
MDArrayView< Scalar > getDataNonConst(bool includePadding=true)
Get a non-const view of the data as an MDArrayView.
Definition: Domi_MDVector.hpp:2386
int getCommDim(int axis) const
Get the communicator size along the given axis.
Definition: Domi_MDVector.hpp:1686
Slice getGlobalRankBounds(int axis, bool withBndryPad=false) const
Get the global loop bounds on this processor along the specified axis.
Definition: Domi_MDVector.hpp:1756
MDVector< Scalar > & operator=(const MDVector< Scalar > &source)
Assignment operator.
Definition: Domi_MDVector.hpp:1619
int getLowerBndryPad(int axis) const
Get the size of the lower boundary padding along the given axis.
Definition: Domi_MDVector.hpp:1846
Invalid argument exception type.
Definition: Domi_Exceptions.hpp:53
MDMap Error exception type.
Definition: Domi_Exceptions.hpp:102