|
Tpetra Matrix/Vector Services Version of the Day
|
00001 // @HEADER 00002 // *********************************************************************** 00003 // 00004 // Tpetra: Templated Linear Algebra Services Package 00005 // Copyright (2008) Sandia Corporation 00006 // 00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, 00008 // the U.S. Government retains certain rights in this software. 00009 // 00010 // Redistribution and use in source and binary forms, with or without 00011 // modification, are permitted provided that the following conditions are 00012 // met: 00013 // 00014 // 1. Redistributions of source code must retain the above copyright 00015 // notice, this list of conditions and the following disclaimer. 00016 // 00017 // 2. Redistributions in binary form must reproduce the above copyright 00018 // notice, this list of conditions and the following disclaimer in the 00019 // documentation and/or other materials provided with the distribution. 00020 // 00021 // 3. Neither the name of the Corporation nor the names of the 00022 // contributors may be used to endorse or promote products derived from 00023 // this software without specific prior written permission. 00024 // 00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00036 // 00037 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00038 // 00039 // ************************************************************************ 00040 // @HEADER 00041 00042 #ifndef TPETRA_MULTIVECTOR_DEF_HPP 00043 #define TPETRA_MULTIVECTOR_DEF_HPP 00044 00045 #include <Kokkos_NodeTrace.hpp> 00046 00047 #include <Teuchos_Assert.hpp> 00048 #include <Teuchos_as.hpp> 00049 00050 #include "Tpetra_Vector.hpp" 00051 00052 #ifdef DOXYGEN_USE_ONLY 00053 #include "Tpetra_MultiVector_decl.hpp" 00054 #endif 00055 00056 namespace Tpetra { 00057 00058 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00059 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00060 MultiVector (const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& map, 00061 size_t NumVectors, 00062 bool zeroOut) : /* default is true */ 00063 DistObject<Scalar,LocalOrdinal,GlobalOrdinal,Node> (map), 00064 lclMV_ (map->getNode ()) 00065 { 00066 using Teuchos::ArrayRCP; 00067 using Teuchos::RCP; 00068 00069 TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument, 00070 "Tpetra::MultiVector::MultiVector(): NumVectors must be strictly positive."); 00071 const size_t myLen = getLocalLength(); 00072 if (myLen > 0) { 00073 RCP<Node> node = map->getNode(); 00074 ArrayRCP<Scalar> data = node->template allocBuffer<Scalar>(myLen*NumVectors); 00075 MVT::initializeValues(lclMV_,myLen,NumVectors,data,myLen); 00076 if (zeroOut) { 00077 MVT::Init(lclMV_, Teuchos::ScalarTraits<Scalar>::zero()); 00078 } 00079 } 00080 else { 00081 MVT::initializeValues(lclMV_,0,NumVectors,Teuchos::null,0); 00082 } 00083 } 00084 00085 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00086 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00087 MultiVector (const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& source) : 00088 DistObject<Scalar,LocalOrdinal,GlobalOrdinal,Node> (source), 00089 lclMV_ (MVT::getNode (source.lclMV_)) 00090 { 00091 using Teuchos::ArrayRCP; 00092 using Teuchos::RCP; 00093 00094 // copy data from the source MultiVector into this multivector 00095 RCP<Node> node = MVT::getNode(source.lclMV_); 00096 const LocalOrdinal myLen = getLocalLength(); 00097 const size_t numVecs = source.getNumVectors(); 00098 if (myLen > 0) { 00099 // allocate data 00100 ArrayRCP<Scalar> data = node->template allocBuffer<Scalar>(myLen*numVecs); 00101 MVT::initializeValues(lclMV_,myLen,numVecs,data,myLen); 00102 // copy data 00103 { 00104 ArrayRCP<Scalar> dstdata = data; 00105 ArrayRCP<const Scalar> srcdata = MVT::getValues(source.lclMV_); 00106 for (size_t j = 0; j < numVecs; ++j) { 00107 ArrayRCP<const Scalar> srcj = source.getSubArrayRCP(srcdata,j); 00108 KOKKOS_NODE_TRACE("MultiVector::MultiVector(MV)") 00109 node->template copyBuffers<Scalar>(myLen,srcj,dstdata); 00110 dstdata += myLen; 00111 } 00112 } 00113 } 00114 else { 00115 MVT::initializeValues(lclMV_,0,numVecs,Teuchos::null,0); 00116 } 00117 } 00118 00119 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00120 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00121 MultiVector (const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& map, 00122 const Teuchos::ArrayRCP<Scalar>& view, 00123 size_t LDA, 00124 size_t NumVectors, 00125 EPrivateHostViewConstructor /* dummy */) : 00126 DistObject<Scalar,LocalOrdinal,GlobalOrdinal,Node> (map), 00127 lclMV_ (map->getNode ()) 00128 { 00129 using Teuchos::as; 00130 using Teuchos::ArrayRCP; 00131 using Teuchos::RCP; 00132 00133 const std::string tfecfFuncName("MultiVector(view,LDA,NumVector)"); 00134 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(NumVectors < 1, std::invalid_argument, 00135 ": NumVectors must be strictly positive, but you specified NumVectors = " 00136 << NumVectors << "."); 00137 const size_t myLen = getLocalLength(); 00138 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(LDA < myLen, std::invalid_argument, 00139 ": LDA must be large enough to accomodate the local entries."); 00140 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( 00141 as<size_t> (view.size ()) < LDA*(NumVectors - 1) + myLen, 00142 std::invalid_argument, 00143 ": view does not contain enough data to specify the entries in this."); 00144 MVT::initializeValues (lclMV_, myLen, NumVectors, view, myLen); 00145 } 00146 00147 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00148 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00149 MultiVector (const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > &map, 00150 Teuchos::ArrayRCP<Scalar> data, 00151 size_t LDA, 00152 size_t NumVectors, 00153 EPrivateComputeViewConstructor /* dummy */) : 00154 DistObject<Scalar,LocalOrdinal,GlobalOrdinal,Node> (map), 00155 lclMV_ (map->getNode ()) 00156 { 00157 using Teuchos::as; 00158 using Teuchos::ArrayRCP; 00159 using Teuchos::RCP; 00160 00161 const std::string tfecfFuncName("MultiVector(data,LDA,NumVector)"); 00162 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(NumVectors < 1, std::invalid_argument, 00163 ": NumVectors must be strictly positive."); 00164 const size_t myLen = getLocalLength(); 00165 #ifdef HAVE_TPETRA_DEBUG 00166 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( 00167 as<size_t> (data.size ()) < LDA*(NumVectors - 1) + myLen, 00168 std::runtime_error, 00169 ": data does not contain enough data to specify the entries in this."); 00170 #endif 00171 MVT::initializeValues(lclMV_,myLen,NumVectors,data,LDA); 00172 } 00173 00174 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00175 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00176 MultiVector (const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > &map, 00177 Teuchos::ArrayRCP<Scalar> data, 00178 size_t LDA, 00179 Teuchos::ArrayView<const size_t> WhichVectors, 00180 EPrivateComputeViewConstructor /* dummy */) : 00181 DistObject<Scalar,LocalOrdinal,GlobalOrdinal,Node> (map), 00182 lclMV_ (map->getNode ()), 00183 whichVectors_ (WhichVectors) 00184 { 00185 using Teuchos::as; 00186 using Teuchos::ArrayRCP; 00187 using Teuchos::RCP; 00188 00189 const std::string tfecfFuncName("MultiVector(data,LDA,WhichVectors)"); 00190 const size_t myLen = getLocalLength(); 00191 size_t maxVector = *std::max_element(WhichVectors.begin(), WhichVectors.end()); 00192 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(LDA < myLen, std::runtime_error, 00193 ": LDA must be large enough to accomodate the local entries."); 00194 #ifdef HAVE_TPETRA_DEBUG 00195 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( 00196 as<size_t> (data.size ()) < LDA * maxVector + myLen, std::runtime_error, 00197 ": data does not contain enough data to specify the entries in this."); 00198 #endif 00199 if (WhichVectors.size () == 1) { 00200 // shift data so that desired vector is vector 0 00201 maxVector = 0; 00202 data += LDA * WhichVectors[0]; 00203 // kill whichVectors_; we are constant stride 00204 whichVectors_.clear (); 00205 } 00206 if (myLen > 0) { 00207 MVT::initializeValues(lclMV_,myLen,maxVector+1,data,LDA); 00208 } 00209 else { 00210 MVT::initializeValues(lclMV_,0,WhichVectors.size(),Teuchos::null,0); 00211 } 00212 } 00213 00214 00215 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00216 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00217 MultiVector (const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& map, 00218 const Teuchos::ArrayView<const Scalar>& A, 00219 size_t LDA, 00220 size_t NumVectors) : 00221 DistObject<Scalar,LocalOrdinal,GlobalOrdinal,Node> (map), 00222 lclMV_ (map->getNode ()) 00223 { 00224 using Teuchos::as; 00225 using Teuchos::ArrayRCP; 00226 using Teuchos::ArrayView; 00227 using Teuchos::RCP; 00228 00229 const std::string tfecfFuncName("MultiVector(A,LDA)"); 00230 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(NumVectors < 1, std::invalid_argument, 00231 ": NumVectors must be strictly positive."); 00232 const size_t myLen = getLocalLength(); 00233 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(LDA < myLen, std::runtime_error, 00234 ": LDA must be large enough to accomodate the local entries."); 00235 #ifdef HAVE_TPETRA_DEBUG 00236 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( 00237 as<size_t> (A.size ()) < LDA*(NumVectors - 1) + myLen, 00238 std::runtime_error, 00239 ": A does not contain enough data to specify the entries in this."); 00240 #endif 00241 if (myLen > 0) { 00242 RCP<Node> node = MVT::getNode(lclMV_); 00243 ArrayRCP<Scalar> mydata = node->template allocBuffer<Scalar>(myLen*NumVectors); 00244 MVT::initializeValues(lclMV_,myLen,NumVectors,mydata,myLen); 00245 KOKKOS_NODE_TRACE("MultiVector::MultiVector(1D)") 00246 ArrayRCP<Scalar> myview = node->template viewBufferNonConst<Scalar>(Kokkos::WriteOnly,myLen*NumVectors,mydata); 00247 typename ArrayView<const Scalar>::iterator srcit = A.begin(); 00248 for (size_t j = 0; j < NumVectors; ++j) { 00249 std::copy(srcit,srcit+myLen,myview); 00250 srcit += LDA; 00251 myview += myLen; 00252 } 00253 mydata = Teuchos::null; 00254 myview = Teuchos::null; 00255 } 00256 else { 00257 MVT::initializeValues(lclMV_,0,NumVectors,Teuchos::null,0); 00258 } 00259 } 00260 00261 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00262 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00263 MultiVector (const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& map, 00264 const Teuchos::ArrayView<const ArrayView<const Scalar> >& ArrayOfPtrs, 00265 size_t NumVectors) : 00266 DistObject<Scalar,LocalOrdinal,GlobalOrdinal,Node> (map), 00267 lclMV_ (map->getNode ()) 00268 { 00269 using Teuchos::as; 00270 using Teuchos::ArrayRCP; 00271 using Teuchos::ArrayView; 00272 using Teuchos::RCP; 00273 00274 const std::string tfecfFuncName("MultiVector(ArrayOfPtrs)"); 00275 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( 00276 NumVectors < 1 || NumVectors != as<size_t>(ArrayOfPtrs.size()), 00277 std::runtime_error, 00278 ": ArrayOfPtrs.size() must be strictly positive and as large as ArrayOfPtrs."); 00279 const size_t myLen = getLocalLength(); 00280 if (myLen > 0) { 00281 RCP<Node> node = MVT::getNode(lclMV_); 00282 ArrayRCP<Scalar> mydata = node->template allocBuffer<Scalar>(myLen*NumVectors); 00283 MVT::initializeValues(lclMV_,myLen,NumVectors,mydata,myLen); 00284 KOKKOS_NODE_TRACE("MultiVector::MultiVector(2D)") 00285 ArrayRCP<Scalar> myview = node->template viewBufferNonConst<Scalar>(Kokkos::WriteOnly,myLen*NumVectors,mydata); 00286 for (size_t j = 0; j < NumVectors; ++j) { 00287 #ifdef HAVE_TPETRA_DEBUG 00288 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( 00289 as<size_t>(ArrayOfPtrs[j].size()) != getLocalLength(), 00290 std::runtime_error, 00291 ": ArrayOfPtrs[" << j << "].size() (== " << ArrayOfPtrs[j].size() 00292 << ") is not equal to getLocalLength() (== " << getLocalLength()); 00293 #endif 00294 typename ArrayView<const Scalar>::iterator src = ArrayOfPtrs[j].begin(); 00295 std::copy(src,src+myLen,myview); 00296 myview += myLen; 00297 } 00298 myview = Teuchos::null; 00299 mydata = Teuchos::null; 00300 } 00301 else { 00302 MVT::initializeValues(lclMV_,0,NumVectors,Teuchos::null,0); 00303 } 00304 } 00305 00306 00307 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00308 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::~MultiVector() { 00309 } 00310 00311 00312 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00313 bool MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::isConstantStride() const { 00314 return whichVectors_.empty(); 00315 } 00316 00317 00318 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00319 size_t MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalLength() const { 00320 return this->getMap()->getNodeNumElements(); 00321 } 00322 00323 00324 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00325 global_size_t MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getGlobalLength() const { 00326 return this->getMap()->getGlobalNumElements(); 00327 } 00328 00329 00330 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00331 size_t MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getStride() const { 00332 if (isConstantStride()) { 00333 return MVT::getStride(lclMV_); 00334 } 00335 return 0; 00336 } 00337 00338 00339 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00340 bool 00341 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00342 checkSizes(const DistObject<Scalar,LocalOrdinal,GlobalOrdinal,Node> &sourceObj) 00343 { 00344 using Teuchos::RCP; 00345 using Teuchos::rcp_dynamic_cast; 00346 using Teuchos::rcpFromRef; 00347 00348 typedef MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> MV; 00349 // rcp_dynamic_cast gives us superior cast failure output to dynamic_cast. 00350 RCP<const MV> A = rcp_dynamic_cast<const MV> (rcpFromRef (sourceObj), true); 00351 00352 // This method is called in DistObject::doTransfer(). By that 00353 // point, we've already constructed an Import or Export object 00354 // using the two multivectors' Maps, which means that (hopefully) 00355 // we've already checked other attributes of the multivectors. 00356 // Thus, all we need to do here is check the number of columns. 00357 return (A->getNumVectors() == this->getNumVectors()); 00358 } 00359 00360 00361 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00362 void MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::copyAndPermute( 00363 const DistObject<Scalar,LocalOrdinal,GlobalOrdinal,Node> & sourceObj, 00364 size_t numSameIDs, 00365 const Teuchos::ArrayView<const LocalOrdinal> &permuteToLIDs, 00366 const Teuchos::ArrayView<const LocalOrdinal> &permuteFromLIDs) 00367 { 00368 using Teuchos::ArrayRCP; 00369 using Teuchos::ArrayView; 00370 using Teuchos::RCP; 00371 00372 const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &sourceMV = dynamic_cast<const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &>(sourceObj); 00373 typename ArrayView<const LocalOrdinal>::iterator pTo, pFrom; 00374 // any other error will be caught by Teuchos 00375 TEUCHOS_TEST_FOR_EXCEPTION(permuteToLIDs.size() != permuteFromLIDs.size(), std::runtime_error, 00376 "Tpetra::MultiVector::copyAndPermute(): permuteToLIDs and permuteFromLIDs must have the same size."); 00377 RCP<Node> node = MVT::getNode(lclMV_); 00378 const size_t numCols = getNumVectors(); 00379 00380 // Copy the vectors one at a time. This works whether or not the 00381 // multivectors have constant stride. 00382 for (size_t j = 0; j < numCols; ++j) { 00383 // Get host views of the current source and destination column. 00384 // 00385 // FIXME (mfh 10 Mar 2012) Copying should really be done on the 00386 // device. There's no need to bring everything back to the 00387 // host, copy there, and then copy back. 00388 ArrayRCP<const Scalar> srcptr = sourceMV.getSubArrayRCP(sourceMV.cview_,j); 00389 ArrayRCP< Scalar> dstptr = getSubArrayRCP(ncview_,j); 00390 00391 // The first numImportIDs IDs are the same between source and 00392 // target, so we can just copy the data. (This favors faster 00393 // contiguous access whenever we can do it.) 00394 std::copy(srcptr,srcptr+numSameIDs,dstptr); 00395 00396 // For the remaining GIDs, execute the permutations. This may 00397 // involve noncontiguous access of both source and destination 00398 // vectors, depending on the LID lists. 00399 // 00400 // FIXME (mfh 09 Apr 2012) Permutation should be done on the 00401 // device, not on the host. 00402 for (pTo = permuteToLIDs.begin(), pFrom = permuteFromLIDs.begin(); 00403 pTo != permuteToLIDs.end(); ++pTo, ++pFrom) { 00404 dstptr[*pTo] = srcptr[*pFrom]; 00405 } 00406 } 00407 } 00408 00409 00410 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00411 void MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::packAndPrepare( 00412 const DistObject<Scalar,LocalOrdinal,GlobalOrdinal,Node> & sourceObj, 00413 const ArrayView<const LocalOrdinal> &exportLIDs, 00414 Array<Scalar> &exports, 00415 const ArrayView<size_t> &numExportPacketsPerLID, 00416 size_t& constantNumPackets, 00417 Distributor & /* distor */ ) 00418 { 00419 using Teuchos::Array; 00420 using Teuchos::ArrayView; 00421 using Teuchos::as; 00422 typedef MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> MV; 00423 00424 const MV& sourceMV = dynamic_cast<const MV&> (sourceObj); 00425 00426 /* The layout in the export for MultiVectors is as follows: 00427 exports = { all of the data from row exportLIDs.front() ; 00428 .... 00429 all of the data from row exportLIDs.back() } 00430 This doesn't have the best locality, but is necessary because 00431 the data for a Packet (all data associated with an LID) is 00432 required to be contiguous. */ 00433 TEUCHOS_TEST_FOR_EXCEPTION(as<int> (numExportPacketsPerLID.size ()) != exportLIDs.size (), 00434 std::runtime_error, "Tpetra::MultiVector::packAndPrepare(): size of num" 00435 "ExportPacketsPerLID buffer should be the same as exportLIDs. numExport" 00436 "PacketsPerLID.size() = " << numExportPacketsPerLID.size() << ", but " 00437 "exportLIDs.size() = " << exportLIDs.size() << "."); 00438 const KMV& srcData = sourceMV.lclMV_; 00439 const size_t numCols = sourceMV.getNumVectors (); 00440 const size_t stride = MVT::getStride (srcData); 00441 constantNumPackets = numCols; 00442 exports.resize (numCols * exportLIDs.size ()); 00443 typename ArrayView<const LocalOrdinal>::iterator idptr; 00444 typename Array<Scalar>::iterator expptr; 00445 expptr = exports.begin(); 00446 00447 if (sourceMV.isConstantStride ()) { 00448 size_t i = 0; 00449 for (idptr = exportLIDs.begin(); idptr != exportLIDs.end(); ++idptr, ++i) { 00450 for (size_t j = 0; j < numCols; ++j) { 00451 *expptr++ = sourceMV.cview_[j*stride + (*idptr)]; 00452 } 00453 //we shouldn't need to set numExportPacketsPerLID[i] since we have set 00454 //constantNumPackets to a nonzero value. But we'll set it anyway, since 00455 //I'm not sure if the API will remain the way it is. 00456 numExportPacketsPerLID[i] = numCols; 00457 } 00458 } 00459 else { 00460 size_t i = 0; 00461 for (idptr = exportLIDs.begin(); idptr != exportLIDs.end(); ++idptr, ++i) { 00462 for (size_t j = 0; j < numCols; ++j) { 00463 *expptr++ = sourceMV.cview_[sourceMV.whichVectors_[j]*stride + (*idptr)]; 00464 } 00465 numExportPacketsPerLID[i] = numCols; 00466 } 00467 } 00468 } 00469 00470 00471 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00472 void 00473 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00474 unpackAndCombine (const ArrayView<const LocalOrdinal> &importLIDs, 00475 const ArrayView<const Scalar> &imports, 00476 const ArrayView<size_t> &numPacketsPerLID, 00477 size_t constantNumPackets, 00478 Distributor & /* distor */, 00479 CombineMode CM) 00480 { 00481 using Teuchos::ArrayView; 00482 using Teuchos::as; 00483 typedef Teuchos::ScalarTraits<Scalar> SCT; 00484 00485 const std::string tfecfFuncName("unpackAndCombine()"); 00486 /* The layout in the export for MultiVectors is as follows: 00487 imports = { all of the data from row exportLIDs.front() ; 00488 .... 00489 all of the data from row exportLIDs.back() } 00490 This doesn't have the best locality, but is necessary because 00491 the data for a Packet (all data associated with an LID) is 00492 required to be contiguous. */ 00493 #ifdef HAVE_TPETRA_DEBUG 00494 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(as<size_t>(imports.size()) != getNumVectors()*importLIDs.size(), std::runtime_error, 00495 ": Imports buffer size must be consistent with the amount of data to be exported."); 00496 #endif 00497 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(as<size_t>(constantNumPackets) == as<size_t>(0), std::runtime_error, 00498 ": constantNumPackets input argument must be nonzero."); 00499 00500 const size_t myStride = MVT::getStride(lclMV_); 00501 const size_t numVecs = getNumVectors(); 00502 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(as<size_t>(numPacketsPerLID.size()) != as<size_t>(importLIDs.size()), std::runtime_error, 00503 ": numPacketsPerLID must have the same length as importLIDs."); 00504 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(as<size_t>(numVecs) != as<size_t>(constantNumPackets), std::runtime_error, 00505 ": constantNumPackets must equal numVecs."); 00506 00507 if (numVecs > 0 && importLIDs.size()) { 00508 typename ArrayView<const Scalar>::iterator impptr; 00509 typename ArrayView<const LocalOrdinal>::iterator idptr; 00510 impptr = imports.begin(); 00511 // NOTE (mfh 10 Mar 2012) If you want to implement custom 00512 // combine modes, start editing here. Also, if you trust 00513 // inlining, it would be nice to condense this code by using a 00514 // binary function object f: 00515 // 00516 // ncview_[...] = f (ncview_[...], *impptr++); 00517 if (CM == INSERT || CM == REPLACE) { 00518 if (isConstantStride()) { 00519 for (idptr = importLIDs.begin(); idptr != importLIDs.end(); ++idptr) { 00520 for (size_t j = 0; j < numVecs; ++j) { 00521 ncview_[myStride*j + *idptr] = *impptr++; 00522 } 00523 } 00524 } 00525 else { 00526 for (idptr = importLIDs.begin(); idptr != importLIDs.end(); ++idptr) { 00527 for (size_t j = 0; j < numVecs; ++j) { 00528 ncview_[myStride*whichVectors_[j] + *idptr] = *impptr++; 00529 } 00530 } 00531 } 00532 } 00533 else if (CM == ADD) { 00534 if (isConstantStride()) { 00535 for (idptr = importLIDs.begin(); idptr != importLIDs.end(); ++idptr) { 00536 for (size_t j = 0; j < numVecs; ++j) { 00537 ncview_[myStride*j + *idptr] += *impptr++; 00538 } 00539 } 00540 } 00541 else { 00542 for (idptr = importLIDs.begin(); idptr != importLIDs.end(); ++idptr) { 00543 for (size_t j = 0; j < numVecs; ++j) { 00544 ncview_[myStride*whichVectors_[j] + *idptr] += *impptr++; 00545 } 00546 } 00547 } 00548 } 00549 else if (CM == ABSMAX) { 00550 if (isConstantStride()) { 00551 for (idptr = importLIDs.begin(); idptr != importLIDs.end(); ++idptr) { 00552 for (size_t j = 0; j < numVecs; ++j) { 00553 Scalar &curval = ncview_[myStride*j + *idptr]; 00554 const Scalar &newval = *impptr++; 00555 curval = std::max( SCT::magnitude(curval), SCT::magnitude(newval) ); 00556 } 00557 } 00558 } 00559 else { 00560 for (idptr = importLIDs.begin(); idptr != importLIDs.end(); ++idptr) { 00561 for (size_t j = 0; j < numVecs; ++j) { 00562 Scalar &curval = ncview_[myStride*whichVectors_[j] + *idptr]; 00563 const Scalar &newval = *impptr++; 00564 curval = std::max( SCT::magnitude(curval), SCT::magnitude(newval) ); 00565 } 00566 } 00567 } 00568 } 00569 else { 00570 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(CM != ADD && CM != REPLACE && CM != INSERT && CM != ABSMAX, std::invalid_argument, 00571 ": Invalid CombineMode: " << CM); 00572 } 00573 } 00574 } 00575 00576 00577 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00578 inline size_t MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getNumVectors() const { 00579 if (isConstantStride()) { 00580 return MVT::getNumCols(lclMV_); 00581 } 00582 else { 00583 return whichVectors_.size(); 00584 } 00585 } 00586 00587 00588 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00589 void 00590 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00591 dot (const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A, 00592 const Teuchos::ArrayView<Scalar> &dots) const 00593 { 00594 using Teuchos::Array; 00595 using Teuchos::ArrayRCP; 00596 using Teuchos::as; 00597 using Teuchos::arcp_const_cast; 00598 00599 const std::string tfecfFuncName("dot()"); 00600 const size_t myLen = getLocalLength(); 00601 const size_t numVecs = getNumVectors(); 00602 #ifdef HAVE_TPETRA_DEBUG 00603 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( !this->getMap()->isCompatible(*A.getMap()), std::runtime_error, 00604 ": MultiVectors do not have compatible Maps:" << std::endl 00605 << "this->getMap(): " << std::endl << *this->getMap() 00606 << "A.getMap(): " << std::endl << *A.getMap() << std::endl); 00607 #else 00608 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( getLocalLength() != A.getLocalLength(), std::runtime_error, 00609 ": MultiVectors do not have the same local length."); 00610 #endif 00611 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(A.getNumVectors() != numVecs, std::runtime_error, 00612 ": MultiVectors must have the same number of vectors."); 00613 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(Teuchos::as<size_t>(dots.size()) != numVecs, std::runtime_error, 00614 ": dots.size() must be as large as the number of vectors in *this and A."); 00615 if (isConstantStride() && A.isConstantStride()) { 00616 MVT::Dot(lclMV_,A.lclMV_,dots); 00617 } 00618 else { 00619 KMV v(MVT::getNode(lclMV_)), a(MVT::getNode(lclMV_)); 00620 ArrayRCP<Scalar> vptr = arcp_const_cast<Scalar>(MVT::getValues(lclMV_)); 00621 ArrayRCP<Scalar> avptr = arcp_const_cast<Scalar>(MVT::getValues(A.lclMV_)); 00622 for (size_t j=0; j < numVecs; ++j) { 00623 ArrayRCP<Scalar> vj = getSubArrayRCP( vptr,j); 00624 ArrayRCP<Scalar> avj = A.getSubArrayRCP(avptr,j); 00625 MVT::initializeValues(a,myLen, 1, avj, myLen); 00626 MVT::initializeValues(v,myLen, 1, vj, myLen); 00627 dots[j] = MVT::Dot((const KMV&)v,(const KMV &)a); 00628 } 00629 } 00630 if (this->isDistributed()) { 00631 Array<Scalar> ldots(dots); 00632 Teuchos::reduceAll(*this->getMap()->getComm(),Teuchos::REDUCE_SUM,as<int>(numVecs),ldots.getRawPtr(),dots.getRawPtr()); 00633 } 00634 } 00635 00636 00637 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00638 void 00639 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00640 norm2 (const Teuchos::ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>& norms) const 00641 { 00642 using Teuchos::arcp_const_cast; 00643 using Teuchos::Array; 00644 using Teuchos::ArrayRCP; 00645 using Teuchos::ArrayView; 00646 using Teuchos::as; 00647 using Teuchos::reduceAll; 00648 using Teuchos::REDUCE_SUM; 00649 using Teuchos::ScalarTraits; 00650 typedef typename ScalarTraits<Scalar>::magnitudeType MT; 00651 typedef ScalarTraits<MT> STM; 00652 00653 const size_t numVecs = this->getNumVectors(); 00654 TEUCHOS_TEST_FOR_EXCEPTION(as<size_t>(norms.size()) != numVecs, 00655 std::runtime_error, 00656 "Tpetra::MultiVector::norm2(norms): norms.size() must be as large as the number of vectors in *this."); 00657 if (isConstantStride ()) { 00658 MVT::Norm2Squared (lclMV_,norms); 00659 } 00660 else { 00661 KMV v (MVT::getNode (lclMV_)); 00662 ArrayRCP<Scalar> vi; 00663 for (size_t i=0; i < numVecs; ++i) { 00664 vi = arcp_const_cast<Scalar> (MVT::getValues (lclMV_, whichVectors_[i])); 00665 MVT::initializeValues (v, MVT::getNumRows (lclMV_), 1, vi, MVT::getStride (lclMV_)); 00666 norms[i] = MVT::Norm2Squared (v); 00667 } 00668 } 00669 if (this->isDistributed ()) { 00670 Array<MT> lnorms (norms); 00671 // FIXME (mfh 25 Apr 2012) Somebody explain to me why we're 00672 // calling Teuchos::reduceAll when MultiVector has a perfectly 00673 // good reduce() function. 00674 reduceAll (*this->getMap ()->getComm (), REDUCE_SUM, as<int> (numVecs), 00675 lnorms.getRawPtr (), norms.getRawPtr ()); 00676 } 00677 for (typename ArrayView<MT>::iterator n = norms.begin(); n != norms.begin()+numVecs; ++n) { 00678 (*n) = STM::squareroot (*n); 00679 } 00680 } 00681 00682 00683 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00684 void 00685 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00686 normWeighted (const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& weights, 00687 const Teuchos::ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType> &norms) const 00688 { 00689 using Teuchos::arcp_const_cast; 00690 using Teuchos::Array; 00691 using Teuchos::ArrayRCP; 00692 using Teuchos::ArrayView; 00693 using Teuchos::as; 00694 using Teuchos::reduceAll; 00695 using Teuchos::REDUCE_SUM; 00696 using Teuchos::ScalarTraits; 00697 typedef ScalarTraits<Scalar> SCT; 00698 typedef typename SCT::magnitudeType Mag; 00699 00700 const std::string tfecfFuncName("normWeighted()"); 00701 00702 const Mag OneOverN = ScalarTraits<Mag>::one() / getGlobalLength(); 00703 bool OneW = false; 00704 const size_t numVecs = this->getNumVectors(); 00705 if (weights.getNumVectors() == 1) { 00706 OneW = true; 00707 } 00708 else { 00709 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(weights.getNumVectors() != numVecs, std::runtime_error, 00710 ": MultiVector of weights must contain either one vector or the same number of vectors as this."); 00711 } 00712 #ifdef HAVE_TPETRA_DEBUG 00713 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( !this->getMap()->isCompatible(*weights.getMap()), std::runtime_error, 00714 ": MultiVectors do not have compatible Maps:" << std::endl 00715 << "this->getMap(): " << std::endl << *this->getMap() 00716 << "weights.getMap(): " << std::endl << *weights.getMap() << std::endl); 00717 #else 00718 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( getLocalLength() != weights.getLocalLength(), std::runtime_error, 00719 ": MultiVectors do not have the same local length."); 00720 #endif 00721 const size_t myLen = getLocalLength(); 00722 00723 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(as<size_t>(norms.size()) != numVecs, std::runtime_error, 00724 ": norms.size() must be as large as the number of vectors in *this."); 00725 if (isConstantStride() && weights.isConstantStride()) { 00726 MVT::WeightedNorm(lclMV_,weights.lclMV_,norms); 00727 } 00728 else { 00729 KMV v(MVT::getNode(lclMV_)), w(MVT::getNode(lclMV_)); 00730 ArrayRCP<Scalar> vptr = arcp_const_cast<Scalar> (MVT::getValues (lclMV_)); 00731 ArrayRCP<Scalar> wptr = arcp_const_cast<Scalar> (MVT::getValues (weights.lclMV_)); 00732 ArrayRCP<Scalar> wj = wptr.persistingView (0,myLen); 00733 MVT::initializeValues (w,myLen, 1, wj, myLen); 00734 for (size_t j=0; j < numVecs; ++j) { 00735 ArrayRCP<Scalar> vj = getSubArrayRCP (vptr, j); 00736 MVT::initializeValues(v,myLen, 1, vj, myLen); 00737 if (! OneW) { 00738 wj = weights.getSubArrayRCP(wptr,j); 00739 MVT::initializeValues (w, myLen, 1, wj, myLen); 00740 } 00741 norms[j] = MVT::WeightedNorm ((const KMV&)v, (const KMV &)w); 00742 } 00743 } 00744 if (this->isDistributed ()) { 00745 Array<Mag> lnorms(norms); 00746 reduceAll (*this->getMap ()->getComm (), REDUCE_SUM, as<int>(numVecs), 00747 lnorms.getRawPtr (), norms.getRawPtr ()); 00748 } 00749 for (typename ArrayView<Mag>::iterator n = norms.begin(); n != norms.begin()+numVecs; ++n) { 00750 *n = ScalarTraits<Mag>::squareroot(*n * OneOverN); 00751 } 00752 } 00753 00754 00755 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00756 void 00757 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00758 norm1 (const ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType> &norms) const 00759 { 00760 using Teuchos::Array; 00761 using Teuchos::ArrayRCP; 00762 using Teuchos::arcp_const_cast; 00763 using Teuchos::as; 00764 using Teuchos::reduceAll; 00765 using Teuchos::REDUCE_SUM; 00766 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Mag; 00767 00768 const size_t numVecs = this->getNumVectors(); 00769 TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::as<size_t>(norms.size()) != numVecs, std::runtime_error, 00770 "Tpetra::MultiVector::norm1(norms): norms.size() must be as large as the number of vectors in *this."); 00771 if (isConstantStride()) { 00772 MVT::Norm1(lclMV_,norms); 00773 } 00774 else { 00775 KMV v(MVT::getNode(lclMV_)); 00776 ArrayRCP<Scalar> vj; 00777 for (size_t j=0; j < numVecs; ++j) { 00778 vj = arcp_const_cast<Scalar> (MVT::getValues(lclMV_,whichVectors_[j]) ); 00779 MVT::initializeValues (v, MVT::getNumRows(lclMV_), 1, vj, MVT::getStride (lclMV_)); 00780 norms[j] = MVT::Norm1 ((const KMV&)v); 00781 } 00782 } 00783 if (this->isDistributed ()) { 00784 Array<Mag> lnorms (norms); 00785 reduceAll (*this->getMap ()->getComm (), REDUCE_SUM, as<int> (numVecs), lnorms.getRawPtr (), norms.getRawPtr ()); 00786 } 00787 } 00788 00789 00790 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00791 void 00792 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00793 normInf (const Teuchos::ArrayView<typename Teuchos::ScalarTraits<Scalar>::magnitudeType> &norms) const 00794 { 00795 using Teuchos::Array; 00796 using Teuchos::ArrayRCP; 00797 using Teuchos::arcp_const_cast; 00798 using Teuchos::as; 00799 using Teuchos::reduceAll; 00800 using Teuchos::REDUCE_MAX; 00801 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Mag; 00802 00803 const size_t numVecs = this->getNumVectors(); 00804 TEUCHOS_TEST_FOR_EXCEPTION(as<size_t>(norms.size()) != numVecs, std::runtime_error, 00805 "Tpetra::MultiVector::normInf(norms): norms.size() must be as large as the number of vectors in *this."); 00806 if (isConstantStride ()) { 00807 MVT::NormInf(lclMV_,norms); 00808 } 00809 else { 00810 KMV v(MVT::getNode(lclMV_)); 00811 ArrayRCP<Scalar> vj; 00812 for (size_t j=0; j < numVecs; ++j) { 00813 vj = Teuchos::arcp_const_cast<Scalar>( MVT::getValues(lclMV_,whichVectors_[j]) ); 00814 MVT::initializeValues(v,MVT::getNumRows(lclMV_), 1, vj, MVT::getStride(lclMV_)); 00815 norms[j] = MVT::NormInf((const KMV&)v); 00816 } 00817 } 00818 if (this->isDistributed()) { 00819 Array<Mag> lnorms(norms); 00820 reduceAll (*this->getMap ()->getComm (), REDUCE_MAX, as<int> (numVecs), lnorms.getRawPtr (), norms.getRawPtr ()); 00821 } 00822 } 00823 00824 00825 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00826 void 00827 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00828 meanValue (const Teuchos::ArrayView<Scalar> &means) const 00829 { 00830 using Teuchos::Array; 00831 using Teuchos::ArrayRCP; 00832 using Teuchos::arcp_const_cast; 00833 using Teuchos::as; 00834 using Teuchos::arcp_const_cast; 00835 using Teuchos::reduceAll; 00836 using Teuchos::REDUCE_SUM; 00837 typedef Teuchos::ScalarTraits<Scalar> SCT; 00838 00839 const size_t numVecs = getNumVectors(); 00840 const size_t myLen = getLocalLength(); 00841 TEUCHOS_TEST_FOR_EXCEPTION(as<size_t>(means.size()) != numVecs, std::runtime_error, 00842 "Tpetra::MultiVector::meanValue(): means.size() must be as large as the number of vectors in *this."); 00843 // compute local components of the means 00844 // sum these across all nodes 00845 if (isConstantStride()) { 00846 MVT::Sum(lclMV_,means); 00847 } 00848 else { 00849 KMV v(MVT::getNode(lclMV_)); 00850 ArrayRCP<Scalar> vptr = arcp_const_cast<Scalar>(MVT::getValues(lclMV_)); 00851 for (size_t j=0; j < numVecs; ++j) { 00852 ArrayRCP<Scalar> vj = getSubArrayRCP( vptr,j); 00853 MVT::initializeValues(v,myLen, 1, vj, myLen); 00854 means[j] = MVT::Sum((const KMV &)v); 00855 } 00856 } 00857 if (this->isDistributed()) { 00858 Array<Scalar> lmeans(means); 00859 // only combine if we are a distributed MV 00860 reduceAll (*this->getMap ()->getComm (), REDUCE_SUM, as<int> (numVecs), 00861 lmeans.getRawPtr (), means.getRawPtr ()); 00862 } 00863 // mfh 12 Apr 2012: Don't take out the cast from the ordinal type 00864 // to the magnitude type, since operator/ (std::complex<T>, int) 00865 // isn't necessarily defined. 00866 const Scalar OneOverN = 00867 SCT::one() / as<typename SCT::magnitudeType> (getGlobalLength ()); 00868 for (typename ArrayView<Scalar>::iterator i = means.begin(); 00869 i != means.begin() + numVecs; 00870 ++i) { 00871 (*i) = (*i) * OneOverN; 00872 } 00873 } 00874 00875 00876 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00877 void 00878 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00879 randomize() 00880 { 00881 if (isConstantStride ()) { 00882 MVT::Random (lclMV_); 00883 } 00884 else { 00885 const size_t numVecs = this->getNumVectors (); 00886 KMV v (MVT::getNode (lclMV_)); 00887 Teuchos::ArrayRCP<Scalar> vj; 00888 for (size_t j = 0; j < numVecs; ++j) { 00889 vj = MVT::getValuesNonConst (lclMV_, whichVectors_[j]); 00890 MVT::initializeValues (v, MVT::getNumRows (lclMV_), 1, vj, MVT::getStride (lclMV_)); 00891 MVT::Random (v); 00892 } 00893 } 00894 } 00895 00896 00897 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00898 void 00899 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00900 putScalar (const Scalar &alpha) 00901 { 00902 const size_t numVecs = getNumVectors(); 00903 if (isConstantStride()) { 00904 MVT::Init(lclMV_,alpha); 00905 } 00906 else { 00907 KMV v(MVT::getNode(lclMV_)); 00908 Teuchos::ArrayRCP<Scalar> vj; 00909 for (size_t j=0; j < numVecs; ++j) { 00910 vj = MVT::getValuesNonConst(lclMV_,whichVectors_[j]); 00911 MVT::initializeValues(v,MVT::getNumRows(lclMV_), 1, vj, MVT::getStride(lclMV_)); 00912 MVT::Init(v,alpha); 00913 } 00914 } 00915 } 00916 00917 00918 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00919 void 00920 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00921 replaceMap (const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > &map) 00922 { 00923 #ifdef HAVE_TEUCHOS_DEBUG 00924 TEUCHOS_TEST_FOR_EXCEPTION(! this->getMap ()->isCompatible (*map), 00925 std::invalid_argument, "Tpetra::MultiVector::replaceMap(): The input map " 00926 "is not compatible with this multivector's current map. The replaceMap() " 00927 "method is not for data redistribution; use Import or Export for that."); 00928 #endif // HAVE_TEUCHOS_DEBUG 00929 this->map_ = map; 00930 } 00931 00932 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00933 void 00934 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00935 scale (const Scalar &alpha) 00936 { 00937 using Teuchos::arcp_const_cast; 00938 using Teuchos::ArrayRCP; 00939 00940 // NOTE: can't substitute putScalar(0.0) for scale(0.0), because 00941 // the former will overwrite NaNs present in the MultiVector, while the 00942 // semantics of this call require multiplying them by 0, which IEEE requires to be NaN 00943 const size_t numVecs = getNumVectors(); 00944 if (alpha == Teuchos::ScalarTraits<Scalar>::one()) { 00945 // do nothing 00946 } 00947 else if (isConstantStride()) { 00948 MVT::Scale(lclMV_,alpha); 00949 } 00950 else { 00951 KMV v(MVT::getNode(lclMV_)); 00952 ArrayRCP<Scalar> vj; 00953 for (size_t j=0; j < numVecs; ++j) { 00954 vj = arcp_const_cast<Scalar> (MVT::getValues(lclMV_, whichVectors_[j])); 00955 MVT::initializeValues (v, MVT::getNumRows (lclMV_), 1, vj, MVT::getStride (lclMV_)); 00956 MVT::Scale (v, alpha); 00957 } 00958 } 00959 } 00960 00961 00962 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00963 void 00964 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00965 scale (Teuchos::ArrayView<const Scalar> alphas) 00966 { 00967 using Teuchos::arcp_const_cast; 00968 using Teuchos::ArrayRCP; 00969 using Teuchos::as; 00970 00971 const size_t numVecs = this->getNumVectors(); 00972 TEUCHOS_TEST_FOR_EXCEPTION(as<size_t>(alphas.size()) != numVecs, std::runtime_error, 00973 "Tpetra::MultiVector::scale(alphas): alphas.size() must be as large as the number of vectors in *this."); 00974 KMV vec(MVT::getNode(lclMV_)); 00975 const size_t myLen = MVT::getNumRows(lclMV_); 00976 if (myLen == 0) return; 00977 ArrayRCP<Scalar> mybuf = MVT::getValuesNonConst(lclMV_); 00978 for (size_t j = 0; j < numVecs; ++j) { 00979 if (alphas[j] == Teuchos::ScalarTraits<Scalar>::one()) { 00980 // do nothing: NaN * 1.0 == NaN, Number*1.0 == Number 00981 } 00982 else { 00983 ArrayRCP<Scalar> mybufj = getSubArrayRCP(mybuf,j); 00984 MVT::initializeValues(vec,myLen,1,mybufj,myLen); 00985 MVT::Scale(vec,alphas[j]); 00986 } 00987 } 00988 } 00989 00990 00991 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 00992 void 00993 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 00994 scale (const Scalar &alpha, const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A) 00995 { 00996 using Teuchos::arcp_const_cast; 00997 using Teuchos::ArrayRCP; 00998 using Teuchos::as; 00999 01000 const std::string tfecfFuncName("scale(alpha,A)"); 01001 01002 const size_t numVecs = getNumVectors(), 01003 myLen = getLocalLength(); 01004 #ifdef HAVE_TPETRA_DEBUG 01005 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(! this->getMap ()->isCompatible (*A.getMap ()), 01006 std::runtime_error, ": MultiVectors do not have compatible Maps:" << std::endl 01007 << "this->getMap(): " << std::endl << *(this->getMap ()) 01008 << "A.getMap(): " << std::endl << *(A.getMap ()) << std::endl); 01009 #else 01010 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( getLocalLength() != A.getLocalLength(), std::runtime_error, 01011 ": MultiVectors do not have the same local length."); 01012 #endif 01013 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(A.getNumVectors() != numVecs, std::runtime_error, 01014 ": MultiVectors must have the same number of vectors."); 01015 if (isConstantStride() && A.isConstantStride()) { 01016 // set me == alpha*A 01017 MVT::Scale(lclMV_,alpha,(const KMV&)A.lclMV_); 01018 } 01019 else { 01020 KMV v(MVT::getNode(lclMV_)), a(MVT::getNode(lclMV_)); 01021 ArrayRCP<Scalar> vptr = MVT::getValuesNonConst (lclMV_); 01022 ArrayRCP<Scalar> avptr = arcp_const_cast<Scalar> (MVT::getValues (A.lclMV_)); 01023 for (size_t j=0; j < numVecs; ++j) { 01024 ArrayRCP<Scalar> vj = getSubArrayRCP (vptr, j); 01025 ArrayRCP<Scalar> avj = A.getSubArrayRCP (avptr, j); 01026 MVT::initializeValues (a,myLen, 1, avj, myLen); 01027 MVT::initializeValues (v,myLen, 1, vj, myLen); 01028 MVT::Scale (v, alpha, (const KMV &)a); 01029 } 01030 } 01031 } 01032 01033 01034 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01035 void 01036 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01037 reciprocal (const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A) 01038 { 01039 using Teuchos::arcp_const_cast; 01040 using Teuchos::ArrayRCP; 01041 using Teuchos::as; 01042 01043 const std::string tfecfFuncName("reciprocal()"); 01044 01045 #ifdef HAVE_TPETRA_DEBUG 01046 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( !this->getMap()->isCompatible(*A.getMap()), std::runtime_error, 01047 ": MultiVectors do not have compatible Maps:" << std::endl 01048 << "this->getMap(): " << std::endl << *this->getMap() 01049 << "A.getMap(): " << std::endl << *A.getMap() << std::endl); 01050 #else 01051 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( getLocalLength() != A.getLocalLength(), std::runtime_error, 01052 ": MultiVectors do not have the same local length."); 01053 #endif 01054 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(A.getNumVectors() != this->getNumVectors(), std::runtime_error, 01055 ": MultiVectors must have the same number of vectors."); 01056 const size_t numVecs = getNumVectors(); 01057 const size_t myLen = getLocalLength(); 01058 try { 01059 if (isConstantStride() && A.isConstantStride()) { 01060 MVT::Recip(lclMV_,(const KMV&)A.lclMV_); 01061 } 01062 else { 01063 KMV v(MVT::getNode(lclMV_)), a(MVT::getNode(lclMV_)); 01064 ArrayRCP<Scalar> vptr = MVT::getValuesNonConst(lclMV_), 01065 avptr = arcp_const_cast<Scalar>(MVT::getValues(A.lclMV_)); 01066 for (size_t j=0; j < numVecs; ++j) { 01067 ArrayRCP<Scalar> vj = getSubArrayRCP( vptr,j), 01068 avj = A.getSubArrayRCP(avptr,j); 01069 MVT::initializeValues(a,myLen, 1, avj, myLen); 01070 MVT::initializeValues(v,myLen, 1, vj, myLen); 01071 MVT::Recip(v,(const KMV &)a); 01072 } 01073 } 01074 } 01075 catch (std::runtime_error &e) { 01076 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(true,std::runtime_error, 01077 ": caught exception from Kokkos:" << std::endl 01078 << e.what() << std::endl); 01079 } 01080 } 01081 01082 01083 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01084 void MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::abs(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A) { 01085 using Teuchos::arcp_const_cast; 01086 using Teuchos::ArrayRCP; 01087 using Teuchos::as; 01088 01089 const std::string tfecfFuncName("abs()"); 01090 #ifdef HAVE_TPETRA_DEBUG 01091 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( !this->getMap()->isCompatible(*A.getMap()), std::runtime_error, 01092 ": MultiVectors do not have compatible Maps:" << std::endl 01093 << "this->getMap(): " << std::endl << *this->getMap() 01094 << "A.getMap(): " << std::endl << *A.getMap() << std::endl); 01095 #else 01096 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( getLocalLength() != A.getLocalLength(), std::runtime_error, 01097 ": MultiVectors do not have the same local length."); 01098 #endif 01099 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(A.getNumVectors() != this->getNumVectors(), std::runtime_error, 01100 ": MultiVectors must have the same number of vectors."); 01101 const size_t myLen = getLocalLength(); 01102 const size_t numVecs = getNumVectors(); 01103 if (isConstantStride() && A.isConstantStride()) { 01104 MVT::Abs(lclMV_,(const KMV&)A.lclMV_); 01105 } 01106 else { 01107 KMV v(MVT::getNode(lclMV_)), a(MVT::getNode(lclMV_)); 01108 ArrayRCP<Scalar> vptr = MVT::getValuesNonConst(lclMV_), 01109 avptr = arcp_const_cast<Scalar>(MVT::getValues(A.lclMV_)); 01110 for (size_t j=0; j < numVecs; ++j) { 01111 ArrayRCP<Scalar> vj = getSubArrayRCP( vptr,j), 01112 avj = A.getSubArrayRCP(avptr,j); 01113 MVT::initializeValues(a,myLen, 1, avj, myLen); 01114 MVT::initializeValues(v,myLen, 1, vj, myLen); 01115 MVT::Abs(v,(const KMV &)a); 01116 } 01117 } 01118 } 01119 01120 01121 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01122 void MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::update( 01123 const Scalar &alpha, const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A, 01124 const Scalar &beta) 01125 { 01126 using Teuchos::arcp_const_cast; 01127 using Teuchos::ArrayRCP; 01128 using Teuchos::as; 01129 01130 const std::string tfecfFuncName("update()"); 01131 // this = beta*this + alpha*A 01132 // must support case where &this == &A 01133 // can't short circuit on alpha==0.0 or beta==0.0, because 0.0*NaN != 0.0 01134 #ifdef HAVE_TPETRA_DEBUG 01135 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( !this->getMap()->isCompatible(*A.getMap()), std::runtime_error, 01136 ": MultiVectors do not have compatible Maps:" << std::endl 01137 << "this->getMap(): " << std::endl << this->getMap() 01138 << "A.getMap(): " << std::endl << *A.getMap() << std::endl); 01139 #else 01140 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( getLocalLength() != A.getLocalLength(), std::runtime_error, 01141 ": MultiVectors do not have the same local length."); 01142 #endif 01143 const size_t myLen = getLocalLength(); 01144 const size_t numVecs = getNumVectors(); 01145 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(A.getNumVectors() != this->getNumVectors(), std::runtime_error, 01146 ": MultiVectors must have the same number of vectors."); 01147 if (isConstantStride() && A.isConstantStride()) { 01148 MVT::GESUM(lclMV_,alpha,(const KMV&)A.lclMV_,beta); 01149 } 01150 else { 01151 KMV v(MVT::getNode(lclMV_)), a(MVT::getNode(lclMV_)); 01152 ArrayRCP<Scalar> vptr = MVT::getValuesNonConst(lclMV_), 01153 avptr = arcp_const_cast<Scalar>(MVT::getValues(A.lclMV_)); 01154 for (size_t j=0; j < numVecs; ++j) { 01155 ArrayRCP<Scalar> vj = getSubArrayRCP( vptr,j), 01156 avj = A.getSubArrayRCP(avptr,j); 01157 MVT::initializeValues(a,myLen, 1, avj, myLen); 01158 MVT::initializeValues(v,myLen, 1, vj, myLen); 01159 MVT::GESUM(v,alpha,(const KMV &)a,beta); 01160 } 01161 } 01162 } 01163 01164 01165 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01166 void MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::update( 01167 const Scalar &alpha, const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &A, 01168 const Scalar &beta, const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &B, 01169 const Scalar &gamma) 01170 { 01171 using Teuchos::arcp_const_cast; 01172 using Teuchos::ArrayRCP; 01173 using Teuchos::as; 01174 01175 const std::string tfecfFuncName("update()"); 01176 // this = alpha*A + beta*B + gamma*this 01177 // must support case where &this == &A or &this == &B 01178 // can't short circuit on alpha==0.0 or beta==0.0 or gamma==0.0, because 0.0*NaN != 0.0 01179 #ifdef HAVE_TPETRA_DEBUG 01180 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( !this->getMap()->isCompatible(*A.getMap()) || !this->getMap()->isCompatible(*B.getMap()), 01181 std::runtime_error, 01182 ": MultiVectors do not have compatible Maps:" << std::endl 01183 << "this->getMap(): " << std::endl << *this->getMap() 01184 << "A.getMap(): " << std::endl << *A.getMap() << std::endl 01185 << "B.getMap(): " << std::endl << *B.getMap() << std::endl); 01186 #else 01187 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( getLocalLength() != A.getLocalLength() || getLocalLength() != B.getLocalLength(), std::runtime_error, 01188 ": MultiVectors do not have the same local length."); 01189 #endif 01190 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(A.getNumVectors() != this->getNumVectors() || B.getNumVectors() != this->getNumVectors(), std::runtime_error, 01191 ": MultiVectors must have the same number of vectors."); 01192 const size_t myLen = getLocalLength(); 01193 const size_t numVecs = getNumVectors(); 01194 if (isConstantStride() && A.isConstantStride() && B.isConstantStride()) { 01195 MVT::GESUM(lclMV_,alpha,(const KMV&)A.lclMV_,beta,(const KMV&)B.lclMV_,gamma); 01196 } 01197 else { 01198 KMV v(MVT::getNode(lclMV_)), a(MVT::getNode(lclMV_)), b(MVT::getNode(lclMV_)); 01199 ArrayRCP<Scalar> vptr = MVT::getValuesNonConst(lclMV_), 01200 avptr = arcp_const_cast<Scalar>(MVT::getValues(A.lclMV_)), 01201 bvptr = arcp_const_cast<Scalar>(MVT::getValues(B.lclMV_)); 01202 for (size_t j=0; j < numVecs; ++j) { 01203 ArrayRCP<Scalar> vj = getSubArrayRCP( vptr,j), 01204 avj = A.getSubArrayRCP(avptr,j), 01205 bvj = B.getSubArrayRCP(bvptr,j); 01206 MVT::initializeValues(b,myLen, 1, bvj, myLen); 01207 MVT::initializeValues(a,myLen, 1, avj, myLen); 01208 MVT::initializeValues(v,myLen, 1, vj, myLen); 01209 MVT::GESUM(v,alpha,(const KMV&)a,beta,(const KMV&)b,gamma); 01210 } 01211 } 01212 } 01213 01214 01215 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01216 Teuchos::ArrayRCP<const Scalar> 01217 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01218 getData (size_t j) const 01219 { 01220 Teuchos::RCP<Node> node = MVT::getNode(lclMV_); 01221 KOKKOS_NODE_TRACE("MultiVector::getData()") 01222 return node->template viewBuffer<Scalar> (getLocalLength (), getSubArrayRCP (MVT::getValues(lclMV_), j)); 01223 } 01224 01225 01226 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01227 Teuchos::ArrayRCP<Scalar> 01228 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01229 getDataNonConst(size_t j) 01230 { 01231 Teuchos::RCP<Node> node = MVT::getNode(lclMV_); 01232 KOKKOS_NODE_TRACE("MultiVector::getDataNonConst()") 01233 return node->template viewBufferNonConst<Scalar> (Kokkos::ReadWrite, getLocalLength (), getSubArrayRCP (MVT::getValuesNonConst (lclMV_), j)); 01234 } 01235 01236 01237 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01238 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& 01239 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01240 operator= (const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &source) 01241 { 01242 const std::string tfecfFuncName("operator=()"); 01243 // Check for special case of this=Source, in which case we do nothing. 01244 if (this != &source) { 01245 #ifdef HAVE_TPETRA_DEBUG 01246 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( !this->getMap()->isCompatible(*source.getMap()), std::invalid_argument, 01247 ": MultiVectors do not have compatible Maps:" << std::endl 01248 << "this->getMap(): " << std::endl << *this->getMap() 01249 << "source.getMap(): " << std::endl << *source.getMap() << std::endl); 01250 #else 01251 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( getLocalLength() != source.getLocalLength(), std::invalid_argument, 01252 ": MultiVectors do not have the same local length."); 01253 #endif 01254 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(source.getNumVectors() != getNumVectors(), std::invalid_argument, 01255 ": MultiVectors must have the same number of vectors."); 01256 Teuchos::RCP<Node> node = MVT::getNode (lclMV_); 01257 const size_t numVecs = getNumVectors(); 01258 if (isConstantStride() && source.isConstantStride() && getLocalLength()==getStride() && source.getLocalLength()==source.getStride()) { 01259 // Both multivectors' data are stored contiguously, so we can 01260 // copy in one call. 01261 KOKKOS_NODE_TRACE("MultiVector::operator=()") 01262 node->template copyBuffers<Scalar>(getLocalLength()*numVecs, MVT::getValues(source.lclMV_), MVT::getValuesNonConst(lclMV_) ); 01263 } 01264 else { 01265 // We have to copy the columns one at a time. 01266 for (size_t j=0; j < numVecs; ++j) { 01267 KOKKOS_NODE_TRACE("MultiVector::operator=()") 01268 node->template copyBuffers<Scalar>(getLocalLength(), source.getSubArrayRCP(MVT::getValues(source.lclMV_),j), 01269 getSubArrayRCP(MVT::getValuesNonConst(lclMV_),j) ); 01270 } 01271 } 01272 } 01273 return(*this); 01274 } 01275 01276 01277 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01278 Teuchos::RCP<MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > 01279 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01280 subCopy (const Teuchos::ArrayView<const size_t> &cols) const 01281 { 01282 using Teuchos::RCP; 01283 using Teuchos::rcp; 01284 01285 TEUCHOS_TEST_FOR_EXCEPTION(cols.size() < 1, std::runtime_error, 01286 "Tpetra::MultiVector::subCopy(cols): cols must contain at least one column."); 01287 size_t numCopyVecs = cols.size(); 01288 const bool zeroData = false; 01289 RCP<Node> node = MVT::getNode(lclMV_); 01290 RCP<MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > mv; 01291 // mv is allocated with constant stride 01292 mv = rcp (new MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> (this->getMap (), numCopyVecs, zeroData)); 01293 // copy data from *this into mv 01294 for (size_t j=0; j<numCopyVecs; ++j) { 01295 KOKKOS_NODE_TRACE("MultiVector::subCopy()") 01296 node->template copyBuffers<Scalar> (getLocalLength (), 01297 getSubArrayRCP (MVT::getValues (lclMV_), cols[j]), 01298 MVT::getValuesNonConst (mv->lclMV_, j)); 01299 } 01300 return mv; 01301 } 01302 01303 01304 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01305 Teuchos::RCP<MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > 01306 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01307 subCopy (const Teuchos::Range1D &colRng) const 01308 { 01309 using Teuchos::RCP; 01310 using Teuchos::rcp; 01311 01312 TEUCHOS_TEST_FOR_EXCEPTION(colRng.size() == 0, std::runtime_error, 01313 "Tpetra::MultiVector::subCopy(Range1D): range must include at least one vector."); 01314 size_t numCopyVecs = colRng.size(); 01315 const bool zeroData = false; 01316 RCP<Node> node = MVT::getNode(lclMV_); 01317 RCP<MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > mv; 01318 // mv is allocated with constant stride 01319 mv = rcp (new MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> (this->getMap (), numCopyVecs, zeroData)); 01320 // copy data from *this into mv 01321 for (size_t js=colRng.lbound(), jd=0; jd<numCopyVecs; ++jd, ++js) { 01322 KOKKOS_NODE_TRACE("MultiVector::subCopy()") 01323 node->template copyBuffers<Scalar> (getLocalLength (), 01324 getSubArrayRCP (MVT::getValues (lclMV_), js), 01325 MVT::getValuesNonConst (mv->lclMV_, jd)); 01326 } 01327 return mv; 01328 } 01329 01330 01331 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01332 Teuchos::RCP<const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > 01333 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01334 offsetView (const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& subMap, 01335 size_t offset) const 01336 { 01337 using Teuchos::arcp_const_cast; 01338 using Teuchos::ArrayRCP; 01339 using Teuchos::RCP; 01340 using Teuchos::rcp; 01341 01342 typedef const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> CMV; 01343 TEUCHOS_TEST_FOR_EXCEPTION( subMap->getNodeNumElements() + offset > this->getLocalLength(), std::runtime_error, 01344 "Tpetra::MultiVector::offsetView(subMap,offset): sizes are not sane.\noffset == " << offset << "\nsubMap: " << subMap->description() << "\nthis->rowMap: " << this->getMap()->description()); 01345 const size_t numVecs = this->getNumVectors(), 01346 myStride = MVT::getStride(lclMV_), 01347 newLen = subMap->getNodeNumElements(); 01348 // this is const, so the lclMV_ is const, so that we can only get const buffers 01349 // we will cast away the const; this is okay, because 01350 // a) the constructor doesn't modify the data, and 01351 // b) we are encapsulating in a const MV before returning 01352 ArrayRCP<const Scalar> cbuf = MVT::getValues(lclMV_); 01353 ArrayRCP<Scalar> ncbuf = arcp_const_cast<Scalar> (cbuf); 01354 RCP<CMV> constViewMV; 01355 if (isConstantStride()) { 01356 // view goes from first entry of first vector to last entry of last vector 01357 ArrayRCP<Scalar> subdata = ncbuf.persistingView( offset, myStride * (numVecs-1) + newLen ); 01358 constViewMV = rcp (new MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> (subMap, subdata, myStride, numVecs, COMPUTE_VIEW_CONSTRUCTOR)); 01359 } 01360 else { 01361 // use same which index, but with an offset start pointer 01362 size_t maxSubVecIndex = *std::max_element(whichVectors_.begin(), whichVectors_.end()); 01363 ArrayRCP<Scalar> subdata = ncbuf.persistingView( offset, myStride * maxSubVecIndex + newLen ); 01364 constViewMV = rcp (new MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> (subMap, subdata, myStride, whichVectors_, COMPUTE_VIEW_CONSTRUCTOR)); 01365 } 01366 return constViewMV; 01367 } 01368 01369 01370 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01371 Teuchos::RCP<MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > 01372 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01373 offsetViewNonConst (const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& subMap, 01374 size_t offset) 01375 { 01376 using Teuchos::arcp_const_cast; 01377 using Teuchos::ArrayRCP; 01378 using Teuchos::RCP; 01379 using Teuchos::rcp; 01380 01381 typedef MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> MV; 01382 TEUCHOS_TEST_FOR_EXCEPTION( subMap->getNodeNumElements() + offset > this->getLocalLength(), std::runtime_error, 01383 "Tpetra::MultiVector::offsetView(subMap,offset): sizes are not sane.\noffset == " 01384 << offset << "\nsubMap: " << subMap->description() << "\nthis->rowMap: " 01385 << this->getMap()->description()); 01386 const size_t numVecs = this->getNumVectors(), 01387 myStride = MVT::getStride(lclMV_), 01388 newLen = subMap->getNodeNumElements(); 01389 // this is const, so the lclMV_ is const, so that we can only get const buffers 01390 // we will cast away the const; this is okay, because 01391 // a) the constructor doesn't modify the data, and 01392 // b) we are encapsulating in a const MV before returning 01393 ArrayRCP<Scalar> buf = MVT::getValuesNonConst(lclMV_); 01394 RCP<MV> subViewMV; 01395 if (isConstantStride()) { 01396 // view goes from first entry of first vector to last entry of last vector 01397 ArrayRCP<Scalar> subdata = buf.persistingView( offset, myStride * (numVecs-1) + newLen ); 01398 subViewMV = rcp (new MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> (subMap, subdata, myStride, numVecs, COMPUTE_VIEW_CONSTRUCTOR)); 01399 } 01400 else { 01401 // use same which index, but with an offset start pointer 01402 size_t maxSubVecIndex = *std::max_element(whichVectors_.begin(), whichVectors_.end()); 01403 ArrayRCP<Scalar> subdata = buf.persistingView( offset, myStride * maxSubVecIndex + newLen ); 01404 subViewMV = rcp (new MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> (subMap, subdata, myStride, whichVectors_, COMPUTE_VIEW_CONSTRUCTOR)); 01405 } 01406 return subViewMV; 01407 } 01408 01409 01410 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01411 Teuchos::RCP<const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > 01412 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01413 subView (const ArrayView<const size_t> &cols) const 01414 { 01415 using Teuchos::arcp_const_cast; 01416 using Teuchos::Array; 01417 using Teuchos::ArrayRCP; 01418 using Teuchos::RCP; 01419 using Teuchos::rcp; 01420 using Teuchos::rcp_const_cast; 01421 typedef MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> MV; 01422 01423 TEUCHOS_TEST_FOR_EXCEPTION(cols.size() == 0, std::runtime_error, 01424 "Tpetra::MultiVector::subView(ArrayView): range must include at least one vector."); 01425 // this is const, so the lclMV_ is const, so that we can only get const buffers 01426 // we will cast away the const; this is okay, because 01427 // a) the constructor doesn't modify the data, and 01428 // b) we are encapsulating in a const MV before returning 01429 const size_t myStride = MVT::getStride(lclMV_), 01430 myLen = MVT::getNumRows(lclMV_), 01431 numViewCols = cols.size(); 01432 // use the smallest view possible of the buffer: from the first element of the minInd vector to the last element of the maxInd vector 01433 // this minimizes overlap between views, and keeps view of the minimum amount necessary, in order to allow node to achieve maximum efficiency. 01434 // adjust the indices appropriately; shift so that smallest index is 0 01435 ArrayRCP<const Scalar> cbuf = MVT::getValues(lclMV_); 01436 ArrayRCP<Scalar> ncbuf = arcp_const_cast<Scalar> (cbuf); 01437 Array<size_t> newCols (numViewCols); 01438 size_t minInd = Teuchos::OrdinalTraits<size_t>::max(); 01439 size_t maxInd = Teuchos::OrdinalTraits<size_t>::zero(); 01440 if (isConstantStride()) { 01441 for (size_t j=0; j < numViewCols; ++j) { 01442 newCols[j] = cols[j]; 01443 if (newCols[j] < minInd) minInd = newCols[j]; 01444 if (maxInd < newCols[j]) maxInd = newCols[j]; 01445 } 01446 } 01447 else { 01448 for (size_t j=0; j < numViewCols; ++j) { 01449 newCols[j] = whichVectors_[cols[j]]; 01450 if (newCols[j] < minInd) minInd = newCols[j]; 01451 if (maxInd < newCols[j]) maxInd = newCols[j]; 01452 } 01453 } 01454 ArrayRCP<Scalar> minbuf = ncbuf.persistingView(minInd * myStride, myStride * (maxInd - minInd) + myLen); 01455 for (size_t j=0; j < numViewCols; ++j) { 01456 newCols[j] -= minInd; 01457 } 01458 RCP<const MV> constViewMV; 01459 return rcp_const_cast<const MV> (rcp (new MV (this->getMap (), minbuf, myStride, newCols (), COMPUTE_VIEW_CONSTRUCTOR))); 01460 } 01461 01462 01463 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01464 Teuchos::RCP<const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > 01465 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01466 subView (const Teuchos::Range1D &colRng) const 01467 { 01468 using Teuchos::arcp_const_cast; 01469 using Teuchos::Array; 01470 using Teuchos::ArrayRCP; 01471 using Teuchos::RCP; 01472 using Teuchos::rcp; 01473 typedef MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> MV; 01474 01475 TEUCHOS_TEST_FOR_EXCEPTION(colRng.size() == 0, std::runtime_error, 01476 "Tpetra::MultiVector::subView(Range1D): range must include at least one vector."); 01477 size_t numViewVecs = colRng.size(); 01478 // this is const, so the lclMV_ is const, so that we can only get const buffers 01479 // we will cast away the const; this is okay, because 01480 // a) the constructor doesn't modify the data, and 01481 // b) we are encapsulating in a const MV before returning 01482 ArrayRCP<const Scalar> cbuf = MVT::getValues(lclMV_); 01483 ArrayRCP<Scalar> ncbuf = arcp_const_cast<Scalar>(cbuf); 01484 // resulting MultiVector is constant stride only if *this is 01485 if (isConstantStride()) { 01486 // view goes from first entry of first vector to last entry of last vector 01487 ArrayRCP<Scalar> subdata = ncbuf.persistingView( MVT::getStride(lclMV_) * colRng.lbound(), 01488 MVT::getStride(lclMV_) * (numViewVecs-1) + getLocalLength() ); 01489 return rcp (new MV (this->getMap (), subdata, MVT::getStride (lclMV_), numViewVecs, COMPUTE_VIEW_CONSTRUCTOR)); 01490 } 01491 // otherwise, use a subset of this whichVectors_ to construct new multivector 01492 Array<size_t> whchvecs( whichVectors_.begin()+colRng.lbound(), whichVectors_.begin()+colRng.ubound()+1 ); 01493 return rcp (new MV (this->getMap (), ncbuf, MVT::getStride (lclMV_), whchvecs, COMPUTE_VIEW_CONSTRUCTOR)); 01494 } 01495 01496 01497 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01498 Teuchos::RCP<MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > 01499 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01500 subViewNonConst (const ArrayView<const size_t> &cols) 01501 { 01502 using Teuchos::as; 01503 using Teuchos::arcp_const_cast; 01504 using Teuchos::Array; 01505 using Teuchos::ArrayRCP; 01506 using Teuchos::RCP; 01507 using Teuchos::rcp; 01508 typedef MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> MV; 01509 01510 TEUCHOS_TEST_FOR_EXCEPTION(cols.size() == 0, std::runtime_error, 01511 "Tpetra::MultiVector::subViewNonConst(ArrayView): range must include at least one vector."); 01512 if (isConstantStride()) { 01513 return rcp (new MV (this->getMap (), MVT::getValuesNonConst (lclMV_), MVT::getStride (lclMV_), cols, COMPUTE_VIEW_CONSTRUCTOR)); 01514 } 01515 // else, lookup current whichVectors_ using cols 01516 Array<size_t> newcols(cols.size()); 01517 for (size_t j = 0; j < as<size_t> (cols.size ()); ++j) { 01518 newcols[j] = whichVectors_[cols[j]]; 01519 } 01520 return rcp (new MV (this->getMap (), MVT::getValuesNonConst (lclMV_), MVT::getStride (lclMV_), newcols (), COMPUTE_VIEW_CONSTRUCTOR)); 01521 } 01522 01523 01524 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01525 Teuchos::RCP<MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > 01526 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01527 subViewNonConst (const Teuchos::Range1D &colRng) 01528 { 01529 using Teuchos::as; 01530 using Teuchos::arcp_const_cast; 01531 using Teuchos::Array; 01532 using Teuchos::ArrayRCP; 01533 using Teuchos::RCP; 01534 using Teuchos::rcp; 01535 typedef MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> MV; 01536 01537 TEUCHOS_TEST_FOR_EXCEPTION(colRng.size() == 0, std::runtime_error, 01538 "Tpetra::MultiVector::subViewNonConst(Range1D): range must include at least one vector."); 01539 size_t numViewVecs = colRng.size(); 01540 // resulting MultiVector is constant stride only if *this is 01541 if (isConstantStride ()) { 01542 // view goes from first entry of first vector to last entry of last vector 01543 const size_t stride = MVT::getStride(lclMV_); 01544 ArrayRCP<Scalar> data = MVT::getValuesNonConst(lclMV_); 01545 ArrayRCP<Scalar> subdata = data.persistingView( stride * colRng.lbound(), 01546 stride * (numViewVecs-1) + getLocalLength() ); 01547 return rcp (new MV (this->getMap (), subdata, stride, numViewVecs, COMPUTE_VIEW_CONSTRUCTOR)); 01548 } 01549 // otherwise, use a subset of this whichVectors_ to construct new multivector 01550 Array<size_t> whchvecs( whichVectors_.begin()+colRng.lbound(), whichVectors_.begin()+colRng.ubound()+1 ); 01551 const size_t stride = MVT::getStride(lclMV_); 01552 ArrayRCP<Scalar> data = MVT::getValuesNonConst(lclMV_); 01553 return rcp (new MV (this->getMap (), data, stride, whchvecs, COMPUTE_VIEW_CONSTRUCTOR)); 01554 } 01555 01556 01557 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01558 Teuchos::RCP<const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > 01559 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01560 getVector (size_t j) const 01561 { 01562 using Teuchos::arcp_const_cast; 01563 using Teuchos::ArrayRCP; 01564 using Teuchos::rcp; 01565 using Teuchos::rcp_const_cast; 01566 typedef Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> V; 01567 01568 #ifdef HAVE_TPETRA_DEBUG 01569 TEUCHOS_TEST_FOR_EXCEPTION( vectorIndexOutOfRange(j), std::runtime_error, 01570 "Tpetra::MultiVector::getVector(j): index j (== " << j << ") exceeds valid column range for this multivector."); 01571 #endif 01572 // this is const, so lclMV_ is const, so we get const buff 01573 // it is safe to cast away the const because we will wrap it in a const Vector below 01574 ArrayRCP<Scalar> ncbuff; 01575 if (getLocalLength() > 0) { 01576 ArrayRCP<const Scalar> cbuff = getSubArrayRCP(MVT::getValues(lclMV_),j); 01577 ncbuff = arcp_const_cast<Scalar>(cbuff); 01578 } 01579 return rcp_const_cast<const V> (rcp (new V (this->getMap (), ncbuff, COMPUTE_VIEW_CONSTRUCTOR))); 01580 } 01581 01582 01583 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01584 Teuchos::RCP<Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > 01585 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getVectorNonConst(size_t j) 01586 { 01587 using Teuchos::ArrayRCP; 01588 using Teuchos::rcp; 01589 typedef Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> V; 01590 01591 #ifdef HAVE_TPETRA_DEBUG 01592 TEUCHOS_TEST_FOR_EXCEPTION( vectorIndexOutOfRange(j), std::runtime_error, 01593 "Tpetra::MultiVector::getVectorNonConst(j): index j (== " << j << ") exceeds valid column range for this multivector."); 01594 #endif 01595 ArrayRCP<Scalar> ncbuff; 01596 if (getLocalLength() > 0) { 01597 ncbuff = getSubArrayRCP (MVT::getValuesNonConst (lclMV_), j); 01598 } 01599 return rcp (new V (this->getMap (), ncbuff, COMPUTE_VIEW_CONSTRUCTOR)); 01600 } 01601 01602 01603 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01604 void 01605 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01606 get1dCopy (Teuchos::ArrayView<Scalar> A, size_t LDA) const 01607 { 01608 using Teuchos::ArrayRCP; 01609 using Teuchos::ArrayView; 01610 using Teuchos::RCP; 01611 using Teuchos::rcp; 01612 01613 const std::string tfecfFuncName("get1dCopy(A,LDA)"); 01614 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(LDA < getLocalLength(), std::runtime_error, 01615 ": specified stride is not large enough for local vector length."); 01616 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(Teuchos::as<size_t>(A.size()) < LDA*(getNumVectors()-1)+getLocalLength(), std::runtime_error, 01617 ": specified stride/storage is not large enough for the number of vectors."); 01618 RCP<Node> node = MVT::getNode(lclMV_); 01619 const size_t myStride = MVT::getStride(lclMV_), 01620 numCols = getNumVectors(), 01621 myLen = getLocalLength(); 01622 if (myLen > 0) { 01623 ArrayRCP<const Scalar> mydata = MVT::getValues(lclMV_); 01624 KOKKOS_NODE_TRACE("MultiVector::getVectorNonConst()") 01625 ArrayRCP<const Scalar> myview = node->template viewBuffer<Scalar>(myStride*(numCols-1)+myLen,mydata); 01626 typename ArrayView<Scalar>::iterator Aptr = A.begin(); 01627 for (size_t j=0; j<numCols; j++) { 01628 ArrayRCP<const Scalar> myviewj = getSubArrayRCP(myview,j); 01629 std::copy(myviewj,myviewj+myLen,Aptr); 01630 Aptr += LDA; 01631 } 01632 myview = Teuchos::null; 01633 mydata = Teuchos::null; 01634 } 01635 } 01636 01637 01638 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01639 void 01640 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01641 get2dCopy (Teuchos::ArrayView<const Teuchos::ArrayView<Scalar> > ArrayOfPtrs) const 01642 { 01643 using Teuchos::as; 01644 using Teuchos::ArrayRCP; 01645 using Teuchos::ArrayView; 01646 using Teuchos::RCP; 01647 using Teuchos::rcp; 01648 01649 const std::string tfecfFuncName("get2dCopy(ArrayOfPtrs)"); 01650 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(Teuchos::as<size_t>(ArrayOfPtrs.size()) != getNumVectors(), std::runtime_error, 01651 ": Array of pointers must contain as many pointers as the MultiVector has rows."); 01652 RCP<Node> node = MVT::getNode(lclMV_); 01653 const size_t numCols = getNumVectors(), 01654 myLen = getLocalLength(); 01655 if (myLen > 0) { 01656 ArrayRCP<const Scalar> mybuff = MVT::getValues(lclMV_); 01657 KOKKOS_NODE_TRACE("MultiVector::get2dCopy()") 01658 ArrayRCP<const Scalar> myview = node->template viewBuffer<Scalar>(mybuff.size(), mybuff); 01659 for (size_t j=0; j<numCols; ++j) { 01660 #ifdef HAVE_TPETRA_DEBUG 01661 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(as<size_t>(ArrayOfPtrs[j].size()) != getLocalLength(), std::runtime_error, 01662 ": The ArrayView provided in ArrayOfPtrs[" << j << "] was not large enough to contain the local entries."); 01663 #endif 01664 ArrayRCP<const Scalar> myviewj = getSubArrayRCP(myview,j); 01665 std::copy(myviewj,myviewj+myLen,ArrayOfPtrs[j].begin()); 01666 } 01667 myview = Teuchos::null; 01668 mybuff = Teuchos::null; 01669 } 01670 } 01671 01672 01673 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01674 Teuchos::ArrayRCP<const Scalar> 01675 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::get1dView () const 01676 { 01677 TEUCHOS_TEST_FOR_EXCEPTION(!isConstantStride(), std::runtime_error, 01678 "Tpetra::MultiVector::get1dView() requires that this MultiVector have constant stride."); 01679 Teuchos::RCP<Node> node = MVT::getNode(lclMV_); 01680 KOKKOS_NODE_TRACE("MultiVector::get1dView()") 01681 return node->template viewBuffer<Scalar>( getStride()*(getNumVectors()-1)+getLocalLength(), MVT::getValues(lclMV_) ); 01682 } 01683 01684 01685 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01686 Teuchos::ArrayRCP<Scalar> 01687 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::get1dViewNonConst () 01688 { 01689 TEUCHOS_TEST_FOR_EXCEPTION(!isConstantStride(), std::runtime_error, 01690 "Tpetra::MultiVector::get1dViewNonConst(): requires that this MultiVector have constant stride."); 01691 Teuchos::RCP<Node> node = MVT::getNode(lclMV_); 01692 KOKKOS_NODE_TRACE("MultiVector::get1dViewNonConst()") 01693 return node->template viewBufferNonConst<Scalar>( Kokkos::ReadWrite, getStride()*(getNumVectors()-1)+getLocalLength(), MVT::getValuesNonConst(lclMV_) ); 01694 } 01695 01696 01697 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01698 Teuchos::ArrayRCP<Teuchos::ArrayRCP<Scalar> > 01699 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::get2dViewNonConst() 01700 { 01701 using Teuchos::arcp; 01702 using Teuchos::ArrayRCP; 01703 using Teuchos::RCP; 01704 01705 RCP<Node> node = MVT::getNode(lclMV_); 01706 ArrayRCP<ArrayRCP<Scalar> > views = arcp<ArrayRCP<Scalar> > (getNumVectors ()); 01707 if (isConstantStride ()) { 01708 const size_t myStride = MVT::getStride(lclMV_), 01709 numCols = getNumVectors(), 01710 myLen = getLocalLength(); 01711 if (myLen > 0) { 01712 KOKKOS_NODE_TRACE("MultiVector::get2dViewNonConst()") 01713 ArrayRCP<Scalar> myview = node->template viewBufferNonConst<Scalar>(Kokkos::ReadWrite,myStride*(numCols-1)+myLen,MVT::getValuesNonConst(lclMV_)); 01714 for (size_t j=0; j<numCols; ++j) { 01715 views[j] = myview.persistingView(0,myLen); 01716 myview += myStride; 01717 } 01718 } 01719 } 01720 else { 01721 const size_t myStride = MVT::getStride(lclMV_), 01722 numCols = MVT::getNumCols(lclMV_), 01723 myCols = getNumVectors(), 01724 myLen = MVT::getNumRows(lclMV_); 01725 if (myLen > 0) { 01726 KOKKOS_NODE_TRACE("MultiVector::get2dViewNonConst()") 01727 ArrayRCP<Scalar> myview = node->template viewBufferNonConst<Scalar>(Kokkos::ReadWrite,myStride*(numCols-1)+myLen,MVT::getValuesNonConst(lclMV_)); 01728 for (size_t j=0; j<myCols; ++j) { 01729 views[j] = myview.persistingView(whichVectors_[j]*myStride,myLen); 01730 } 01731 } 01732 } 01733 return views; 01734 } 01735 01736 01737 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01738 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const Scalar> > 01739 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::get2dView() const 01740 { 01741 using Teuchos::arcp; 01742 using Teuchos::ArrayRCP; 01743 using Teuchos::RCP; 01744 01745 RCP<Node> node = MVT::getNode(lclMV_); 01746 ArrayRCP< ArrayRCP<const Scalar> > views = arcp<ArrayRCP<const Scalar> >(getNumVectors()); 01747 if (isConstantStride()) { 01748 const size_t myStride = MVT::getStride(lclMV_), 01749 numCols = getNumVectors(), 01750 myLen = getLocalLength(); 01751 if (myLen > 0) { 01752 KOKKOS_NODE_TRACE("MultiVector::get2dView()") 01753 ArrayRCP<const Scalar> myview = node->template viewBuffer<Scalar>(myStride*(numCols-1)+myLen,MVT::getValues(lclMV_)); 01754 for (size_t j=0; j<numCols; ++j) { 01755 views[j] = myview.persistingView(0,myLen); 01756 myview += myStride; 01757 } 01758 } 01759 } 01760 else { 01761 const size_t myStride = MVT::getStride(lclMV_), 01762 numCols = MVT::getNumCols(lclMV_), 01763 myCols = getNumVectors(), 01764 myLen = MVT::getNumRows(lclMV_); 01765 if (myLen > 0) { 01766 KOKKOS_NODE_TRACE("MultiVector::get2dView()") 01767 ArrayRCP<const Scalar> myview = node->template viewBuffer<Scalar>(myStride*(numCols-1)+myLen,MVT::getValues(lclMV_)); 01768 for (size_t j=0; j<myCols; ++j) { 01769 views[j] = myview.persistingView(whichVectors_[j]*myStride,myLen); 01770 } 01771 } 01772 } 01773 return views; 01774 } 01775 01776 01777 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01778 void 01779 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01780 multiply (Teuchos::ETransp transA, 01781 Teuchos::ETransp transB, 01782 const Scalar &alpha, 01783 const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A, 01784 const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& B, 01785 const Scalar &beta) 01786 { 01787 using Teuchos::NO_TRANS; // enums 01788 using Teuchos::TRANS; 01789 using Teuchos::CONJ_TRANS; 01790 using Teuchos::null; 01791 using Teuchos::ScalarTraits; // traits 01792 using Teuchos::as; 01793 using Teuchos::RCP; 01794 using Teuchos::rcp; 01795 typedef MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> MV; 01796 01797 // This routine performs a variety of matrix-matrix multiply operations, interpreting 01798 // the MultiVector (this-aka C , A and B) as 2D matrices. Variations are due to 01799 // the fact that A, B and C can be local replicated or global distributed 01800 // MultiVectors and that we may or may not operate with the transpose of 01801 // A and B. Possible cases are: 01802 // Num 01803 // OPERATIONS cases Notes 01804 // 1) C(local) = A^X(local) * B^X(local) 4 (X=Trans or Not, No comm needed) 01805 // 2) C(local) = A^T(distr) * B (distr) 1 (2D dot product, replicate C) 01806 // 3) C(distr) = A (distr) * B^X(local) 2 (2D vector update, no comm needed) 01807 // 01808 // The following operations are not meaningful for 1D distributions: 01809 // 01810 // u1) C(local) = A^T(distr) * B^T(distr) 1 01811 // u2) C(local) = A (distr) * B^X(distr) 2 01812 // u3) C(distr) = A^X(local) * B^X(local) 4 01813 // u4) C(distr) = A^X(local) * B^X(distr) 4 01814 // u5) C(distr) = A^T(distr) * B^X(local) 2 01815 // u6) C(local) = A^X(distr) * B^X(local) 4 01816 // u7) C(distr) = A^X(distr) * B^X(local) 4 01817 // u8) C(local) = A^X(local) * B^X(distr) 4 01818 // 01819 // Total of 32 case (2^5). 01820 01821 const std::string errPrefix("Tpetra::MultiVector::multiply(transOpA,transOpB,alpha,A,B,beta): "); 01822 01823 TEUCHOS_TEST_FOR_EXCEPTION( ScalarTraits<Scalar>::isComplex && (transA == TRANS || transB == TRANS), std::invalid_argument, 01824 errPrefix << "non-conjugate transpose not supported for complex types."); 01825 transA = (transA == NO_TRANS ? NO_TRANS : CONJ_TRANS); 01826 transB = (transB == NO_TRANS ? NO_TRANS : CONJ_TRANS); 01827 01828 // Compute effective dimensions, w.r.t. transpose operations on 01829 size_t A_nrows = (transA==CONJ_TRANS) ? A.getNumVectors() : A.getLocalLength(); 01830 size_t A_ncols = (transA==CONJ_TRANS) ? A.getLocalLength() : A.getNumVectors(); 01831 size_t B_nrows = (transB==CONJ_TRANS) ? B.getNumVectors() : B.getLocalLength(); 01832 size_t B_ncols = (transB==CONJ_TRANS) ? B.getLocalLength() : B.getNumVectors(); 01833 01834 Scalar beta_local = beta; // local copy of beta; might be reassigned below 01835 01836 TEUCHOS_TEST_FOR_EXCEPTION( getLocalLength() != A_nrows || getNumVectors() != B_ncols || A_ncols != B_nrows, std::runtime_error, 01837 errPrefix << "dimension of *this, op(A) and op(B) must be consistent."); 01838 01839 bool A_is_local = !A.isDistributed(); 01840 bool B_is_local = !B.isDistributed(); 01841 bool C_is_local = !this->isDistributed(); 01842 bool Case1 = ( C_is_local && A_is_local && B_is_local); // Case 1: C(local) = A^X(local) * B^X(local) 01843 bool Case2 = ( C_is_local && !A_is_local && !B_is_local && transA==CONJ_TRANS && transB==NO_TRANS); // Case 2: C(local) = A^T(distr) * B (distr) 01844 bool Case3 = (!C_is_local && !A_is_local && B_is_local && transA==NO_TRANS ); // Case 3: C(distr) = A (distr) * B^X(local) 01845 01846 // Test that we are considering a meaningful cases 01847 TEUCHOS_TEST_FOR_EXCEPTION( !Case1 && !Case2 && !Case3, std::runtime_error, 01848 errPrefix << "multiplication of op(A) and op(B) into *this is not a supported use case."); 01849 01850 if (beta != ScalarTraits<Scalar>::zero() && Case2) 01851 { 01852 // if Case2, then C is local and contributions must be summed across all nodes 01853 // however, if beta != 0, then accumulate beta*C into the sum 01854 // when summing across all nodes, we only want to accumulate this once, so 01855 // set beta == 0 on all nodes except node 0 01856 int MyPID = this->getMap()->getComm()->getRank(); 01857 if (MyPID!=0) beta_local = ScalarTraits<Scalar>::zero(); 01858 } 01859 01860 // Check if A, B, C have constant stride, if not then make temp copy (strided) 01861 RCP<const MV> Atmp, Btmp; 01862 RCP<MV> Ctmp; 01863 if (isConstantStride() == false) Ctmp = rcp (new MV (*this)); 01864 else Ctmp = rcp(this,false); 01865 01866 if (A.isConstantStride() == false) Atmp = rcp (new MV (A)); 01867 else Atmp = rcp(&A,false); 01868 01869 if (B.isConstantStride() == false) Btmp = rcp (new MV (B)); 01870 else Btmp = rcp(&B,false); 01871 01872 #ifdef HAVE_TEUCHOS_DEBUG 01873 TEUCHOS_TEST_FOR_EXCEPTION(!Ctmp->isConstantStride() || !Btmp->isConstantStride() || !Atmp->isConstantStride(), std::logic_error, 01874 errPrefix << "failed making temporary strided copies of input multivectors."); 01875 #endif 01876 01877 KMV &C_mv = Ctmp->lclMV_; 01878 { 01879 // get local multivectors 01880 const KMV &A_mv = Atmp->lclMV_; 01881 const KMV &B_mv = Btmp->lclMV_; 01882 // do the multiply (GEMM) 01883 MVT::GEMM(C_mv,transA,transB,alpha,A_mv,B_mv,beta_local); 01884 } 01885 01886 // Dispose of (possibly) extra copies of A, B 01887 Atmp = null; 01888 Btmp = null; 01889 01890 RCP<Node> node = MVT::getNode(lclMV_); 01891 // If *this was not strided, copy the data from the strided version and then delete it 01892 if (! isConstantStride ()) { 01893 // *this is not strided, we must put data from Ctmp into *this 01894 TEUCHOS_TEST_FOR_EXCEPT(&C_mv != &lclMV_); 01895 const size_t numVecs = MVT::getNumCols(lclMV_); 01896 for (size_t j=0; j < numVecs; ++j) { 01897 KOKKOS_NODE_TRACE("MultiVector::multiply()") 01898 node->template copyBuffers<Scalar>(getLocalLength(),MVT::getValues(C_mv,j),MVT::getValuesNonConst(lclMV_,whichVectors_[j])); 01899 } 01900 } 01901 01902 // If Case 2 then sum up *this and distribute it to all processors. 01903 if (Case2) { 01904 this->reduce(); 01905 } 01906 } 01907 01908 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01909 void 01910 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 01911 elementWiseMultiply (Scalar scalarAB, 01912 const Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& A, 01913 const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& B, 01914 Scalar scalarThis) 01915 { 01916 using Teuchos::arcp_const_cast; 01917 01918 const std::string tfecfFuncName("elementWiseMultiply()"); 01919 01920 #ifdef HAVE_TPETRA_DEBUG 01921 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( !this->getMap()->isCompatible(*A.getMap()) || !this->getMap()->isCompatible(*B.getMap()), std::runtime_error, 01922 ": MultiVectors do not have compatible Maps:" << std::endl 01923 << "this->getMap(): " << std::endl << *this->getMap() 01924 << "A.getMap(): " << std::endl << *A.getMap() << std::endl 01925 << "B.getMap(): " << std::endl << *B.getMap() << std::endl); 01926 #else 01927 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( getLocalLength() != A.getLocalLength() || getLocalLength() != B.getLocalLength(), std::runtime_error, 01928 ": MultiVectors do not have the same local length."); 01929 #endif 01930 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(B.getNumVectors() != this->getNumVectors(), std::runtime_error, 01931 ": MultiVectors 'this' and B must have the same number of vectors."); 01932 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(A.getNumVectors() != 1, std::runtime_error, 01933 ": MultiVectors A must have just 1 vector."); 01934 try { 01935 MVT::ElemMult(lclMV_,scalarThis,scalarAB,(const KMV &)A.lclMV_,(const KMV &)B.lclMV_); 01936 } 01937 catch (std::runtime_error &e) { 01938 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(true,std::runtime_error, 01939 ": caught exception from Kokkos:" << std::endl 01940 << e.what() << std::endl); 01941 } 01942 } 01943 01944 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 01945 void 01946 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::reduce() 01947 { 01948 using Teuchos::Array; 01949 using Teuchos::ArrayView; 01950 using Teuchos::ArrayRCP; 01951 using Teuchos::as; 01952 using Teuchos::RCP; 01953 using Teuchos::reduceAll; 01954 using Teuchos::REDUCE_SUM; 01955 01956 // This method should be called only for locally replicated 01957 // MultiVectors (! isDistributed ()). 01958 TEUCHOS_TEST_FOR_EXCEPTION(this->isDistributed (), std::runtime_error, 01959 "Tpetra::MultiVector::reduce() should only be called with locally " 01960 "replicated or otherwise not distributed MultiVector objects."); 01961 RCP<const Teuchos::Comm<int> > comm = this->getMap()->getComm(); 01962 if (comm->getSize() == 1) return; 01963 RCP<Node> node = MVT::getNode(lclMV_); 01964 // sum the data across all multivectors 01965 // need to have separate packed buffers for send and receive 01966 // if we're packed, we'll just set the receive buffer as our data, the send as a copy 01967 // if we're not packed, we'll use allocated buffers for both. 01968 ArrayView<Scalar> target; 01969 const size_t myStride = MVT::getStride(lclMV_), 01970 numCols = MVT::getNumCols(lclMV_), 01971 myLen = MVT::getNumRows(lclMV_); 01972 Array<Scalar> sourceBuffer(numCols*myLen), tmparr(0); 01973 bool packed = isConstantStride() && (myStride == myLen); 01974 KOKKOS_NODE_TRACE("MultiVector::reduce()") 01975 ArrayRCP<Scalar> bufView = node->template viewBufferNonConst<Scalar>( 01976 Kokkos::ReadWrite,myStride*(numCols-1)+myLen, 01977 MVT::getValuesNonConst(lclMV_) ); 01978 if (packed) { 01979 // copy data from view to sourceBuffer, reduce into view below 01980 target = bufView(0,myLen*numCols); 01981 std::copy(target.begin(),target.end(),sourceBuffer.begin()); 01982 } 01983 else { 01984 // copy data into sourceBuffer, reduce from sourceBuffer into tmparr below, copy back to view after that 01985 tmparr.resize(myLen*numCols); 01986 Scalar *sptr = sourceBuffer.getRawPtr(); 01987 ArrayRCP<const Scalar> vptr = bufView; 01988 for (size_t j=0; j<numCols; ++j) { 01989 std::copy(vptr,vptr+myLen,sptr); 01990 sptr += myLen; 01991 vptr += myStride; 01992 } 01993 target = tmparr(); 01994 } 01995 // reduce 01996 reduceAll<int,Scalar> (*comm, REDUCE_SUM, as<int> (numCols*myLen), 01997 sourceBuffer.getRawPtr (), target.getRawPtr ()); 01998 if (! packed) { 01999 // copy tmparr back into view 02000 const Scalar *sptr = tmparr.getRawPtr(); 02001 ArrayRCP<Scalar> vptr = bufView; 02002 for (size_t j=0; j<numCols; ++j) { 02003 std::copy (sptr, sptr+myLen, vptr); 02004 sptr += myLen; 02005 vptr += myStride; 02006 } 02007 } 02008 bufView = Teuchos::null; 02009 } 02010 02011 02012 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02013 void 02014 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 02015 replaceLocalValue (LocalOrdinal MyRow, 02016 size_t VectorIndex, 02017 const Scalar &ScalarValue) 02018 { 02019 #ifdef HAVE_TPETRA_DEBUG 02020 const std::string tfecfFuncName("replaceLocalValue()"); 02021 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(MyRow < this->getMap()->getMinLocalIndex() || MyRow > this->getMap()->getMaxLocalIndex(), 02022 std::runtime_error, ": row index is invalid."); 02023 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( vectorIndexOutOfRange(VectorIndex), 02024 std::runtime_error, ": vector index is invalid."); 02025 #endif 02026 getDataNonConst(VectorIndex)[MyRow] = ScalarValue; 02027 } 02028 02029 02030 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02031 void 02032 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 02033 sumIntoLocalValue (LocalOrdinal MyRow, 02034 size_t VectorIndex, 02035 const Scalar &ScalarValue) 02036 { 02037 #ifdef HAVE_TPETRA_DEBUG 02038 const std::string tfecfFuncName("sumIntoLocalValue"); 02039 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(MyRow < this->getMap()->getMinLocalIndex() || MyRow > this->getMap()->getMaxLocalIndex(), std::runtime_error, 02040 ": row index is invalid."); 02041 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( vectorIndexOutOfRange(VectorIndex), std::runtime_error, 02042 ": vector index is invalid."); 02043 #endif 02044 getDataNonConst(VectorIndex)[MyRow] += ScalarValue; 02045 } 02046 02047 02048 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02049 void 02050 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 02051 replaceGlobalValue (GlobalOrdinal GlobalRow, size_t VectorIndex, const Scalar &ScalarValue) 02052 { 02053 LocalOrdinal MyRow = this->getMap()->getLocalElement(GlobalRow); 02054 #ifdef HAVE_TPETRA_DEBUG 02055 const std::string tfecfFuncName("replaceGlobalValue()"); 02056 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(MyRow == Teuchos::OrdinalTraits<LocalOrdinal>::invalid(), std::runtime_error, 02057 ": row index is not present on this processor."); 02058 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( vectorIndexOutOfRange(VectorIndex), std::runtime_error, 02059 ": vector index is invalid."); 02060 #endif 02061 getDataNonConst(VectorIndex)[MyRow] = ScalarValue; 02062 } 02063 02064 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02065 void MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::sumIntoGlobalValue(GlobalOrdinal GlobalRow, size_t VectorIndex, const Scalar &ScalarValue) 02066 { 02067 LocalOrdinal MyRow = this->getMap()->getLocalElement(GlobalRow); 02068 #ifdef HAVE_TEUCHOS_DEBUG 02069 const std::string tfecfFuncName("sumIntoGlobalValue()"); 02070 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(MyRow == Teuchos::OrdinalTraits<LocalOrdinal>::invalid(), std::runtime_error, 02071 ": row index is not present on this processor."); 02072 TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC( vectorIndexOutOfRange(VectorIndex), std::runtime_error, 02073 ": vector index is invalid."); 02074 #endif 02075 getDataNonConst(VectorIndex)[MyRow] += ScalarValue; 02076 } 02077 02078 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02079 template <class T> 02080 Teuchos::ArrayRCP<T> 02081 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 02082 getSubArrayRCP (Teuchos::ArrayRCP<T> arr, 02083 size_t j) const 02084 { 02085 const size_t stride = MVT::getStride (lclMV_); 02086 const size_t myLen = getLocalLength (); 02087 Teuchos::ArrayRCP<T> ret; 02088 if (isConstantStride()) { 02089 ret = arr.persistingView (j*stride, myLen); 02090 } 02091 else { 02092 ret = arr.persistingView (whichVectors_[j]*stride, myLen); 02093 } 02094 return ret; 02095 } 02096 02097 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02098 const Kokkos::MultiVector<Scalar,Node>& 02099 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalMV() const { 02100 return lclMV_; 02101 } 02102 02103 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02104 Kokkos::MultiVector<Scalar,Node>& 02105 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalMVNonConst() { 02106 return lclMV_; 02107 } 02108 02109 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02110 std::string 02111 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::description() const 02112 { 02113 std::ostringstream oss; 02114 oss << Teuchos::Describable::description(); 02115 oss << "{length = "<<getGlobalLength() 02116 << ", getNumVectors = "<<getNumVectors() 02117 << ", isConstantStride() = "<<isConstantStride() 02118 << "}"; 02119 return oss.str(); 02120 } 02121 02122 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02123 void 02124 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 02125 describe (Teuchos::FancyOStream &out, 02126 const Teuchos::EVerbosityLevel verbLevel) const 02127 { 02128 using Teuchos::ArrayRCP; 02129 using Teuchos::RCP; 02130 using Teuchos::VERB_DEFAULT; 02131 using Teuchos::VERB_NONE; 02132 using Teuchos::VERB_LOW; 02133 using Teuchos::VERB_MEDIUM; 02134 using Teuchos::VERB_HIGH; 02135 using Teuchos::VERB_EXTREME; 02136 using std::endl; 02137 using std::setw; 02138 02139 // Set default verbosity if applicable. 02140 const Teuchos::EVerbosityLevel vl = 02141 (verbLevel == VERB_DEFAULT) ? VERB_LOW : verbLevel; 02142 02143 RCP<const Teuchos::Comm<int> > comm = this->getMap()->getComm(); 02144 const int myImageID = comm->getRank(); 02145 const int numImages = comm->getSize(); 02146 02147 if (vl != VERB_NONE) { 02148 // Don't set the tab level unless we're printing something. 02149 Teuchos::OSTab tab (out); 02150 02151 if (myImageID == 0) { // >= VERB_LOW prints description() 02152 out << this->description() << endl; 02153 } 02154 for (int imageCtr = 0; imageCtr < numImages; ++imageCtr) { 02155 if (myImageID == imageCtr) { 02156 if (vl != VERB_LOW) { 02157 // At verbosity > VERB_LOW, each process prints something. 02158 out << "Process " << myImageID << ":" << endl; 02159 02160 Teuchos::OSTab procTab (out); 02161 // >= VERB_MEDIUM: print the local vector length. 02162 out << "local length=" << getLocalLength(); 02163 if (vl != VERB_MEDIUM) { 02164 // >= VERB_HIGH: print isConstantStride() and getStride() 02165 if (isConstantStride()) { 02166 out << ", constant stride=" << getStride() << endl; 02167 } 02168 else { 02169 out << ", not constant stride" << endl; 02170 } 02171 if (vl == VERB_EXTREME) { 02172 // VERB_EXTREME: print all the values in the multivector. 02173 out << "Values:" << endl; 02174 ArrayRCP<ArrayRCP<const Scalar> > X = this->get2dView(); 02175 for (size_t i = 0; i < getLocalLength(); ++i) { 02176 for (size_t j = 0; j < getNumVectors(); ++j) { 02177 out << X[j][i]; 02178 if (j < getNumVectors() - 1) { 02179 out << " "; 02180 } 02181 } // for each column 02182 out << endl; 02183 } // for each row 02184 } // if vl == VERB_EXTREME 02185 } // if (vl != VERB_MEDIUM) 02186 else { // vl == VERB_LOW 02187 out << endl; 02188 } 02189 } // if vl != VERB_LOW 02190 } // if it is my process' turn to print 02191 } // for each process in the communicator 02192 } // if vl != VERB_NONE 02193 } 02194 02195 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02196 void 02197 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 02198 createViews() const 02199 { 02200 TPETRA_EFFICIENCY_WARNING(! cview_.is_null (), std::runtime_error, 02201 "::createViews(): The const view " 02202 "has already been created and is therefore not null. (For" 02203 "Tpetra developers: cview_.total_count() = " << cview_.total_count () 02204 << ". This " 02205 "means that MultiVector is either creating a view unnecessarily, or " 02206 "hanging on to a view beyond its needed scope. This " 02207 "probably does not affect correctness, it but does affect total memory " 02208 "use. Please report this performance bug to the Tpetra developers."); 02209 02210 Teuchos::RCP<Node> node = this->getMap ()->getNode (); 02211 if (cview_.is_null () && getLocalLength () > 0) { 02212 Teuchos::ArrayRCP<const Scalar> buff = MVT::getValues (lclMV_); 02213 KOKKOS_NODE_TRACE("MultiVector::createViews()") 02214 cview_ = node->template viewBuffer<Scalar> (buff.size (), buff); 02215 } 02216 } 02217 02218 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02219 void 02220 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>:: 02221 createViewsNonConst (Kokkos::ReadWriteOption rwo) 02222 { 02223 02224 TPETRA_EFFICIENCY_WARNING(! ncview_.is_null (), std::runtime_error, 02225 "::createViewsNonConst(): The nonconst view " 02226 "has already been created and is therefore not null. (For" 02227 "Tpetra developers: ncview_.total_count() = " << ncview_.total_count () 02228 << ". This " 02229 "means that MultiVector is either creating a view unnecessarily, or " 02230 "hanging on to a view beyond its needed scope. This " 02231 "probably does not affect correctness, it but does affect total memory " 02232 "use. Please report this performance bug to the Tpetra developers."); 02233 02234 Teuchos::RCP<Node> node = this->getMap ()->getNode (); 02235 if (ncview_.is_null () && getLocalLength () > 0) { 02236 Teuchos::ArrayRCP<Scalar> buff = MVT::getValuesNonConst (lclMV_); 02237 KOKKOS_NODE_TRACE("MultiVector::createViewsNonConst()") 02238 ncview_ = node->template viewBufferNonConst<Scalar> (rwo, buff.size (), buff); 02239 } 02240 } 02241 02242 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node> 02243 void 02244 MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::releaseViews () const 02245 { 02246 const int constViewCount = cview_.total_count (); 02247 const int nonconstViewCount = ncview_.total_count (); 02248 02249 // This prevents unused variable compiler warnings, in case Tpetra 02250 // efficiency warnings aren't enabled. 02251 (void) constViewCount; 02252 (void) nonconstViewCount; 02253 02254 TPETRA_EFFICIENCY_WARNING(constViewCount > 1 || nonconstViewCount > 1, 02255 std::runtime_error, "::releaseViews(): Either the const view or the " 02256 "nonconst view has a reference count greater than 1. (For Tpetra " 02257 "developers: cview_.total_count() = " << constViewCount << " and ncview_." 02258 "total_count() = " << nonconstViewCount << ". This means that release" 02259 "Views() won't actually free memory. This probably does not affect " 02260 "correctness, it but does affect total memory use. Please report this " 02261 "performance bug to the Tpetra developers."); 02262 02263 cview_ = Teuchos::null; 02264 ncview_ = Teuchos::null; 02265 } 02266 } // namespace Tpetra 02267 02268 // 02269 // Explicit instantiation macro 02270 // 02271 // Must be expanded from within the Tpetra namespace! 02272 // 02273 02274 #define TPETRA_MULTIVECTOR_INSTANT(SCALAR,LO,GO,NODE) \ 02275 \ 02276 template class MultiVector< SCALAR , LO , GO , NODE >; \ 02277 02278 02279 #endif // TPETRA_MULTIVECTOR_DEF_HPP
1.7.4