42 #ifndef THYRA_TPETRA_MULTIVECTOR_HPP
43 #define THYRA_TPETRA_MULTIVECTOR_HPP
57 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
62 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
65 const RCP<
const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
66 const RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraMultiVector
69 initializeImpl(tpetraVectorSpace, domainSpace, tpetraMultiVector);
73 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
76 const RCP<
const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
77 const RCP<
const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraMultiVector
80 initializeImpl(tpetraVectorSpace, domainSpace, tpetraMultiVector);
84 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
88 return tpetraMultiVector_.getNonconstObj();
92 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
96 return tpetraMultiVector_;
103 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
114 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
118 tpetraMultiVector_.getNonconstObj()->putScalar(alpha);
122 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
126 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
131 tpetraMultiVector_.getNonconstObj()->assign(*tmv);
134 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
135 tpetraMultiVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
136 tpetraMultiVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
138 tpetraMultiVector_.getNonconstObj()->sync_host ();
139 tpetraMultiVector_.getNonconstObj()->modify_host ();
141 MultiVectorDefaultBase<Scalar>::assignMultiVecImpl(mv);
146 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
150 tpetraMultiVector_.getNonconstObj()->scale(alpha);
154 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
157 const MultiVectorBase<Scalar>& mv
160 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
166 tpetraMultiVector_.getNonconstObj()->update(alpha, *tmv, ST::one());
169 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
170 tpetraMultiVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
171 tpetraMultiVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
173 tpetraMultiVector_.getNonconstObj()->sync_host ();
174 tpetraMultiVector_.getNonconstObj()->modify_host ();
176 MultiVectorDefaultBase<Scalar>::updateImpl(alpha, mv);
181 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
184 const ArrayView<
const Ptr<
const MultiVectorBase<Scalar> > >& mv,
193 typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
196 bool allCastsSuccessful =
true;
198 auto mvIter = mv.begin();
199 auto tmvIter = tmvs.begin();
200 for (; mvIter != mv.end(); ++mvIter, ++tmvIter) {
201 tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromPtr(*mvIter));
205 allCastsSuccessful =
false;
213 auto len = tmvs.
size();
215 tpetraMultiVector_.getNonconstObj()->scale(beta);
216 }
else if (len == 1 && allCastsSuccessful) {
217 tpetraMultiVector_.getNonconstObj()->update(alpha[0], *tmvs[0], beta);
218 }
else if (len == 2 && allCastsSuccessful) {
219 tpetraMultiVector_.getNonconstObj()->update(alpha[0], *tmvs[0], alpha[1], *tmvs[1], beta);
220 }
else if (allCastsSuccessful) {
222 auto tmvIter = tmvs.begin();
223 auto alphaIter = alpha.
begin();
228 for (; tmvIter != tmvs.end(); ++tmvIter) {
229 if (tmvIter->getRawPtr() == tpetraMultiVector_.getConstObj().getRawPtr()) {
236 tmvIter = tmvs.
begin();
240 if ((tmvs.size() % 2) == 0) {
241 tpetraMultiVector_.getNonconstObj()->scale(beta);
243 tpetraMultiVector_.getNonconstObj()->update(*alphaIter, *(*tmvIter), beta);
247 for (; tmvIter != tmvs.end(); tmvIter+=2, alphaIter+=2) {
248 tpetraMultiVector_.getNonconstObj()->update(
249 *alphaIter, *(*tmvIter), *(alphaIter+1), *(*(tmvIter+1)), ST::one());
253 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
254 tpetraMultiVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
255 tpetraMultiVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
257 tpetraMultiVector_.getNonconstObj()->sync_host ();
258 tpetraMultiVector_.getNonconstObj()->modify_host ();
260 MultiVectorDefaultBase<Scalar>::linearCombinationImpl(alpha, mv, beta);
265 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
267 const MultiVectorBase<Scalar>& mv,
271 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
276 tpetraMultiVector_.getConstObj()->dot(*tmv, prods);
279 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
280 tpetraMultiVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
281 tpetraMultiVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
283 tpetraMultiVector_.getNonconstObj()->sync_host ();
284 tpetraMultiVector_.getNonconstObj()->modify_host ();
286 MultiVectorDefaultBase<Scalar>::dotsImpl(mv, prods);
291 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
293 const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
296 tpetraMultiVector_.getConstObj()->norm1(norms);
300 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
302 const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
305 tpetraMultiVector_.getConstObj()->norm2(norms);
309 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
311 const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
314 tpetraMultiVector_.getConstObj()->normInf(norms);
318 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
325 return constTpetraVector<Scalar>(
327 tpetraMultiVector_->getVector(j)
332 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
339 return tpetraVector<Scalar>(
341 tpetraMultiVector_.getNonconstObj()->getVectorNonConst(j)
346 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
349 const Range1D& col_rng_in
352 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
353 std::cerr <<
"\nTpetraMultiVector::subView(Range1D) const called!\n";
355 const Range1D colRng = this->validateColRange(col_rng_in);
358 this->getConstTpetraMultiVector()->subView(colRng);
361 tpetraVectorSpace<Scalar>(
362 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
363 tpetraView->getNumVectors(),
364 tpetraView->getMap()->getComm()
368 return constTpetraMultiVector(
376 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
379 const Range1D& col_rng_in
382 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
383 std::cerr <<
"\nTpetraMultiVector::subView(Range1D) called!\n";
385 const Range1D colRng = this->validateColRange(col_rng_in);
388 this->getTpetraMultiVector()->subViewNonConst(colRng);
391 tpetraVectorSpace<Scalar>(
392 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
393 tpetraView->getNumVectors(),
394 tpetraView->getMap()->getComm()
398 return tpetraMultiVector(
406 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
412 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
413 std::cerr <<
"\nTpetraMultiVector::subView(ArrayView) const called!\n";
416 Array<std::size_t> cols(cols_in.
size());
417 for (Array<std::size_t>::size_type i = 0; i < cols.size(); ++i)
418 cols[i] = static_cast<std::size_t>(cols_in[i]);
421 this->getConstTpetraMultiVector()->subView(cols());
424 tpetraVectorSpace<Scalar>(
425 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
426 tpetraView->getNumVectors(),
427 tpetraView->getMap()->getComm()
431 return constTpetraMultiVector(
439 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
445 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
446 std::cerr <<
"\nTpetraMultiVector::subView(ArrayView) called!\n";
449 Array<std::size_t> cols(cols_in.
size());
450 for (Array<std::size_t>::size_type i = 0; i < cols.size(); ++i)
451 cols[i] = static_cast<std::size_t>(cols_in[i]);
454 this->getTpetraMultiVector()->subViewNonConst(cols());
457 tpetraVectorSpace<Scalar>(
458 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
459 tpetraView->getNumVectors(),
460 tpetraView->getMap()->getComm()
464 return tpetraMultiVector(
472 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
476 const ArrayView<
const Ptr<
const MultiVectorBase<Scalar> > > &multi_vecs,
477 const ArrayView<
const Ptr<MultiVectorBase<Scalar> > > &targ_multi_vecs,
478 const ArrayView<
const Ptr<RTOpPack::ReductTarget> > &reduct_objs,
479 const Ordinal primary_global_offset
485 for (
auto itr = multi_vecs.begin(); itr != multi_vecs.end(); ++itr) {
486 Ptr<const TMV> tmv = Teuchos::ptr_dynamic_cast<
const TMV>(*itr);
488 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
489 Teuchos::rcp_const_cast<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(
490 tmv->getConstTpetraMultiVector())->
template sync<Kokkos::HostSpace>();
492 Teuchos::rcp_const_cast<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(
493 tmv->getConstTpetraMultiVector())-> sync_host ();
499 for (
auto itr = targ_multi_vecs.begin(); itr != targ_multi_vecs.end(); ++itr) {
500 Ptr<TMV> tmv = Teuchos::ptr_dynamic_cast<TMV>(*itr);
502 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
503 tmv->getTpetraMultiVector()->template sync<Kokkos::HostSpace>();
504 tmv->getTpetraMultiVector()->template modify<Kokkos::HostSpace>();
506 tmv->getTpetraMultiVector()->sync_host ();
507 tmv->getTpetraMultiVector()->modify_host ();
512 MultiVectorAdapterBase<Scalar>::mvMultiReductApplyOpImpl(
513 primary_op, multi_vecs, targ_multi_vecs, reduct_objs, primary_global_offset);
517 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
520 const Range1D &rowRng,
521 const Range1D &colRng,
526 typedef typename Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
527 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
528 Teuchos::rcp_const_cast<TMV>(
529 tpetraMultiVector_.getConstObj())->
template sync<Kokkos::HostSpace>();
531 Teuchos::rcp_const_cast<TMV>(
532 tpetraMultiVector_.getConstObj())->sync_host ();
535 SpmdMultiVectorDefaultBase<Scalar>::
536 acquireDetachedMultiVectorViewImpl(rowRng, colRng, sub_mv);
540 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
543 const Range1D &rowRng,
544 const Range1D &colRng,
549 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
550 tpetraMultiVector_.getNonconstObj()->template sync<Kokkos::HostSpace>();
551 tpetraMultiVector_.getNonconstObj()->template modify<Kokkos::HostSpace>();
553 tpetraMultiVector_.getNonconstObj()->sync_host ();
554 tpetraMultiVector_.getNonconstObj()->modify_host ();
557 SpmdMultiVectorDefaultBase<Scalar>::
558 acquireNonconstDetachedMultiVectorViewImpl(rowRng, colRng, sub_mv);
562 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
568 SpmdMultiVectorDefaultBase<Scalar>::
569 commitNonconstDetachedMultiVectorViewImpl(sub_mv);
572 typedef typename Tpetra::MultiVector<
573 Scalar,LocalOrdinal,GlobalOrdinal,Node>::execution_space execution_space;
574 tpetraMultiVector_.getNonconstObj()->template sync<execution_space>();
625 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
629 return tpetraVectorSpace_;
633 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
635 const Ptr<ArrayRCP<Scalar> > &localValues,
const Ptr<Ordinal> &leadingDim
638 *localValues = tpetraMultiVector_.getNonconstObj()->get1dViewNonConst();
639 *leadingDim = tpetraMultiVector_->getStride();
643 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
645 const Ptr<ArrayRCP<const Scalar> > &localValues,
const Ptr<Ordinal> &leadingDim
648 *localValues = tpetraMultiVector_->get1dView();
649 *leadingDim = tpetraMultiVector_->getStride();
653 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
655 const EOpTransp M_trans,
656 const MultiVectorBase<Scalar> &X,
657 const Ptr<MultiVectorBase<Scalar> > &Y,
663 typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
671 typedef typename TMV::execution_space execution_space;
672 Teuchos::rcp_const_cast<TMV>(X_tpetra)->
template sync<execution_space>();
673 Y_tpetra->template sync<execution_space>();
674 Teuchos::rcp_const_cast<TMV>(
675 tpetraMultiVector_.getConstObj())->
template sync<execution_space>();
680 "Error, conjugation without transposition is not allowed for complex scalar types!");
698 Y_tpetra->template modify<execution_space>();
699 Y_tpetra->multiply(trans,
Teuchos::NO_TRANS, alpha, *tpetraMultiVector_.getConstObj(), *X_tpetra, beta);
701 #ifdef KOKKOS_ENABLE_DEPRECATED_CODE
702 Teuchos::rcp_const_cast<TMV>(
703 tpetraMultiVector_.getConstObj())->
template sync<Kokkos::HostSpace>();
705 Teuchos::rcp_const_cast<TMV>(
706 tpetraMultiVector_.getConstObj())->sync_host ();
708 SpmdMultiVectorDefaultBase<Scalar>::euclideanApply(M_trans, X, Y, alpha, beta);
716 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
717 template<
class TpetraMultiVector_t>
720 const RCP<
const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
731 tpetraVectorSpace_ = tpetraVectorSpace;
732 domainSpace_ = domainSpace;
733 tpetraMultiVector_.initialize(tpetraMultiVector);
734 this->updateSpmdSpace();
738 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
743 using Teuchos::rcp_dynamic_cast;
747 RCP<TMV> tmv = rcp_dynamic_cast<TMV>(mv);
749 return tmv->getTpetraMultiVector();
752 RCP<TV> tv = rcp_dynamic_cast<TV>(mv);
754 return tv->getTpetraVector();
760 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
765 using Teuchos::rcp_dynamic_cast;
771 return tmv->getConstTpetraMultiVector();
776 return tv->getConstTpetraVector();
786 #endif // THYRA_TPETRA_MULTIVECTOR_HPP
virtual void updateImpl(Scalar alpha, const MultiVectorBase< Scalar > &mv)
Concrete implementation of Thyra::MultiVector in terms of Tpetra::MultiVector.
TpetraMultiVector()
Construct to uninitialized.
RCP< const VectorBase< Scalar > > colImpl(Ordinal j) const
void getNonconstLocalMultiVectorDataImpl(const Ptr< ArrayRCP< Scalar > > &localValues, const Ptr< Ordinal > &leadingDim)
Concrete implementation of an SPMD vector space for Tpetra.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetraMultiVector()
Extract the underlying non-const Tpetra::MultiVector object.
virtual void assignImpl(Scalar alpha)
RCP< MultiVectorBase< Scalar > > nonconstContigSubViewImpl(const Range1D &colRng)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void euclideanApply(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
void initializeImpl(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< const ScalarProdVectorSpaceBase< Scalar > > &domainSpace, const RCP< TpetraMultiVector_t > &tpetraMultiVector)
RCP< const MultiVectorBase< Scalar > > nonContigSubViewImpl(const ArrayView< const int > &cols_in) const
void initialize(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< const ScalarProdVectorSpaceBase< Scalar > > &domainSpace, const RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraMultiVector)
Initialize.
#define TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE(index, lower_inclusive, upper_exclusive)
virtual void scaleImpl(Scalar alpha)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
RCP< const MultiVectorBase< Scalar > > contigSubViewImpl(const Range1D &colRng) const
virtual void mvMultiReductApplyOpImpl(const RTOpPack::RTOpT< Scalar > &primary_op, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &multi_vecs, const ArrayView< const Ptr< MultiVectorBase< Scalar > > > &targ_multi_vecs, const ArrayView< const Ptr< RTOpPack::ReductTarget > > &reduct_objs, const Ordinal primary_global_offset) const
RCP< MultiVectorBase< Scalar > > nonconstNonContigSubViewImpl(const ArrayView< const int > &cols_in)
void constInitialize(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< const ScalarProdVectorSpaceBase< Scalar > > &domainSpace, const RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraMultiVector)
Initialize.
virtual void norms2Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
virtual void linearCombinationImpl(const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta)
Concrete Thyra::SpmdVectorBase using Tpetra::Vector.
virtual void dotsImpl(const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
RCP< VectorBase< Scalar > > nonconstColImpl(Ordinal j)
bool nonnull(const boost::shared_ptr< T > &p)
void getLocalMultiVectorDataImpl(const Ptr< ArrayRCP< const Scalar > > &localValues, const Ptr< Ordinal > &leadingDim) const
virtual void assignMultiVecImpl(const MultiVectorBase< Scalar > &mv)
virtual void norms1Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraMultiVector() const
Extract the underlying const Tpetra::MultiVector object.
RCP< const SpmdVectorSpaceBase< Scalar > > spmdSpaceImpl() const
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
virtual void normsInfImpl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
RCP< const ScalarProdVectorSpaceBase< Scalar > > domainScalarProdVecSpc() const
void acquireNonconstDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
void acquireDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
void commitNonconstDetachedMultiVectorViewImpl(RTOpPack::SubMultiVectorView< Scalar > *sub_mv)