10#ifndef THYRA_TPETRA_MULTIVECTOR_HPP
11#define THYRA_TPETRA_MULTIVECTOR_HPP
13#include "Thyra_TpetraMultiVector_decl.hpp"
14#include "Thyra_TpetraVectorSpace.hpp"
15#include "Thyra_TpetraVector.hpp"
16#include "Teuchos_Assert.hpp"
25template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
30template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
41template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
52template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
56 return tpetraMultiVector_.getNonconstObj();
60template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
64 return tpetraMultiVector_;
71template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
82template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
86 tpetraMultiVector_.getNonconstObj()->putScalar(alpha);
90template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
99 tpetraMultiVector_.getNonconstObj()->assign(*tmv);
106template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
110 tpetraMultiVector_.getNonconstObj()->scale(alpha);
114template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
126 tpetraMultiVector_.getNonconstObj()->update(alpha, *tmv, ST::one());
133template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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) {
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());
209template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
220 tpetraMultiVector_.getConstObj()->dot(*tmv, prods);
227template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
232 tpetraMultiVector_.getConstObj()->norm1(
norms);
236template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
241 tpetraMultiVector_.getConstObj()->norm2(
norms);
245template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
250 tpetraMultiVector_.getConstObj()->normInf(
norms);
254template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
263 tpetraMultiVector_->getVector(j)
268template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
277 tpetraMultiVector_.getNonconstObj()->getVectorNonConst(j)
282template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
288#ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
289 std::cerr <<
"\nTpetraMultiVector::subView(Range1D) const called!\n";
298 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
299 tpetraView->getNumVectors(),
300 tpetraView->getMap()->getComm()
312template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
318#ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
319 std::cerr <<
"\nTpetraMultiVector::subView(Range1D) called!\n";
328 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
329 tpetraView->getNumVectors(),
330 tpetraView->getMap()->getComm()
342template <
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";
354 cols[i] =
static_cast<std::size_t
>(cols_in[i]);
361 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
362 tpetraView->getNumVectors(),
363 tpetraView->getMap()->getComm()
375template <
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";
387 cols[i] =
static_cast<std::size_t
>(cols_in[i]);
394 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
395 tpetraView->getNumVectors(),
396 tpetraView->getMap()->getComm()
408template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
415 const Ordinal primary_global_offset
420 primary_op, multi_vecs, targ_multi_vecs, reduct_objs, primary_global_offset);
424template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
437template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
450template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
509template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
513 return tpetraVectorSpace_;
517template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
522 *localValues = tpetraMultiVector_.getNonconstObj()->get1dViewNonConst();
523 *leadingDim = tpetraMultiVector_->getStride();
527template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
532 *localValues = tpetraMultiVector_->get1dView();
533 *leadingDim = tpetraMultiVector_->getStride();
537template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
539 const EOpTransp M_trans,
547 typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
553 if (nonnull(X_tpetra) && nonnull(Y_tpetra)) {
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);
586template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
587template<
class TpetraMultiVector_t>
588void TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::initializeImpl(
602 domainSpace_ = domainSpace;
604 this->updateSpmdSpace();
608template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
609RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,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();
627 return Teuchos::null;
630template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
635 using Teuchos::rcp_dynamic_cast;
636 typedef Thyra::TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
637 typedef Thyra::TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
641 return tmv->getConstTpetraMultiVector();
646 return tv->getConstTpetraVector();
649 return Teuchos::null;
RCP< const VectorSpaceBase< Scalar > > domain() const
Returns this->domainScalarProdVecSpc().
Interface for a collection of column vectors called a multi-vector.
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 =0
Apply a reduction/transformation operator column by column and return an array of the reduction objec...
void norms(const MultiVectorBase< Scalar > &V, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms)
Column-wise multi-vector natural norm.
virtual void dotsImpl(const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
Default implementation of dots using RTOps.
virtual void linearCombinationImpl(const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta)
Default implementation of linear_combination using RTOps.
virtual void assignMultiVecImpl(const MultiVectorBase< Scalar > &mv)
Default implementation of assign(MV) using RTOps.
virtual void updateImpl(Scalar alpha, const MultiVectorBase< Scalar > &mv)
Default implementation of update using RTOps.
Range1D validateColRange(const Range1D &rowCol) const
Validate and resize the column range.
void acquireDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
void euclideanApply(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
Uses GEMM() and Teuchos::reduceAll() to implement.
void acquireNonconstDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
void commitNonconstDetachedMultiVectorViewImpl(RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
Concrete implementation of Thyra::MultiVector in terms of Tpetra::MultiVector.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetraMultiVector()
Extract the underlying non-const Tpetra::MultiVector object.
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.
RCP< const SpmdVectorSpaceBase< Scalar > > spmdSpaceImpl() const
RCP< const ScalarProdVectorSpaceBase< Scalar > > domainScalarProdVecSpc() const
RCP< const TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > constTpetraMultiVector(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)
Nonmember constructor for const TpetraMultiVector.
RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraMultiVector() const
Extract the underlying const Tpetra::MultiVector object.
virtual void norms2Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
void commitNonconstDetachedMultiVectorViewImpl(RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
virtual void updateImpl(Scalar alpha, const MultiVectorBase< Scalar > &mv)
void acquireNonconstDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
RCP< const VectorBase< Scalar > > colImpl(Ordinal j) const
void acquireDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
void getLocalMultiVectorDataImpl(const Ptr< ArrayRCP< const Scalar > > &localValues, const Ptr< Ordinal > &leadingDim) const
RCP< MultiVectorBase< Scalar > > nonconstNonContigSubViewImpl(const ArrayView< const int > &cols_in)
RCP< const MultiVectorBase< Scalar > > nonContigSubViewImpl(const ArrayView< const int > &cols_in) const
RCP< MultiVectorBase< Scalar > > nonconstContigSubViewImpl(const Range1D &colRng)
RCP< const MultiVectorBase< Scalar > > contigSubViewImpl(const Range1D &colRng) const
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 norms1Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
TpetraMultiVector()
Construct to uninitialized.
RCP< VectorBase< Scalar > > nonconstColImpl(Ordinal j)
virtual void dotsImpl(const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
virtual void euclideanApply(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
virtual void assignMultiVecImpl(const MultiVectorBase< Scalar > &mv)
virtual void scaleImpl(Scalar alpha)
void getNonconstLocalMultiVectorDataImpl(const Ptr< ArrayRCP< Scalar > > &localValues, const Ptr< Ordinal > &leadingDim)
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< TpetraMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > tpetraMultiVector(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< const ScalarProdVectorSpaceBase< Scalar > > &domainSpace, const RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraMultiVector)
Nonmember constructor for non-const TpetraMultiVector.
virtual void assignImpl(Scalar alpha)
virtual void linearCombinationImpl(const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta)
virtual void normsInfImpl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
Concrete implementation of an SPMD vector space for Tpetra.
RCP< TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > tpetraVectorSpace(const RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &tpetraMap)
Nonmember constructor that creats a serial vector space.
Concrete Thyra::SpmdVectorBase using Tpetra::Vector.
RCP< const TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > constTpetraVector(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector)
Nonmember constructor for TpetraVector.
RCP< TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > tpetraVector(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector)
Nonmember constructor for TpetraVector.
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE(index, lower_inclusive, upper_exclusive)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
bool nonnull(const std::shared_ptr< T > &p)
NOTRANS
Type for the dimension of a vector space. `**/ typedef Teuchos::Ordinal Ordinal;.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)