Thyra Version of the Day
Loading...
Searching...
No Matches
Thyra_TpetraVector_decl.hpp
1// @HEADER
2// *****************************************************************************
3// Thyra: Interfaces and Support for Abstract Numerical Algorithms
4//
5// Copyright 2004 NTESS and the Thyra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef THYRA_TPETRA_VECTOR_DECL_HPP
11#define THYRA_TPETRA_VECTOR_DECL_HPP
12
13
14#include "Thyra_SpmdVectorDefaultBase.hpp"
15#include "Thyra_TpetraVectorSpace_decl.hpp"
16#include "Tpetra_Vector.hpp"
17#include "Teuchos_ConstNonconstObjectContainer.hpp"
18
19
20namespace Thyra {
21
22
27template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
29 : virtual public SpmdVectorDefaultBase<Scalar>
30{
31public:
32
35
38
40 void initialize(
42 const RCP<Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVector
43 );
44
46 void constInitialize(
48 const RCP<const Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVector
49 );
50
54
58
60
66
67 // Should these Impl functions should alsp be protected???
68//protected:
69
75
79 void getNonconstLocalVectorDataImpl(const Ptr<ArrayRCP<Scalar> > &localValues);
81 void getLocalVectorDataImpl(const Ptr<ArrayRCP<const Scalar> > &localValues) const;
83
84protected:
85
88
90 virtual void randomizeImpl(Scalar l, Scalar u);
91
93 virtual void absImpl(const VectorBase<Scalar>& x);
94
96 virtual void reciprocalImpl(const VectorBase<Scalar>& x);
97
99 virtual void eleWiseScaleImpl(const VectorBase<Scalar>& x);
100
104
106 virtual void applyOpImpl(
107 const RTOpPack::RTOpT<Scalar> &op,
108 const ArrayView<const Ptr<const VectorBase<Scalar> > > &vecs,
109 const ArrayView<const Ptr<VectorBase<Scalar> > > &targ_vecs,
110 const Ptr<RTOpPack::ReductTarget> &reduct_obj,
111 const Ordinal global_offset
112 ) const;
113
116 const Range1D& rng,
118 ) const;
119
122 const Range1D& rng,
124 );
125
129 );
130
132
136 virtual void assignImpl(Scalar alpha);
137
139 virtual void assignMultiVecImpl(const MultiVectorBase<Scalar>& mv);
140
142 virtual void scaleImpl(Scalar alpha);
143
145 virtual void updateImpl(
146 Scalar alpha,
148 );
149
151 virtual void linearCombinationImpl(
152 const ArrayView<const Scalar>& alpha,
153 const ArrayView<const Ptr<const MultiVectorBase<Scalar> > >& mv,
154 const Scalar& beta
155 );
156
158 virtual void dotsImpl(
159 const MultiVectorBase<Scalar>& mv,
160 const ArrayView<Scalar>& prods
161 ) const;
162
164 virtual void norms1Impl(
166 ) const;
167
169 virtual void norms2Impl(
171 ) const;
172
174 virtual void normsInfImpl(
176 ) const;
177
179
182
184 void applyImpl(
185 const EOpTransp M_trans,
187 const Ptr<MultiVectorBase<Scalar> > &Y,
188 const Scalar alpha,
189 const Scalar beta
190 ) const;
191
193
194private:
195
196 // ///////////////////////////////////////
197 // Private data members
198
200 tpetraVectorSpace_;
201
203 domainSpace_;
204
206 tpetraVector_;
207
208 typedef Tpetra::MultiVector<Scalar, LocalOrdinal,GlobalOrdinal,Node> TpetraMultiVector_t;
209
210 // ////////////////////////////////////
211 // Private member functions
212
213 template<class TpetraVector_t>
214 void initializeImpl(
215 const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
217 );
218
219 // Non-throwing Tpetra Vector/MultiVector extraction methods.
220 // Return null if casting failed.
222 getTpetraMultiVector(const RCP<MultiVectorBase<Scalar> >& mv) const;
223
225 getConstTpetraMultiVector(const RCP<const MultiVectorBase<Scalar> >& mv) const;
226
228 getTpetraVector(const RCP<VectorBase<Scalar> >& v) const;
229
231 getConstTpetraVector(const RCP<const VectorBase<Scalar> >& v) const;
232
233};
234
235
240template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
241inline
244 const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
245 const RCP<Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVector
246 )
247{
250 v->initialize(tpetraVectorSpace, tpetraVector);
251 return v;
252}
253
254
259template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
260inline
263 const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
264 const RCP<const Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVector
265 )
266{
269 v->constInitialize(tpetraVectorSpace, tpetraVector);
270 return v;
271}
272
273
274} // end namespace Thyra
275
276
277#endif // THYRA_TPETRA_VECTOR_DECL_HPP
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.
Base class for SPMD vectors that can provide views of contiguous elements in a process.
Concrete implementation of an SPMD vector space for Tpetra.
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)
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.
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.
EOpTransp
Enumeration for determining how a linear operator is applied. `*.
Teuchos::Ordinal Ordinal
Type for the dimension of a vector space. `*.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)