Anasazi Version of the Day
Loading...
Searching...
No Matches
AnasaziThyraAdapter.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Anasazi: Block Eigensolvers Package
4//
5// Copyright 2004 NTESS and the Anasazi contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
15#ifndef ANASAZI_THYRA_ADAPTER_HPP
16#define ANASAZI_THYRA_ADAPTER_HPP
17
20#include "AnasaziConfigDefs.hpp"
21
22#include <Thyra_DetachedMultiVectorView.hpp>
23#include <Thyra_MultiVectorBase.hpp>
24#include <Thyra_MultiVectorStdOps.hpp>
25#include <Thyra_VectorStdOps.hpp>
26
27#include <Teuchos_Ptr.hpp>
28#include <Teuchos_ArrayRCP.hpp>
29#include <Teuchos_ArrayView.hpp>
30
31namespace Anasazi {
32
34 //
35 // Implementation of the Anasazi::MultiVecTraits for Thyra::MultiVectorBase
36 //
38
47 template<class ScalarType>
48 class MultiVecTraits< ScalarType, Thyra::MultiVectorBase<ScalarType> >
49 {
50 private:
51 typedef Thyra::MultiVectorBase<ScalarType> TMVB;
52 typedef Teuchos::ScalarTraits<ScalarType> ST;
53 typedef typename ST::magnitudeType magType;
54
55 public:
56
59
64 static Teuchos::RCP<TMVB> Clone( const TMVB & mv, const int numvecs )
65 {
66 Teuchos::RCP<TMVB> c = Thyra::createMembers( mv.range(), numvecs );
67 return c;
68 }
69
74 static Teuchos::RCP<TMVB>
75 CloneCopy (const TMVB& mv)
76 {
77 const int numvecs = mv.domain()->dim();
78 // create the new multivector
79 Teuchos::RCP< TMVB > cc = Thyra::createMembers (mv.range(), numvecs);
80 // copy the data from the source multivector to the new multivector
81 Thyra::assign (Teuchos::outArg (*cc), mv);
82 return cc;
83 }
84
90 static Teuchos::RCP< TMVB > CloneCopy( const TMVB & mv, const std::vector<int>& index )
91 {
92 const int numvecs = index.size();
93 // create the new multivector
94 Teuchos::RCP<TMVB > cc = Thyra::createMembers (mv.range(), numvecs);
95 // create a view to the relevant part of the source multivector
96 Teuchos::RCP<const TMVB > view = mv.subView ( Teuchos::arrayViewFromVector( index ) );
97 // copy the data from the relevant view to the new multivector
98 Thyra::assign (Teuchos::outArg (*cc), *view);
99 return cc;
100 }
101
102 static Teuchos::RCP<TMVB>
103 CloneCopy (const TMVB& mv, const Teuchos::Range1D& index)
104 {
105 const int numVecs = index.size();
106 // Create the new multivector
107 Teuchos::RCP<TMVB> cc = Thyra::createMembers (mv.range(), numVecs);
108 // Create a view to the relevant part of the source multivector
109 Teuchos::RCP<const TMVB> view = mv.subView (index);
110 // Copy the data from the view to the new multivector.
111 Thyra::assign (Teuchos::outArg (*cc), *view);
112 return cc;
113 }
114
120 static Teuchos::RCP< TMVB > CloneViewNonConst( TMVB & mv, const std::vector<int>& index )
121 {
122 int numvecs = index.size();
123
124 // We do not assume that the indices are sorted, nor do we check that
125 // index.size() > 0. This code is fail-safe, in the sense that a zero
126 // length index vector will pass the error on the Thyra.
127
128 // Thyra has two ways to create an indexed View:
129 // * contiguous (via a range of columns)
130 // * indexed (via a vector of column indices)
131 // The former is significantly more efficient than the latter, in terms of
132 // computations performed with/against the created view.
133 // We will therefore check to see if the given indices are contiguous, and
134 // if so, we will use the contiguous view creation method.
135
136 int lb = index[0];
137 bool contig = true;
138 for (int i=0; i<numvecs; i++) {
139 if (lb+i != index[i]) contig = false;
140 }
141
142 Teuchos::RCP< TMVB > cc;
143 if (contig) {
144 const Thyra::Range1D rng(lb,lb+numvecs-1);
145 // create a contiguous view to the relevant part of the source multivector
146 cc = mv.subView(rng);
147 }
148 else {
149 // create an indexed view to the relevant part of the source multivector
150 cc = mv.subView( Teuchos::arrayViewFromVector( index ) );
151 }
152 return cc;
153 }
154
155 static Teuchos::RCP<TMVB>
156 CloneViewNonConst (TMVB& mv, const Teuchos::Range1D& index)
157 {
158 // We let Thyra be responsible for checking that the index range
159 // is nonempty.
160 //
161 // Create and return a contiguous view to the relevant part of
162 // the source multivector.
163 return mv.subView (index);
164 }
165
171 static Teuchos::RCP<const TMVB > CloneView( const TMVB & mv, const std::vector<int>& index )
172 {
173 int numvecs = index.size();
174
175 // We do not assume that the indices are sorted, nor do we check that
176 // index.size() > 0. This code is fail-safe, in the sense that a zero
177 // length index vector will pass the error on the Thyra.
178
179 // Thyra has two ways to create an indexed View:
180 // * contiguous (via a range of columns)
181 // * indexed (via a vector of column indices)
182 // The former is significantly more efficient than the latter, in terms of
183 // computations performed with/against the created view.
184 // We will therefore check to see if the given indices are contiguous, and
185 // if so, we will use the contiguous view creation method.
186
187 int lb = index[0];
188 bool contig = true;
189 for (int i=0; i<numvecs; i++) {
190 if (lb+i != index[i]) contig = false;
191 }
192
193 Teuchos::RCP< const TMVB > cc;
194 if (contig) {
195 const Thyra::Range1D rng(lb,lb+numvecs-1);
196 // create a contiguous view to the relevant part of the source multivector
197 cc = mv.subView(rng);
198 }
199 else {
200 // create an indexed view to the relevant part of the source multivector
201 cc = mv.subView(Teuchos::arrayViewFromVector( index ) );
202 }
203 return cc;
204 }
205
206 static Teuchos::RCP<const TMVB>
207 CloneView (const TMVB& mv, const Teuchos::Range1D& index)
208 {
209 // We let Thyra be responsible for checking that the index range
210 // is nonempty.
211 //
212 // Create and return a contiguous view to the relevant part of
213 // the source multivector.
214 return mv.subView (index);
215 }
216
218
221
223 static ptrdiff_t GetGlobalLength( const TMVB & mv )
224 { return mv.range()->dim(); }
225
227 static int GetNumberVecs( const TMVB & mv )
228 { return mv.domain()->dim(); }
229
231
234
237 static void MvTimesMatAddMv( const ScalarType alpha, const TMVB & A,
238 const Teuchos::SerialDenseMatrix<int,ScalarType>& B,
239 const ScalarType beta, TMVB & mv )
240 {
241 int m = B.numRows();
242 int n = B.numCols();
243 // Create a view of the B object!
244 Teuchos::RCP< const TMVB >
245 B_thyra = Thyra::createMembersView(
246 A.domain(),
247 RTOpPack::ConstSubMultiVectorView<ScalarType>(
248 0, m, 0, n,
249 Teuchos::arcpFromArrayView(Teuchos::arrayView(&B(0,0), B.stride()*B.numCols())), B.stride()
250 )
251 );
252 // perform the operation via A: mv <- alpha*A*B_thyra + beta*mv
253 A.apply(Thyra::NOTRANS,*B_thyra,Teuchos::outArg (mv),alpha,beta);
254 }
255
258 static void MvAddMv( const ScalarType alpha, const TMVB & A,
259 const ScalarType beta, const TMVB & B, TMVB & mv )
260 {
261 using Teuchos::tuple; using Teuchos::ptrInArg; using Teuchos::inoutArg;
262
263 Thyra::linear_combination<ScalarType>(
264 tuple(alpha, beta)(), tuple(ptrInArg(A), ptrInArg(B))(), ST::zero(), inoutArg(mv));
265 }
266
269 static void MvTransMv( const ScalarType alpha, const TMVB & A, const TMVB & mv,
270 Teuchos::SerialDenseMatrix<int,ScalarType>& B )
271 {
272 // Create a multivector to hold the result (m by n)
273 int m = A.domain()->dim();
274 int n = mv.domain()->dim();
275 // Create a view of the B object!
276 Teuchos::RCP< TMVB >
277 B_thyra = Thyra::createMembersView(
278 A.domain(),
279 RTOpPack::SubMultiVectorView<ScalarType>(
280 0, m, 0, n,
281 Teuchos::arcpFromArrayView(Teuchos::arrayView(&B(0,0), B.stride()*B.numCols())), B.stride()
282 )
283 );
284 A.apply(Thyra::CONJTRANS,mv,B_thyra.ptr(),alpha,Teuchos::ScalarTraits<ScalarType>::zero());
285 }
286
290 static void MvDot( const TMVB & mv, const TMVB & A, std::vector<ScalarType> &b )
291 { Thyra::dots(mv,A,Teuchos::arrayViewFromVector( b )); }
292
295 static void
296 MvScale (TMVB& mv,
297 const ScalarType alpha)
298 {
299 Thyra::scale (alpha, Teuchos::inOutArg (mv));
300 }
301
304 static void
305 MvScale (TMVB& mv,
306 const std::vector<ScalarType>& alpha)
307 {
308 for (unsigned int i=0; i<alpha.size(); i++) {
309 Thyra::scale (alpha[i], mv.col(i).ptr());
310 }
311 }
312
314
317
321 static void MvNorm( const TMVB & mv, std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> &normvec )
322 { Thyra::norms_2(mv,Teuchos::arrayViewFromVector( normvec )); }
323
325
328
331 static void SetBlock( const TMVB & A, const std::vector<int>& index, TMVB & mv )
332 {
333 // Extract the "numvecs" columns of mv indicated by the index vector.
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++) {
338 indexA[i] = i;
339 }
340 // Thyra::assign requires that both arguments have the same number of
341 // vectors. Enforce this, by shrinking one to match the other.
342 if ( numAcols < numvecs ) {
343 // A does not have enough columns to satisfy index_plus. Shrink
344 // index_plus.
345 numvecs = numAcols;
346 }
347 else if ( numAcols > numvecs ) {
348 numAcols = numvecs;
349 indexA.resize( numAcols );
350 }
351 // create a view to the relevant part of the source multivector
352 Teuchos::RCP< const TMVB > relsource = A.subView( Teuchos::arrayViewFromVector( indexA ) );
353 // create a view to the relevant part of the destination multivector
354 Teuchos::RCP< TMVB > reldest = mv.subView( Teuchos::arrayViewFromVector( index ) );
355 // copy the data to the destination multivector subview
356 Thyra::assign (Teuchos::outArg (*reldest), *relsource);
357 }
358
359 static void
360 SetBlock (const TMVB& A, const Teuchos::Range1D& index, TMVB& mv)
361 {
362 const int numColsA = A.domain()->dim();
363 const int numColsMv = mv.domain()->dim();
364 // 'index' indexes into mv; it's the index set of the target.
365 const bool validIndex = index.lbound() >= 0 && index.ubound() < numColsMv;
366 // We can't take more columns out of A than A has.
367 const bool validSource = index.size() <= numColsA;
368
369 if (! validIndex || ! validSource)
370 {
371 std::ostringstream os;
372 os << "Anasazi::MultiVecTraits<Scalar, Thyra::MultiVectorBase<Scalar> "
373 ">::SetBlock(A, [" << index.lbound() << ", " << index.ubound()
374 << "], mv): ";
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!");
386 }
387
388 // View of the relevant column(s) of the target multivector mv.
389 // We avoid view creation overhead by only creating a view if
390 // the index range is different than [0, (# columns in mv) - 1].
391 Teuchos::RCP<TMVB> mv_view;
392 if (index.lbound() == 0 && index.ubound()+1 == numColsMv)
393 mv_view = Teuchos::rcpFromRef (mv); // Non-const, non-owning RCP
394 else
395 mv_view = mv.subView (index);
396
397 // View of the relevant column(s) of the source multivector A.
398 // If A has fewer columns than mv_view, then create a view of
399 // the first index.size() columns of A.
400 Teuchos::RCP<const TMVB> A_view;
401 if (index.size() == numColsA)
402 A_view = Teuchos::rcpFromRef (A); // Const, non-owning RCP
403 else
404 A_view = A.subView (Teuchos::Range1D(0, index.size()-1));
405
406 // Copy the data to the destination multivector.
407 Thyra::assign (Teuchos::outArg (*mv_view), *A_view);
408 }
409
410 static void
411 Assign (const TMVB& A, TMVB& mv)
412 {
413 using Teuchos::ptr;
414 using Teuchos::Range1D;
415 using Teuchos::RCP;
416
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!");
428 }
429 // Copy the data to the destination multivector.
430 if (numColsA == numColsMv) {
431 Thyra::assign (outArg (mv), A);
432 }
433 else {
434 RCP<TMVB> mv_view = CloneViewNonConst (mv, Range1D(0, numColsA-1));
435 Thyra::assign (outArg (*mv_view), A);
436 }
437 }
438
441 static void MvRandom( TMVB & mv )
442 {
443 // Thyra::randomize generates via a uniform distribution on [l,u]
444 // We will use this to generate on [-1,1]
445 Thyra::randomize(-Teuchos::ScalarTraits<ScalarType>::one(),
446 Teuchos::ScalarTraits<ScalarType>::one(),
447 Teuchos::outArg (mv));
448 }
449
452 static void
453 MvInit (TMVB& mv,
454 ScalarType alpha = Teuchos::ScalarTraits<ScalarType>::zero())
455 {
456 Thyra::assign (Teuchos::outArg (mv), alpha);
457 }
458
460
463
466 static void MvPrint( const TMVB & mv, std::ostream& os )
467 {
468 Teuchos::RCP<Teuchos::FancyOStream> out = Teuchos::getFancyOStream(Teuchos::rcp(&os,false));
469 out->setf(std::ios_base::scientific);
470 out->precision(16);
471 mv.describe(*out,Teuchos::VERB_EXTREME);
472 }
473
475
476 };
477
478
480 //
481 // Implementation of the Anasazi::OperatorTraits for Thyra::LinearOpBase
482 //
484
494 template <class ScalarType>
495 class OperatorTraits < ScalarType, Thyra::MultiVectorBase<ScalarType>, Thyra::LinearOpBase<ScalarType> >
496 {
497 public:
498
502 static void Apply ( const Thyra::LinearOpBase< ScalarType >& Op, const Thyra::MultiVectorBase< ScalarType > & x, Thyra::MultiVectorBase< ScalarType > & y )
503 {
504 Op.apply(Thyra::NOTRANS,x,Teuchos::outArg (y), Teuchos::ScalarTraits<ScalarType>::one(),Teuchos::ScalarTraits<ScalarType>::zero());
505 }
506
507 };
508
509} // end of Anasazi namespace
510
511#endif
512// end of file ANASAZI_THYRA_ADAPTER_HPP
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.