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     tpetraMultiVector_.getNonconstObj()->sync_host ();
 
  135     tpetraMultiVector_.getNonconstObj()->modify_host ();
 
  136     MultiVectorDefaultBase<Scalar>::assignMultiVecImpl(mv);
 
  141 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  145   tpetraMultiVector_.getNonconstObj()->scale(alpha);
 
  149 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  152   const MultiVectorBase<Scalar>& mv
 
  155   auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
 
  161     tpetraMultiVector_.getNonconstObj()->update(alpha, *tmv, ST::one());
 
  164     tpetraMultiVector_.getNonconstObj()->sync_host ();
 
  165     tpetraMultiVector_.getNonconstObj()->modify_host ();
 
  166     MultiVectorDefaultBase<Scalar>::updateImpl(alpha, mv);
 
  171 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  174   const ArrayView<
const Ptr<
const MultiVectorBase<Scalar> > >& mv,
 
  183   typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
 
  186   bool allCastsSuccessful = 
true;
 
  188     auto mvIter = mv.begin();
 
  189     auto tmvIter = tmvs.begin();
 
  190     for (; mvIter != mv.end(); ++mvIter, ++tmvIter) {
 
  191       tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromPtr(*mvIter));
 
  195         allCastsSuccessful = 
false;
 
  203   auto len = tmvs.
size();
 
  205     tpetraMultiVector_.getNonconstObj()->scale(beta);
 
  206   } 
else if (len == 1 && allCastsSuccessful) {
 
  207     tpetraMultiVector_.getNonconstObj()->update(alpha[0], *tmvs[0], beta);
 
  208   } 
else if (len == 2 && allCastsSuccessful) {
 
  209     tpetraMultiVector_.getNonconstObj()->update(alpha[0], *tmvs[0], alpha[1], *tmvs[1], beta);
 
  210   } 
