42 #ifndef THYRA_TPETRA_VECTOR_HPP
43 #define THYRA_TPETRA_VECTOR_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,Node>(
105 tpetraVector_.getConstObj()->getMap()->getComm()
116 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
120 return tpetraVectorSpace_;
127 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
129 const Ptr<ArrayRCP<Scalar> > &localValues )
131 *localValues = tpetraVector_.getNonconstObj()->get1dViewNonConst();
135 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
137 const Ptr<ArrayRCP<const Scalar> > &localValues )
const
139 *localValues = tpetraVector_->get1dView();
146 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
154 if (!tpetraVector_.getNonconstObj()->isDistributed()) {
155 auto comm = tpetraVector_.getNonconstObj()->getMap()->getComm();
156 if (tpetraVector_.getConstObj()->getMap()->getComm()->getRank() == 0)
157 tpetraVector_.getNonconstObj()->randomize(l, u);
160 tpetraVector_.getNonconstObj()->reduce();
162 tpetraVector_.getNonconstObj()->randomize(l, u);
167 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
169 const VectorBase<Scalar>& x
172 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
177 tpetraVector_.getNonconstObj()->abs(*tx);
180 tpetraVector_.getNonconstObj()->sync_host ();
181 tpetraVector_.getNonconstObj()->modify_host ();
182 VectorDefaultBase<Scalar>::absImpl(x);
187 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
189 const VectorBase<Scalar>& x
192 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
197 tpetraVector_.getNonconstObj()->reciprocal(*tx);
200 tpetraVector_.getNonconstObj()->sync_host ();
201 tpetraVector_.getNonconstObj()->modify_host ();
202 VectorDefaultBase<Scalar>::reciprocalImpl(x);
207 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
209 const VectorBase<Scalar>& x
212 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
218 tpetraVector_.getNonconstObj()->elementWiseMultiply(
219 ST::one(), *tx, *tpetraVector_.getConstObj(), ST::zero());
222 tpetraVector_.getNonconstObj()->sync_host ();
223 tpetraVector_.getNonconstObj()->modify_host ();
224 VectorDefaultBase<Scalar>::eleWiseScaleImpl(x);
229 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
232 const VectorBase<Scalar>& x
235 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
243 = Tpetra::createVector<Scalar>(tx->getMap());
244 temp->elementWiseMultiply(
245 ST::one(), *tx, *tpetraVector_.getConstObj(), ST::zero());
246 return ST::magnitude(ST::squareroot(tpetraVector_.getConstObj()->dot(*temp)));
249 tpetraVector_.getNonconstObj()->sync_host ();
250 return VectorDefaultBase<Scalar>::norm2WeightedImpl(x);
255 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
258 const ArrayView<
const Ptr<
const VectorBase<Scalar> > > &vecs,
259 const ArrayView<
const Ptr<VectorBase<Scalar> > > &targ_vecs,
260 const Ptr<RTOpPack::ReductTarget> &reduct_obj,
265 for (
auto itr = vecs.begin(); itr != vecs.end(); ++itr) {
266 auto tv = this->getConstTpetraVector(Teuchos::rcpFromPtr(*itr));
268 typedef Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
269 Teuchos::rcp_const_cast<TV>(tv)->sync_host ();
274 for (
auto itr = targ_vecs.begin(); itr != targ_vecs.end(); ++itr) {
275 auto tv = this->getTpetraVector(Teuchos::rcpFromPtr(*itr));
282 SpmdVectorDefaultBase<Scalar>::applyOpImpl(op, vecs, targ_vecs, reduct_obj, global_offset);
286 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
294 typedef typename Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
295 Teuchos::rcp_const_cast<TV>(
296 tpetraVector_.getConstObj())->sync_host ();
298 SpmdVectorDefaultBase<Scalar>::acquireDetachedVectorViewImpl(rng, sub_vec);
302 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
310 tpetraVector_.getNonconstObj()->sync_host ();
311 tpetraVector_.getNonconstObj()->modify_host ();
313 SpmdVectorDefaultBase<Scalar>::acquireNonconstDetachedVectorViewImpl(rng, sub_vec);
317 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
323 SpmdVectorDefaultBase<Scalar>::commitNonconstDetachedVectorViewImpl(sub_vec);
326 typedef typename Tpetra::Vector<
327 Scalar,LocalOrdinal,GlobalOrdinal,Node>::execution_space execution_space;
328 tpetraVector_.getNonconstObj()->template sync<execution_space>();
335 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
338 tpetraVector_.getNonconstObj()->putScalar(alpha);
342 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
346 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
351 tpetraVector_.getNonconstObj()->assign(*tmv);
354 tpetraVector_.getNonconstObj()->sync_host ();
355 tpetraVector_.getNonconstObj()->modify_host ();
356 MultiVectorDefaultBase<Scalar>::assignMultiVecImpl(mv);
361 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
364 tpetraVector_.getNonconstObj()->scale(alpha);
368 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
371 const MultiVectorBase<Scalar>& mv
374 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
380 tpetraVector_.getNonconstObj()->update(alpha, *tmv, ST::one());
383 tpetraVector_.getNonconstObj()->sync_host();
384 tpetraVector_.getNonconstObj()->modify_host();
385 MultiVectorDefaultBase<Scalar>::updateImpl(alpha, mv);
390 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
393 const ArrayView<
const Ptr<
const MultiVectorBase<Scalar> > >& mv,
404 bool allCastsSuccessful =
true;
406 auto mvIter = mv.begin();
407 auto tmvIter = tmvs.begin();
408 for (; mvIter != mv.end(); ++mvIter, ++tmvIter) {
409 tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromPtr(*mvIter));
413 allCastsSuccessful =
false;
421 auto len = mv.size();
423 tpetraVector_.getNonconstObj()->scale(beta);
424 }
else if (len == 1 && allCastsSuccessful) {
425 tpetraVector_.getNonconstObj()->update(alpha[0], *tmvs[0], beta);
426 }
else if (len == 2 && allCastsSuccessful) {
427 tpetraVector_.getNonconstObj()->update(alpha[0], *tmvs[0], alpha[1], *tmvs[1], beta);
428 }
else if (allCastsSuccessful) {
430 auto tmvIter = tmvs.begin();
431 auto alphaIter = alpha.
begin();
436 for (; tmvIter != tmvs.end(); ++tmvIter) {
437 if (tmvIter->getRawPtr() == tpetraVector_.getConstObj().getRawPtr()) {
445 tmvIter = tmvs.begin();
449 if ((tmvs.size() % 2) == 0) {
450 tpetraVector_.getNonconstObj()->scale(beta);
452 tpetraVector_.getNonconstObj()->update(*alphaIter, *(*tmvIter), beta);
456 for (; tmvIter != tmvs.end(); tmvIter+=2, alphaIter+=2) {
457 tpetraVector_.getNonconstObj()->update(
458 *alphaIter, *(*tmvIter), *(alphaIter+1), *(*(tmvIter+1)), ST::one());
462 tpetraVector_.getNonconstObj()->sync_host ();
463 tpetraVector_.getNonconstObj()->modify_host ();
464 MultiVectorDefaultBase<Scalar>::linearCombinationImpl(alpha, mv, beta);
469 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
471 const MultiVectorBase<Scalar>& mv,
475 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
480 tpetraVector_.getConstObj()->dot(*tmv, prods);
483 tpetraVector_.getNonconstObj()->sync_host ();
484 tpetraVector_.getNonconstObj()->modify_host ();
485 MultiVectorDefaultBase<Scalar>::dotsImpl(mv, prods);
490 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
492 const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
495 tpetraVector_.getConstObj()->norm1(norms);
499 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
501 const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
504 tpetraVector_.getConstObj()->norm2(norms);
508 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
510 const ArrayView<
typename ScalarTraits<Scalar>::magnitudeType>& norms
513 tpetraVector_.getConstObj()->normInf(norms);
517 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
519 const EOpTransp M_trans,
520 const MultiVectorBase<Scalar> &X,
521 const Ptr<MultiVectorBase<Scalar> > &Y,
527 typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
528 typedef Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
536 typedef typename TMV::execution_space execution_space;
537 Teuchos::rcp_const_cast<TMV>(X_tpetra)->
template sync<execution_space>();
538 Y_tpetra->template sync<execution_space>();
539 Teuchos::rcp_const_cast<TV>(tpetraVector_.getConstObj())->
template sync<execution_space>();
544 "Error, conjugation without transposition is not allowed for complex scalar types!");
562 Y_tpetra->template modify<execution_space>();
563 Y_tpetra->multiply(trans,
Teuchos::NO_TRANS, alpha, *tpetraVector_.getConstObj(), *X_tpetra, beta);
565 Teuchos::rcp_const_cast<TV>(tpetraVector_.getConstObj())->sync_host ();
566 VectorDefaultBase<Scalar>::applyImpl(M_trans, X, Y, alpha, beta);
575 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
576 template<
class TpetraVector_t>
586 tpetraVectorSpace_ = tpetraVectorSpace;
587 tpetraVector_.initialize(tpetraVector);
588 this->updateSpmdSpace();
592 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
597 using Teuchos::rcp_dynamic_cast;
601 RCP<TMV> tmv = rcp_dynamic_cast<TMV>(mv);
603 return tmv->getTpetraMultiVector();
606 RCP<TV> tv = rcp_dynamic_cast<TV>(mv);
608 return tv->getTpetraVector();
615 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
620 using Teuchos::rcp_dynamic_cast;
626 return tmv->getConstTpetraMultiVector();
631 return tv->getConstTpetraVector();
638 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
644 RCP<TV> tv = Teuchos::rcp_dynamic_cast<TV>(v);
646 return tv->getTpetraVector();
653 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
661 return tv->getConstTpetraVector();
671 #endif // THYRA_TPETRA_VECTOR_HPP
void constInitialize(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector)
Initialize.
Concrete implementation of Thyra::MultiVector in terms of Tpetra::MultiVector.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetraMultiVector(const RCP< MultiVectorBase< Scalar > > &mv) const
void applyImpl(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
Concrete implementation of an SPMD vector space for Tpetra.
void commitNonconstDetachedVectorViewImpl(RTOpPack::SubVectorView< Scalar > *sub_vec)
void initialize(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector)
Initialize.
RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraMultiVector(const RCP< const MultiVectorBase< Scalar > > &mv) const
TpetraVector()
Construct to uninitialized.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void norms1Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
virtual void applyOpImpl(const RTOpPack::RTOpT< Scalar > &op, const ArrayView< const Ptr< const VectorBase< Scalar > > > &vecs, const ArrayView< const Ptr< VectorBase< Scalar > > > &targ_vecs, const Ptr< RTOpPack::ReductTarget > &reduct_obj, const Ordinal global_offset) const
virtual void norms2Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
void getLocalVectorDataImpl(const Ptr< ArrayRCP< const Scalar > > &localValues) const
virtual void assignMultiVecImpl(const MultiVectorBase< Scalar > &mv)
void acquireNonconstDetachedVectorViewImpl(const Range1D &rng, RTOpPack::SubVectorView< Scalar > *sub_vec)
RCP< const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraVector() const
Get the embedded non-const Tpetra::Vector.
virtual void reciprocalImpl(const VectorBase< Scalar > &x)
void acquireDetachedVectorViewImpl(const Range1D &rng, RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm2WeightedImpl(const VectorBase< Scalar > &x) const
virtual void linearCombinationImpl(const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta)
virtual void absImpl(const VectorBase< Scalar > &x)
void getNonconstLocalVectorDataImpl(const Ptr< ArrayRCP< Scalar > > &localValues)
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVector_t
RCP< const SpmdVectorSpaceBase< Scalar > > spmdSpaceImpl() const
void initializeImpl(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< TpetraVector_t > &tpetraVector)
Concrete Thyra::SpmdVectorBase using Tpetra::Vector.
virtual void dotsImpl(const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetraVector()
Get the embedded non-const Tpetra::Vector.
bool nonnull(const boost::shared_ptr< T > &p)
virtual void assignImpl(Scalar alpha)
virtual void eleWiseScaleImpl(const VectorBase< Scalar > &x)
virtual void normsInfImpl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
RCP< const VectorSpaceBase< Scalar > > domain() const
virtual void randomizeImpl(Scalar l, Scalar u)
virtual void updateImpl(Scalar alpha, const MultiVectorBase< Scalar > &mv)
virtual void scaleImpl(Scalar alpha)
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)