Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_TpetraMultiVector_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Xpetra: A linear algebra interface package
4//
5// Copyright 2012 NTESS and the Xpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef XPETRA_TPETRAMULTIVECTOR_DEF_HPP
11#define XPETRA_TPETRAMULTIVECTOR_DEF_HPP
13
14#include "Xpetra_TpetraMap.hpp" //TMP
15#include "Xpetra_Utils.hpp"
16#include "Xpetra_TpetraImport.hpp"
17#include "Xpetra_TpetraExport.hpp"
18
20#include "Tpetra_MultiVector.hpp"
21#include "Tpetra_Vector.hpp"
22#include "Tpetra_Details_Random.hpp"
23
24namespace Xpetra {
25
27template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
29 : vec_(Teuchos::rcp(new Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(map), NumVectors, zeroOut))) {
30 // TAW 1/30/2016: even though Tpetra allows numVecs == 0, Epetra does not. Introduce exception to keep behavior of Epetra and Tpetra consistent.
31 TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument, "Xpetra::TpetraMultiVector(map,numVecs,zeroOut): numVecs = " << NumVectors << " < 1.");
32}
33
35template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
38 : vec_(Teuchos::rcp(new Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(source), copyOrView))) {}
39
41template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
43 TpetraMultiVector(const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &map, const Teuchos::ArrayView<const Scalar> &A, size_t LDA, size_t NumVectors)
44 : vec_(Teuchos::rcp(new Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(map), A, LDA, NumVectors))) {
45 // TAW 1/30/2016: even though Tpetra allows numVecs == 0, Epetra does not. Introduce exception to keep behavior of Epetra and Tpetra consistent.
46 TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument, "Xpetra::TpetraMultiVector(map,A,LDA,numVecs): numVecs = " << NumVectors << " < 1.");
47}
48
50template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
53 : vec_(Teuchos::rcp(new Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(map), ArrayOfPtrs, NumVectors))) {
54 // TAW 1/30/2016: even though Tpetra allows numVecs == 0, Epetra does not. Introduce exception to keep behavior of Epetra and Tpetra consistent.
55 TEUCHOS_TEST_FOR_EXCEPTION(NumVectors < 1, std::invalid_argument, "Xpetra::TpetraMultiVector(map,ArrayOfPtrs,numVecs): numVecs = " << NumVectors << " < 1.");
56}
57
59template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
64template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
66 replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) {
67 XPETRA_MONITOR("TpetraMultiVector::replaceGlobalValue");
68 vec_->replaceGlobalValue(globalRow, vectorIndex, value);
69}
72template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
74 sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value) {
75 XPETRA_MONITOR("TpetraMultiVector::sumIntoGlobalValue");
76 vec_->sumIntoGlobalValue(globalRow, vectorIndex, value);
77}
78
80template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
82 replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) {
83 XPETRA_MONITOR("TpetraMultiVector::replaceLocalValue");
84 vec_->replaceLocalValue(myRow, vectorIndex, value);
86
88template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
90 sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value) {
91 XPETRA_MONITOR("TpetraMultiVector::sumIntoLocalValue");
92 vec_->sumIntoLocalValue(myRow, vectorIndex, value);
93}
94
96template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
98 putScalar(const Scalar &value) {
99 XPETRA_MONITOR("TpetraMultiVector::putScalar");
100 vec_->putScalar(value);
101}
104template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
106 reduce() {
107 XPETRA_MONITOR("TpetraMultiVector::reduce");
108 vec_->reduce();
109}
110
111template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
114 getVector(size_t j) const {
115 XPETRA_MONITOR("TpetraMultiVector::getVector");
116 return toXpetra(vec_->getVector(j));
118
120template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
124 XPETRA_MONITOR("TpetraMultiVector::getVectorNonConst");
125 return toXpetra(vec_->getVectorNonConst(j));
127
129template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
132 getData(size_t j) const {
133 XPETRA_MONITOR("TpetraMultiVector::getData");
134 return vec_->getData(j);
135}
136
138template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
141 getDataNonConst(size_t j) {
142 XPETRA_MONITOR("TpetraMultiVector::getDataNonConst");
143 return vec_->getDataNonConst(j);
144}
145
147template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
150 XPETRA_MONITOR("TpetraMultiVector::get1dCopy");
151 vec_->get1dCopy(A, LDA);
153
155template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
158 XPETRA_MONITOR("TpetraMultiVector::get2dCopy");
159 vec_->get2dCopy(ArrayOfPtrs);
160}
163template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
166 get1dView() const {
167 XPETRA_MONITOR("TpetraMultiVector::get1dView");
168 return vec_->get1dView();
169}
172template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
175 get2dView() const {
176 XPETRA_MONITOR("TpetraMultiVector::get2dView");
177 return vec_->get2dView();
178}
179
181template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
185 XPETRA_MONITOR("TpetraMultiVector::get1dViewNonConst");
186 return vec_->get1dViewNonConst();
187}
188
190template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
194 XPETRA_MONITOR("TpetraMultiVector::get2dViewNonConst");
195 return vec_->get2dViewNonConst();
197
199template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
202 XPETRA_MONITOR("TpetraMultiVector::dot");
203 vec_->dot(toTpetra(A), dots);
204}
207template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
213
215template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
218 XPETRA_MONITOR("TpetraMultiVector::reciprocal");
219 vec_->reciprocal(toTpetra(A));
221
223template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
225 scale(const Scalar &alpha) {
226 XPETRA_MONITOR("TpetraMultiVector::scale");
227 vec_->scale(alpha);
229
231template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
234 XPETRA_MONITOR("TpetraMultiVector::scale");
235 vec_->scale(alpha);
237
239template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
242 XPETRA_MONITOR("TpetraMultiVector::scale");
243 vec_->scale(alpha, toTpetra(A));
244}
245
247template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
249 update(const Scalar &alpha, const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &A, const Scalar &beta) {
250 XPETRA_MONITOR("TpetraMultiVector::update");
251 vec_->update(alpha, toTpetra(A), beta);
253
255template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
257 update(const Scalar &alpha, const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &A, const Scalar &beta, const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &B, const Scalar &gamma) {
258 XPETRA_MONITOR("TpetraMultiVector::update");
259 vec_->update(alpha, toTpetra(A), beta, toTpetra(B), gamma);
261
263template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
266 XPETRA_MONITOR("TpetraMultiVector::norm1");
267 vec_->norm1(norms);
268}
269
271template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
274 XPETRA_MONITOR("TpetraMultiVector::norm2");
275 vec_->norm2(norms);
276}
277
279template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
282 XPETRA_MONITOR("TpetraMultiVector::normInf");
283 vec_->normInf(norms);
284}
285
287template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
289 meanValue(const Teuchos::ArrayView<Scalar> &means) const {
290 XPETRA_MONITOR("TpetraMultiVector::meanValue");
291 vec_->meanValue(means);
292}
293
295template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
298 XPETRA_MONITOR("TpetraMultiVector::multiply");
299 vec_->multiply(transA, transB, alpha, toTpetra(A), toTpetra(B), beta);
300}
301
303template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
305 getNumVectors() const {
306 XPETRA_MONITOR("TpetraMultiVector::getNumVectors");
307 return vec_->getNumVectors();
308}
309
311template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
313 getLocalLength() const {
314 XPETRA_MONITOR("TpetraMultiVector::getLocalLength");
315 return vec_->getLocalLength();
316}
317
319template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
321 getGlobalLength() const {
322 XPETRA_MONITOR("TpetraMultiVector::getGlobalLength");
323 return vec_->getGlobalLength();
324}
325
326// \brief Checks to see if the local length, number of vectors and size of Scalar type match
327template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
330 XPETRA_MONITOR("TpetraMultiVector::isSameSize");
331 return vec_->isSameSize(toTpetra(vec));
332}
333
335template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
337 description() const {
338 XPETRA_MONITOR("TpetraMultiVector::description");
339 return vec_->description();
340}
341
343template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
345 describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const {
346 XPETRA_MONITOR("TpetraMultiVector::describe");
347 vec_->describe(out, verbLevel);
348}
349
351template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
353 randomize(bool bUseXpetraImplementation) {
354 XPETRA_MONITOR("TpetraMultiVector::randomize");
355
356 if (bUseXpetraImplementation)
358 else
359 vec_->randomize();
360}
361
363template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
365 randomize(const Scalar &minVal, const Scalar &maxVal, bool bUseXpetraImplementation) {
366 XPETRA_MONITOR("TpetraMultiVector::randomize");
367
368 if (bUseXpetraImplementation)
370 else
371 vec_->randomize(minVal, maxVal);
372}
373
374template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
377 getMap() const {
378 XPETRA_MONITOR("TpetraMultiVector::getMap");
379 return toXpetra(vec_->getMap());
380}
381
382template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
385 XPETRA_MONITOR("TpetraMultiVector::doImport");
386
387 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
389 this->getTpetra_MultiVector()->doImport(*v, toTpetra(importer), toTpetra(CM));
390}
391
392template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
395 XPETRA_MONITOR("TpetraMultiVector::beginImport");
396
397 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
399 this->getTpetra_MultiVector()->beginImport(*v, toTpetra(importer), toTpetra(CM));
400}
401
402template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
405 XPETRA_MONITOR("TpetraMultiVector::endImport");
406
407 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
409 this->getTpetra_MultiVector()->endImport(*v, toTpetra(importer), toTpetra(CM));
410}
411
412template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
415 XPETRA_MONITOR("TpetraMultiVector::doExport");
416
417 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
419 this->getTpetra_MultiVector()->doExport(*v, toTpetra(importer), toTpetra(CM));
420}
421
422template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
425 XPETRA_MONITOR("TpetraMultiVector::beginExport");
426
427 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
429 this->getTpetra_MultiVector()->beginExport(*v, toTpetra(importer), toTpetra(CM));
430}
431
432template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
435 XPETRA_MONITOR("TpetraMultiVector::endExport");
436
437 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
439 this->getTpetra_MultiVector()->endExport(*v, toTpetra(importer), toTpetra(CM));
440}
441
442template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
445 XPETRA_MONITOR("TpetraMultiVector::doImport");
446
447 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
449 this->getTpetra_MultiVector()->doImport(*v, toTpetra(exporter), toTpetra(CM));
450}
451
452template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
455 XPETRA_MONITOR("TpetraMultiVector::beginImport");
456
457 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
459 this->getTpetra_MultiVector()->beginImport(*v, toTpetra(exporter), toTpetra(CM));
460}
461
462template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
465 XPETRA_MONITOR("TpetraMultiVector::endImport");
466
467 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, source, tSource, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
469 this->getTpetra_MultiVector()->endImport(*v, toTpetra(exporter), toTpetra(CM));
470}
471
472template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
475 XPETRA_MONITOR("TpetraMultiVector::doExport");
476
477 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
479 this->getTpetra_MultiVector()->doExport(*v, toTpetra(exporter), toTpetra(CM));
480}
481
482template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
485 XPETRA_MONITOR("TpetraMultiVector::beginExport");
486
487 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
489 this->getTpetra_MultiVector()->beginExport(*v, toTpetra(exporter), toTpetra(CM));
490}
491
492template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
495 XPETRA_MONITOR("TpetraMultiVector::endExport");
496
497 XPETRA_DYNAMIC_CAST(const TpetraMultiVectorClass, dest, tDest, "Xpetra::TpetraMultiVector::doImport only accept Xpetra::TpetraMultiVector as input arguments."); // TODO: remove and use toTpetra()
499 this->getTpetra_MultiVector()->endExport(*v, toTpetra(exporter), toTpetra(CM));
500}
501
502template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
505 XPETRA_MONITOR("TpetraMultiVector::replaceMap");
506 this->getTpetra_MultiVector()->replaceMap(toTpetra(map));
507}
508
510template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
514
516template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
520
522template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
524 setSeed(unsigned int seed) {
525 XPETRA_MONITOR("TpetraMultiVector::seedrandom");
527 // Tell Tpetra to update its RNG pool for the new random seed
528 Tpetra::Details::Static_Random_XorShift64_Pool<typename Node::device_type::execution_space>::resetPool(getMap()->getComm()->getRank());
529}
530
531template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
534 getLocalViewHost(Tpetra::Access::ReadOnlyStruct) const {
535 return subview(vec_->getLocalViewHost(Tpetra::Access::ReadOnly), Kokkos::ALL(), Kokkos::ALL());
536}
537
538template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
541 getLocalViewDevice(Tpetra::Access::ReadOnlyStruct) const {
542 return subview(vec_->getLocalViewDevice(Tpetra::Access::ReadOnly), Kokkos::ALL(), Kokkos::ALL());
543}
544
545template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
548 getLocalViewHost(Tpetra::Access::OverwriteAllStruct) const {
549 return subview(vec_->getLocalViewHost(Tpetra::Access::OverwriteAll), Kokkos::ALL(), Kokkos::ALL());
550}
551
552template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
555 getLocalViewDevice(Tpetra::Access::OverwriteAllStruct) const {
556 return subview(vec_->getLocalViewDevice(Tpetra::Access::OverwriteAll), Kokkos::ALL(), Kokkos::ALL());
557}
558
559template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
562 getLocalViewHost(Tpetra::Access::ReadWriteStruct) const {
563 return subview(vec_->getLocalViewHost(Tpetra::Access::ReadWrite), Kokkos::ALL(), Kokkos::ALL());
564}
565
566template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
569 getLocalViewDevice(Tpetra::Access::ReadWriteStruct) const {
570 return subview(vec_->getLocalViewDevice(Tpetra::Access::ReadWrite), Kokkos::ALL(), Kokkos::ALL());
571}
572
575template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
579 const this_type *rhsPtr = dynamic_cast<const this_type *>(&rhs);
581 rhsPtr == NULL, std::invalid_argument,
582 "Xpetra::MultiVector::operator=:"
583 " The left-hand side (LHS) of the assignment has a different type than "
584 "the right-hand side (RHS). The LHS has type Xpetra::TpetraMultiVector"
585 " (which means it wraps a Tpetra::MultiVector), but the RHS has some "
586 "other type. This probably means that the RHS wraps an "
587 "Epetra_MultiVector. Xpetra::MultiVector does not currently implement "
588 "assignment from an Epetra object to a Tpetra object, though this could"
589 " be added with sufficient interest.");
590
592 RCP<const TMV> rhsImpl = rhsPtr->getTpetra_MultiVector();
593 RCP<TMV> lhsImpl = this->getTpetra_MultiVector();
594
596 rhsImpl.is_null(), std::logic_error,
597 "Xpetra::MultiVector::operator= "
598 "(in Xpetra::TpetraMultiVector::assign): *this (the right-hand side of "
599 "the assignment) has a null RCP<Tpetra::MultiVector> inside. Please "
600 "report this bug to the Xpetra developers.");
602 lhsImpl.is_null(), std::logic_error,
603 "Xpetra::MultiVector::operator= "
604 "(in Xpetra::TpetraMultiVector::assign): The left-hand side of the "
605 "assignment has a null RCP<Tpetra::MultiVector> inside. Please report "
606 "this bug to the Xpetra developers.");
607
608 Tpetra::deep_copy(*lhsImpl, *rhsImpl);
609}
610
611} // namespace Xpetra
612
613// Following header file inculsion is needed for the dynamic_cast to TpetraVector in
614// elementWiseMultiply (because we cannot dynamic_cast if target is not a complete type)
615// It is included here to avoid circular dependency between Vector and MultiVector
616// TODO: there is certainly a more elegant solution...
617#include "Xpetra_TpetraVector.hpp"
618
619namespace Xpetra {
620
621template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
623 elementWiseMultiply(Scalar scalarAB,
626 Scalar scalarThis) {
627 XPETRA_MONITOR("TpetraMultiVector::elementWiseMultiply");
628
629 // XPETRA_DYNAMIC_CAST won't take TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>
630 // as an argument, hence the following typedef.
632 XPETRA_DYNAMIC_CAST(const tpv, A, tA, "Xpetra::TpetraMultiVectorMatrix->multiply() only accept Xpetra::TpetraMultiVector as input arguments.");
633 XPETRA_DYNAMIC_CAST(const TpetraMultiVector, B, tB, "Xpetra::TpetraMultiVectorMatrix->multiply() only accept Xpetra::TpetraMultiVector as input arguments.");
634 vec_->elementWiseMultiply(scalarAB, *tA.getTpetra_Vector(), *tB.getTpetra_MultiVector(), scalarThis);
635}
636
637} // namespace Xpetra
638
639#define XPETRA_TPETRAMULTIVECTOR_SHORT
640#endif // XPETRA_TPETRAMULTIVECTOR_HPP
#define XPETRA_MONITOR(funcName)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
bool is_null() const
virtual void Xpetra_randomize()
Set multi-vector values to random numbers. XPetra implementation.
void replaceMap(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
void replaceGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Replace value, using global (row) index.
void doExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
bool isSameSize(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &vec) const
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.
virtual dual_view_type::t_dev_const_um getLocalViewDevice(Tpetra::Access::ReadOnlyStruct) 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::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(size_t j)
Return a Vector which is a nonconst view of column j.
void doImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
virtual void assign(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &rhs)
Implementation of the assignment operator (operator=); does a deep copy.
void beginExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
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 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 putScalar(const Scalar &value)
Set all values in the multivector with the given value.
Teuchos::ArrayRCP< Scalar > get1dViewNonConst()
Nonconst persisting (1-D) view of this multivector's local values.
size_t getLocalLength() const
Local number of rows on the calling process.
size_t getNumVectors() const
Number of columns in the multivector.
void sumIntoGlobalValue(GlobalOrdinal globalRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using global (row) index.
void setSeed(unsigned int seed)
Set seed for Random function.
void reduce()
Sum values of a locally replicated multivector across all processes.
void replaceLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Replace value, using local (row) index.
void norm2(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
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 endExport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
std::string description() const
A simple one-line description of this object.
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 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.
virtual dual_view_type::t_host_const_um getLocalViewHost(Tpetra::Access::ReadOnlyStruct) const
global_size_t getGlobalLength() const
Global number of rows in the multivector.
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,...
virtual ~TpetraMultiVector()
Destructor (virtual for memory safety of derived classes).
void endImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
void abs(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A)
Put element-wise absolute values of input Multi-vector in target: A = abs(this).
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 sumIntoLocalValue(LocalOrdinal myRow, size_t vectorIndex, const Scalar &value)
Add value to existing value, using local (row) index.
Teuchos::ArrayRCP< const Scalar > get1dView() const
Const persisting (1-D) view of this multivector's local values.
Teuchos::ArrayRCP< const Scalar > getData(size_t j) const
Const view of the local values in a particular vector of this multivector.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_MultiVector() const
Get the underlying Tpetra multivector.
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's local values.
void beginImport(const DistObject< Scalar, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Teuchos::ArrayRCP< Teuchos::ArrayRCP< Scalar > > get2dViewNonConst()
Return non-const persisting pointers to values.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
The Map describing the parallel distribution of this object.
void norm1(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute 1-norm of each vector in multi-vector.
void normInf(const Teuchos::ArrayView< typename Teuchos::ScalarTraits< Scalar >::magnitudeType > &norms) const
Compute Inf-norm of each vector in multi-vector.
void randomize(bool bUseXpetraImplementation=false)
Set multi-vector values to random numbers.
void get1dCopy(Teuchos::ArrayView< Scalar > A, size_t LDA) const
Fill the given array with a copy of this multivector's local values.
Teuchos::ArrayRCP< Teuchos::ArrayRCP< const Scalar > > get2dView() const
Return const persisting pointers to values.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > graph)
size_t global_size_t
Global size_t object.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
CombineMode
Xpetra::Combine Mode enumerable type.
static void seedrandom(unsigned int s)