Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_EpetraMultiVector.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_EPETRAMULTIVECTOR_HPP
47 #define XPETRA_EPETRAMULTIVECTOR_HPP
48 
49 /* this file is automatically generated - do not edit (see script/epetra.py) */
50 
51 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
52 #include <Kokkos_Core.hpp>
53 #include <Kokkos_DualView.hpp>
54 #endif
55 
56 
58 
59 #include "Xpetra_MultiVector.hpp"
60 #include "Xpetra_Vector.hpp"
61 
62 #include "Xpetra_EpetraMap.hpp"
63 #include "Xpetra_EpetraExport.hpp"
64 #include "Xpetra_Utils.hpp"
65 #include "Xpetra_EpetraUtils.hpp"
66 #include "Xpetra_EpetraImport.hpp"
67 #include "Xpetra_Exceptions.hpp"
68 #include "Epetra_SerialComm.h"
69 
70 #include <Epetra_MultiVector.h>
71 #include <Epetra_Vector.h>
72 
73 namespace Xpetra {
74 
75  // TODO: move that elsewhere
76  template<class GlobalOrdinal, class Node>
77  const Epetra_MultiVector & toEpetra(const MultiVector<double,int,GlobalOrdinal,Node> &);
78  template<class GlobalOrdinal, class Node>
79  Epetra_MultiVector & toEpetra(MultiVector<double, int,GlobalOrdinal,Node> &);
80  template<class GlobalOrdinal, class Node>
81  RCP<MultiVector<double, int, GlobalOrdinal, Node> > toXpetra(RCP<Epetra_MultiVector> vec);
82 
83  // we need this forward declaration
84 #ifndef DOXYGEN_SHOULD_SKIP_THIS
85  template<class GlobalOrdinal, class Node> class EpetraVectorT;
86 #endif
87 
88  template<class EpetraGlobalOrdinal, class Node>
90  : public virtual MultiVector<double, int, EpetraGlobalOrdinal, Node>
91  {
92  typedef double Scalar;
93  typedef int LocalOrdinal;
94  typedef EpetraGlobalOrdinal GlobalOrdinal;
95 
96  public:
97 
99 
100 
102  EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true) {
103  TEUCHOS_TEST_FOR_EXCEPTION(true, Xpetra::Exceptions::RuntimeError,
104  "Xpetra::EpetraMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
105  }
106 
109  TEUCHOS_TEST_FOR_EXCEPTION(true, Xpetra::Exceptions::RuntimeError,
110  "Xpetra::EpetraMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
111  }
112 
114  EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors) {
115  TEUCHOS_TEST_FOR_EXCEPTION(true, Xpetra::Exceptions::RuntimeError,
116  "Xpetra::EpetraMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
117  }
118 
120  virtual ~EpetraMultiVectorT() {}
121 
123 
125 
126 
128  void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { }
129 
131  void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { }
132 
134  void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { }
135 
137  void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { }
138 
140  void putScalar(const Scalar &value) { }
141 
143 
145 
146 
148  Teuchos::RCP< const Vector< double, int, GlobalOrdinal, Node > > getVector(size_t j) const {
149  return Teuchos::null;
150  }
151 
153  Teuchos::RCP< Vector< double, int, GlobalOrdinal, Node > > getVectorNonConst(size_t j) {
154  return Teuchos::null;
155  }
156 
158  Teuchos::ArrayRCP< const Scalar > getData(size_t j) const {
159  return ArrayRCP<const Scalar>();
160  }
161 
163  Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j) {
164  return ArrayRCP<Scalar>();
165  }
166 
168 
170 
171 
173  void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const { }
174 
177 
180 
182  void scale(const Scalar &alpha) { }
183 
185  void scale (Teuchos::ArrayView< const Scalar > alpha) { }
186 
188  void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta) { }
189 
192 
194  void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { }
195 
197  void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { }
198 
200  void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { }
201 
203  void meanValue(const Teuchos::ArrayView< Scalar > &means) const { }
204 
206  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) { }
207 
210 
212 
214 
215 
217  size_t getNumVectors() const { return 0; }
218 
220  size_t getLocalLength() const { return 0; }
221 
223  global_size_t getGlobalLength() const { return 0; }
224 
225  // \brief Checks to see if the local length, number of vectors and size of Scalar type match
226  bool isSameSize(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & vec) const { return false; }
227 
229 
231 
232 
234  std::string description() const { return std::string(""); }
235 
237  void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const { }
238 
240 
242  void randomize(bool bUseXpetraImplementation = false) { }
243 
245  //{@
246 
248  Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const { return Teuchos::null; }
249 
252 
255 
258 
261 
263  void replaceMap(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map) { }
264 
266 
268 
269 
271  EpetraMultiVectorT(const RCP<Epetra_MultiVector> &vec) { //TODO removed const
272  TEUCHOS_TEST_FOR_EXCEPTION(true, Xpetra::Exceptions::RuntimeError,
273  "Xpetra::EpetraMultiVector only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)");
274  }
275 
277  RCP<Epetra_MultiVector> getEpetra_MultiVector() const { return Teuchos::null; }
278 
280  void setSeed(unsigned int seed) { }
281 
282 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
283 
285 
296  template<class TargetDeviceType>
297  typename Kokkos::Impl::if_c<
298  std::is_same<
299  typename dual_view_type::t_dev_um::execution_space::memory_space,
300  typename TargetDeviceType::memory_space>::value,
301  typename dual_view_type::t_dev_um,
302  typename dual_view_type::t_host_um>::type
303  getLocalView () const {
304  typename Kokkos::Impl::if_c<
305  std::is_same<
306  typename dual_view_type::t_dev_um::execution_space::memory_space,
307  typename TargetDeviceType::memory_space>::value,
308  typename dual_view_type::t_dev_um,
309  typename dual_view_type::t_host_um>::type dummy;
310  return dummy;
311  }
312 
313  typename dual_view_type::t_host_um getHostLocalView () const {
314  return typename dual_view_type::t_host_um();
315  }
316 
317  typename dual_view_type::t_dev_um getDeviceLocalView() const {
318  throw std::runtime_error("Epetra does not support device views!");
319 #ifndef __NVCC__ //prevent nvcc warning
320  typename dual_view_type::t_dev_um ret;
321 #endif
322  TEUCHOS_UNREACHABLE_RETURN(ret);
323  }
324 
325 #endif
326 
328 
329  protected:
332  virtual void
334 
335  }; // EpetraMultiVectorT class
336 
337  // specialization on GO=int and Node=EpetraNode
338 #ifndef XPETRA_EPETRA_NO_32BIT_GLOBAL_INDICES
339  template<>
341  : public virtual MultiVector<double, int, int, EpetraNode>
342  {
343  typedef double Scalar;
344  typedef int LocalOrdinal;
345  typedef int GlobalOrdinal;
346  typedef EpetraNode Node;
347 
348  public:
349 
351 
352 
354  EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
355  : vec_(Teuchos::rcp(new Epetra_MultiVector(toEpetra<GlobalOrdinal,Node>(map), Teuchos::as<int>(NumVectors), zeroOut))) { }
356 
359  : vec_(Teuchos::rcp(new Epetra_MultiVector(toEpetra<GlobalOrdinal,Node>(source)))) { }
360 
362  EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors) {
363  //TODO: input argument 'NumVectors' is not necessary in both Xpetra and Tpetra interface. Should it be removed?
364 
365  const std::string tfecfFuncName("MultiVector(ArrayOfPtrs)");
366  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(NumVectors < 1 || NumVectors != Teuchos::as<size_t>(ArrayOfPtrs.size()), std::runtime_error,
367  ": ArrayOfPtrs.size() must be strictly positive and as large as ArrayOfPtrs.");
368 
369  #ifdef HAVE_XPETRA_DEBUG
370  // This cannot be tested by Epetra itself
371  {
372  size_t localLength = map->getNodeNumElements();
373  for(int j=0; j<ArrayOfPtrs.size(); j++) {
374  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(Teuchos::as<size_t>(ArrayOfPtrs[j].size()) != localLength, std::runtime_error,
375  ": ArrayOfPtrs[" << j << "].size() (== " << ArrayOfPtrs[j].size() <<
376  ") is not equal to getLocalLength() (== " << localLength);
377 
378  }
379  }
380  #endif
381 
382  // Convert Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > to double**
383  Array<const double*> arrayOfRawPtrs(ArrayOfPtrs.size());
384  for(int i=0; i<ArrayOfPtrs.size(); i++) {
385  arrayOfRawPtrs[i] = ArrayOfPtrs[i].getRawPtr();
386  }
387  double** rawArrayOfRawPtrs = const_cast<double**>(arrayOfRawPtrs.getRawPtr()); // This const_cast should be fine, because Epetra_DataAccess=Copy.
388 
389  vec_ = Teuchos::rcp(new Epetra_MultiVector(Copy, toEpetra<GlobalOrdinal,Node>(map), rawArrayOfRawPtrs, static_cast<int>(NumVectors)));
390  }
391 
393  virtual ~EpetraMultiVectorT() {}
394 
396 
398 
399 
401  void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("EpetraMultiVectorT::replaceGlobalValue"); vec_->ReplaceGlobalValue(globalRow, Teuchos::as<int>(vectorIndex), value); }
402 
404  void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("EpetraMultiVectorT::sumIntoGlobalValue"); vec_->SumIntoGlobalValue(globalRow, Teuchos::as<int>(vectorIndex), value); }
405 
407  void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("EpetraMultiVectorT::replaceLocalValue"); vec_->ReplaceMyValue(myRow, Teuchos::as<int>(vectorIndex), value); }
408 
410  void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("EpetraMultiVectorT::sumIntoLocalValue"); vec_->SumIntoMyValue(myRow, Teuchos::as<int>(vectorIndex), value); }
411 
413  void putScalar(const Scalar &value) { XPETRA_MONITOR("EpetraMultiVectorT::putScalar"); vec_->PutScalar(value); }
414 
416 
418 
419 
421  Teuchos::RCP< const Vector< double, int, int, EpetraNode > > getVector(size_t j) const;
422 
424  Teuchos::RCP< Vector< double, int, int, EpetraNode > > getVectorNonConst(size_t j);
425 
427  Teuchos::ArrayRCP< const Scalar > getData(size_t j) const {
428  XPETRA_MONITOR("EpetraMultiVectorT::getData");
429 
430  double ** arrayOfPointers;
431 
432  vec_->ExtractView(&arrayOfPointers);
433 
434  double * data = arrayOfPointers[j];
435  int localLength = vec_->MyLength();
436 
437  return ArrayRCP<double>(data, 0, localLength, false); // no ownership
438  }
439 
441  Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j) {
442  XPETRA_MONITOR("EpetraMultiVectorT::getDataNonConst");
443 
444  double ** arrayOfPointers;
445 
446  vec_->ExtractView(&arrayOfPointers);
447 
448  double * data = arrayOfPointers[j];
449  int localLength = vec_->MyLength();
450 
451  return ArrayRCP<double>(data, 0, localLength, false); // no ownership
452  }
453 
455 
457 
458 
460  void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const {
461  XPETRA_MONITOR("EpetraMultiVectorT::dot");
462 
463  XPETRA_DYNAMIC_CAST(const EpetraMultiVectorT<GlobalOrdinal XPETRA_COMMA Node>, A, eA, "This Xpetra::EpetraMultiVectorT method only accept Xpetra::EpetraMultiVectorT as input arguments.");
464  vec_->Dot(*eA.getEpetra_MultiVector(), dots.getRawPtr());
465  }
466 
468  void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("EpetraMultiVectorT::abs"); vec_->Abs(toEpetra<GlobalOrdinal,Node>(A)); }
469 
471  void reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("EpetraMultiVectorT::reciprocal"); vec_->Reciprocal(toEpetra<GlobalOrdinal,Node>(A)); }
472 
474  void scale(const Scalar &alpha) { XPETRA_MONITOR("EpetraMultiVectorT::scale"); vec_->Scale(alpha); }
475 
477  void scale (Teuchos::ArrayView< const Scalar > alpha) {
478  XPETRA_MONITOR("EpetraMultiVectorT::scale");
479  // Epetra, unlike Tpetra, doesn't implement this version of
480  // scale(). Deal with this by scaling one column at a time.
481  const size_t numVecs = this->getNumVectors ();
482  for (size_t j = 0; j < numVecs; ++j) {
483  Epetra_Vector *v = (*vec_)(j);
484  v->Scale (alpha[j]);
485  }
486  }
487 
489  void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta) { XPETRA_MONITOR("EpetraMultiVectorT::update"); vec_->Update(alpha, toEpetra<GlobalOrdinal,Node>(A), beta); }
490 
492  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) { XPETRA_MONITOR("EpetraMultiVectorT::update"); vec_->Update(alpha, toEpetra<GlobalOrdinal,Node>(A), beta, toEpetra<GlobalOrdinal,Node>(B), gamma); }
493 
495  void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("EpetraMultiVectorT::norm1"); vec_->Norm1(norms.getRawPtr()); }
496 
498  void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("EpetraMultiVectorT::norm2"); vec_->Norm2(norms.getRawPtr()); }
499 
501  void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("EpetraMultiVectorT::normInf"); vec_->NormInf(norms.getRawPtr()); }
502 
504  void meanValue(const Teuchos::ArrayView< Scalar > &means) const { XPETRA_MONITOR("EpetraMultiVectorT::meanValue"); vec_->MeanValue(means.getRawPtr()); } //TODO: modify ArrayView size ??
505 
507  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) { XPETRA_MONITOR("EpetraMultiVectorT::multiply"); vec_->Multiply(toEpetra(transA), toEpetra(transB), alpha, toEpetra(A), toEpetra(B), beta); }
508 
510  void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis) { XPETRA_MONITOR("EpetraMultiVectorT::elementWiseMultiply"); vec_->Multiply(scalarAB, toEpetra<GlobalOrdinal,Node>(A), toEpetra<GlobalOrdinal,Node>(B), scalarThis); }
511 
513 
515 
516 
518  size_t getNumVectors() const { XPETRA_MONITOR("EpetraMultiVectorT::getNumVectors"); return vec_->NumVectors(); }
519 
521  size_t getLocalLength() const { XPETRA_MONITOR("EpetraMultiVectorT::getLocalLength"); return vec_->MyLength(); }
522 
524  global_size_t getGlobalLength() const { XPETRA_MONITOR("EpetraMultiVectorT::getGlobalLength"); return vec_->GlobalLength64(); }
525 
528  XPETRA_MONITOR("EpetraMultiVectorT::isSameSize");
529  auto vv = toEpetra<GlobalOrdinal,Node>(vec);
530  return ( (vec_->MyLength() == vv.MyLength()) && (vec_->NumVectors() == vv.NumVectors()));
531  }
532 
534 
536 
537 
539  std::string description() const {
540  XPETRA_MONITOR("EpetraMultiVectorT::description");
541  TEUCHOS_TEST_FOR_EXCEPTION(1, Xpetra::Exceptions::NotImplemented, "TODO");
542  TEUCHOS_UNREACHABLE_RETURN("TODO");
543  }
544 
546  void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel /* verbLevel */=Teuchos::Describable::verbLevel_default) const {
547  XPETRA_MONITOR("EpetraMultiVectorT::describe");
548  vec_->Print(out);
549  }
550 
552 
554  void randomize(bool bUseXpetraImplementation = false) {
555  XPETRA_MONITOR("EpetraMultiVectorT::randomize");
556 
557  if (bUseXpetraImplementation)
559  else
560  vec_->Random();
561  }
562 
564  //{@
565 
567  Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const { XPETRA_MONITOR("EpetraMultiVectorT::getMap"); return toXpetra<GlobalOrdinal,Node>(vec_->Map()); }
568 
571  XPETRA_MONITOR("EpetraMultiVectorT::doImport");
572 
573  XPETRA_DYNAMIC_CAST(const EpetraMultiVectorT<GlobalOrdinal XPETRA_COMMA Node>, source, tSource, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraMultiVectorT as input arguments.");
574  XPETRA_DYNAMIC_CAST(const EpetraImportT<GlobalOrdinal XPETRA_COMMA Node>, importer, tImporter, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraImportT as input arguments.");
575 
576  RCP<Epetra_MultiVector> v = tSource.getEpetra_MultiVector();
577  int err = this->getEpetra_MultiVector()->Import(*v, *tImporter.getEpetra_Import(), toEpetra(CM));
578  TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra is " << err);
579  }
580 
583  XPETRA_MONITOR("EpetraMultiVectorT::doExport");
584 
585  XPETRA_DYNAMIC_CAST(const EpetraMultiVectorT<GlobalOrdinal XPETRA_COMMA Node>, dest, tDest, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraMultiVectorT as input arguments.");
586  XPETRA_DYNAMIC_CAST(const EpetraImportT<GlobalOrdinal XPETRA_COMMA Node>, importer, tImporter, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraImportT as input arguments.");
587 
588  RCP<Epetra_MultiVector> v = tDest.getEpetra_MultiVector();
589  int err = this->getEpetra_MultiVector()->Export(*v, *tImporter.getEpetra_Import(), toEpetra(CM));
590  TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
591  }
592 
595  XPETRA_MONITOR("EpetraMultiVectorT::doImport");
596 
597  XPETRA_DYNAMIC_CAST(const EpetraMultiVectorT<GlobalOrdinal XPETRA_COMMA Node>, source, tSource, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraMultiVectorT as input arguments.");
598  XPETRA_DYNAMIC_CAST(const EpetraExportT<GlobalOrdinal XPETRA_COMMA Node>, exporter, tExporter, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraImportT as input arguments.");
599 
600  RCP<Epetra_MultiVector> v = tSource.getEpetra_MultiVector();
601  int err = this->getEpetra_MultiVector()->Import(*v, *tExporter.getEpetra_Export(), toEpetra(CM));
602  TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
603  }
604 
607  XPETRA_MONITOR("EpetraMultiVectorT::doExport");
608 
609  XPETRA_DYNAMIC_CAST(const EpetraMultiVectorT<GlobalOrdinal XPETRA_COMMA Node>, dest, tDest, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraMultiVectorT as input arguments.");
610  XPETRA_DYNAMIC_CAST(const EpetraExportT<GlobalOrdinal XPETRA_COMMA Node>, exporter, tExporter, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraImportT as input arguments.");
611 
612  RCP<Epetra_MultiVector> v = tDest.getEpetra_MultiVector();
613  int err = this->getEpetra_MultiVector()->Export(*v, *tExporter.getEpetra_Export(), toEpetra(CM));
614  TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
615  }
616 
618  void replaceMap(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map) {
619  XPETRA_MONITOR("EpetraMultiVectorT::replaceMap");
620  int err = 0;
621  if (!map.is_null()) {
622  err = this->getEpetra_MultiVector()->ReplaceMap(toEpetra<GlobalOrdinal,Node>(map));
623 
624  } else {
625  // Replace map with a dummy map to avoid potential hangs later
626  Epetra_SerialComm SComm;
627  Epetra_Map NewMap((GlobalOrdinal) vec_->MyLength(), (GlobalOrdinal) vec_->Map().IndexBase64(), SComm);
628  err = this->getEpetra_MultiVector()->ReplaceMap(NewMap);
629  }
630  TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
631  }
632 
634 
636 
637 
639  EpetraMultiVectorT(const RCP<Epetra_MultiVector> &vec) : vec_(vec) { } //TODO removed const
640 
642  RCP<Epetra_MultiVector> getEpetra_MultiVector() const { return vec_; }
643 
645  void setSeed(unsigned int seed) {
646  XPETRA_MONITOR("EpetraMultiVectorT::seedrandom");
647 
648  Teuchos::ScalarTraits< Scalar >::seedrandom(seed);
649  vec_->SetSeed(seed);
650  }
651 
652 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
653 
655 
666  template<class TargetDeviceType>
667  typename Kokkos::Impl::if_c<
668  std::is_same<
669  typename dual_view_type::t_dev_um::execution_space::memory_space,
670  typename TargetDeviceType::memory_space>::value,
671  typename dual_view_type::t_dev_um,
672  typename dual_view_type::t_host_um>::type
673  getLocalView () const {
674  return this->MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::template getLocalView<TargetDeviceType>();
675  }
676 
677  typename dual_view_type::t_host_um getHostLocalView () const {
678  typedef Kokkos::View< typename dual_view_type::t_host::data_type ,
679  Kokkos::LayoutLeft,
680  typename dual_view_type::t_host::device_type ,
681  Kokkos::MemoryUnmanaged> epetra_view_type;
682 
683  // access Epetra multivector data
684  double* data = NULL;
685  int myLDA;
686  vec_->ExtractView(&data, &myLDA);
687  int localLength = vec_->MyLength();
688  int numVectors = getNumVectors();
689 
690  // create view
691  epetra_view_type test = epetra_view_type(data, localLength, numVectors);
692  typename dual_view_type::t_host_um ret = subview(test, Kokkos::ALL(), Kokkos::ALL());
693 
694  return ret;
695  }
696 
697  typename dual_view_type::t_dev_um getDeviceLocalView() const {
698  throw std::runtime_error("Epetra does not support device views!");
699 #ifndef __NVCC__ //prevent nvcc warning
700  typename dual_view_type::t_dev_um ret;
701 #endif
702  TEUCHOS_UNREACHABLE_RETURN(ret);
703  }
704 
705 #endif
706 
708 
709  protected:
712  virtual void
714  typedef EpetraMultiVectorT this_type;
715  const this_type* rhsPtr = dynamic_cast<const this_type*> (&rhs);
716  TEUCHOS_TEST_FOR_EXCEPTION(
717  rhsPtr == NULL, std::invalid_argument, "Xpetra::MultiVector::operator=: "
718  "The left-hand side (LHS) of the assignment has a different type than "
719  "the right-hand side (RHS). The LHS has type Xpetra::EpetraMultiVectorT "
720  "(which means it wraps an Epetra_MultiVector), but the RHS has some "
721  "other type. This probably means that the RHS wraps a Tpetra::Multi"
722  "Vector. Xpetra::MultiVector does not currently implement assignment "
723  "from a Tpetra object to an Epetra object, though this could be added "
724  "with sufficient interest.");
725 
726  RCP<const Epetra_MultiVector> rhsImpl = rhsPtr->getEpetra_MultiVector ();
727  RCP<Epetra_MultiVector> lhsImpl = this->getEpetra_MultiVector ();
728 
729  TEUCHOS_TEST_FOR_EXCEPTION(
730  rhsImpl.is_null (), std::logic_error, "Xpetra::MultiVector::operator= "
731  "(in Xpetra::EpetraMultiVectorT::assign): *this (the right-hand side of "
732  "the assignment) has a null RCP<Epetra_MultiVector> inside. Please "
733  "report this bug to the Xpetra developers.");
734  TEUCHOS_TEST_FOR_EXCEPTION(
735  lhsImpl.is_null (), std::logic_error, "Xpetra::MultiVector::operator= "
736  "(in Xpetra::EpetraMultiVectorT::assign): The left-hand side of the "
737  "assignment has a null RCP<Epetra_MultiVector> inside. Please report "
738  "this bug to the Xpetra developers.");
739 
740  // Epetra_MultiVector's assignment operator does a deep copy.
741  *lhsImpl = *rhsImpl;
742  }
743 
744  private:
746  RCP< Epetra_MultiVector > vec_;
747 
748  }; // EpetraMultiVectorT class (specialization on GO=int, NO=EpetraNode
749 #endif
750 
751  // specialization on GO=long long and EpetraNode
752 #ifndef XPETRA_EPETRA_NO_64BIT_GLOBAL_INDICES
753  template<>
754  class EpetraMultiVectorT<long long, EpetraNode>
755  : public virtual MultiVector<double, int, long long, EpetraNode>
756  {
757  typedef double Scalar;
758  typedef int LocalOrdinal;
759  typedef long long GlobalOrdinal;
760  typedef EpetraNode Node;
761 
762  public:
763 
765 
766 
768  EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
769  : vec_(Teuchos::rcp(new Epetra_MultiVector(toEpetra<GlobalOrdinal,Node>(map), Teuchos::as<int>(NumVectors), zeroOut))) { }
770 
773  : vec_(Teuchos::rcp(new Epetra_MultiVector(toEpetra<GlobalOrdinal,Node>(source)))) { }
774 
776  EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors) {
777  //TODO: input argument 'NumVectors' is not necessary in both Xpetra and Tpetra interface. Should it be removed?
778 
779  const std::string tfecfFuncName("MultiVector(ArrayOfPtrs)");
780  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(NumVectors < 1 || NumVectors != Teuchos::as<size_t>(ArrayOfPtrs.size()), std::runtime_error,
781  ": ArrayOfPtrs.size() must be strictly positive and as large as ArrayOfPtrs.");
782 
783  #ifdef HAVE_XPETRA_DEBUG
784  // This cannot be tested by Epetra itself
785  {
786  size_t localLength = map->getNodeNumElements();
787  for(int j=0; j<ArrayOfPtrs.size(); j++) {
788  TEUCHOS_TEST_FOR_EXCEPTION_CLASS_FUNC(Teuchos::as<size_t>(ArrayOfPtrs[j].size()) != localLength, std::runtime_error,
789  ": ArrayOfPtrs[" << j << "].size() (== " << ArrayOfPtrs[j].size() <<
790  ") is not equal to getLocalLength() (== " << localLength);
791 
792  }
793  }
794  #endif
795 
796  // Convert Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > to double**
797  Array<const double*> arrayOfRawPtrs(ArrayOfPtrs.size());
798  for(int i=0; i<ArrayOfPtrs.size(); i++) {
799  arrayOfRawPtrs[i] = ArrayOfPtrs[i].getRawPtr();
800  }
801  double** rawArrayOfRawPtrs = const_cast<double**>(arrayOfRawPtrs.getRawPtr()); // This const_cast should be fine, because Epetra_DataAccess=Copy.
802 
803  vec_ = Teuchos::rcp(new Epetra_MultiVector(Copy, toEpetra<GlobalOrdinal,Node>(map), rawArrayOfRawPtrs, NumVectors));
804  }
805 
807  virtual ~EpetraMultiVectorT() {}
808 
810 
812 
813 
815  void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("EpetraMultiVectorT::replaceGlobalValue"); vec_->ReplaceGlobalValue(globalRow, Teuchos::as<int>(vectorIndex), value); }
816 
818  void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("EpetraMultiVectorT::sumIntoGlobalValue"); vec_->SumIntoGlobalValue(globalRow, Teuchos::as<int>(vectorIndex), value); }
819 
821  void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("EpetraMultiVectorT::replaceLocalValue"); vec_->ReplaceMyValue(myRow, Teuchos::as<int>(vectorIndex), value); }
822 
824  void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) { XPETRA_MONITOR("EpetraMultiVectorT::sumIntoLocalValue"); vec_->SumIntoMyValue(myRow, Teuchos::as<int>(vectorIndex), value); }
825 
827  void putScalar(const Scalar &value) { XPETRA_MONITOR("EpetraMultiVectorT::putScalar"); vec_->PutScalar(value); }
828 
830 
832 
833 
835  Teuchos::RCP< const Vector< double, int, long long, EpetraNode > > getVector(size_t j) const;
836 
838  Teuchos::RCP< Vector< double, int, long long, EpetraNode > > getVectorNonConst(size_t j);
839 
841  Teuchos::ArrayRCP< const Scalar > getData(size_t j) const {
842  XPETRA_MONITOR("EpetraMultiVectorT::getData");
843 
844  double ** arrayOfPointers;
845 
846  vec_->ExtractView(&arrayOfPointers);
847 
848  double * data = arrayOfPointers[j];
849  int localLength = vec_->MyLength();
850 
851  return ArrayRCP<double>(data, 0, localLength, false); // no ownership
852  }
853 
855  Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j) {
856  XPETRA_MONITOR("EpetraMultiVectorT::getDataNonConst");
857 
858  double ** arrayOfPointers;
859 
860  vec_->ExtractView(&arrayOfPointers);
861 
862  double * data = arrayOfPointers[j];
863  int localLength = vec_->MyLength();
864 
865  return ArrayRCP<double>(data, 0, localLength, false); // no ownership
866  }
867 
869 
871 
872 
874  void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const {
875  XPETRA_MONITOR("EpetraMultiVectorT::dot");
876 
877  XPETRA_DYNAMIC_CAST(const EpetraMultiVectorT<GlobalOrdinal XPETRA_COMMA Node>, A, eA, "This Xpetra::EpetraMultiVectorT method only accept Xpetra::EpetraMultiVectorT as input arguments.");
878  vec_->Dot(*eA.getEpetra_MultiVector(), dots.getRawPtr());
879  }
880 
882  void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("EpetraMultiVectorT::abs"); vec_->Abs(toEpetra<GlobalOrdinal,Node>(A)); }
883 
885  void reciprocal(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A) { XPETRA_MONITOR("EpetraMultiVectorT::reciprocal"); vec_->Reciprocal(toEpetra<GlobalOrdinal,Node>(A)); }
886 
888  void scale(const Scalar &alpha) { XPETRA_MONITOR("EpetraMultiVectorT::scale"); vec_->Scale(alpha); }
889 
891  void scale (Teuchos::ArrayView< const Scalar > alpha) {
892  XPETRA_MONITOR("EpetraMultiVectorT::scale");
893  // Epetra, unlike Tpetra, doesn't implement this version of
894  // scale(). Deal with this by scaling one column at a time.
895  const size_t numVecs = this->getNumVectors ();
896  for (size_t j = 0; j < numVecs; ++j) {
897  Epetra_Vector *v = (*vec_)(j);
898  v->Scale (alpha[j]);
899  }
900  }
901 
903  void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta) { XPETRA_MONITOR("EpetraMultiVectorT::update"); vec_->Update(alpha, toEpetra<GlobalOrdinal,Node>(A), beta); }
904 
906  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) { XPETRA_MONITOR("EpetraMultiVectorT::update"); vec_->Update(alpha, toEpetra<GlobalOrdinal,Node>(A), beta, toEpetra<GlobalOrdinal,Node>(B), gamma); }
907 
909  void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("EpetraMultiVectorT::norm1"); vec_->Norm1(norms.getRawPtr()); }
910 
912  void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("EpetraMultiVectorT::norm2"); vec_->Norm2(norms.getRawPtr()); }
913 
915  void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const { XPETRA_MONITOR("EpetraMultiVectorT::normInf"); vec_->NormInf(norms.getRawPtr()); }
916 
918  void meanValue(const Teuchos::ArrayView< Scalar > &means) const { XPETRA_MONITOR("EpetraMultiVectorT::meanValue"); vec_->MeanValue(means.getRawPtr()); } //TODO: modify ArrayView size ??
919 
921  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) { XPETRA_MONITOR("EpetraMultiVectorT::multiply"); vec_->Multiply(toEpetra(transA), toEpetra(transB), alpha, toEpetra(A), toEpetra(B), beta); }
922 
924  void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis) { XPETRA_MONITOR("EpetraMultiVectorT::elementWiseMultiply"); vec_->Multiply(scalarAB, toEpetra<GlobalOrdinal,Node>(A), toEpetra<GlobalOrdinal,Node>(B), scalarThis); }
925 
927 
929 
930 
932  size_t getNumVectors() const { XPETRA_MONITOR("EpetraMultiVectorT::getNumVectors"); return vec_->NumVectors(); }
933 
935  size_t getLocalLength() const { XPETRA_MONITOR("EpetraMultiVectorT::getLocalLength"); return vec_->MyLength(); }
936 
938  global_size_t getGlobalLength() const { XPETRA_MONITOR("EpetraMultiVectorT::getGlobalLength"); return vec_->GlobalLength64(); }
939 
942  XPETRA_MONITOR("EpetraMultiVectorT::isSameSize");
943  auto vv = toEpetra<GlobalOrdinal,Node>(vec);
944  return ( (vec_->MyLength() == vv.MyLength()) && (vec_->NumVectors() == vv.NumVectors()));
945  }
946 
948 
950 
951 
953  std::string description() const {
954  XPETRA_MONITOR("EpetraMultiVectorT::description");
955  TEUCHOS_TEST_FOR_EXCEPTION(true, Xpetra::Exceptions::NotImplemented, "TODO");
956  //return "TODO"; // unreachable
957  }
958 
960  void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel /* verbLevel */=Teuchos::Describable::verbLevel_default) const {
961  XPETRA_MONITOR("EpetraMultiVectorT::describe");
962  vec_->Print(out);
963  }
964 
966 
968  void randomize(bool bUseXpetraImplementation = false) {
969  XPETRA_MONITOR("EpetraMultiVectorT::randomize");
970 
971  if (bUseXpetraImplementation)
973  else
974  vec_->Random();
975  }
976 
978  //{@
979 
981  Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const { XPETRA_MONITOR("EpetraMultiVectorT::getMap"); return toXpetra<GlobalOrdinal,Node>(vec_->Map()); }
982 
985  XPETRA_MONITOR("EpetraMultiVectorT::doImport");
986 
987  XPETRA_DYNAMIC_CAST(const EpetraMultiVectorT<GlobalOrdinal XPETRA_COMMA Node>, source, tSource, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraMultiVectorT as input arguments.");
988  XPETRA_DYNAMIC_CAST(const EpetraImportT<GlobalOrdinal XPETRA_COMMA Node>, importer, tImporter, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraImportT as input arguments.");
989 
990  RCP<Epetra_MultiVector> v = tSource.getEpetra_MultiVector();
991  int err = this->getEpetra_MultiVector()->Import(*v, *tImporter.getEpetra_Import(), toEpetra(CM));
992  TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra is " << err);
993  }
994 
997  XPETRA_MONITOR("EpetraMultiVectorT::doExport");
998 
999  XPETRA_DYNAMIC_CAST(const EpetraMultiVectorT<GlobalOrdinal XPETRA_COMMA Node>, dest, tDest, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraMultiVectorT as input arguments.");
1000  XPETRA_DYNAMIC_CAST(const EpetraImportT<GlobalOrdinal XPETRA_COMMA Node>, importer, tImporter, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraImportT as input arguments.");
1001 
1002  RCP<Epetra_MultiVector> v = tDest.getEpetra_MultiVector();
1003  int err = this->getEpetra_MultiVector()->Export(*v, *tImporter.getEpetra_Import(), toEpetra(CM));
1004  TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
1005  }
1006 
1009  XPETRA_MONITOR("EpetraMultiVectorT::doImport");
1010 
1011  XPETRA_DYNAMIC_CAST(const EpetraMultiVectorT<GlobalOrdinal XPETRA_COMMA Node>, source, tSource, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraMultiVectorT as input arguments.");
1012  XPETRA_DYNAMIC_CAST(const EpetraExportT<GlobalOrdinal XPETRA_COMMA Node>, exporter, tExporter, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraImportT as input arguments.");
1013 
1014  RCP<Epetra_MultiVector> v = tSource.getEpetra_MultiVector();
1015  int err = this->getEpetra_MultiVector()->Import(*v, *tExporter.getEpetra_Export(), toEpetra(CM));
1016  TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
1017  }
1018 
1021  XPETRA_MONITOR("EpetraMultiVectorT::doExport");
1022 
1023  XPETRA_DYNAMIC_CAST(const EpetraMultiVectorT<GlobalOrdinal XPETRA_COMMA Node>, dest, tDest, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraMultiVectorT as input arguments.");
1024  XPETRA_DYNAMIC_CAST(const EpetraExportT<GlobalOrdinal XPETRA_COMMA Node>, exporter, tExporter, "Xpetra::EpetraMultiVectorT::doImport only accept Xpetra::EpetraImportT as input arguments.");
1025 
1026  RCP<Epetra_MultiVector> v = tDest.getEpetra_MultiVector();
1027  int err = this->getEpetra_MultiVector()->Export(*v, *tExporter.getEpetra_Export(), toEpetra(CM));
1028  TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
1029  }
1030 
1032  void replaceMap(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map) {
1033  XPETRA_MONITOR("EpetraMultiVectorT::replaceMap");
1034  int err = 0;
1035  if (!map.is_null()) {
1036  err = this->getEpetra_MultiVector()->ReplaceMap(toEpetra<GlobalOrdinal,Node>(map));
1037 
1038  } else {
1039  // Replace map with a dummy map to avoid potential hangs later
1040  Epetra_SerialComm SComm;
1041  Epetra_Map NewMap((GlobalOrdinal) vec_->MyLength(), (GlobalOrdinal) vec_->Map().IndexBase64(), SComm);
1042  err = this->getEpetra_MultiVector()->ReplaceMap(NewMap);
1043  }
1044  TEUCHOS_TEST_FOR_EXCEPTION(err != 0, std::runtime_error, "Catch error code returned by Epetra.");
1045  }
1046 
1048 
1050 
1051 
1053  EpetraMultiVectorT(const RCP<Epetra_MultiVector> &vec) : vec_(vec) { } //TODO removed const
1054 
1056  RCP<Epetra_MultiVector> getEpetra_MultiVector() const { return vec_; }
1057 
1059  void setSeed(unsigned int seed) {
1060  XPETRA_MONITOR("EpetraMultiVectorT::seedrandom");
1061 
1062  Teuchos::ScalarTraits< Scalar >::seedrandom(seed);
1063  vec_->SetSeed(seed);
1064  }
1065 
1066 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
1067 
1069 
1080  template<class TargetDeviceType>
1081  typename Kokkos::Impl::if_c<
1082  std::is_same<
1083  typename dual_view_type::t_dev_um::execution_space::memory_space,
1084  typename TargetDeviceType::memory_space>::value,
1085  typename dual_view_type::t_dev_um,
1086  typename dual_view_type::t_host_um>::type
1087  getLocalView () const {
1088  return this->MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >::template getLocalView<TargetDeviceType>();
1089  }
1090 
1091  typename dual_view_type::t_host_um getHostLocalView () const {
1092  typedef Kokkos::View< typename dual_view_type::t_host::data_type ,
1093  Kokkos::LayoutLeft,
1094  typename dual_view_type::t_host::device_type ,
1095  Kokkos::MemoryUnmanaged> epetra_view_type;
1096 
1097  // access Epetra multivector data
1098  double* data = NULL;
1099  int myLDA;
1100  vec_->ExtractView(&data, &myLDA);
1101  int localLength = vec_->MyLength();
1102  int numVectors = getNumVectors();
1103 
1104  // create view
1105  epetra_view_type test = epetra_view_type(data, localLength, numVectors);
1106  typename dual_view_type::t_host_um ret = subview(test, Kokkos::ALL(), Kokkos::ALL());
1107 
1108  return ret;
1109  }
1110 
1111  typename dual_view_type::t_dev_um getDeviceLocalView() const {
1112  throw std::runtime_error("Epetra does not support device views!");
1113 #ifndef __NVCC__ //prevent nvcc warning
1114  typename dual_view_type::t_dev_um ret;
1115 #endif
1116  TEUCHOS_UNREACHABLE_RETURN(ret);
1117  }
1118 
1119 #endif
1120 
1122 
1123  protected:
1126  virtual void
1128  typedef EpetraMultiVectorT this_type;
1129  const this_type* rhsPtr = dynamic_cast<const this_type*> (&rhs);
1130  TEUCHOS_TEST_FOR_EXCEPTION(
1131  rhsPtr == NULL, std::invalid_argument, "Xpetra::MultiVector::operator=: "
1132  "The left-hand side (LHS) of the assignment has a different type than "
1133  "the right-hand side (RHS). The LHS has type Xpetra::EpetraMultiVectorT "
1134  "(which means it wraps an Epetra_MultiVector), but the RHS has some "
1135  "other type. This probably means that the RHS wraps a Tpetra::Multi"
1136  "Vector. Xpetra::MultiVector does not currently implement assignment "
1137  "from a Tpetra object to an Epetra object, though this could be added "
1138  "with sufficient interest.");
1139 
1140  RCP<const Epetra_MultiVector> rhsImpl = rhsPtr->getEpetra_MultiVector ();
1141  RCP<Epetra_MultiVector> lhsImpl = this->getEpetra_MultiVector ();
1142 
1143  TEUCHOS_TEST_FOR_EXCEPTION(
1144  rhsImpl.is_null (), std::logic_error, "Xpetra::MultiVector::operator= "
1145  "(in Xpetra::EpetraMultiVectorT::assign): *this (the right-hand side of "
1146  "the assignment) has a null RCP<Epetra_MultiVector> inside. Please "
1147  "report this bug to the Xpetra developers.");
1148  TEUCHOS_TEST_FOR_EXCEPTION(
1149  lhsImpl.is_null (), std::logic_error, "Xpetra::MultiVector::operator= "
1150  "(in Xpetra::EpetraMultiVectorT::assign): The left-hand side of the "
1151  "assignment has a null RCP<Epetra_MultiVector> inside. Please report "
1152  "this bug to the Xpetra developers.");
1153 
1154  // Epetra_MultiVector's assignment operator does a deep copy.
1155  *lhsImpl = *rhsImpl;
1156  }
1157 
1158  private:
1160  RCP< Epetra_MultiVector > vec_;
1161 
1162  }; // EpetraMultiVectorT class (specialization on GO=long long, NO=EpetraNode
1163 #endif
1164 
1165 } // Xpetra namespace
1166 
1167 #include "Xpetra_EpetraVector.hpp" // to avoid incomplete type instantiated above in out-of-body functions.
1168 
1169 #endif // XPETRA_EPETRAMULTIVECTOR_HPP
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const
Compute the dot product of each corresponding pair of vectors (columns) in A and B.
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: this = gamma*this + alpha*A + beta*B.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
Set multi-vector values from array of pointers using Teuchos memory management classes. (copy).
size_t getNumVectors() const
Number of columns in the multivector.
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: this = beta*this + alpha*A.
std::string description() const
A simple one-line description of this object.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const
Compute the dot product of each corresponding pair of vectors (columns) in A and B.
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: this = gamma*this + alpha*A + beta*B.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
Set multi-vector values from array of pointers using Teuchos memory management classes. (copy).
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: this = gamma*this + alpha*A + beta*B.
void scale(const Scalar &alpha)
Scale in place: this = alpha*this.
void scale(const Scalar &alpha)
Scale in place: this = alpha*this.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
Update: this = beta*this + alpha*A.
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
RCP< Epetra_MultiVector > vec_
The Epetra_MultiVector which this class wraps.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
size_t getLocalLength() const
Local number of rows on the calling process.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
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...
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic MultiVector constuctor.
EpetraMultiVectorT(const RCP< Epetra_MultiVector > &vec)
EpetraMultiVectorT constructor to wrap a Epetra_MultiVector object.
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 doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
Exception throws to report errors in the internal logical of the program.
virtual void Xpetra_randomize()
Set multi-vector values to random numbers. XPetra implementation.
const Epetra_CrsGraph & toEpetra(const RCP< const CrsGraph< int, GlobalOrdinal, Node > > &graph)
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
Checks to see if the local length, number of vectors and size of Scalar type match.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic MultiVector constuctor.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
RCP< Epetra_MultiVector > getEpetra_MultiVector() const
Get the underlying Epetra multivector.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
Multiply a Vector A elementwise by a MultiVector B.
global_size_t getGlobalLength() const
Global number of rows in the multivector.
void norm1(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
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 replaceMap(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
Replace the underlying Map in place.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
EpetraMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
MultiVector copy constructor.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
Multiply a Vector A elementwise by a MultiVector B.
RCP< Epetra_MultiVector > getEpetra_MultiVector() const
Get the underlying Epetra multivector.
global_size_t getGlobalLength() const
Global number of rows in the multivector.
void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
std::string description() const
A simple one-line description of this object.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
size_t getNumVectors() const
Number of columns in the multivector.
void setSeed(unsigned int seed)
Set seed for Random function.
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
EpetraMultiVectorT(const RCP< Epetra_MultiVector > &vec)
EpetraMultiVectorT constructor to wrap a Epetra_MultiVector object.
virtual ~EpetraMultiVectorT()
MultiVector destructor.
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const Teuchos::ArrayView< const Teuchos::ArrayView< const Scalar > > &ArrayOfPtrs, size_t NumVectors)
Set multi-vector values from array of pointers using Teuchos memory management classes. (copy).
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
EpetraMultiVectorT(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Basic MultiVector constuctor.
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...
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
size_t getNumVectors() const
Number of columns in the multivector.
void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
void replaceMap(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
Replace the underlying Map in place.
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).
global_size_t getGlobalLength() const
Global number of rows in the multivector.
Exception throws when you call an unimplemented method of Xpetra.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_t global_size_t
Global size_t object.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)
Update: this = beta*this + alpha*A.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
size_t getLocalLength() const
Local number of rows on the calling process.
void sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
Checks to see if the local length, number of vectors and size of Scalar type match.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)
View of the local values in a particular vector of this multivector.
void setSeed(unsigned int seed)
Set seed for Random function.
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).
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
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...
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).
EpetraMultiVectorT(const RCP< Epetra_MultiVector > &vec)
EpetraMultiVectorT constructor to wrap a Epetra_MultiVector object.
void putScalar(const Scalar &value)
Set all values in the multivector with the given value.
void scale(Teuchos::ArrayView< const Scalar > alpha)
Scale the current values of a multi-vector, this[j] = alpha[j]*this[j].
void setSeed(unsigned int seed)
Set seed for Random function.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
void scale(const Scalar &alpha)
Scale in place: this = alpha*this.
RCP< Epetra_MultiVector > vec_
The Epetra_MultiVector which this class wraps.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel=Teuchos::Describable::verbLevel_default) const
Print the object with the given verbosity level to a FancyOStream.
void replaceMap(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
Replace the underlying Map in place.
EpetraMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
MultiVector copy constructor.
RCP< Epetra_MultiVector > getEpetra_MultiVector() const
Get the underlying Epetra multivector.
CombineMode
Xpetra::Combine Mode enumerable type.
#define XPETRA_MONITOR(funcName)
void normInf(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
void elementWiseMultiply(Scalar scalarAB, const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, Scalar scalarThis)
Multiply a Vector A elementwise by a MultiVector B.
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).
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
EpetraMultiVectorT(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source)
MultiVector copy constructor.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
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).
std::string description() const
A simple one-line description of this object.
Teuchos::RCP< const Vector< double, int, GlobalOrdinal, Node > > getVector(size_t j) const
Return a Vector which is a const view of column j.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
Teuchos::RCP< Vector< double, int, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
void dot(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Teuchos::ArrayView< Scalar > &dots) const
Compute the dot product of each corresponding pair of vectors (columns) in A and B.
void norm2(const Teuchos::ArrayView< Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const