Thyra Version of the Day
Loading...
Searching...
No Matches
Thyra_VectorBase.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_VECTOR_BASE_DECL_HPP
11#define THYRA_VECTOR_BASE_DECL_HPP
12
13
14#include "Thyra_OperatorVectorTypes.hpp"
15#include "Thyra_MultiVectorBase_decl.hpp"
16#include "RTOpPack_RTOpT.hpp"
17#include "RTOpPack_SparseSubVectorT.hpp"
18
19
20namespace Thyra {
21
22
113template<class Scalar>
114class VectorBase : virtual public MultiVectorBase<Scalar>
115{
116public:
117
118#ifdef THYRA_INJECT_USING_DECLARATIONS
119 using MultiVectorBase<Scalar>::apply;
120#endif
121
124
125 // Overloading assign for VectorBase argument
126 using MultiVectorBase<Scalar>::assign;
127
135 { assignVecImpl(x); }
136
146 void randomize(Scalar l, Scalar u)
147 { randomizeImpl(l,u); }
148
149 // Overloading update for VectorBase argument.
150 using MultiVectorBase<Scalar>::update;
151
158 void update(
159 Scalar alpha,
160 const VectorBase<Scalar>& x)
161 { updateVecImpl(alpha, x); }
162
163 // Overloading linear_combination for VectorBase arguments.
165
190 const ArrayView<const Scalar>& alpha,
191 const ArrayView<const Ptr<const VectorBase<Scalar> > >& x,
192 const Scalar& beta
193 )
194 { linearCombinationVecImpl(alpha, x, beta); }
195
199 Scalar dot(const VectorBase<Scalar>& x) const
200 { return dotImpl(x); }
201
206 norm_1() const
207 { return norm1Impl(); }
208
214 norm_2() const
215 { return norm2Impl(); }
216
223 { return norm2WeightedImpl(x); }
224
230 norm_inf() const
231 { return normInfImpl(); }
232
234
237
255
257
260
266 const RTOpPack::RTOpT<Scalar> &op,
267 const ArrayView<const Ptr<const VectorBase<Scalar> > > &vecs,
268 const ArrayView<const Ptr<VectorBase<Scalar> > > &targ_vecs,
269 const Ptr<RTOpPack::ReductTarget> &reduct_obj,
270 const Ordinal global_offset
271 ) const
272 {
273 applyOpImpl(op, vecs, targ_vecs, reduct_obj, global_offset);
274 }
275
277
280
306 virtual RCP<VectorBase<Scalar> > clone_v() const = 0;
307
309
312
319 ) const
320 { acquireDetachedVectorViewImpl(rng,sub_vec); }
321
330
339
348
355 )
356 { setSubVectorImpl(sub_vec); }
357
359
360protected:
361
364
368 virtual void assignVecImpl(const VectorBase<Scalar>& x) = 0;
369
373 virtual void randomizeImpl(Scalar l, Scalar u) = 0;
374
378 virtual void absImpl(const VectorBase<Scalar>& x) = 0;
379
383 virtual void reciprocalImpl(const VectorBase<Scalar>& x) = 0;
384
388 virtual void eleWiseScaleImpl(const VectorBase<Scalar>& x) = 0;
389
393 virtual void updateVecImpl(
394 Scalar alpha,
395 const VectorBase<Scalar>& x) = 0;
396
401 const ArrayView<const Scalar>& alpha,
402 const ArrayView<const Ptr<const VectorBase<Scalar> > >& x,
403 const Scalar& beta
404 ) = 0;
405
409 virtual Scalar dotImpl(const VectorBase<Scalar>& x) const = 0;
410
415 norm1Impl() const = 0;
416
421 norm2Impl() const = 0;
422
428
433 normInfImpl() const = 0;
434
451 virtual void applyOpImpl(
452 const RTOpPack::RTOpT<Scalar> &op,
453 const ArrayView<const Ptr<const VectorBase<Scalar> > > &vecs,
454 const ArrayView<const Ptr<VectorBase<Scalar> > > &targ_vecs,
455 const Ptr<RTOpPack::ReductTarget> &reduct_obj,
456 const Ordinal global_offset
457 ) const = 0;
458
502 ) const = 0;
503
525 ) const = 0;
526
576 const Range1D& rng, RTOpPack::SubVectorView<Scalar>* sub_vec
577 ) = 0;
578
603 ) = 0;
604
628 virtual void setSubVectorImpl(
630 ) = 0;
631
633
634private:
635
636 // Not defined and not to be called
638 operator=(const VectorBase<Scalar>&);
639
640public:
641
642 // These are functions that may be removed in the future and should not be
643 // called by client code.
644
645 // Don't call this directly. Use non-member Thyra::abs(). This is because
646 // this member function may disappear in the future. (see Trilinos GitHub
647 // issue #330)
648 void abs(const VectorBase<Scalar>& x)
649 { absImpl(x); }
650
651 // Don't call this directly. Use non-member Thyra::reciprocal(). This is because
652 // this member function may disappear in the future. (see Trilinos GitHub
653 // issue #330)
654 void reciprocal(const VectorBase<Scalar>& x)
655 { reciprocalImpl(x); }
656
657 // Don't call this directly. Use non-member Thyra::ele_wise_scale(). This
658 // is because this member function may disappear in the future. (see
659 // Trilinos GitHub issue #330)
660 void ele_wise_scale(const VectorBase<Scalar>& x)
661 { eleWiseScaleImpl(x); }
662
663};
664
665
736template<class Scalar>
737inline
739 const RTOpPack::RTOpT<Scalar> &op,
740 const ArrayView<const Ptr<const VectorBase<Scalar> > > &vecs,
741 const ArrayView<const Ptr<VectorBase<Scalar> > > &targ_vecs,
742 const Ptr<RTOpPack::ReductTarget> &reduct_obj,
743 const Ordinal global_offset = 0
744 )
745{
746 if (vecs.size())
747 vecs[0]->applyOp(op, vecs, targ_vecs, reduct_obj, global_offset);
748 else if (targ_vecs.size())
749 targ_vecs[0]->applyOp(op, vecs, targ_vecs, reduct_obj, global_offset);
750}
751
752
753} // end namespace Thyra
754
755
756#endif // THYRA_VECTOR_BASE_DECL_HPP
void apply(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
Apply the linear operator to a multi-vector : Y = alpha*op(M)*X + beta*Y.
Interface for a collection of column vectors called a multi-vector.
Abstract interface for finite-dimensional dense vectors.
virtual void acquireDetachedVectorViewImpl(const Range1D &rng, RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const =0
Get a non-mutable explicit view of a sub-vector.
virtual void commitNonconstDetachedVectorViewImpl(RTOpPack::SubVectorView< Scalar > *sub_vec)=0
Commit changes for a mutable explicit view of a sub-vector.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm2WeightedImpl(const VectorBase< Scalar > &x) const =0
Virtual implementation for NVI norm_2 (weighted).
virtual void assignVecImpl(const VectorBase< Scalar > &x)=0
Virtual implementation for NVI assign.
Teuchos::ScalarTraits< Scalar >::magnitudeType norm_1() const
One (1) norm: result = ||v||1.
void acquireDetachedView(const Range1D &rng, RTOpPack::SubVectorView< Scalar > *sub_vec)
Calls acquireNonconstDetachedVectorViewImpl().
Scalar dot(const VectorBase< Scalar > &x) const
Euclidean dot product: result = x^H * this.
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 =0
Apply a reduction/transformation operator over a set of vectors: op(op(v[0]...v[nv-1],...
void applyOp(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=0)
Apply a reduction/transformation operator over a set of vectors: op(op(v[0]...v[nv-1],...
virtual RCP< VectorBase< Scalar > > clone_v() const =0
Returns a seprate cloned copy of *this vector with the same values but different storage.
virtual void randomizeImpl(Scalar l, Scalar u)=0
Virtual implementation for NVI randomize.
void randomize(Scalar l, Scalar u)
Random vector generation:
virtual void acquireNonconstDetachedVectorViewImpl(const Range1D &rng, RTOpPack::SubVectorView< Scalar > *sub_vec)=0
Get a mutable explicit view of a sub-vector.
void linear_combination(const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const VectorBase< Scalar > > > &x, const Scalar &beta)
Linear combination:
virtual void updateVecImpl(Scalar alpha, const VectorBase< Scalar > &x)=0
Virtual implementation for NVI update.
void setSubVector(const RTOpPack::SparseSubVectorT< Scalar > &sub_vec)
Calls setSubVectorImpl().
virtual void releaseDetachedVectorViewImpl(RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const =0
Free an explicit view of a sub-vector.
void commitDetachedView(RTOpPack::SubVectorView< Scalar > *sub_vec)
Calls commitDetachedView().
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm2Impl() const =0
Virtual implementation for NVI norm_2.
virtual void setSubVectorImpl(const RTOpPack::SparseSubVectorT< Scalar > &sub_vec)=0
Set a specific sub-vector.
void applyOp(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 applyOpImpl().
Teuchos::ScalarTraits< Scalar >::magnitudeType norm_2(const VectorBase< Scalar > &x) const
Weighted Euclidean (2) norm: result = ||v||2.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType normInfImpl() const =0
Virtual implementation for NVI norm_inf.
void assign(const VectorBase< Scalar > &x)
Vector assignment:
virtual void linearCombinationVecImpl(const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const VectorBase< Scalar > > > &x, const Scalar &beta)=0
Virtual implementation for NVI linear_combination.
virtual Scalar dotImpl(const VectorBase< Scalar > &x) const =0
Virtual implementation for NVI dot.
void acquireDetachedView(const Range1D &rng, RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const
Calls acquireDetachedVectorViewImpl().
Teuchos::ScalarTraits< Scalar >::magnitudeType norm_2() const
Euclidean (2) norm: result = ||v||2.
virtual RCP< const VectorSpaceBase< Scalar > > space() const =0
Return a smart pointer to the vector space that this vector belongs to.
void releaseDetachedView(RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const
Calls releaseDetachedVectorViewImpl().
Teuchos::ScalarTraits< Scalar >::magnitudeType norm_inf() const
Infinity norm: result = ||v||inf.
virtual void reciprocalImpl(const VectorBase< Scalar > &x)=0
Virtual implementation for NVI reciprocal.
virtual void eleWiseScaleImpl(const VectorBase< Scalar > &x)=0
Virtual implementation for NVI ele_wise_scale.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm1Impl() const =0
Virtual implementation for NVI norm_1.
virtual void absImpl(const VectorBase< Scalar > &x)=0
Virtual implementation for NVI abs.
void update(Scalar alpha, const VectorBase< Scalar > &x)
AXPY:
Teuchos::Ordinal Ordinal
Type for the dimension of a vector space. `*.