15#ifndef ANASAZI_THYRA_ADAPTER_HPP
16#define ANASAZI_THYRA_ADAPTER_HPP
22#include <Thyra_DetachedMultiVectorView.hpp>
23#include <Thyra_MultiVectorBase.hpp>
24#include <Thyra_MultiVectorStdOps.hpp>
25#include <Thyra_VectorStdOps.hpp>
27#include <Teuchos_Ptr.hpp>
28#include <Teuchos_ArrayRCP.hpp>
29#include <Teuchos_ArrayView.hpp>
47 template<
class ScalarType>
51 typedef Thyra::MultiVectorBase<ScalarType> TMVB;
52 typedef Teuchos::ScalarTraits<ScalarType> ST;
53 typedef typename ST::magnitudeType magType;
64 static Teuchos::RCP<TMVB>
Clone(
const TMVB & mv,
const int numvecs )
66 Teuchos::RCP<TMVB> c = Thyra::createMembers( mv.range(), numvecs );
74 static Teuchos::RCP<TMVB>
77 const int numvecs = mv.domain()->dim();
79 Teuchos::RCP< TMVB > cc = Thyra::createMembers (mv.range(), numvecs);
81 Thyra::assign (Teuchos::outArg (*cc), mv);
90 static Teuchos::RCP< TMVB >
CloneCopy(
const TMVB & mv,
const std::vector<int>& index )
92 const int numvecs = index.size();
94 Teuchos::RCP<TMVB > cc = Thyra::createMembers (mv.range(), numvecs);
96 Teuchos::RCP<const TMVB > view = mv.subView ( Teuchos::arrayViewFromVector( index ) );
98 Thyra::assign (Teuchos::outArg (*cc), *view);
102 static Teuchos::RCP<TMVB>
103 CloneCopy (
const TMVB& mv,
const Teuchos::Range1D& index)
105 const int numVecs = index.size();
107 Teuchos::RCP<TMVB> cc = Thyra::createMembers (mv.range(), numVecs);
109 Teuchos::RCP<const TMVB> view = mv.subView (index);
111 Thyra::assign (Teuchos::outArg (*cc), *view);
122 int numvecs = index.size();
138 for (
int i=0; i<numvecs; i++) {
139 if (lb+i != index[i]) contig =
false;
142 Teuchos::RCP< TMVB > cc;
144 const Thyra::Range1D rng(lb,lb+numvecs-1);
146 cc = mv.subView(rng);
150 cc = mv.subView( Teuchos::arrayViewFromVector( index ) );
155 static Teuchos::RCP<TMVB>
163 return mv.subView (index);
171 static Teuchos::RCP<const TMVB >
CloneView(
const TMVB & mv,
const std::vector<int>& index )
173 int numvecs = index.size();
189 for (
int i=0; i<numvecs; i++) {
190 if (lb+i != index[i]) contig =
false;
193 Teuchos::RCP< const TMVB > cc;
195 const Thyra::Range1D rng(lb,lb+numvecs-1);
197 cc = mv.subView(rng);
201 cc = mv.subView(Teuchos::arrayViewFromVector( index ) );
206 static Teuchos::RCP<const TMVB>
207 CloneView (
const TMVB& mv,
const Teuchos::Range1D& index)
214 return mv.subView (index);
224 {
return mv.range()->dim(); }
228 {
return mv.domain()->dim(); }
238 const Teuchos::SerialDenseMatrix<int,ScalarType>& B,
239 const ScalarType beta, TMVB & mv )
244 Teuchos::RCP< const TMVB >
245 B_thyra = Thyra::createMembersView(
247 RTOpPack::ConstSubMultiVectorView<ScalarType>(
249 Teuchos::arcpFromArrayView(Teuchos::arrayView(&B(0,0), B.stride()*B.numCols())), B.stride()
253 A.apply(Thyra::NOTRANS,*B_thyra,Teuchos::outArg (mv),alpha,beta);
258 static void MvAddMv(
const ScalarType alpha,
const TMVB & A,
259 const ScalarType beta,
const TMVB & B, TMVB & mv )
261 using Teuchos::tuple;
using Teuchos::ptrInArg;
using Teuchos::inoutArg;
263 Thyra::linear_combination<ScalarType>(
264 tuple(alpha, beta)(), tuple(ptrInArg(A), ptrInArg(B))(), ST::zero(), inoutArg(mv));
269 static void MvTransMv(
const ScalarType alpha,
const TMVB & A,
const TMVB & mv,
270 Teuchos::SerialDenseMatrix<int,ScalarType>& B )
273 int m = A.domain()->dim();
274 int n = mv.domain()->dim();
277 B_thyra = Thyra::createMembersView(
279 RTOpPack::SubMultiVectorView<ScalarType>(
281 Teuchos::arcpFromArrayView(Teuchos::arrayView(&B(0,0), B.stride()*B.numCols())), B.stride()
284 A.apply(Thyra::CONJTRANS,mv,B_thyra.ptr(),alpha,Teuchos::ScalarTraits<ScalarType>::zero());
290 static void MvDot(
const TMVB & mv,
const TMVB & A, std::vector<ScalarType> &b )
291 { Thyra::dots(mv,A,Teuchos::arrayViewFromVector( b )); }
297 const ScalarType alpha)
299 Thyra::scale (alpha, Teuchos::inOutArg (mv));
306 const std::vector<ScalarType>& alpha)
308 for (
unsigned int i=0; i<alpha.size(); i++) {
309 Thyra::scale (alpha[i], mv.col(i).ptr());
321 static void MvNorm(
const TMVB & mv, std::vector<
typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> &normvec )
322 { Thyra::norms_2(mv,Teuchos::arrayViewFromVector( normvec )); }
331 static void SetBlock(
const TMVB & A,
const std::vector<int>& index, TMVB & mv )
334 int numvecs = index.size();
335 std::vector<int> indexA(numvecs);
336 int numAcols = A.domain()->dim();
337 for (
int i=0; i<numvecs; i++) {
342 if ( numAcols < numvecs ) {
347 else if ( numAcols > numvecs ) {
349 indexA.resize( numAcols );
352 Teuchos::RCP< const TMVB > relsource = A.subView( Teuchos::arrayViewFromVector( indexA ) );
354 Teuchos::RCP< TMVB > reldest = mv.subView( Teuchos::arrayViewFromVector( index ) );
356 Thyra::assign (Teuchos::outArg (*reldest), *relsource);
360 SetBlock (
const TMVB& A,
const Teuchos::Range1D& index, TMVB& mv)
362 const int numColsA = A.domain()->dim();
363 const int numColsMv = mv.domain()->dim();
365 const bool validIndex = index.lbound() >= 0 && index.ubound() < numColsMv;
367 const bool validSource = index.size() <= numColsA;
369 if (! validIndex || ! validSource)
371 std::ostringstream os;
372 os <<
"Anasazi::MultiVecTraits<Scalar, Thyra::MultiVectorBase<Scalar> "
373 ">::SetBlock(A, [" << index.lbound() <<
", " << index.ubound()
375 TEUCHOS_TEST_FOR_EXCEPTION(index.lbound() < 0, std::invalid_argument,
376 os.str() <<
"Range lower bound must be nonnegative.");
377 TEUCHOS_TEST_FOR_EXCEPTION(index.ubound() >= numColsMv, std::invalid_argument,
378 os.str() <<
"Range upper bound must be less than "
379 "the number of columns " << numColsA <<
" in the "
380 "'mv' output argument.");
381 TEUCHOS_TEST_FOR_EXCEPTION(index.size() > numColsA, std::invalid_argument,
382 os.str() <<
"Range must have no more elements than"
383 " the number of columns " << numColsA <<
" in the "
384 "'A' input argument.");
385 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Should never get here!");
391 Teuchos::RCP<TMVB> mv_view;
392 if (index.lbound() == 0 && index.ubound()+1 == numColsMv)
393 mv_view = Teuchos::rcpFromRef (mv);
395 mv_view = mv.subView (index);
400 Teuchos::RCP<const TMVB> A_view;
401 if (index.size() == numColsA)
402 A_view = Teuchos::rcpFromRef (A);
404 A_view = A.subView (Teuchos::Range1D(0, index.size()-1));
407 Thyra::assign (Teuchos::outArg (*mv_view), *A_view);
411 Assign (
const TMVB& A, TMVB& mv)
414 using Teuchos::Range1D;
417 const int numColsA = A.domain()->dim();
418 const int numColsMv = mv.domain()->dim();
419 if (numColsA > numColsMv) {
420 const std::string prefix (
"Anasazi::MultiVecTraits<Scalar, "
421 "Thyra::MultiVectorBase<Scalar>"
422 " >::Assign(A, mv): ");
423 TEUCHOS_TEST_FOR_EXCEPTION(numColsA > numColsMv, std::invalid_argument,
424 prefix <<
"Input multivector 'A' has "
425 << numColsA <<
" columns, but output multivector "
426 "'mv' has only " << numColsMv <<
" columns.");
427 TEUCHOS_TEST_FOR_EXCEPTION(
true, std::logic_error,
"Should never get here!");
430 if (numColsA == numColsMv) {
431 Thyra::assign (outArg (mv), A);
435 Thyra::assign (outArg (*mv_view), A);
445 Thyra::randomize(-Teuchos::ScalarTraits<ScalarType>::one(),
446 Teuchos::ScalarTraits<ScalarType>::one(),
447 Teuchos::outArg (mv));
454 ScalarType alpha = Teuchos::ScalarTraits<ScalarType>::zero())
456 Thyra::assign (Teuchos::outArg (mv), alpha);
466 static void MvPrint(
const TMVB & mv, std::ostream& os )
468 Teuchos::RCP<Teuchos::FancyOStream> out = Teuchos::getFancyOStream(Teuchos::rcp(&os,
false));
469 out->setf(std::ios_base::scientific);
471 mv.describe(*out,Teuchos::VERB_EXTREME);
494 template <
class ScalarType>
495 class OperatorTraits < ScalarType, Thyra::MultiVectorBase<ScalarType>, Thyra::LinearOpBase<ScalarType> >
502 static void Apply (
const Thyra::LinearOpBase< ScalarType >& Op,
const Thyra::MultiVectorBase< ScalarType > & x, Thyra::MultiVectorBase< ScalarType > & y )
504 Op.apply(Thyra::NOTRANS,x,Teuchos::outArg (y), Teuchos::ScalarTraits<ScalarType>::one(),Teuchos::ScalarTraits<ScalarType>::zero());
Anasazi header file which uses auto-configuration information to include necessary C++ headers.
Declaration of basic traits for the multivector type.
Virtual base class which defines basic traits for the operator type.
static void MvTimesMatAddMv(const ScalarType alpha, const TMVB &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, TMVB &mv)
Update mv with .
static void MvRandom(TMVB &mv)
Replace the vectors in mv with random vectors.
static Teuchos::RCP< const TMVB > CloneView(const TMVB &mv, const std::vector< int > &index)
Creates a new const MultiVectorBase that shares the selected contents of mv (shallow copy).
static Teuchos::RCP< TMVB > CloneViewNonConst(TMVB &mv, const std::vector< int > &index)
Creates a new MultiVectorBase that shares the selected contents of mv (shallow copy).
static void SetBlock(const TMVB &A, const std::vector< int > &index, TMVB &mv)
Copy the vectors in A to a set of vectors in mv indicated by the indices given in index.
static void MvTransMv(const ScalarType alpha, const TMVB &A, const TMVB &mv, Teuchos::SerialDenseMatrix< int, ScalarType > &B)
Compute a dense matrix B through the matrix-matrix multiply .
static Teuchos::RCP< TMVB > CloneCopy(const TMVB &mv, const std::vector< int > &index)
Creates a new MultiVectorBase and copies the selected contents of mv into the new vector (deep copy).
static Teuchos::RCP< TMVB > Clone(const TMVB &mv, const int numvecs)
Creates a new empty MultiVectorBase containing numvecs columns.
static void MvInit(TMVB &mv, ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
Replace each element of the vectors in mv with alpha.
static ptrdiff_t GetGlobalLength(const TMVB &mv)
Obtain the vector length of mv.
static void MvAddMv(const ScalarType alpha, const TMVB &A, const ScalarType beta, const TMVB &B, TMVB &mv)
Replace mv with .
static void MvScale(TMVB &mv, const ScalarType alpha)
Scale each element of the vectors in *this with alpha.
static Teuchos::RCP< TMVB > CloneCopy(const TMVB &mv)
Creates a new MultiVectorBase and copies contents of mv into the new vector (deep copy).
static void MvDot(const TMVB &mv, const TMVB &A, std::vector< ScalarType > &b)
Compute a vector b where the components are the individual dot-products of the i-th columns of A and ...
static void MvScale(TMVB &mv, const std::vector< ScalarType > &alpha)
Scale each element of the i-th vector in *this with alpha[i].
static int GetNumberVecs(const TMVB &mv)
Obtain the number of vectors in mv.
static void MvPrint(const TMVB &mv, std::ostream &os)
Print the mv multi-vector to the os output stream.
static void MvNorm(const TMVB &mv, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec)
Compute the 2-norm of each individual vector of mv. Upon return, normvec[i] holds the value of ,...
Traits class which defines basic operations on multivectors.
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy).
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
Creates a new const MV that shares the selected contents of mv (shallow copy).
static void Assign(const MV &A, MV &mv)
mv := A
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
static void SetBlock(const MV &A, const std::vector< int > &index, MV &mv)
Copy the vectors in A to a set of vectors in mv indicated by the indices given in index.
static void Apply(const Thyra::LinearOpBase< ScalarType > &Op, const Thyra::MultiVectorBase< ScalarType > &x, Thyra::MultiVectorBase< ScalarType > &y)
This method takes the MultiVectorBase x and applies the LinearOpBase Op to it resulting in the MultiV...
Virtual base class which defines basic traits for the operator type.
Namespace Anasazi contains the classes, structs, enums and utilities used by the Anasazi package.