42 #ifndef THYRA_TPETRA_VECTOR_HPP
43 #define THYRA_TPETRA_VECTOR_HPP
46 #include "Thyra_TpetraVector_decl.hpp"
47 #include "Thyra_TpetraMultiVector.hpp"
55 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
60 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
63 const RCP<Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVector
66 initializeImpl(tpetraVectorSpace, tpetraVector);
70 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
73 const RCP<
const Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVector
76 initializeImpl(tpetraVectorSpace, tpetraVector);
80 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
84 return tpetraVector_.getNonconstObj();
88 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
97 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
101 if (domainSpace_.is_null()) {
102 domainSpace_ = tpetraVectorSpace<Scalar>(
103 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal>(
105 tpetraVector_.getConstObj()->getMap()->getComm(),
106 tpetraVector_.getConstObj()->getMap()->getNode()
117 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
121 return tpetraVectorSpace_;
128 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
132 *localValues = tpetraVector_.getNonconstObj()->get1dViewNonConst();
136 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
140 *localValues = tpetraVector_->get1dView();
147 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
155 if (!tpetraVector_.getNonconstObj()->isDistributed()) {
156 auto comm = tpetraVector_.getNonconstObj()->getMap()->getComm();
157 if (tpetraVector_.getConstObj()->getMap()->getComm()->getRank() == 0)
158 tpetraVector_.getNonconstObj()->randomize(l, u);
161 tpetraVector_.getNonconstObj()->reduce();
163 tpetraVector_.getNonconstObj()->randomize(l, u);
168 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
173 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
178 tpetraVector_.getNonconstObj()->abs(*tx);
181 tpetraVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
182 tpetraVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
188 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
193 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
198 tpetraVector_.getNonconstObj()->reciprocal(*tx);
201 tpetraVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
202 tpetraVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
208 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
213 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
219 tpetraVector_.getNonconstObj()->elementWiseMultiply(
220 ST::one(), *tx, *tpetraVector_.getConstObj(), ST::zero());
223 tpetraVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
224 tpetraVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
230 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
236 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
244 = Tpetra::createVector<Scalar>(tx->getMap());
245 temp->elementWiseMultiply(
246 ST::one(), *tx, *tpetraVector_.getConstObj(), ST::zero());
247 return ST::magnitude(ST::squareroot(tpetraVector_.getConstObj()->dot(*temp)));
250 tpetraVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
256 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
266 for (
auto itr = vecs.begin(); itr != vecs.end(); ++itr) {
267 auto tv = this->getConstTpetraVector(Teuchos::rcpFromPtr(*itr));
269 typedef Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
270 Teuchos::rcp_const_cast<TV>(tv)->template sync<Kokkos::HostSpace>();
275 for (
auto itr = targ_vecs.begin(); itr != targ_vecs.end(); ++itr) {
276 auto tv = this->getTpetraVector(Teuchos::rcpFromPtr(*itr));
278 tv->template sync<Kokkos::HostSpace>();
279 tv->template modify<Kokkos::HostSpace>();
287 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
295 typedef typename Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
296 Teuchos::rcp_const_cast<TV>(
297 tpetraVector_.getConstObj())->
template sync<Kokkos::HostSpace>();
303 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
311 tpetraVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
312 tpetraVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
318 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
327 typedef typename Tpetra::Vector<
328 Scalar,LocalOrdinal,GlobalOrdinal,Node>::execution_space execution_space;
329 tpetraVector_.getNonconstObj()->template sync<execution_space>();
336 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
339 tpetraVector_.getNonconstObj()->putScalar(alpha);
343 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
347 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
352 tpetraVector_.getNonconstObj()->assign(*tmv);
355 tpetraVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
356 tpetraVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
362 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
365 tpetraVector_.getNonconstObj()->scale(alpha);
369 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
375 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
381 tpetraVector_.getNonconstObj()->update(alpha, *tmv, ST::one());
384 tpetraVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
385 tpetraVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
391 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
405 bool allCastsSuccessful =
true;
407 auto mvIter = mv.begin();
408 auto tmvIter = tmvs.begin();
409 for (; mvIter != mv.end(); ++mvIter, ++tmvIter) {
410 tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromPtr(*mvIter));
414 allCastsSuccessful =
false;
422 auto len = mv.size();
424 tpetraVector_.getNonconstObj()->scale(beta);
425 }
else if (len == 1 && allCastsSuccessful) {
426 tpetraVector_.getNonconstObj()->update(alpha[0], *tmvs[0], beta);
427 }
else if (len == 2 && allCastsSuccessful) {
428 tpetraVector_.getNonconstObj()->update(alpha[0], *tmvs[0], alpha[1], *tmvs[1], beta);
429 }
else if (allCastsSuccessful) {
431 auto tmvIter = tmvs.begin();
432 auto alphaIter = alpha.
begin();
437 for (; tmvIter != tmvs.end(); ++tmvIter) {
438 if (tmvIter->getRawPtr() == tpetraVector_.getConstObj().getRawPtr()) {
440 tmv = Teuchos::rcp(
new TpetraMultiVector_t(
441 *tpetraVector_.getConstObj(), Teuchos::Copy));
446 tmvIter = tmvs.begin();
450 if ((tmvs.size() % 2) == 0) {
451 tpetraVector_.getNonconstObj()->scale(beta);
453 tpetraVector_.getNonconstObj()->update(*alphaIter, *(*tmvIter), beta);
457 for (; tmvIter != tmvs.end(); tmvIter+=2, alphaIter+=2) {
458 tpetraVector_.getNonconstObj()->update(
459 *alphaIter, *(*tmvIter), *(alphaIter+1), *(*(tmvIter+1)), ST::one());
463 tpetraVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
464 tpetraVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
470 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
476 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
481 tpetraVector_.getConstObj()->dot(*tmv, prods);
484 tpetraVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
485 tpetraVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
491 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
496 tpetraVector_.getConstObj()->norm1(norms);
500 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
505 tpetraVector_.getConstObj()->norm2(norms);
509 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
514 tpetraVector_.getConstObj()->normInf(norms);
518 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
528 typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
529 typedef Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
537 typedef typename TMV::execution_space execution_space;
538 Teuchos::rcp_const_cast<TMV>(X_tpetra)->template sync<execution_space>();
539 Y_tpetra->template sync<execution_space>();
540 Teuchos::rcp_const_cast<TV>(tpetraVector_.getConstObj())->
template sync<execution_space>();
545 "Error, conjugation without transposition is not allowed for complex scalar types!");
550 trans = Teuchos::NO_TRANS;
553 trans = Teuchos::NO_TRANS;
559 trans = Teuchos::CONJ_TRANS;
563 Y_tpetra->template modify<execution_space>();
564 Y_tpetra->multiply(trans, Teuchos::NO_TRANS, alpha, *tpetraVector_.getConstObj(), *X_tpetra, beta);
566 Teuchos::rcp_const_cast<TV>(tpetraVector_.getConstObj())->
template sync<Kokkos::HostSpace>();
576 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
577 template<
class TpetraVector_t>
587 tpetraVectorSpace_ = tpetraVectorSpace;
588 tpetraVector_.initialize(tpetraVector);
589 this->updateSpmdSpace();
593 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
594 RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
595 TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
596 getTpetraMultiVector(
const RCP<MultiVectorBase<Scalar> >& mv)
const
598 using Teuchos::rcp_dynamic_cast;
599 typedef TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
600 typedef TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
602 RCP<TMV> tmv = rcp_dynamic_cast<TMV>(mv);
604 return tmv->getTpetraMultiVector();
607 RCP<TV> tv = rcp_dynamic_cast<TV>(mv);
609 return tv->getTpetraVector();
612 return Teuchos::null;
616 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
617 RCP<const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
618 TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
619 getConstTpetraMultiVector(
const RCP<
const MultiVectorBase<Scalar> >& mv)
const
621 using Teuchos::rcp_dynamic_cast;
622 typedef TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
623 typedef TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
625 RCP<const TMV> tmv = rcp_dynamic_cast<const TMV>(mv);
627 return tmv->getConstTpetraMultiVector();
630 RCP<const TV> tv = rcp_dynamic_cast<const TV>(mv);
632 return tv->getConstTpetraVector();
635 return Teuchos::null;
639 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
640 RCP<Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
644 typedef TpetraVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> TV;
645 RCP<TV> tv = Teuchos::rcp_dynamic_cast<TV>(v);
647 return tv->getTpetraVector();
649 return Teuchos::null;
654 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
655 RCP<const Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
659 typedef TpetraVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> TV;
660 RCP<const TV> tv = Teuchos::rcp_dynamic_cast<const TV>(v);
662 return tv->getConstTpetraVector();
664 return Teuchos::null;
672 #endif // THYRA_TPETRA_VECTOR_HPP