Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_TpetraMultiVector_decl.hpp
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Thyra: Interfaces and Support for Abstract Numerical Algorithms
5 // Copyright (2004) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (bartlettra@ornl.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef THYRA_TPETRA_MULTIVECTOR_DECL_HPP
43 #define THYRA_TPETRA_MULTIVECTOR_DECL_HPP
44 
45 #include "Thyra_SpmdMultiVectorDefaultBase.hpp"
46 #include "Thyra_TpetraVectorSpace_decl.hpp"
47 #include "Tpetra_MultiVector.hpp"
48 #include "Teuchos_ConstNonconstObjectContainer.hpp"
49 
50 
51 namespace Thyra {
52 
53 
61 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
63  : virtual public SpmdMultiVectorDefaultBase<Scalar>
64 {
65 public:
66 
69 
72 
75  void initialize(
76  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
77  const RCP<const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
78  const RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraMultiVector
79  );
80 
83  void constInitialize(
84  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
85  const RCP<const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
86  const RCP<const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraMultiVector
87  );
88 
92 
96 
98 
103  domainScalarProdVecSpc() const;
105 
106 protected:
107 
111  virtual void assignImpl(Scalar alpha);
112 
114  virtual void assignMultiVecImpl(const MultiVectorBase<Scalar>& mv);
115 
117  virtual void scaleImpl(Scalar alpha);
118 
120  virtual void updateImpl(
121  Scalar alpha,
122  const MultiVectorBase<Scalar>& mv
123  );
124 
126  virtual void linearCombinationImpl(
127  const ArrayView<const Scalar>& alpha,
128  const ArrayView<const Ptr<const MultiVectorBase<Scalar> > >& mv,
129  const Scalar& beta
130  );
131 
133  virtual void dotsImpl(
134  const MultiVectorBase<Scalar>& mv,
135  const ArrayView<Scalar>& prods
136  ) const;
137 
139  virtual void norms1Impl(
141  ) const;
142 
144  virtual void norms2Impl(
146  ) const;
147 
149  virtual void normsInfImpl(
151  ) const;
152 
157 
160  contigSubViewImpl(const Range1D& colRng) const;
163  nonconstContigSubViewImpl(const Range1D& colRng);
166  nonContigSubViewImpl(const ArrayView<const int>& cols_in) const;
170 
172  virtual void mvMultiReductApplyOpImpl(
173  const RTOpPack::RTOpT<Scalar> &primary_op,
174  const ArrayView<const Ptr<const MultiVectorBase<Scalar> > > &multi_vecs,
175  const ArrayView<const Ptr<MultiVectorBase<Scalar> > > &targ_multi_vecs,
176  const ArrayView<const Ptr<RTOpPack::ReductTarget> > &reduct_objs,
177  const Ordinal primary_global_offset
178  ) const;
179 
182  const Range1D &rowRng,
183  const Range1D &colRng,
185  ) const;
186 
189  const Range1D &rowRng,
190  const Range1D &colRng,
192  );
193 
197  );
198 
199 // /** \brief . */
200 // RCP<const MultiVectorBase<Scalar> >
201 // nonContigSubViewImpl( const ArrayView<const int> &cols ) const;
202 // /** \brief . */
203 // RCP<MultiVectorBase<Scalar> >
204 // nonconstNonContigSubViewImpl( const ArrayView<const int> &cols );
206 
213  const Ptr<ArrayRCP<Scalar> > &localValues, const Ptr<Ordinal> &leadingDim
214  );
217  const Ptr<ArrayRCP<const Scalar> > &localValues, const Ptr<Ordinal> &leadingDim
218  ) const;
219 
221 
225  virtual void euclideanApply(
226  const EOpTransp M_trans,
227  const MultiVectorBase<Scalar> &X,
228  const Ptr<MultiVectorBase<Scalar> > &Y,
229  const Scalar alpha,
230  const Scalar beta
231  ) const;
232 
234 
235 private:
236 
237  // ///////////////////////////////////////
238  // Private data members
239 
243  tpetraMultiVector_;
244 
245  // ////////////////////////////////////
246  // Private member functions
247 
248  template<class TpetraMultiVector_t>
249  void initializeImpl(
250  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
251  const RCP<const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
253  );
254 
255  // Non-throwing Tpetra MultiVector extraction methods.
256  // Return null if casting failed.
259 
261  getConstTpetraMultiVector(const RCP<const MultiVectorBase<Scalar> >& mv) const;
262 
263 };
264 
265 
270 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
273  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
274  const RCP<const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
275  const RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraMultiVector
276  )
277 {
280  tmv->initialize(tpetraVectorSpace, domainSpace, tpetraMultiVector);
281  return tmv;
282 }
283 
284 
289 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
292  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
293  const RCP<const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
294  const RCP<const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraMultiVector
295  )
296 {
299  tmv->constInitialize(tpetraVectorSpace, domainSpace, tpetraMultiVector);
300  return tmv;
301 }
302 
303 
304 } // end namespace Thyra
305 
306 
307 #endif // THYRA_TPETRA_MULTIVECTOR_DECL_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)
EOpTransp
Enumeration for determining how a linear operator is applied. `*.
RCP< MultiVectorBase< Scalar > > nonconstContigSubViewImpl(const Range1D &colRng)
virtual void euclideanApply(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
Base node implementation class for SPMD multi-vectors.
RCP< const MultiVectorBase< Scalar > > nonContigSubViewImpl(const ArrayView< const int > &cols_in) 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 TpetraMultiVector.
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.
virtual void scaleImpl(Scalar alpha)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::Ordinal Ordinal
Type for the dimension of a vector space. `*.
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
Interface for a collection of column vectors called a multi-vector.
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 TpetraMultiVector.
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.
void norms(const MultiVectorBase< Scalar > &V, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms)
Column-wise multi-vector natural norm.
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)
virtual void dotsImpl(const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
RCP< VectorBase< Scalar > > nonconstColImpl(Ordinal j)
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
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)