10 #ifndef THYRA_TPETRA_MULTIVECTOR_HPP
11 #define THYRA_TPETRA_MULTIVECTOR_HPP
25 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
30 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
33 const RCP<
const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
34 const RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraMultiVector
37 initializeImpl(tpetraVectorSpace, domainSpace, tpetraMultiVector);
41 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
44 const RCP<
const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
45 const RCP<
const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraMultiVector
48 initializeImpl(tpetraVectorSpace, domainSpace, tpetraMultiVector);
52 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
56 return tpetraMultiVector_.getNonconstObj();
60 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
64 return tpetraMultiVector_;
71 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
82 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
86 tpetraMultiVector_.getNonconstObj()->putScalar(alpha);
90 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
94 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
99 tpetraMultiVector_.getNonconstObj()->assign(*tmv);
101 MultiVectorDefaultBase<Scalar>::assignMultiVecImpl(mv);
106 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
110 tpetraMultiVector_.getNonconstObj()->scale(alpha);
114 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
117 const MultiVectorBase<Scalar>& mv
120 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
126 tpetraMultiVector_.getNonconstObj()->update(alpha, *tmv, ST::one());
128 MultiVectorDefaultBase<Scalar>::updateImpl(alpha, mv);
133 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
136 const ArrayView<
const Ptr<
const MultiVectorBase<Scalar> > >& mv,
145 typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
148 bool allCastsSuccessful =
true;
150 auto mvIter = mv.begin();
151 auto tmvIter = tmvs.begin();
152 for (; mvIter != mv.end(); ++mvIter, ++tmvIter) {
153 tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromPtr(*mvIter));
157 allCastsSuccessful =
false;
165 auto len = tmvs.
size();
167 tpetraMultiVector_.getNonconstObj()->scale(beta);
168 }
else if (len == 1 && allCastsSuccessful) {
169 tpetraMultiVector_.getNonconstObj()->update(alpha[0], *tmvs[0], beta);
170 }
else if (len == 2 && allCastsSuccessful) {
171 tpetraMultiVector_.getNonconstObj()->update(alpha[0], *tmvs[0], alpha[1], *tmvs[1], beta);
172 }
else if (allCastsSuccessful) {
174 auto tmvIter = tmvs.begin();
175 auto alphaIter = alpha.
begin();
180 for (; tmvIter != tmvs.end(); ++tmvIter) {
181 if (tmvIter->getRawPtr() == tpetraMultiVector_.getConstObj().getRawPtr()) {
188 tmvIter = tmvs.
begin();
192 if ((tmvs.size() % 2) == 0) {
193 tpetraMultiVector_.getNonconstObj()->scale(beta);
195 tpetraMultiVector_.getNonconstObj()->update(*alphaIter, *(*tmvIter), beta);
199 for (; tmvIter != tmvs.end(); tmvIter+=2, alphaIter+=2) {
200 tpetraMultiVector_.getNonconstObj()->update(
201 *alphaIter, *(*tmvIter), *(alphaIter+1), *(*(tmvIter+1)), ST::one());
204 MultiVectorDefaultBase<Scalar>::linearCombinationImpl(alpha, mv, beta);
209 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
211 const MultiVectorBase<Scalar>& mv,
215 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
220 tpetraMultiVector_.getConstObj()->dot(*tmv, prods);
222 MultiVectorDefaultBase<Scalar>::dotsImpl(mv, prods);
227 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
229 const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
232 tpetraMultiVector_.getConstObj()->norm1(norms);
236 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
238 const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
241 tpetraMultiVector_.getConstObj()->norm2(norms);
245 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
247 const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
250 tpetraMultiVector_.getConstObj()->normInf(norms);
254 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
261 return constTpetraVector<Scalar>(
263 tpetraMultiVector_->getVector(j)
268 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
275 return tpetraVector<Scalar>(
277 tpetraMultiVector_.getNonconstObj()->getVectorNonConst(j)
282 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
285 const Range1D& col_rng_in
288 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
289 std::cerr <<
"\nTpetraMultiVector::subView(Range1D) const called!\n";
291 const Range1D colRng = this->validateColRange(col_rng_in);
294 this->getConstTpetraMultiVector()->subView(colRng);
297 tpetraVectorSpace<Scalar>(
298 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
299 tpetraView->getNumVectors(),
300 tpetraView->getMap()->getComm()
304 return constTpetraMultiVector(
312 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
315 const Range1D& col_rng_in
318 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
319 std::cerr <<
"\nTpetraMultiVector::subView(Range1D) called!\n";
321 const Range1D colRng = this->validateColRange(col_rng_in);
324 this->getTpetraMultiVector()->subViewNonConst(colRng);
327 tpetraVectorSpace<Scalar>(
328 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
329 tpetraView->getNumVectors(),
330 tpetraView->getMap()->getComm()
334 return tpetraMultiVector(
342 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
348 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
349 std::cerr <<
"\nTpetraMultiVector::subView(ArrayView) const called!\n";
352 Array<std::size_t> cols(cols_in.
size());
353 for (Array<std::size_t>::size_type i = 0; i < cols.size(); ++i)
354 cols[i] = static_cast<std::size_t>(cols_in[i]);
357 this->getConstTpetraMultiVector()->subView(cols());
360 tpetraVectorSpace<Scalar>(
361 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
362 tpetraView->getNumVectors(),
363 tpetraView->getMap()->getComm()
367 return constTpetraMultiVector(
375 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
381 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
382 std::cerr <<
"\nTpetraMultiVector::subView(ArrayView) called!\n";
385 Array<std::size_t> cols(cols_in.
size());
386 for (Array<std::size_t>::size_type i = 0; i < cols.size(); ++i)
387 cols[i] = static_cast<std::size_t>(cols_in[i]);
390 this->getTpetraMultiVector()->subViewNonConst(cols());
393 tpetraVectorSpace<Scalar>(
394 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
395 tpetraView->getNumVectors(),
396 tpetraView->getMap()->getComm()
400 return tpetraMultiVector(
408 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
412 const ArrayView<
const Ptr<
const MultiVectorBase<Scalar> > > &multi_vecs,
413 const ArrayView<
const Ptr<MultiVectorBase<Scalar> > > &targ_multi_vecs,
414 const ArrayView<
const Ptr<RTOpPack::ReductTarget> > &reduct_objs,
415 const Ordinal primary_global_offset
419 MultiVectorAdapterBase<Scalar>::mvMultiReductApplyOpImpl(
420 primary_op, multi_vecs, targ_multi_vecs, reduct_objs, primary_global_offset);
424 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
427 const Range1D &rowRng,
428 const Range1D &colRng,
432 SpmdMultiVectorDefaultBase<Scalar>::
433 acquireDetachedMultiVectorViewImpl(rowRng, colRng, sub_mv);
437 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
440 const Range1D &rowRng,
441 const Range1D &colRng,
445 SpmdMultiVectorDefaultBase<Scalar>::
446 acquireNonconstDetachedMultiVectorViewImpl(rowRng, colRng, sub_mv);
450 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
456 SpmdMultiVectorDefaultBase<Scalar>::
457 commitNonconstDetachedMultiVectorViewImpl(sub_mv);
509 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
513 return tpetraVectorSpace_;
517 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
519 const Ptr<ArrayRCP<Scalar> > &localValues,
const Ptr<Ordinal> &leadingDim
522 *localValues = tpetraMultiVector_.getNonconstObj()->get1dViewNonConst();
523 *leadingDim = tpetraMultiVector_->getStride();
527 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
529 const Ptr<ArrayRCP<const Scalar> > &localValues,
const Ptr<Ordinal> &leadingDim
532 *localValues = tpetraMultiVector_->get1dView();
533 *leadingDim = tpetraMultiVector_->getStride();
537 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
539 const EOpTransp M_trans,
540 const MultiVectorBase<Scalar> &X,
541 const Ptr<MultiVectorBase<Scalar> > &Y,
547 typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
557 "Error, conjugation without transposition is not allowed for complex scalar types!");
575 Y_tpetra->multiply(trans,
Teuchos::NO_TRANS, alpha, *tpetraMultiVector_.getConstObj(), *X_tpetra, beta);
578 SpmdMultiVectorDefaultBase<Scalar>::euclideanApply(M_trans, X, Y, alpha, beta);
586 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
587 template<
class TpetraMultiVector_t>
590 const RCP<
const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
601 tpetraVectorSpace_ = tpetraVectorSpace;
602 domainSpace_ = domainSpace;
603 tpetraMultiVector_.initialize(tpetraMultiVector);
604 this->updateSpmdSpace();
608 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
613 using Teuchos::rcp_dynamic_cast;
617 RCP<TMV> tmv = rcp_dynamic_cast<TMV>(mv);
619 return tmv->getTpetraMultiVector();
622 RCP<TV> tv = rcp_dynamic_cast<TV>(mv);
624 return tv->getTpetraVector();
630 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
635 using Teuchos::rcp_dynamic_cast;
641 return tmv->getConstTpetraMultiVector();
646 return tv->getConstTpetraVector();
656 #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)