10#ifndef THYRA_TPETRA_VECTOR_HPP
11#define THYRA_TPETRA_VECTOR_HPP
14#include "Thyra_TpetraVector_decl.hpp"
15#include "Thyra_TpetraMultiVector.hpp"
23template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
28template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
31 const RCP<Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &
tpetraVector
38template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
41 const RCP<
const Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &
tpetraVector
48template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
52 return tpetraVector_.getNonconstObj();
56template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
65template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
69 if (domainSpace_.is_null()) {
71 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
73 tpetraVector_.getConstObj()->getMap()->getComm()
84template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
88 return tpetraVectorSpace_;
95template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
99 *localValues = tpetraVector_.getNonconstObj()->get1dViewNonConst();
103template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
107 *localValues = tpetraVector_->get1dView();
114template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
122 if (!tpetraVector_.getNonconstObj()->isDistributed()) {
123 auto comm = tpetraVector_.getNonconstObj()->getMap()->getComm();
124 if (tpetraVector_.getConstObj()->getMap()->getComm()->getRank() == 0)
125 tpetraVector_.getNonconstObj()->randomize(l, u);
128 tpetraVector_.getNonconstObj()->reduce();
130 tpetraVector_.getNonconstObj()->randomize(l, u);
135template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
145 tpetraVector_.getNonconstObj()->abs(*tx);
152template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
162 tpetraVector_.getNonconstObj()->reciprocal(*tx);
169template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
180 tpetraVector_.getNonconstObj()->elementWiseMultiply(
181 ST::one(), *tx, *tpetraVector_.getConstObj(), ST::zero());
188template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
202 = Tpetra::createVector<Scalar>(tx->getMap());
203 temp->elementWiseMultiply(
204 ST::one(), *tx, *tpetraVector_.getConstObj(), ST::zero());
205 return ST::magnitude(ST::squareroot(tpetraVector_.getConstObj()->dot(*temp)));
212template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
218 const Ordinal global_offset
225template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
236template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
248template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
261template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
264 tpetraVector_.getNonconstObj()->putScalar(alpha);
268template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
272 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
277 tpetraVector_.getNonconstObj()->assign(*tmv);
284template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
287 tpetraVector_.getNonconstObj()->scale(alpha);
291template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
297 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
303 tpetraVector_.getNonconstObj()->update(alpha, *tmv, ST::one());
310template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
324 bool allCastsSuccessful =
true;
326 auto mvIter = mv.begin();
327 auto tmvIter = tmvs.
begin();
328 for (; mvIter != mv.end(); ++mvIter, ++tmvIter) {
329 tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromPtr(*mvIter));
333 allCastsSuccessful =
false;
341 auto len = mv.size();
343 tpetraVector_.getNonconstObj()->scale(beta);
344 }
else if (len == 1 && allCastsSuccessful) {
345 tpetraVector_.getNonconstObj()->update(alpha[0], *tmvs[0], beta);
346 }
else if (len == 2 && allCastsSuccessful) {
347 tpetraVector_.getNonconstObj()->update(alpha[0], *tmvs[0], alpha[1], *tmvs[1], beta);
348 }
else if (allCastsSuccessful) {
350 auto tmvIter = tmvs.
begin();
351 auto alphaIter = alpha.
begin();
356 for (; tmvIter != tmvs.
end(); ++tmvIter) {
357 if (tmvIter->getRawPtr() == tpetraVector_.getConstObj().getRawPtr()) {
365 tmvIter = tmvs.
begin();
369 if ((tmvs.
size() % 2) == 0) {
370 tpetraVector_.getNonconstObj()->scale(beta);
372 tpetraVector_.getNonconstObj()->update(*alphaIter, *(*tmvIter), beta);
376 for (; tmvIter != tmvs.
end(); tmvIter+=2, alphaIter+=2) {
377 tpetraVector_.getNonconstObj()->update(
378 *alphaIter, *(*tmvIter), *(alphaIter+1), *(*(tmvIter+1)), ST::one());
386template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
392 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
397 tpetraVector_.getConstObj()->dot(*tmv, prods);
404template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
409 tpetraVector_.getConstObj()->norm1(
norms);
413template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
418 tpetraVector_.getConstObj()->norm2(
norms);
422template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
427 tpetraVector_.getConstObj()->normInf(
norms);
431template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
433 const EOpTransp M_trans,
441 typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
447 if (nonnull(X_tpetra) && nonnull(Y_tpetra)) {
451 "Error, conjugation without transposition is not allowed for complex scalar types!");
469 Y_tpetra->multiply(trans,
Teuchos::NO_TRANS, alpha, *tpetraVector_.getConstObj(), *X_tpetra, beta);
481template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
482template<
class TpetraVector_t>
483void TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::initializeImpl(
494 this->updateSpmdSpace();
498template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
499RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
500TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
501getTpetraMultiVector(
const RCP<MultiVectorBase<Scalar> >& mv)
const
503 using Teuchos::rcp_dynamic_cast;
504 typedef TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
505 typedef TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
507 RCP<TMV> tmv = rcp_dynamic_cast<TMV>(mv);
509 return tmv->getTpetraMultiVector();
512 RCP<TV> tv = rcp_dynamic_cast<TV>(mv);
514 return tv->getTpetraVector();
517 return Teuchos::null;
521template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
526 using Teuchos::rcp_dynamic_cast;
532 return tmv->getConstTpetraMultiVector();
537 return tv->getConstTpetraVector();
540 return Teuchos::null;
544template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
550 RCP<TV> tv = Teuchos::rcp_dynamic_cast<TV>(v);
552 return tv->getTpetraVector();
554 return Teuchos::null;
559template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
567 return tv->getConstTpetraVector();
569 return Teuchos::null;
Interface for a collection of column vectors called a multi-vector.
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.
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
Calls applyOpImplWithComm(null,op,...).
void acquireNonconstDetachedVectorViewImpl(const Range1D &rng, RTOpPack::SubVectorView< Scalar > *sub_vec)
Implemented through this->getLocalData().
void acquireDetachedVectorViewImpl(const Range1D &rng, RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const
Implemented through this->getLocalData().
void commitNonconstDetachedVectorViewImpl(RTOpPack::SubVectorView< Scalar > *sub_vec)
Implemented through this->commitLocalData().
Concrete implementation of Thyra::MultiVector in terms of Tpetra::MultiVector.
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.
void acquireNonconstDetachedVectorViewImpl(const Range1D &rng, RTOpPack::SubVectorView< Scalar > *sub_vec)
void getLocalVectorDataImpl(const Ptr< ArrayRCP< const Scalar > > &localValues) const
virtual void assignImpl(Scalar alpha)
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm2WeightedImpl(const VectorBase< Scalar > &x) const
virtual void randomizeImpl(Scalar l, Scalar u)
virtual void eleWiseScaleImpl(const VectorBase< Scalar > &x)
RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetraVector()
Get the embedded non-const Tpetra::Vector.
RCP< const VectorSpaceBase< Scalar > > domain() const
virtual void linearCombinationImpl(const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta)
RCP< const SpmdVectorSpaceBase< Scalar > > spmdSpaceImpl() const
virtual void dotsImpl(const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
void getNonconstLocalVectorDataImpl(const Ptr< ArrayRCP< Scalar > > &localValues)
virtual void scaleImpl(Scalar alpha)
TpetraVector()
Construct to uninitialized.
virtual void absImpl(const VectorBase< Scalar > &x)
void constInitialize(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector)
Initialize.
virtual void assignMultiVecImpl(const MultiVectorBase< Scalar > &mv)
virtual void updateImpl(Scalar alpha, const MultiVectorBase< Scalar > &mv)
virtual void normsInfImpl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
void commitNonconstDetachedVectorViewImpl(RTOpPack::SubVectorView< Scalar > *sub_vec)
virtual void norms1Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
virtual void reciprocalImpl(const VectorBase< Scalar > &x)
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
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.
void acquireDetachedVectorViewImpl(const Range1D &rng, RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const
void applyImpl(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
virtual void norms2Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
RCP< const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraVector() const
Get the embedded non-const Tpetra::Vector.
void initialize(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector)
Initialize.
Abstract interface for finite-dimensional dense vectors.
void applyImpl(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
. Applies vector or its adjoint (transpose) as a linear operator.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm2WeightedImpl(const VectorBase< Scalar > &x) const
Default implementation of norm_2 (weighted) using RTOps.
virtual void reciprocalImpl(const VectorBase< Scalar > &x)
Default implementation of reciprocal using RTOps.
virtual void absImpl(const VectorBase< Scalar > &x)
Default implementation of abs using RTOps.
virtual void eleWiseScaleImpl(const VectorBase< Scalar > &x)
Default implementation of ele_wise_scale using RTOps.
#define TEUCHOS_ASSERT(assertion_test)
#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)