All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Xpetra_TpetraMultiVector_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Xpetra: A linear algebra interface package
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef XPETRA_TPETRAMULTIVECTOR_DEF_HPP
47 #define XPETRA_TPETRAMULTIVECTOR_DEF_HPP
49 
50 #include "Xpetra_TpetraMap.hpp" //TMP
51 #include "Xpetra_Utils.hpp"
52 #include "Xpetra_TpetraImport.hpp"
53 #include "Xpetra_TpetraExport.hpp"
54 
56 #include "Tpetra_MultiVector.hpp"
57 #include "Tpetra_Vector.hpp"
58 
59 namespace Xpetra {
60 
61 
63  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
65  : vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(map), NumVectors, zeroOut))) {
66  // TAW 1/30/2016: even though Tpetra allows numVecs == 0, Epetra does not. Introduce exception to keep behavior of Epetra and Tpetra consistent.
67  TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument, "Xpetra::TpetraMultiVector(map,numVecs,zeroOut): numVecs = " << NumVectors << " < 1.");
68  }
69 
71  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
74  : vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(Tpetra::createCopy(toTpetra(source))))) { }
75 
77  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
80  : vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(map), A, LDA, NumVectors))) {
81  // TAW 1/30/2016: even though Tpetra allows numVecs == 0, Epetra does not. Introduce exception to keep behavior of Epetra and Tpetra consistent.
82  TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument, "Xpetra::TpetraMultiVector(map,A,LDA,numVecs): numVecs = " << NumVectors << " < 1.");
83  }
84 
86  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
89  : vec_(Teuchos::rcp(new Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(map), ArrayOfPtrs, NumVectors))) {
90  // TAW 1/30/2016: even though Tpetra allows numVecs == 0, Epetra does not. Introduce exception to keep behavior of Epetra and Tpetra consistent.
91  TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument, "Xpetra::TpetraMultiVector(map,ArrayOfPtrs,numVecs): numVecs = " << NumVectors << " < 1.");
92  }
93 
94 
96  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
99 
101  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
103  replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::replaceGlobalValue"); vec_->replaceGlobalValue(globalRow, vectorIndex, value); }
104 
106  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
108  sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::sumIntoGlobalValue"); vec_->sumIntoGlobalValue(globalRow, vectorIndex, value); }
109 
111  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
113  replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::replaceLocalValue"); vec_->replaceLocalValue(myRow, vectorIndex, value); }
114 
116  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
118  sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::sumIntoLocalValue"); vec_->sumIntoLocalValue(myRow, vectorIndex, value); }
119 
121  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
123  putScalar(const Scalar &value) { XPETRA_MONITOR("TpetraMultiVector::putScalar"); vec_->putScalar(value); }
124 
126  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
128  reduce() { XPETRA_MONITOR("TpetraMultiVector::reduce"); vec_->reduce(); }
129 
130 
131  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
134  getVector(size_t j) const { XPETRA_MONITOR("TpetraMultiVector::getVector"); return toXpetra(vec_->getVector(j)); }
135 
137  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
140  getVectorNonConst(size_t j) { XPETRA_MONITOR("TpetraMultiVector::getVectorNonConst"); return toXpetra(vec_->getVectorNonConst(j)); }
141 
143  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
146  getData(size_t j) const { XPETRA_MONITOR("TpetraMultiVector::getData"); return vec_->getData(j); }
147 
149  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
152  getDataNonConst(size_t j) { XPETRA_MONITOR("TpetraMultiVector::getDataNonConst"); return vec_->getDataNonConst(j); }
153 
155  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
156  void
158  get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const { XPETRA_MONITOR("TpetraMultiVector::get1dCopy"); vec_->get1dCopy(A, LDA); }
159 
161  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
162  void
164  get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const { XPETRA_MONITOR("TpetraMultiVector::get2dCopy"); vec_->get2dCopy(ArrayOfPtrs); }
165 
167  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
170  get1dView() const { XPETRA_MONITOR("TpetraMultiVector::get1dView"); return vec_->get1dView(); }
171 
173  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
176  get2dView() const { XPETRA_MONITOR("TpetraMultiVector::get2dView"); return vec_->get2dView(); }
177 
179  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
182  get1dViewNonConst() { XPETRA_MONITOR("TpetraMultiVector::get1dViewNonConst"); return vec_->get1dViewNonConst(); }
183 
185  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
188  get2dViewNonConst() { XPETRA_MONITOR("TpetraMultiVector::get2dViewNonConst"); return vec_->get2dViewNonConst(); }
189 
191  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
193  dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const { XPETRA_MONITOR("TpetraMultiVector::dot"); vec_->dot(toTpetra(A), dots); }
194 
196  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
198  abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("TpetraMultiVector::abs"); vec_->abs(toTpetra(A)); }
199 
201  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
203  reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("TpetraMultiVector::reciprocal"); vec_->reciprocal(toTpetra(A)); }
204 
206  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
208  scale(const Scalar &alpha) { XPETRA_MONITOR("TpetraMultiVector::scale"); vec_->scale(alpha); }
209 
211  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
213  scale(Teuchos::ArrayView< const Scalar > alpha) { XPETRA_MONITOR("TpetraMultiVector::scale"); vec_->scale(alpha); }
214 
216  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
218  scale(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("TpetraMultiVector::scale"); vec_->scale(alpha, toTpetra(A)); }
219 
221  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
223  update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta) { XPETRA_MONITOR("TpetraMultiVector::update"); vec_->update(alpha, toTpetra(A), beta); }
224 
226  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
228  update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &gamma) { XPETRA_MONITOR("TpetraMultiVector::update"); vec_->update(alpha, toTpetra(A), beta, toTpetra(B), gamma); }
229 
231  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
233  norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("TpetraMultiVector::norm1"); vec_->norm1(norms); }
234 
236  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
238  norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("TpetraMultiVector::norm2"); vec_->norm2(norms); }
239 
241  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
243  normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("TpetraMultiVector::normInf"); vec_->normInf(norms); }
244 
246  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
248  meanValue(const Teuchos::ArrayView< Scalar > &means) const { XPETRA_MONITOR("TpetraMultiVector::meanValue"); vec_->meanValue(means); }
249 
251  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
253  multiply(Teuchos::ETransp transA, Teuchos::ETransp transB, const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &beta) { XPETRA_MONITOR("TpetraMultiVector::multiply"); vec_->multiply(transA, transB, alpha, toTpetra(A), toTpetra(B), beta); }
254 
255 
257  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
259  getNumVectors() const { XPETRA_MONITOR("TpetraMultiVector::getNumVectors"); return vec_->getNumVectors(); }
260 
262  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
264  getLocalLength() const { XPETRA_MONITOR("TpetraMultiVector::getLocalLength"); return vec_->getLocalLength(); }
265 
267  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
269  getGlobalLength() const { XPETRA_MONITOR("TpetraMultiVector::getGlobalLength"); return vec_->getGlobalLength(); }
270 
271  // \brief Checks to see if the local length, number of vectors and size of Scalar type match
272  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
274  isSameSize(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & vec) const { XPETRA_MONITOR("TpetraMultiVector::isSameSize"); return vec_->isSameSize(toTpetra(vec));}
275 
277  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
279  description() const { XPETRA_MONITOR("TpetraMultiVector::description"); return vec_->description(); }
280 
282  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
284  describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const { XPETRA_MONITOR("TpetraMultiVector::describe"); vec_->describe(out, verbLevel); }
285 
287  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
289  randomize(bool bUseXpetraImplementation) {
290  XPETRA_MONITOR("TpetraMultiVector::randomize");
291 
292  if(bUseXpetraImplementation)
294  else
295  vec_->randomize();
296  }
297 
298  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
301  getMap() const { XPETRA_MONITOR("TpetraMultiVector::getMap"); return toXpetra(vec_->getMap()); }
302 
303  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
306  XPETRA_MONITOR("TpetraMultiVector::doImport");
307 
308  XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); //TODO: remove and use toTpetra()
310  this->getTpetra_MultiVector()->doImport(*v, toTpetra(importer), toTpetra(CM));
311  }
312 
313  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
316  XPETRA_MONITOR("TpetraMultiVector::doExport");
317 
318  XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); //TODO: remove and use toTpetra()
320  this->getTpetra_MultiVector()->doExport(*v, toTpetra(importer), toTpetra(CM));
321 
322  }
323 
324  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
327  XPETRA_MONITOR("TpetraMultiVector::doImport");
328 
329  XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); //TODO: remove and use toTpetra()
331  this->getTpetra_MultiVector()->doImport(*v, toTpetra(exporter), toTpetra(CM));
332 
333  }
334 
335  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
338  XPETRA_MONITOR("TpetraMultiVector::doExport");
339 
340  XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); //TODO: remove and use toTpetra()
342  this->getTpetra_MultiVector()->doExport(*v, toTpetra(exporter), toTpetra(CM));
343 
344  }
345 
346  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
349  XPETRA_MONITOR("TpetraMultiVector::replaceMap");
350  this->getTpetra_MultiVector()->replaceMap(toTpetra(map));
351  }
352 
353 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
354  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
355  template<class Node2>
358  clone(const RCP<Node2> &node2) const {
359  XPETRA_MONITOR("TpetraMultiVector::clone");
361  }
362 #endif
363 
364 
366  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
368  TpetraMultiVector(const Teuchos::RCP<Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > &vec) : vec_(vec) { } //TODO removed const
369 
371  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
374  getTpetra_MultiVector() const { return vec_; }
375 
377  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
379  setSeed(unsigned int seed) { XPETRA_MONITOR("TpetraMultiVector::seedrandom"); Teuchos::ScalarTraits< Scalar >::seedrandom(seed); }
380 
381 
382 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
383 #if 0
394  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
395  template<class TargetDeviceType>
396  typename Kokkos::Impl::if_c<
397  Kokkos::Impl::is_same<
398  typename dual_view_type::t_dev_um::execution_space::memory_space,
399  typename TargetDeviceType::memory_space>::value,
400  typename dual_view_type::t_dev_um,
401  typename dual_view_type::t_host_um>::type
403  getLocalView () const {
404  return this->MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::template getLocalView<TargetDeviceType>();
405  }
406 #endif
407 
408  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
409  typename TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dual_view_type::t_host_um
410  TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
411  getHostLocalView () const {
412  return subview(vec_->template getLocalView<typename TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dual_view_type::host_mirror_space> (),
413  Kokkos::ALL(), Kokkos::ALL());
414  }
415 
416  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
417  typename TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dual_view_type::t_dev_um
418  TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
419  getDeviceLocalView() const {
420  return subview(vec_->template getLocalView<typename TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::dual_view_type::t_dev_um::execution_space> (),
421  Kokkos::ALL(), Kokkos::ALL());
422  }
423 
424 #endif
425 
428  template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
431  {
433  const this_type* rhsPtr = dynamic_cast<const this_type*> (&rhs);
435  rhsPtr == NULL, std::invalid_argument, "Xpetra::MultiVector::operator=:"
436  " The left-hand side (LHS) of the assignment has a different type than "
437  "the right-hand side (RHS). The LHS has type Xpetra::TpetraMultiVector"
438  " (which means it wraps a Tpetra::MultiVector), but the RHS has some "
439  "other type. This probably means that the RHS wraps an "
440  "Epetra_MultiVector. Xpetra::MultiVector does not currently implement "
441  "assignment from an Epetra object to a Tpetra object, though this could"
442  " be added with sufficient interest.");
443 
444  typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> TMV;
445  RCP<const TMV> rhsImpl = rhsPtr->getTpetra_MultiVector ();
446  RCP<TMV> lhsImpl = this->getTpetra_MultiVector ();
447 
449  rhsImpl.is_null (), std::logic_error, "Xpetra::MultiVector::operator= "
450  "(in Xpetra::TpetraMultiVector::assign): *this (the right-hand side of "
451  "the assignment) has a null RCP<Tpetra::MultiVector> inside. Please "
452  "report this bug to the Xpetra developers.");
454  lhsImpl.is_null (), std::logic_error, "Xpetra::MultiVector::operator= "
455  "(in Xpetra::TpetraMultiVector::assign): The left-hand side of the "
456  "assignment has a null RCP<Tpetra::MultiVector> inside. Please report "
457  "this bug to the Xpetra developers.");
458 
459  Tpetra::deep_copy (*lhsImpl, *rhsImpl);
460  }
461 
462 
463 #ifdef HAVE_XPETRA_EPETRA
464 
465 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
466  (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
467 
468  // specialization for TpetraMultiVector on EpetraNode and GO=int
469  template <class Scalar>
470  class TpetraMultiVector<Scalar,int,int,EpetraNode>
471  : public virtual MultiVector< Scalar, int, int, EpetraNode >
472  {
473  typedef int LocalOrdinal;
474  typedef int GlobalOrdinal;
475  typedef EpetraNode Node;
476 
477  // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
479 
480  public:
481 
483 
484 
488  }
489 
491  TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true) {
493  }
494 
498  }
499 
503  }
504 
508  }
509 
510 
512  virtual ~TpetraMultiVector() { }
513 
515 
517 
518 
520  void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { }
521 
523  void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { }
524 
526  void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { }
527 
529  void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { }
530 
532  void putScalar(const Scalar &value) { }
533 
535  void reduce() { }
536 
538 
540 
541 
544 
547 
550 
553 
555  void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const { }
556 
558  void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const { }
559 
562 
565 
568 
571 
573 
575 
576 
579 
582 
585 
587  void scale(const Scalar &alpha) { }
588 
591 
593  void scale(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { }
594 
596  void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta) { }
597 
599  void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &gamma) { }
600 
603 
606 
609 
611  void meanValue(const Teuchos::ArrayView< Scalar > &means) const { }
612 
615 
617 
619 
620 
622  size_t getNumVectors() const { return 0; }
623 
625  size_t getLocalLength() const { return 0; }
626 
628  global_size_t getGlobalLength() const { return 0; }
629 
630  // \! Checks to see if the local length, number of vectors and size of Scalar type match
631  bool isSameSize(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & vec) const { return false; }
632 
634 
636 
637 
639  std::string description() const { return std::string(""); }
640 
643 
645 
648 
650  void randomize(bool bUseXpetraImplementation = false) { }
651 
652  //{@
653  // Implements DistObject interface
654 
656 
658 
660 
662 
664 
666 
667 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
668  template<class Node2>
669  RCP<MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node2> > clone(const RCP<Node2> &node2) const { return Teuchos::null; }
670 #endif
671 
672 
674 
675 
677  TpetraMultiVector(const Teuchos::RCP<Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > &vec) {
679  }
680 
683 
685  void setSeed(unsigned int seed) { }
686 
687 
688 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
690  //typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::host_execution_space host_execution_space;
691  //typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dev_execution_space dev_execution_space;
692 
703  template<class TargetDeviceType>
704  typename Kokkos::Impl::if_c<
705  Kokkos::Impl::is_same<
706  typename dual_view_type::t_dev_um::execution_space::memory_space,
707  typename TargetDeviceType::memory_space>::value,
708  typename dual_view_type::t_dev_um,
709  typename dual_view_type::t_host_um>::type
710  getLocalView () const {
711  typename Kokkos::Impl::if_c<
712  Kokkos::Impl::is_same<
713  typename dual_view_type::t_dev_um::execution_space::memory_space,
714  typename TargetDeviceType::memory_space>::value,
715  typename dual_view_type::t_dev_um,
716  typename dual_view_type::t_host_um>::type dummy;
717  return dummy;
718  }
719 
720  typename dual_view_type::t_host_um getHostLocalView () const {
721  //return subview(vec_->template getLocalView<typename dual_view_type::host_mirror_space> (),
722  // Kokkos::ALL(), Kokkos::ALL());
723  return typename dual_view_type::t_host_um();
724  }
725 
726  typename dual_view_type::t_dev_um getDeviceLocalView() const {
727  //return subview(vec_->template getLocalView<typename dual_view_type::t_dev_um::execution_space> (),
728  // Kokkos::ALL(), Kokkos::ALL());
729  return typename dual_view_type::t_dev_um();
730  }
731 
732 #endif
733 
735 
736  protected:
739  virtual void
741  { }
742  }; // TpetraMultiVector class (specialization GO=int, NO=EpetraNode)
743 #endif
744 
745 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \
746  (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))))
747 
748  // specialization for TpetraMultiVector on EpetraNode and GO=long long
749  template <class Scalar>
750  class TpetraMultiVector<Scalar,int,long long,EpetraNode>
751  : public virtual MultiVector< Scalar, int, long long, EpetraNode >
752  {
753  typedef int LocalOrdinal;
754  typedef long long GlobalOrdinal;
755  typedef EpetraNode Node;
756 
757  // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
759 
760  public:
761 
763 
764 
766  TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true) {
768  }
769 
773  }
774 
778  }
779 
783  }
784 
785 
787  virtual ~TpetraMultiVector() { }
788 
790 
792 
793 
795  void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { }
796 
798  void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { }
799 
801  void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { }
802 
804  void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { }
805 
807  void putScalar(const Scalar &value) { }
808 
810  void reduce() { }
811 
813 
815 
816 
819 
822 
825 
828 
830  void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const { }
831 
833  void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const { }
834 
837 
840 
843 
846 
848 
850 
851 
854 
857 
860 
862  void scale(const Scalar &alpha) { }
863 
866 
868  void scale(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { }
869 
871  void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta) { }
872 
874  void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &gamma) { }
875 
878 
881 
884 
886  void meanValue(const Teuchos::ArrayView< Scalar > &means) const { }
887 
890 
892 
894 
895 
897  size_t getNumVectors() const { return 0; }
898 
900  size_t getLocalLength() const { return 0; }
901 
903  global_size_t getGlobalLength() const { return 0; }
904 
905  // \! Checks to see if the local length, number of vectors and size of Scalar type match
906  bool isSameSize(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & vec) const { return false; }
907 
909 
911 
912 
914  std::string description() const { return std::string(""); }
915 
918 
920 
923 
925  void randomize(bool bUseXpetraImplementation = false) { }
926 
927  //{@
928  // Implements DistObject interface
929 
931 
933 
935 
937 
939 
941 
942 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
943  template<class Node2>
944  RCP<MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node2> > clone(const RCP<Node2> &node2) const { return Teuchos::null; }
945 #endif
946 
947 
949 
950 
952  TpetraMultiVector(const Teuchos::RCP<Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node> > &vec) {
954  }
955 
958 
960  void setSeed(unsigned int seed) { }
961 
962 
963 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
965  //typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::host_execution_space host_execution_space;
966  //typedef typename Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>::dev_execution_space dev_execution_space;
967 
978  template<class TargetDeviceType>
979  typename Kokkos::Impl::if_c<
980  Kokkos::Impl::is_same<
981  typename dual_view_type::t_dev_um::execution_space::memory_space,
982  typename TargetDeviceType::memory_space>::value,
983  typename dual_view_type::t_dev_um,
984  typename dual_view_type::t_host_um>::type
985  getLocalView () const {
986  typename Kokkos::Impl::if_c<
987  Kokkos::Impl::is_same<
988  typename dual_view_type::t_dev_um::execution_space::memory_space,
989  typename TargetDeviceType::memory_space>::value,
990  typename dual_view_type::t_dev_um,
991  typename dual_view_type::t_host_um>::type dummy;
992  return dummy;
993  }
994 
995  typename dual_view_type::t_host_um getHostLocalView () const {
996  //return subview(vec_->template getLocalView<typename dual_view_type::host_mirror_space> (),
997  // Kokkos::ALL(), Kokkos::ALL());
998  return typename dual_view_type::t_host_um();
999  }
1000 
1001  typename dual_view_type::t_dev_um getDeviceLocalView() const {
1002  //return subview(vec_->template getLocalView<typename dual_view_type::t_dev_um::execution_space> (),
1003  // Kokkos::ALL(), Kokkos::ALL());
1004  return typename dual_view_type::t_dev_um();
1005  }
1006 
1007 #endif
1008 
1010 
1011  protected:
1014  virtual void
1016  { }
1017  }; // TpetraMultiVector class (specialization GO=int, NO=EpetraNode)
1018 
1019 #endif // TpetraMultiVector class (specialization GO=long long, NO=EpetraNode)
1020 
1021 #endif // HAVE_XPETRA_EPETRA
1022 
1023 } // Xpetra namespace
1024 
1025 // Following header file inculsion is needed for the dynamic_cast to TpetraVector in elementWiseMultiply (because we cannot dynamic_cast if target is not a complete type)
1026 // It is included here to avoid circular dependency between Vector and MultiVector
1027 // TODO: there is certainly a more elegant solution...
1028 #include "Xpetra_TpetraVector.hpp"
1029 
1030 namespace Xpetra {
1031  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
1033  XPETRA_MONITOR("TpetraMultiVector::elementWiseMultiply");
1034 
1035  // XPETRA_DYNAMIC_CAST won't take TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
1036  // as an argument, hence the following typedef.
1038  XPETRA_DYNAMIC_CAST(const tpv, A, tA, "Xpetra::TpetraMultiVectorMatrix->multiply() only accept Xpetra::TpetraMultiVector as input arguments.");
1039  XPETRA_DYNAMIC_CAST(const TpetraMultiVector, B, tB, "Xpetra::TpetraMultiVectorMatrix->multiply() only accept Xpetra::TpetraMultiVector as input arguments.");
1040  vec_->elementWiseMultiply(scalarAB, *tA.getTpetra_Vector(), *tB.getTpetra_MultiVector(), scalarThis);
1041  }
1042 
1043 } // Xpetra namespace
1044 
1045 #define XPETRA_TPETRAMULTIVECTOR_SHORT
1046 #endif // XPETRA_TPETRAMULTIVECTOR_HPP
void replaceMap(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
size_t getLocalLength() const
Local number of rows on the calling process.
RCP< Map< LocalOrdinal, GlobalOrdinal, Node2 > > clone(const Map< LocalOrdinal, GlobalOrdinal, Node1 > &map, const RCP< Node2 > &node2)
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
Teuchos::ArrayRCP< const Scalar > get1dView() const
Const persisting (1-D) view of this multivector&#39;s local values.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_MultiVector() const
Get the underlying Tpetra multivector.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const
Return const persisting pointers to values.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const
Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
void scale(const Scalar &alpha)
Scale the current values of a multi-vector, this = alpha*this.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic constuctor.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void setSeed(unsigned int seed)
Set seed for Random function.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst()
Return non-const persisting pointers to values.
void reduce()
Sum values of a locally replicated multivector across all processes.
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const
Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
std::string description() const
A simple one-line description of this object.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
Create multivector by copying array of views of local data.
void reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
virtual ~TpetraMultiVector()
Destructor (virtual for memory safety of derived classes).
void setSeed(unsigned int seed)
Set seed for Random function.
void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
Element-wise multiply of a Vector A with a TpetraMultiVector B.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic constuctor.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
size_t getNumVectors() const
Number of columns in the multivector.
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
Element-wise multiply of a Vector A with a TpetraMultiVector B.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void Xpetra_randomize()
Set multi-vector values to random numbers. XPetra implementation.
std::string description() const
A simple one-line description of this object.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
size_t getNumVectors() const
Number of columns in the multivector.
void setSeed(unsigned int seed)
Set seed for Random function.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void scale(const Scalar &alpha)
Scale the current values of a multi-vector, this = alpha*this.
TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVectorClass
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void scale(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Replace multi-vector values with scaled values of A, this = alpha*A.
void replaceMap(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const
Return const persisting pointers to values.
void reduce()
Sum values of a locally replicated multivector across all processes.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
#define XPETRA_TPETRA_ETI_EXCEPTION(cl, obj, go, node)
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const
Fill the given array with a copy of this multivector&#39;s local values.
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import data into this object using an Import object (&quot;forward mode&quot;).
void meanValue(const Teuchos::ArrayView< Scalar > &means) const
Compute mean (average) value of each vector in multi-vector. The outcome of this routine is undefined...
TpetraMultiVector(const Teuchos::RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &vec)
TpetraMultiVector constructor to wrap a Tpetra::MultiVector object.
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const
Fill the given array with a copy of this multivector&#39;s local values.
void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
Element-wise multiply of a Vector A with a TpetraMultiVector B.
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
void multiply(Teuchos::ETransp transA, Teuchos::ETransp transB, const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &beta)
Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).
Teuchos::ArrayRCP< Scalar > get1dViewNonConst()
Nonconst persisting (1-D) view of this multivector&#39;s local values.
TpetraMultiVector(const Teuchos::RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &vec)
TpetraMultiVector constructor to wrap a Tpetra::MultiVector object.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst()
Return non-const persisting pointers to values.
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &gamma)
Update multi-vector with scaled values of A and B, this = gamma*this + alpha*A + beta*B.
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
size_t getLocalLength() const
Local number of rows on the calling process.
void multiply(Teuchos::ETransp transA, Teuchos::ETransp transB, const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &beta)
Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const
Fill the given array with a copy of this multivector&#39;s local values.
void multiply(Teuchos::ETransp transA, Teuchos::ETransp transB, const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &beta)
Matrix-matrix multiplication: this = beta*this + alpha*op(A)*op(B).
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
std::string description() const
A simple one-line description of this object.
void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const
Compute dot product of each corresponding pair of vectors, dots[i] = this[i].dot(A[i]).
size_t getNumVectors() const
Number of columns in the multivector.
Teuchos::ArrayRCP< Scalar > get1dViewNonConst()
Nonconst persisting (1-D) view of this multivector&#39;s local values.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
void reduce()
Sum values of a locally replicated multivector across all processes.
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void get2dCopy(Teuchos::ArrayView< const Teuchos::ArrayView< Scalar > > ArrayOfPtrs) const
Fill the given array with a copy of this multivector&#39;s local values.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
void scale(const Scalar &alpha)
Scale the current values of a multi-vector, this = alpha*this.
void reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
virtual ~TpetraMultiVector()
Destructor (virtual for memory safety of derived classes).
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
static void seedrandom(unsigned int s)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_t global_size_t
Global size_t object.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic constuctor.
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, const Scalar &gamma)
Update multi-vector with scaled values of A and B, this = gamma*this + alpha*A + beta*B.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
static const EVerbosityLevel verbLevel_default
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const
Fill the given array with a copy of this multivector&#39;s local values.
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const
Fill the given array with a copy of this multivector&#39;s local values.
void reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise reciprocal values of input Multi-vector in target, this(i,j) = 1/A(i,j).
global_size_t getGlobalLength() const
Global number of rows in the multivector.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
virtual ~TpetraMultiVector()
Destructor (virtual for memory safety of derived classes).
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
Create multivector by copying array of views of local data.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
Teuchos::RCP< const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
size_t getLocalLength() const
Local number of rows on the calling process.
void scale(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Replace multi-vector values with scaled values of A, this = alpha*A.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Scalar > &A, size_t LDA, size_t NumVectors)
Create multivector by copying two-dimensional array of local data.
TpetraMultiVector(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
Copy constructor (performs a deep copy).
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
TpetraMultiVector(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Scalar > &A, size_t LDA, size_t NumVectors)
Create multivector by copying two-dimensional array of local data.
void meanValue(const Teuchos::ArrayView< Scalar > &means) const
Compute mean (average) value of each vector in multi-vector. The outcome of this routine is undefined...
Teuchos::ArrayRCP< const Scalar > get1dView() const
Const persisting (1-D) view of this multivector&#39;s local values.
Teuchos::ArrayRCP< const Scalar > get1dView() const
Const persisting (1-D) view of this multivector&#39;s local values.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export data into this object using an Import object (&quot;reverse mode&quot;).
CombineMode
Xpetra::Combine Mode enumerable type.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_MultiVector() const
Get the underlying Tpetra multivector.
TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVectorClass
#define XPETRA_MONITOR(funcName)
void meanValue(const Teuchos::ArrayView< Scalar > &means) const
Compute mean (average) value of each vector in multi-vector. The outcome of this routine is undefined...
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
Teuchos::ArrayRCP< Scalar > get1dViewNonConst()
Nonconst persisting (1-D) view of this multivector&#39;s local values.
global_size_t getGlobalLength() const
Global number of rows in the multivector.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst()
Return non-const persisting pointers to values.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
TpetraMultiVector(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
Copy constructor (performs a deep copy).
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const
Return const persisting pointers to values.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
global_size_t getGlobalLength() const
Global number of rows in the multivector.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
void replaceMap(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
bool is_null() const
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_MultiVector() const
Get the underlying Tpetra multivector.