Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_TpetraMultiVector_def.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_TPETRA_MULTIVECTOR_HPP
11 #define THYRA_TPETRA_MULTIVECTOR_HPP
12 
13 #include "Thyra_TpetraMultiVector_decl.hpp"
14 #include "Thyra_TpetraVectorSpace.hpp"
15 #include "Thyra_TpetraVector.hpp"
16 #include "Teuchos_Assert.hpp"
17 
18 
19 namespace Thyra {
20 
21 
22 // Constructors/initializers/accessors
23 
24 
25 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
27 {}
28 
29 
30 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
32  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
33  const RCP<const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
35  )
36 {
37  initializeImpl(tpetraVectorSpace, domainSpace, tpetraMultiVector);
38 }
39 
40 
41 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
43  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
44  const RCP<const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
46  )
47 {
48  initializeImpl(tpetraVectorSpace, domainSpace, tpetraMultiVector);
49 }
50 
51 
52 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
55 {
56  return tpetraMultiVector_.getNonconstObj();
57 }
58 
59 
60 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
63 {
64  return tpetraMultiVector_;
65 }
66 
67 
68 // Overridden public functions form MultiVectorAdapterBase
69 
70 
71 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
74 {
75  return domainSpace_;
76 }
77 
78 
79 // Overridden protected functions from MultiVectorBase
80 
81 
82 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
83 void
85 {
86  tpetraMultiVector_.getNonconstObj()->putScalar(alpha);
87 }
88 
89 
90 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
93 {
94  auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
95 
96  // If the cast succeeded, call Tpetra directly.
97  // Otherwise, fall back to the RTOp implementation.
98  if (nonnull(tmv)) {
99  tpetraMultiVector_.getNonconstObj()->assign(*tmv);
100  } else {
102  }
103 }
104 
105 
106 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
107 void
109 {
110  tpetraMultiVector_.getNonconstObj()->scale(alpha);
111 }
112 
113 
114 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
116  Scalar alpha,
117  const MultiVectorBase<Scalar>& mv
118  )
119 {
120  auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
121 
122  // If the cast succeeded, call Tpetra directly.
123  // Otherwise, fall back to the RTOp implementation.
124  if (nonnull(tmv)) {
126  tpetraMultiVector_.getNonconstObj()->update(alpha, *tmv, ST::one());
127  } else {
129  }
130 }
131 
132 
133 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
135  const ArrayView<const Scalar>& alpha,
136  const ArrayView<const Ptr<const MultiVectorBase<Scalar> > >& mv,
137  const Scalar& beta
138  )
139 {
140 #ifdef TEUCHOS_DEBUG
141  TEUCHOS_ASSERT_EQUALITY(alpha.size(), mv.size());
142 #endif
143 
144  // Try to cast mv to an array of this type
146  Teuchos::Array<RCP<const TMV> > tmvs(mv.size());
147  RCP<const TMV> tmv;
148  bool allCastsSuccessful = true;
149  {
150  auto mvIter = mv.begin();
151  auto tmvIter = tmvs.begin();
152  for (; mvIter != mv.end(); ++mvIter, ++tmvIter) {
153  tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromPtr(*mvIter));
154  if (nonnull(tmv)) {
155  *tmvIter = tmv;
156  } else {
157  allCastsSuccessful = false;
158  break;
159  }
160  }
161  }
162 
163  // If casts succeeded, or input arrays are size 0, call Tpetra directly.
164  // Otherwise, fall back to the RTOp implementation.
165  auto len = tmvs.size();
166  if (len == 0) {
167  tpetraMultiVector_.getNonconstObj()->scale(beta);
168  } else if (len == 1 && allCastsSuccessful) {
169  tpetraMultiVector_.getNonconstObj()->update(alpha[0], *tmvs[0], beta);
170  } else if (len == 2 && allCastsSuccessful) {
171  tpetraMultiVector_.getNonconstObj()->update(alpha[0], *tmvs[0], alpha[1], *tmvs[1], beta);
172  } else if (allCastsSuccessful) {
174  auto tmvIter = tmvs.begin();
175  auto alphaIter = alpha.begin();
176 
177  // Check if any entry of tmvs aliases this object's wrapped vector.
178  // If so, replace that entry in the array with a copy.
179  tmv = Teuchos::null;
180  for (; tmvIter != tmvs.end(); ++tmvIter) {
181  if (tmvIter->getRawPtr() == tpetraMultiVector_.getConstObj().getRawPtr()) {
182  if (tmv.is_null()) {
183  tmv = Teuchos::rcp(new TMV(*tpetraMultiVector_.getConstObj(), Teuchos::Copy));
184  }
185  *tmvIter = tmv;
186  }
187  }
188  tmvIter = tmvs.begin();
189 
190  // We add two MVs at a time, so only scale if even num MVs,
191  // and additionally do the first addition if odd num MVs.
192  if ((tmvs.size() % 2) == 0) {
193  tpetraMultiVector_.getNonconstObj()->scale(beta);
194  } else {
195  tpetraMultiVector_.getNonconstObj()->update(*alphaIter, *(*tmvIter), beta);
196  ++tmvIter;
197  ++alphaIter;
198  }
199  for (; tmvIter != tmvs.end(); tmvIter+=2, alphaIter+=2) {
200  tpetraMultiVector_.getNonconstObj()->update(
201  *alphaIter, *(*tmvIter), *(alphaIter+1), *(*(tmvIter+1)), ST::one());
202  }
203  } else {
205  }
206 }
207 
208 
209 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
211  const MultiVectorBase<Scalar>& mv,
212  const ArrayView<Scalar>& prods
213  ) const
214 {
215  auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
216 
217  // If the cast succeeded, call Tpetra directly.
218  // Otherwise, fall back to the RTOp implementation.
219  if (nonnull(tmv)) {
220  tpetraMultiVector_.getConstObj()->dot(*tmv, prods);
221  } else {
223  }
224 }
225 
226 
227 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
229  const ArrayView<typename ScalarTraits<Scalar>::magnitudeType>& norms
230  ) const
231 {
232  tpetraMultiVector_.getConstObj()->norm1(norms);
233 }
234 
235 
236 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
238  const ArrayView<typename ScalarTraits<Scalar>::magnitudeType>& norms
239  ) const
240 {
241  tpetraMultiVector_.getConstObj()->norm2(norms);
242 }
243 
244 
245 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
247  const ArrayView<typename ScalarTraits<Scalar>::magnitudeType>& norms
248  ) const
249 {
250  tpetraMultiVector_.getConstObj()->normInf(norms);
251 }
252 
253 
254 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
257 {
258 #ifdef TEUCHOS_DEBUG
259  TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE(j, 0, this->domain()->dim());
260 #endif
261  return constTpetraVector<Scalar>(
262  tpetraVectorSpace_,
263  tpetraMultiVector_->getVector(j)
264  );
265 }
266 
267 
268 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
271 {
272 #ifdef TEUCHOS_DEBUG
273  TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE(j, 0, this->domain()->dim());
274 #endif
275  return tpetraVector<Scalar>(
276  tpetraVectorSpace_,
277  tpetraMultiVector_.getNonconstObj()->getVectorNonConst(j)
278  );
279 }
280 
281 
282 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
285  const Range1D& col_rng_in
286  ) const
287 {
288 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
289  std::cerr << "\nTpetraMultiVector::subView(Range1D) const called!\n";
290 #endif
291  const Range1D colRng = this->validateColRange(col_rng_in);
292 
294  this->getConstTpetraMultiVector()->subView(colRng);
295 
296  const RCP<const ScalarProdVectorSpaceBase<Scalar> > viewDomainSpace =
297  tpetraVectorSpace<Scalar>(
298  Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
299  tpetraView->getNumVectors(),
300  tpetraView->getMap()->getComm()
301  )
302  );
303 
304  return constTpetraMultiVector(
305  tpetraVectorSpace_,
306  viewDomainSpace,
307  tpetraView
308  );
309 }
310 
311 
312 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
315  const Range1D& col_rng_in
316  )
317 {
318 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
319  std::cerr << "\nTpetraMultiVector::subView(Range1D) called!\n";
320 #endif
321  const Range1D colRng = this->validateColRange(col_rng_in);
322 
324  this->getTpetraMultiVector()->subViewNonConst(colRng);
325 
326  const RCP<const ScalarProdVectorSpaceBase<Scalar> > viewDomainSpace =
327  tpetraVectorSpace<Scalar>(
328  Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
329  tpetraView->getNumVectors(),
330  tpetraView->getMap()->getComm()
331  )
332  );
333 
334  return tpetraMultiVector(
335  tpetraVectorSpace_,
336  viewDomainSpace,
337  tpetraView
338  );
339 }
340 
341 
342 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
345  const ArrayView<const int>& cols_in
346  ) const
347 {
348 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
349  std::cerr << "\nTpetraMultiVector::subView(ArrayView) const called!\n";
350 #endif
351  // Tpetra wants col indices as size_t
352  Array<std::size_t> cols(cols_in.size());
353  for (Array<std::size_t>::size_type i = 0; i < cols.size(); ++i)
354  cols[i] = static_cast<std::size_t>(cols_in[i]);
355 
357  this->getConstTpetraMultiVector()->subView(cols());
358 
359  const RCP<const ScalarProdVectorSpaceBase<Scalar> > viewDomainSpace =
360  tpetraVectorSpace<Scalar>(
361  Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
362  tpetraView->getNumVectors(),
363  tpetraView->getMap()->getComm()
364  )
365  );
366 
367  return constTpetraMultiVector(
368  tpetraVectorSpace_,
369  viewDomainSpace,
370  tpetraView
371  );
372 }
373 
374 
375 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
378  const ArrayView<const int>& cols_in
379  )
380 {
381 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
382  std::cerr << "\nTpetraMultiVector::subView(ArrayView) called!\n";
383 #endif
384  // Tpetra wants col indices as size_t
385  Array<std::size_t> cols(cols_in.size());
386  for (Array<std::size_t>::size_type i = 0; i < cols.size(); ++i)
387  cols[i] = static_cast<std::size_t>(cols_in[i]);
388 
390  this->getTpetraMultiVector()->subViewNonConst(cols());
391 
392  const RCP<const ScalarProdVectorSpaceBase<Scalar> > viewDomainSpace =
393  tpetraVectorSpace<Scalar>(
394  Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
395  tpetraView->getNumVectors(),
396  tpetraView->getMap()->getComm()
397  )
398  );
399 
400  return tpetraMultiVector(
401  tpetraVectorSpace_,
402  viewDomainSpace,
403  tpetraView
404  );
405 }
406 
407 
408 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
411  const RTOpPack::RTOpT<Scalar> &primary_op,
412  const ArrayView<const Ptr<const MultiVectorBase<Scalar> > > &multi_vecs,
413  const ArrayView<const Ptr<MultiVectorBase<Scalar> > > &targ_multi_vecs,
414  const ArrayView<const Ptr<RTOpPack::ReductTarget> > &reduct_objs,
415  const Ordinal primary_global_offset
416  ) const
417 {
418 
420  primary_op, multi_vecs, targ_multi_vecs, reduct_objs, primary_global_offset);
421 }
422 
423 
424 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
427  const Range1D &rowRng,
428  const Range1D &colRng,
430  ) const
431 {
433  acquireDetachedMultiVectorViewImpl(rowRng, colRng, sub_mv);
434 }
435 
436 
437 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
440  const Range1D &rowRng,
441  const Range1D &colRng,
443  )
444 {
446  acquireNonconstDetachedMultiVectorViewImpl(rowRng, colRng, sub_mv);
447 }
448 
449 
450 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
454  )
455 {
458 
459 }
460 
461 
462 /* ToDo: Implement these?
463 
464 
465 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
466 RCP<const MultiVectorBase<Scalar> >
467 TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::nonContigSubViewImpl(
468  const ArrayView<const int> &cols
469  ) const
470 {
471  THYRA_DEBUG_ASSERT_MV_COLS("nonContigSubViewImpl(cols)", cols);
472  const int numCols = cols.size();
473  const ArrayRCP<Scalar> localValuesView = createContiguousCopy(cols);
474  return defaultSpmdMultiVector<Scalar>(
475  spmdRangeSpace_,
476  createSmallScalarProdVectorSpaceBase<Scalar>(spmdRangeSpace_, numCols),
477  localValuesView
478  );
479 }
480 
481 
482 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
483 RCP<MultiVectorBase<Scalar> >
484 TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::nonconstNonContigSubViewImpl(
485  const ArrayView<const int> &cols )
486 {
487  THYRA_DEBUG_ASSERT_MV_COLS("nonContigSubViewImpl(cols)", cols);
488  const int numCols = cols.size();
489  const ArrayRCP<Scalar> localValuesView = createContiguousCopy(cols);
490  const Ordinal localSubDim = spmdRangeSpace_->localSubDim();
491  RCP<CopyBackSpmdMultiVectorEntries<Scalar> > copyBackView =
492  copyBackSpmdMultiVectorEntries<Scalar>(cols, localValuesView.getConst(),
493  localSubDim, localValues_.create_weak(), leadingDim_);
494  return Teuchos::rcpWithEmbeddedObjPreDestroy(
495  new TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>(
496  spmdRangeSpace_,
497  createSmallScalarProdVectorSpaceBase<Scalar>(spmdRangeSpace_, numCols),
498  localValuesView),
499  copyBackView
500  );
501 }
502 
503 */
504 
505 
506 // Overridden protected members from SpmdMultiVectorBase
507 
508 
509 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
512 {
513  return tpetraVectorSpace_;
514 }
515 
516 
517 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
519  const Ptr<ArrayRCP<Scalar> > &localValues, const Ptr<Ordinal> &leadingDim
520  )
521 {
522  *localValues = tpetraMultiVector_.getNonconstObj()->get1dViewNonConst();
523  *leadingDim = tpetraMultiVector_->getStride();
524 }
525 
526 
527 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
529  const Ptr<ArrayRCP<const Scalar> > &localValues, const Ptr<Ordinal> &leadingDim
530  ) const
531 {
532  *localValues = tpetraMultiVector_->get1dView();
533  *leadingDim = tpetraMultiVector_->getStride();
534 }
535 
536 
537 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
539  const EOpTransp M_trans,
540  const MultiVectorBase<Scalar> &X,
541  const Ptr<MultiVectorBase<Scalar> > &Y,
542  const Scalar alpha,
543  const Scalar beta
544  ) const
545 {
546  // Try to extract Tpetra objects from X and Y
548  Teuchos::RCP<const TMV> X_tpetra = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(X));
549  Teuchos::RCP<TMV> Y_tpetra = this->getTpetraMultiVector(Teuchos::rcpFromPtr(Y));
550 
551  // If the cast succeeded, call Tpetra directly.
552  // Otherwise, fall back to the default implementation.
553  if (nonnull(X_tpetra) && nonnull(Y_tpetra)) {
555  TEUCHOS_TEST_FOR_EXCEPTION(ST::isComplex && (M_trans == CONJ),
556  std::logic_error,
557  "Error, conjugation without transposition is not allowed for complex scalar types!");
558 
560  switch (M_trans) {
561  case NOTRANS:
562  trans = Teuchos::NO_TRANS;
563  break;
564  case CONJ:
565  trans = Teuchos::NO_TRANS;
566  break;
567  case TRANS:
568  trans = Teuchos::TRANS;
569  break;
570  case CONJTRANS:
571  trans = Teuchos::CONJ_TRANS;
572  break;
573  }
574 
575  Y_tpetra->multiply(trans, Teuchos::NO_TRANS, alpha, *tpetraMultiVector_.getConstObj(), *X_tpetra, beta);
576  Kokkos::fence();
577  } else {
578  SpmdMultiVectorDefaultBase<Scalar>::euclideanApply(M_trans, X, Y, alpha, beta);
579  }
580 
581 }
582 
583 // private
584 
585 
586 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
587 template<class TpetraMultiVector_t>
589  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
590  const RCP<const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
591  const RCP<TpetraMultiVector_t> &tpetraMultiVector
592  )
593 {
594 #ifdef THYRA_DEBUG
595  TEUCHOS_ASSERT(nonnull(tpetraVectorSpace));
596  TEUCHOS_ASSERT(nonnull(domainSpace));
597  TEUCHOS_ASSERT(nonnull(tpetraMultiVector));
598  // ToDo: Check to make sure that tpetraMultiVector is compatible with
599  // tpetraVectorSpace.
600 #endif
601  tpetraVectorSpace_ = tpetraVectorSpace;
602  domainSpace_ = domainSpace;
603  tpetraMultiVector_.initialize(tpetraMultiVector);
604  this->updateSpmdSpace();
605 }
606 
607 
608 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
609 RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
611 getTpetraMultiVector(const RCP<MultiVectorBase<Scalar> >& mv) const
612 {
613  using Teuchos::rcp_dynamic_cast;
616 
617  RCP<TMV> tmv = rcp_dynamic_cast<TMV>(mv);
618  if (nonnull(tmv)) {
619  return tmv->getTpetraMultiVector();
620  }
621 
622  RCP<TV> tv = rcp_dynamic_cast<TV>(mv);
623  if (nonnull(tv)) {
624  return tv->getTpetraVector();
625  }
626 
627  return Teuchos::null;
628 }
629 
630 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
631 RCP<const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
633 getConstTpetraMultiVector(const RCP<const MultiVectorBase<Scalar> >& mv) const
634 {
635  using Teuchos::rcp_dynamic_cast;
638 
639  RCP<const TMV> tmv = rcp_dynamic_cast<const TMV>(mv);
640  if (nonnull(tmv)) {
641  return tmv->getConstTpetraMultiVector();
642  }
643 
644  RCP<const TV> tv = rcp_dynamic_cast<const TV>(mv);
645  if (nonnull(tv)) {
646  return tv->getConstTpetraVector();
647  }
648 
649  return Teuchos::null;
650 }
651 
652 
653 } // end namespace Thyra
654 
655 
656 #endif // THYRA_TPETRA_MULTIVECTOR_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)
void euclideanApply(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
Uses GEMM() and Teuchos::reduceAll() to implement.
iterator begin() const
virtual void assignMultiVecImpl(const MultiVectorBase< Scalar > &mv)
Default implementation of assign(MV) using RTOps.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Use the non-transposed operator.
virtual void euclideanApply(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
size_type size() 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
virtual void linearCombinationImpl(const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta)
Default implementation of linear_combination using RTOps.
Use the transposed operator with complex-conjugate clements (same as TRANS for real scalar types)...
RCP< const MultiVectorBase< Scalar > > nonContigSubViewImpl(const ArrayView< const int > &cols_in) const
Use the non-transposed operator with complex-conjugate elements (same as NOTRANS for real scalar type...
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.
Use the transposed operator.
#define TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE(index, lower_inclusive, upper_exclusive)
virtual void scaleImpl(Scalar alpha)
Ordinal size_type
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< 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.
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)
Concrete Thyra::SpmdVectorBase using Tpetra::Vector.
void acquireDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
virtual void dotsImpl(const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
RCP< VectorBase< Scalar > > nonconstColImpl(Ordinal j)
void commitNonconstDetachedMultiVectorViewImpl(RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
bool nonnull(const boost::shared_ptr< T > &p)
void getLocalMultiVectorDataImpl(const Ptr< ArrayRCP< const Scalar > > &localValues, const Ptr< Ordinal > &leadingDim) const
size_type size() 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
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
void acquireNonconstDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
virtual void normsInfImpl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
iterator begin()
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
virtual void dotsImpl(const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
Default implementation of dots using RTOps.
void commitNonconstDetachedMultiVectorViewImpl(RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
virtual void updateImpl(Scalar alpha, const MultiVectorBase< Scalar > &mv)
Default implementation of update using RTOps.