else if (allCastsSuccessful) {
 
  212     auto tmvIter = tmvs.begin();
 
  213     auto alphaIter = alpha.
begin();
 
  218     for (; tmvIter != tmvs.end(); ++tmvIter) {
 
  219       if (tmvIter->getRawPtr() == tpetraMultiVector_.getConstObj().getRawPtr()) {
 
  226     tmvIter = tmvs.
begin();
 
  230     if ((tmvs.size() % 2) == 0) {
 
  231       tpetraMultiVector_.getNonconstObj()->scale(beta);
 
  233       tpetraMultiVector_.getNonconstObj()->update(*alphaIter, *(*tmvIter), beta);
 
  237     for (; tmvIter != tmvs.end(); tmvIter+=2, alphaIter+=2) {
 
  238       tpetraMultiVector_.getNonconstObj()->update(
 
  239         *alphaIter, *(*tmvIter), *(alphaIter+1), *(*(tmvIter+1)), ST::one());
 
  243     tpetraMultiVector_.getNonconstObj()->sync_host ();
 
  244     tpetraMultiVector_.getNonconstObj()->modify_host ();
 
  245     MultiVectorDefaultBase<Scalar>::linearCombinationImpl(alpha, mv, beta);
 
  250 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  252     const MultiVectorBase<Scalar>& mv,
 
  256   auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
 
  261     tpetraMultiVector_.getConstObj()->dot(*tmv, prods);
 
  264     tpetraMultiVector_.getNonconstObj()->sync_host ();
 
  265     tpetraMultiVector_.getNonconstObj()->modify_host ();
 
  266     MultiVectorDefaultBase<Scalar>::dotsImpl(mv, prods);
 
  271 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  273   const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
 
  276   tpetraMultiVector_.getConstObj()->norm1(norms);
 
  280 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  282     const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
 
  285   tpetraMultiVector_.getConstObj()->norm2(norms);
 
  289 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  291   const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
 
  294   tpetraMultiVector_.getConstObj()->normInf(norms);
 
  298 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  305   return constTpetraVector<Scalar>(
 
  307     tpetraMultiVector_->getVector(j)
 
  312 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  319   return tpetraVector<Scalar>(
 
  321     tpetraMultiVector_.getNonconstObj()->getVectorNonConst(j)
 
  326 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  329   const Range1D& col_rng_in
 
  332 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT 
  333   std::cerr << 
"\nTpetraMultiVector::subView(Range1D) const called!\n";
 
  335   const Range1D colRng = this->validateColRange(col_rng_in);
 
  338     this->getConstTpetraMultiVector()->subView(colRng);
 
  341     tpetraVectorSpace<Scalar>(
 
  342         Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
 
  343           tpetraView->getNumVectors(),
 
  344           tpetraView->getMap()->getComm()
 
  348   return constTpetraMultiVector(
 
  356 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  359   const Range1D& col_rng_in
 
  362 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT 
  363   std::cerr << 
"\nTpetraMultiVector::subView(Range1D) called!\n";
 
  365   const Range1D colRng = this->validateColRange(col_rng_in);
 
  368     this->getTpetraMultiVector()->subViewNonConst(colRng);
 
  371     tpetraVectorSpace<Scalar>(
 
  372         Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
 
  373           tpetraView->getNumVectors(),
 
  374           tpetraView->getMap()->getComm()
 
  378   return tpetraMultiVector(
 
  386 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  392 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT 
  393   std::cerr << 
"\nTpetraMultiVector::subView(ArrayView) const called!\n";
 
  396   Array<std::size_t> cols(cols_in.
size());
 
  397   for (Array<std::size_t>::size_type i = 0; i < cols.size(); ++i)
 
  398     cols[i] = static_cast<std::size_t>(cols_in[i]);
 
  401     this->getConstTpetraMultiVector()->subView(cols());
 
  404     tpetraVectorSpace<Scalar>(
 
  405         Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
 
  406           tpetraView->getNumVectors(),
 
  407           tpetraView->getMap()->getComm()
 
  411   return constTpetraMultiVector(
 
  419 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  425 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT 
  426   std::cerr << 
"\nTpetraMultiVector::subView(ArrayView) called!\n";
 
  429   Array<std::size_t> cols(cols_in.
size());
 
  430   for (Array<std::size_t>::size_type i = 0; i < cols.size(); ++i)
 
  431     cols[i] = static_cast<std::size_t>(cols_in[i]);
 
  434     this->getTpetraMultiVector()->subViewNonConst(cols());
 
  437     tpetraVectorSpace<Scalar>(
 
  438         Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
 
  439           tpetraView->getNumVectors(),
 
  440           tpetraView->getMap()->getComm()
 
  444   return tpetraMultiVector(
 
  452 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  456   const ArrayView<
const Ptr<
const MultiVectorBase<Scalar> > > &multi_vecs,
 
  457   const ArrayView<
const Ptr<MultiVectorBase<Scalar> > > &targ_multi_vecs,
 
  458   const ArrayView<
const Ptr<RTOpPack::ReductTarget> > &reduct_objs,
 
  459   const Ordinal primary_global_offset
 
  465   for (
auto itr = multi_vecs.begin(); itr != multi_vecs.end(); ++itr) {
 
  466     Ptr<const TMV> tmv = Teuchos::ptr_dynamic_cast<
const TMV>(*itr);
 
  468       Teuchos::rcp_const_cast<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(
 
  469       tmv->getConstTpetraMultiVector())-> sync_host ();
 
  474   for (
auto itr = targ_multi_vecs.begin(); itr != targ_multi_vecs.end(); ++itr) {
 
  475     Ptr<TMV> tmv = Teuchos::ptr_dynamic_cast<TMV>(*itr);
 
  477       tmv->getTpetraMultiVector()->sync_host ();
 
  478       tmv->getTpetraMultiVector()->modify_host ();
 
  482   MultiVectorAdapterBase<Scalar>::mvMultiReductApplyOpImpl(
 
  483     primary_op, multi_vecs, targ_multi_vecs, reduct_objs, primary_global_offset);
 
  487 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  490   const Range1D &rowRng,
 
  491   const Range1D &colRng,
 
  496   typedef typename Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
 
  497   Teuchos::rcp_const_cast<TMV>(
 
  498     tpetraMultiVector_.getConstObj())->sync_host ();
 
  500   SpmdMultiVectorDefaultBase<Scalar>::
 
  501     acquireDetachedMultiVectorViewImpl(rowRng, colRng, sub_mv);
 
  505 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  508   const Range1D &rowRng,
 
  509   const Range1D &colRng,
 
  514   tpetraMultiVector_.getNonconstObj()->sync_host ();
 
  515   tpetraMultiVector_.getNonconstObj()->modify_host ();
 
  517   SpmdMultiVectorDefaultBase<Scalar>::
 
  518     acquireNonconstDetachedMultiVectorViewImpl(rowRng, colRng, sub_mv);
 
  522 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  528   SpmdMultiVectorDefaultBase<Scalar>::
 
  529     commitNonconstDetachedMultiVectorViewImpl(sub_mv);
 
  532   typedef typename Tpetra::MultiVector<
 
  533     Scalar,LocalOrdinal,GlobalOrdinal,Node>::execution_space execution_space;
 
  534   tpetraMultiVector_.getNonconstObj()->template sync<execution_space>();
 
  585 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  589   return tpetraVectorSpace_;
 
  593 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  595   const Ptr<ArrayRCP<Scalar> > &localValues, 
const Ptr<Ordinal> &leadingDim
 
  598   *localValues = tpetraMultiVector_.getNonconstObj()->get1dViewNonConst();
 
  599   *leadingDim = tpetraMultiVector_->getStride();
 
  603 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  605   const Ptr<ArrayRCP<const Scalar> > &localValues, 
const Ptr<Ordinal> &leadingDim
 
  608   *localValues = tpetraMultiVector_->get1dView();
 
  609   *leadingDim = tpetraMultiVector_->getStride();
 
  613 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  615   const EOpTransp M_trans,
 
  616   const MultiVectorBase<Scalar> &X,
 
  617   const Ptr<MultiVectorBase<Scalar> > &Y,
 
  623   typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
 
  631     typedef typename TMV::execution_space execution_space;
 
  632     Teuchos::rcp_const_cast<TMV>(X_tpetra)->
template sync<execution_space>();
 
  633     Y_tpetra->template sync<execution_space>();
 
  634     Teuchos::rcp_const_cast<TMV>(
 
  635       tpetraMultiVector_.getConstObj())->
template sync<execution_space>();
 
  640       "Error, conjugation without transposition is not allowed for complex scalar types!");
 
  658     Y_tpetra->template modify<execution_space>();
 
  659     Y_tpetra->multiply(trans, 
Teuchos::NO_TRANS, alpha, *tpetraMultiVector_.getConstObj(), *X_tpetra, beta);
 
  661     Teuchos::rcp_const_cast<TMV>(
 
  662       tpetraMultiVector_.getConstObj())->sync_host ();
 
  663     SpmdMultiVectorDefaultBase<Scalar>::euclideanApply(M_trans, X, Y, alpha, beta);
 
  671 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  672 template<
class TpetraMultiVector_t>
 
  675   const RCP<
const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
 
  686   tpetraVectorSpace_ = tpetraVectorSpace;
 
  687   domainSpace_ = domainSpace;
 
  688   tpetraMultiVector_.initialize(tpetraMultiVector);
 
  689   this->updateSpmdSpace();
 
  693 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  698   using Teuchos::rcp_dynamic_cast;
 
  702   RCP<TMV> tmv = rcp_dynamic_cast<TMV>(mv);
 
  704     return tmv->getTpetraMultiVector();
 
  707   RCP<TV> tv = rcp_dynamic_cast<TV>(mv);
 
  709     return tv->getTpetraVector();
 
  715 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
  720   using Teuchos::rcp_dynamic_cast;
 
  726     return tmv->getConstTpetraMultiVector();
 
  731     return tv->getConstTpetraVector();
 
  741 #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)