Stratimikos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
BelosThyraAdapter.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Stratimikos: Thyra-based strategies for linear solvers
5 // Copyright (2006) 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 (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
52 #ifndef BELOS_THYRA_ADAPTER_HPP
53 #define BELOS_THYRA_ADAPTER_HPP
54 
55 #include "BelosConfigDefs.hpp"
56 #include "BelosMultiVecTraits.hpp"
57 #include "BelosOperatorTraits.hpp"
58 
59 #include <Thyra_DetachedMultiVectorView.hpp>
60 #include <Thyra_MultiVectorBase.hpp>
61 #include <Thyra_MultiVectorStdOps.hpp>
62 #ifdef HAVE_BELOS_TSQR
63 # include <Thyra_TsqrAdaptor.hpp>
64 #endif // HAVE_BELOS_TSQR
65 
66 #include <Teuchos_TimeMonitor.hpp>
67 
68 namespace Belos {
69 
71  //
72  // Implementation of the Belos::MultiVecTraits for Thyra::MultiVectorBase
73  //
75 
82  template<class ScalarType>
83  class MultiVecTraits< ScalarType, Thyra::MultiVectorBase<ScalarType> >
84  {
85  private:
88  typedef typename ST::magnitudeType magType;
89 
90  public:
91 
94 
99  static Teuchos::RCP<TMVB> Clone( const TMVB& mv, const int numvecs )
100  {
101  Teuchos::RCP<TMVB> c = Thyra::createMembers( mv.range(), numvecs );
102  return c;
103  }
104 
109  static Teuchos::RCP<TMVB> CloneCopy( const TMVB& mv )
110  {
111  int numvecs = mv.domain()->dim();
112  // create the new multivector
113  Teuchos::RCP< TMVB > cc = Thyra::createMembers( mv.range(), numvecs );
114  // copy the data from the source multivector to the new multivector
115  Thyra::assign(cc.ptr(), mv);
116  return cc;
117  }
118 
124  static Teuchos::RCP<TMVB> CloneCopy( const TMVB& mv, const std::vector<int>& index )
125  {
126  int numvecs = index.size();
127  // create the new multivector
128  Teuchos::RCP<TMVB> cc = Thyra::createMembers( mv.range(), numvecs );
129  // create a view to the relevant part of the source multivector
130  Teuchos::RCP<const TMVB> view = mv.subView(index);
131  // copy the data from the relevant view to the new multivector
132  Thyra::assign(cc.ptr(), *view);
133  return cc;
134  }
135 
136  static Teuchos::RCP<TMVB>
137  CloneCopy (const TMVB& mv, const Teuchos::Range1D& index)
138  {
139  const int numVecs = index.size();
140  // Create the new multivector
141  Teuchos::RCP<TMVB> cc = Thyra::createMembers (mv.range(), numVecs);
142  // Create a view to the relevant part of the source multivector
143  Teuchos::RCP<const TMVB> view = mv.subView (index);
144  // Copy the data from the view to the new multivector.
145  Thyra::assign (cc.ptr(), *view);
146  return cc;
147  }
148 
154  static Teuchos::RCP<TMVB> CloneViewNonConst( TMVB& mv, const std::vector<int>& index )
155  {
156  int numvecs = index.size();
157 
158  // We do not assume that the indices are sorted, nor do we check that
159  // index.size() > 0. This code is fail-safe, in the sense that a zero
160  // length index std::vector will pass the error on the Thyra.
161 
162  // Thyra has two ways to create an indexed View:
163  // * contiguous (via a range of columns)
164  // * indexed (via a std::vector of column indices)
165  // The former is significantly more efficient than the latter, in terms of
166  // computations performed with/against the created view.
167  // We will therefore check to see if the given indices are contiguous, and
168  // if so, we will use the contiguous view creation method.
169 
170  int lb = index[0];
171  bool contig = true;
172  for (int i=0; i<numvecs; i++) {
173  if (lb+i != index[i]) contig = false;
174  }
175 
177  if (contig) {
178  const Thyra::Range1D rng(lb,lb+numvecs-1);
179  // create a contiguous view to the relevant part of the source multivector
180  cc = mv.subView(rng);
181  }
182  else {
183  // create an indexed view to the relevant part of the source multivector
184  cc = mv.subView(index);
185  }
186  return cc;
187  }
188 
189  static Teuchos::RCP<TMVB>
190  CloneViewNonConst (TMVB& mv, const Teuchos::Range1D& index)
191  {
192  // We let Thyra be responsible for checking that the index range
193  // is nonempty.
194  //
195  // Create and return a contiguous view to the relevant part of
196  // the source multivector.
197  return mv.subView (index);
198  }
199 
200 
206  static Teuchos::RCP<const TMVB> CloneView( const TMVB& mv, const std::vector<int>& index )
207  {
208  int numvecs = index.size();
209 
210  // We do not assume that the indices are sorted, nor do we check that
211  // index.size() > 0. This code is fail-safe, in the sense that a zero
212  // length index std::vector will pass the error on the Thyra.
213 
214  // Thyra has two ways to create an indexed View:
215  // * contiguous (via a range of columns)
216  // * indexed (via a std::vector of column indices)
217  // The former is significantly more efficient than the latter, in terms of
218  // computations performed with/against the created view.
219  // We will therefore check to see if the given indices are contiguous, and
220  // if so, we will use the contiguous view creation method.
221 
222  int lb = index[0];
223  bool contig = true;
224  for (int i=0; i<numvecs; i++) {
225  if (lb+i != index[i]) contig = false;
226  }
227 
229  if (contig) {
230  const Thyra::Range1D rng(lb,lb+numvecs-1);
231  // create a contiguous view to the relevant part of the source multivector
232  cc = mv.subView(rng);
233  }
234  else {
235  // create an indexed view to the relevant part of the source multivector
236  cc = mv.subView(index);
237  }
238  return cc;
239  }
240 
242  CloneView (const TMVB& mv, const Teuchos::Range1D& index)
243  {
244  // We let Thyra be responsible for checking that the index range
245  // is nonempty.
246  //
247  // Create and return a contiguous view to the relevant part of
248  // the source multivector.
249  return mv.subView (index);
250  }
251 
253 
256 
258  static ptrdiff_t GetGlobalLength( const TMVB& mv ) {
259  return Teuchos::as<ptrdiff_t>(mv.range()->dim());
260  }
261 
263  static int GetNumberVecs( const TMVB& mv )
264  { return mv.domain()->dim(); }
265 
267 
270 
273  static void MvTimesMatAddMv( const ScalarType alpha, const TMVB& A,
275  const ScalarType beta, TMVB& mv )
276  {
277  using Teuchos::arrayView; using Teuchos::arcpFromArrayView;
278 
279  Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("Belos::MVT::MvTimesMatAddMv")));
280 
281  const int m = B.numRows();
282  const int n = B.numCols();
283  // Create a view of the B object!
285  B_thyra = Thyra::createMembersView(
286  A.domain(),
288  0, m, 0, n,
289  arcpFromArrayView(arrayView(&B(0,0), B.stride()*B.numCols())), B.stride()
290  )
291  );
292  // perform the operation via A: mv <- alpha*A*B_thyra + beta*mv
293  Thyra::apply<ScalarType>(A, Thyra::NOTRANS, *B_thyra, Teuchos::outArg(mv), alpha, beta);
294  }
295 
298  static void MvAddMv( const ScalarType alpha, const TMVB& A,
299  const ScalarType beta, const TMVB& B, TMVB& mv )
300  {
302  using Teuchos::tuple; using Teuchos::ptrInArg; using Teuchos::inoutArg;
303 
304  Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("Belos::MVT::MvAddMv")));
305 
306  Thyra::linear_combination<ScalarType>(
307  tuple(alpha, beta)(), tuple(ptrInArg(A), ptrInArg(B))(), ST::zero(), inoutArg(mv));
308  }
309 
312  static void MvScale ( TMVB& mv, const ScalarType alpha )
313  {
314  Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("Belos::MVT::MvScale")));
315 
316  Thyra::scale(alpha, Teuchos::inoutArg(mv));
317  }
318 
321  static void MvScale (TMVB& mv, const std::vector<ScalarType>& alpha)
322  {
323  Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("Belos::MVT::MvScale")));
324 
325  for (unsigned int i=0; i<alpha.size(); i++) {
326  Thyra::scale<ScalarType> (alpha[i], mv.col(i).ptr());
327  }
328  }
329 
332  static void MvTransMv( const ScalarType alpha, const TMVB& A, const TMVB& mv,
334  {
335  using Teuchos::arrayView; using Teuchos::arcpFromArrayView;
336 
337  Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("Belos::MVT::MvTransMv")));
338 
339  // Create a multivector to hold the result (m by n)
340  int m = A.domain()->dim();
341  int n = mv.domain()->dim();
342  // Create a view of the B object!
344  B_thyra = Thyra::createMembersView(
345  A.domain(),
347  0, m, 0, n,
348  arcpFromArrayView(arrayView(&B(0,0), B.stride()*B.numCols())), B.stride()
349  )
350  );
351  Thyra::apply<ScalarType>(A, Thyra::CONJTRANS, mv, B_thyra.ptr(), alpha);
352  }
353 
357  static void MvDot( const TMVB& mv, const TMVB& A, std::vector<ScalarType>& b )
358  {
359  Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("Belos::MVT::MvDot")));
360 
361  Thyra::dots(mv, A, Teuchos::arrayViewFromVector(b));
362  }
363 
365 
368 
372  static void MvNorm( const TMVB& mv, std::vector<magType>& normvec,
373  NormType type = TwoNorm ) {
374  Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("Belos::MVT::MvNorm")));
375 
376  if(type == TwoNorm)
377  Thyra::norms_2(mv, Teuchos::arrayViewFromVector(normvec));
378  else if(type == OneNorm)
379  Thyra::norms_1(mv, Teuchos::arrayViewFromVector(normvec));
380  else if(type == InfNorm)
381  Thyra::norms_inf(mv, Teuchos::arrayViewFromVector(normvec));
382  else
383  TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument,
384  "Belos::MultiVecTraits::MvNorm (Thyra specialization): "
385  "invalid norm type. Must be either TwoNorm, OneNorm or InfNorm");
386  }
387 
389 
392 
395  static void SetBlock( const TMVB& A, const std::vector<int>& index, TMVB& mv )
396  {
397  // Extract the "numvecs" columns of mv indicated by the index std::vector.
398  int numvecs = index.size();
399  std::vector<int> indexA(numvecs);
400  int numAcols = A.domain()->dim();
401  for (int i=0; i<numvecs; i++) {
402  indexA[i] = i;
403  }
404  // Thyra::assign requires that both arguments have the same number of
405  // vectors. Enforce this, by shrinking one to match the other.
406  if ( numAcols < numvecs ) {
407  // A does not have enough columns to satisfy index_plus. Shrink
408  // index_plus.
409  numvecs = numAcols;
410  }
411  else if ( numAcols > numvecs ) {
412  numAcols = numvecs;
413  indexA.resize( numAcols );
414  }
415  // create a view to the relevant part of the source multivector
416  Teuchos::RCP< const TMVB > relsource = A.subView(indexA);
417  // create a view to the relevant part of the destination multivector
418  Teuchos::RCP< TMVB > reldest = mv.subView(index);
419  // copy the data to the destination multivector subview
420  Thyra::assign(reldest.ptr(), *relsource);
421  }
422 
423  static void
424  SetBlock (const TMVB& A, const Teuchos::Range1D& index, TMVB& mv)
425  {
426  const int numColsA = A.domain()->dim();
427  const int numColsMv = mv.domain()->dim();
428  // 'index' indexes into mv; it's the index set of the target.
429  const bool validIndex = index.lbound() >= 0 && index.ubound() < numColsMv;
430  // We can't take more columns out of A than A has.
431  const bool validSource = index.size() <= numColsA;
432 
433  if (! validIndex || ! validSource)
434  {
435  std::ostringstream os;
436  os << "Belos::MultiVecTraits<Scalar, Thyra::MultiVectorBase<Scalar> "
437  ">::SetBlock(A, [" << index.lbound() << ", " << index.ubound()
438  << "], mv): ";
439  TEUCHOS_TEST_FOR_EXCEPTION(index.lbound() < 0, std::invalid_argument,
440  os.str() << "Range lower bound must be nonnegative.");
441  TEUCHOS_TEST_FOR_EXCEPTION(index.ubound() >= numColsMv, std::invalid_argument,
442  os.str() << "Range upper bound must be less than "
443  "the number of columns " << numColsA << " in the "
444  "'mv' output argument.");
445  TEUCHOS_TEST_FOR_EXCEPTION(index.size() > numColsA, std::invalid_argument,
446  os.str() << "Range must have no more elements than"
447  " the number of columns " << numColsA << " in the "
448  "'A' input argument.");
449  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "Should never get here!");
450  }
451 
452  // View of the relevant column(s) of the target multivector mv.
453  // We avoid view creation overhead by only creating a view if
454  // the index range is different than [0, (# columns in mv) - 1].
455  Teuchos::RCP<TMVB> mv_view;
456  if (index.lbound() == 0 && index.ubound()+1 == numColsMv)
457  mv_view = Teuchos::rcpFromRef (mv); // Non-const, non-owning RCP
458  else
459  mv_view = mv.subView (index);
460 
461  // View of the relevant column(s) of the source multivector A.
462  // If A has fewer columns than mv_view, then create a view of
463  // the first index.size() columns of A.
465  if (index.size() == numColsA)
466  A_view = Teuchos::rcpFromRef (A); // Const, non-owning RCP
467  else
468  A_view = A.subView (Teuchos::Range1D(0, index.size()-1));
469 
470  // Copy the data to the destination multivector.
471  Thyra::assign(mv_view.ptr(), *A_view);
472  }
473 
474  static void
475  Assign (const TMVB& A, TMVB& mv)
476  {
477  Teuchos::TimeMonitor tM(*Teuchos::TimeMonitor::getNewTimer(std::string("Belos::MVT::Assign")));
478 
479  const int numColsA = A.domain()->dim();
480  const int numColsMv = mv.domain()->dim();
481  if (numColsA > numColsMv)
482  {
483  std::ostringstream os;
484  os << "Belos::MultiVecTraits<Scalar, Thyra::MultiVectorBase<Scalar>"
485  " >::Assign(A, mv): ";
486  TEUCHOS_TEST_FOR_EXCEPTION(numColsA > numColsMv, std::invalid_argument,
487  os.str() << "Input multivector 'A' has "
488  << numColsA << " columns, but output multivector "
489  "'mv' has only " << numColsMv << " columns.");
490  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error, "Should never get here!");
491  }
492  // Copy the data to the destination multivector.
493  if (numColsA == numColsMv) {
494  Thyra::assign (Teuchos::outArg (mv), A);
495  } else {
496  Teuchos::RCP<TMVB> mv_view =
497  CloneViewNonConst (mv, Teuchos::Range1D(0, numColsA-1));
498  Thyra::assign (mv_view.ptr(), A);
499  }
500  }
501 
504  static void MvRandom( TMVB& mv )
505  {
506  // Thyra::randomize generates via a uniform distribution on [l,u]
507  // We will use this to generate on [-1,1]
508  Thyra::randomize<ScalarType>(
511  Teuchos::outArg(mv));
512  }
513 
515  static void
516  MvInit (TMVB& mv, ScalarType alpha = Teuchos::ScalarTraits<ScalarType>::zero())
517  {
518  Thyra::assign (Teuchos::outArg (mv), alpha);
519  }
520 
522 
525 
528  static void MvPrint( const TMVB& mv, std::ostream& os )
529  { os << describe(mv,Teuchos::VERB_EXTREME); }
530 
532 
533 #ifdef HAVE_BELOS_TSQR
534  typedef Thyra::TsqrAdaptor< ScalarType > tsqr_adaptor_type;
540 #endif // HAVE_BELOS_TSQR
541  };
542 
544  //
545  // Implementation of the Belos::OperatorTraits for Thyra::LinearOpBase
546  //
548 
556  template<class ScalarType>
557  class OperatorTraits <ScalarType,
558  Thyra::MultiVectorBase<ScalarType>,
559  Thyra::LinearOpBase<ScalarType> >
560  {
561  private:
564 
565  public:
581  static void
582  Apply (const TLOB& Op,
583  const TMVB& x,
584  TMVB& y,
585  ETrans trans = NOTRANS)
586  {
587  Thyra::EOpTransp whichOp;
588 
589  // We don't check here whether the operator implements the
590  // requested operation. Call HasApplyTranspose() to check.
591  // Thyra::LinearOpBase implementations are not required to
592  // implement NOTRANS. However, Belos needs NOTRANS
593  // (obviously!), so we assume that Op implements NOTRANS.
594  if (trans == NOTRANS)
595  whichOp = Thyra::NOTRANS;
596  else if (trans == TRANS)
597  whichOp = Thyra::TRANS;
598  else if (trans == CONJTRANS)
599  whichOp = Thyra::CONJTRANS;
600  else
601  TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument,
602  "Belos::OperatorTraits::Apply (Thyra specialization): "
603  "'trans' argument must be neither NOTRANS=" << NOTRANS
604  << ", TRANS=" << TRANS << ", or CONJTRANS=" << CONJTRANS
605  << ", but instead has an invalid value of " << trans << ".");
606  Thyra::apply<ScalarType>(Op, whichOp, x, Teuchos::outArg(y));
607  }
608 
610  static bool HasApplyTranspose (const TLOB& Op)
611  {
613 
614  // Thyra::LinearOpBase's interface lets you check whether the
615  // operator implements any of all four possible combinations of
616  // conjugation and transpose. Belos only needs transpose
617  // (TRANS) if the operator is real; in that case, Apply() does
618  // the same thing with trans = CONJTRANS or TRANS. If the
619  // operator is complex, Belos needs both transpose and conjugate
620  // transpose (CONJTRANS) if the operator is complex.
621  return Op.opSupported (Thyra::TRANS) &&
622  (! STS::isComplex || Op.opSupported (Thyra::CONJTRANS));
623  }
624  };
625 
626 } // end of Belos namespace
627 
628 #endif
629 // end of file BELOS_THYRA_ADAPTER_HPP
static int GetNumberVecs(const TMVB &mv)
Obtain the number of vectors in mv.
virtual RCP< const VectorSpaceBase< Scalar > > range() const =0
EOpTransp
static ptrdiff_t GetGlobalLength(const TMVB &mv)
Obtain the std::vector length of mv.
static void MvPrint(const TMVB &mv, std::ostream &os)
Print the mv multi-std::vector to the os output stream.
Stub adaptor from Thyra::MultiVectorBase to TSQR.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static Teuchos::RCP< const MV > CloneView(const MV &mv, const std::vector< int > &index)
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 std::vector (deep c...
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 void MvInit(MV &mv, const ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
static void MvScale(TMVB &mv, const std::vector< ScalarType > &alpha)
Scale each element of the i-th vector in *this with alpha[i].
static void MvDot(const TMVB &mv, const TMVB &A, std::vector< ScalarType > &b)
Compute a std::vector b where the components are the individual dot-products of the i-th columns of A...
Ordinal ubound() const
static void Assign(const MV &A, MV &mv)
static Teuchos::RCP< TMVB > CloneCopy(const TMVB &mv)
Creates a new MultiVectorBase and copies contents of mv into the new std::vector (deep copy)...
static RCP< Time > getNewTimer(const std::string &name)
static void MvScale(TMVB &mv, const ScalarType alpha)
Scale each element of the vectors in *this with alpha.
static void SetBlock(const MV &A, const std::vector< int > &index, MV &mv)
RCP< const MultiVectorBase< Scalar > > subView(const Range1D &colRng) const
static Teuchos::RCP< MV > CloneViewNonConst(MV &mv, const std::vector< int > &index)
static void MvTimesMatAddMv(const ScalarType alpha, const TMVB &A, const Teuchos::SerialDenseMatrix< int, ScalarType > &B, const ScalarType beta, TMVB &mv)
Update mv with .
Ptr< T > ptr() const
RCP< const VectorBase< Scalar > > col(Ordinal j) const
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 void Apply(const TLOB &Op, const TMVB &x, TMVB &y, ETrans trans=NOTRANS)
Apply Op to x, storing the result in y.
static void MvRandom(TMVB &mv)
Replace the vectors in mv with random vectors.
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Ordinal lbound() const
OrdinalType numCols() const
virtual RCP< const VectorSpaceBase< Scalar > > domain() const =0
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...
Ordinal size() const
static void MvNorm(const TMVB &mv, std::vector< magType > &normvec, NormType type=TwoNorm)
Compute the 2-norm of each individual std::vector of mv. Upon return, normvec[i] holds the value of ...
static bool HasApplyTranspose(const TLOB &Op)
Whether the operator implements applying the transpose.
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).
OrdinalType stride() const
static Teuchos::RCP< TMVB > Clone(const TMVB &mv, const int numvecs)
Creates a new empty MultiVectorBase containing numvecs columns.
static void MvAddMv(const ScalarType alpha, const TMVB &A, const ScalarType beta, const TMVB &B, TMVB &mv)
Replace mv with .
OrdinalType numRows() const
Teuchos::Range1D Range1D

Generated on Fri Jun 5 2020 10:13:28 for Stratimikos by doxygen 1.8.5