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 //
4 // Thyra: Interfaces and Support for Abstract Numerical Algorithms
5 // Copyright (2004) 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 (bartlettra@ornl.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef THYRA_TPETRA_MULTIVECTOR_HPP
43 #define THYRA_TPETRA_MULTIVECTOR_HPP
44 
45 #include "Thyra_TpetraMultiVector_decl.hpp"
46 #include "Thyra_TpetraVectorSpace.hpp"
47 #include "Thyra_TpetraVector.hpp"
48 #include "Teuchos_Assert.hpp"
49 
50 
51 namespace Thyra {
52 
53 
54 // Constructors/initializers/accessors
55 
56 
57 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
59 {}
60 
61 
62 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
64  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
65  const RCP<const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
67  )
68 {
69  initializeImpl(tpetraVectorSpace, domainSpace, tpetraMultiVector);
70 }
71 
72 
73 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
75  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
76  const RCP<const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
78  )
79 {
80  initializeImpl(tpetraVectorSpace, domainSpace, tpetraMultiVector);
81 }
82 
83 
84 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
87 {
88  return tpetraMultiVector_.getNonconstObj();
89 }
90 
91 
92 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
95 {
96  return tpetraMultiVector_;
97 }
98 
99 
100 // Overridden public functions form MultiVectorAdapterBase
101 
102 
103 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
106 {
107  return domainSpace_;
108 }
109 
110 
111 // Overridden protected functions from MultiVectorBase
112 
113 
114 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
115 void
117 {
118  tpetraMultiVector_.getNonconstObj()->putScalar(alpha);
119 }
120 
121 
122 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
125 {
126  auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
127 
128  // If the cast succeeded, call Tpetra directly.
129  // Otherwise, fall back to the RTOp implementation.
130  if (nonnull(tmv)) {
131  tpetraMultiVector_.getNonconstObj()->assign(*tmv);
132  } else {
134  }
135 }
136 
137 
138 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
139 void
141 {
142  tpetraMultiVector_.getNonconstObj()->scale(alpha);
143 }
144 
145 
146 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
148  Scalar alpha,
149  const MultiVectorBase<Scalar>& mv
150  )
151 {
152  auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
153 
154  // If the cast succeeded, call Tpetra directly.
155  // Otherwise, fall back to the RTOp implementation.
156  if (nonnull(tmv)) {
158  tpetraMultiVector_.getNonconstObj()->update(alpha, *tmv, ST::one());
159  } else {
161  }
162 }
163 
164 
165 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
167  const ArrayView<const Scalar>& alpha,
168  const ArrayView<const Ptr<const MultiVectorBase<Scalar> > >& mv,
169  const Scalar& beta
170  )
171 {
172 #ifdef TEUCHOS_DEBUG
173  TEUCHOS_ASSERT_EQUALITY(alpha.size(), mv.size());
174 #endif
175 
176  // Try to cast mv to an array of this type
178  Teuchos::Array<RCP<const TMV> > tmvs(mv.size());
179  RCP<const TMV> tmv;
180  bool allCastsSuccessful = true;
181  {
182  auto mvIter = mv.begin();
183  auto tmvIter = tmvs.begin();
184  for (; mvIter != mv.end(); ++mvIter, ++tmvIter) {
185  tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromPtr(*mvIter));
186  if (nonnull(tmv)) {
187  *tmvIter = tmv;
188  } else {
189  allCastsSuccessful = false;
190  break;
191  }
192  }
193  }
194 
195  // If casts succeeded, or input arrays are size 0, call Tpetra directly.
196  // Otherwise, fall back to the RTOp implementation.
197  auto len = tmvs.size();
198  if (len == 0) {
199  tpetraMultiVector_.getNonconstObj()->scale(beta);
200  } else if (len == 1 && allCastsSuccessful) {
201  tpetraMultiVector_.getNonconstObj()->update(alpha[0], *tmvs[0], beta);
202  } else if (len == 2 && allCastsSuccessful) {
203  tpetraMultiVector_.getNonconstObj()->update(alpha[0], *tmvs[0], alpha[1], *tmvs[1], beta);
204  } else if (allCastsSuccessful) {
206  auto tmvIter = tmvs.begin();
207  auto alphaIter = alpha.begin();
208 
209  // Check if any entry of tmvs aliases this object's wrapped vector.
210  // If so, replace that entry in the array with a copy.
211  tmv = Teuchos::null;
212  for (; tmvIter != tmvs.end(); ++tmvIter) {
213  if (tmvIter->getRawPtr() == tpetraMultiVector_.getConstObj().getRawPtr()) {
214  if (tmv.is_null()) {
215  tmv = Teuchos::rcp(new TMV(*tpetraMultiVector_.getConstObj(), Teuchos::Copy));
216  }
217  *tmvIter = tmv;
218  }
219  }
220  tmvIter = tmvs.begin();
221 
222  // We add two MVs at a time, so only scale if even num MVs,
223  // and additionally do the first addition if odd num MVs.
224  if ((tmvs.size() % 2) == 0) {
225  tpetraMultiVector_.getNonconstObj()->scale(beta);
226  } else {
227  tpetraMultiVector_.getNonconstObj()->update(*alphaIter, *(*tmvIter), beta);
228  ++tmvIter;
229  ++alphaIter;
230  }
231  for (; tmvIter != tmvs.end(); tmvIter+=2, alphaIter+=2) {
232  tpetraMultiVector_.getNonconstObj()->update(
233  *alphaIter, *(*tmvIter), *(alphaIter+1), *(*(tmvIter+1)), ST::one());
234  }
235  } else {
237  }
238 }
239 
240 
241 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
243  const MultiVectorBase<Scalar>& mv,
244  const ArrayView<Scalar>& prods
245  ) const
246 {
247  auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
248 
249  // If the cast succeeded, call Tpetra directly.
250  // Otherwise, fall back to the RTOp implementation.
251  if (nonnull(tmv)) {
252  tpetraMultiVector_.getConstObj()->dot(*tmv, prods);
253  } else {
255  }
256 }
257 
258 
259 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
261  const ArrayView<typename ScalarTraits<Scalar>::magnitudeType>& norms
262  ) const
263 {
264  tpetraMultiVector_.getConstObj()->norm1(norms);
265 }
266 
267 
268 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
270  const ArrayView<typename ScalarTraits<Scalar>::magnitudeType>& norms
271  ) const
272 {
273  tpetraMultiVector_.getConstObj()->norm2(norms);
274 }
275 
276 
277 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
279  const ArrayView<typename ScalarTraits<Scalar>::magnitudeType>& norms
280  ) const
281 {
282  tpetraMultiVector_.getConstObj()->normInf(norms);
283 }
284 
285 
286 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
289 {
290 #ifdef TEUCHOS_DEBUG
291  TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE(j, 0, this->domain()->dim());
292 #endif
293  return constTpetraVector<Scalar>(
294  tpetraVectorSpace_,
295  tpetraMultiVector_->getVector(j)
296  );
297 }
298 
299 
300 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
303 {
304 #ifdef TEUCHOS_DEBUG
305  TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE(j, 0, this->domain()->dim());
306 #endif
307  return tpetraVector<Scalar>(
308  tpetraVectorSpace_,
309  tpetraMultiVector_.getNonconstObj()->getVectorNonConst(j)
310  );
311 }
312 
313 
314 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
317  const Range1D& col_rng_in
318  ) const
319 {
320 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
321  std::cerr << "\nTpetraMultiVector::subView(Range1D) const called!\n";
322 #endif
323  const Range1D colRng = this->validateColRange(col_rng_in);
324 
326  this->getConstTpetraMultiVector()->subView(colRng);
327 
328  const RCP<const ScalarProdVectorSpaceBase<Scalar> > viewDomainSpace =
329  tpetraVectorSpace<Scalar>(
330  Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
331  tpetraView->getNumVectors(),
332  tpetraView->getMap()->getComm()
333  )
334  );
335 
336  return constTpetraMultiVector(
337  tpetraVectorSpace_,
338  viewDomainSpace,
339  tpetraView
340  );
341 }
342 
343 
344 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
347  const Range1D& col_rng_in
348  )
349 {
350 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
351  std::cerr << "\nTpetraMultiVector::subView(Range1D) called!\n";
352 #endif
353  const Range1D colRng = this->validateColRange(col_rng_in);
354 
356  this->getTpetraMultiVector()->subViewNonConst(colRng);
357 
358  const RCP<const ScalarProdVectorSpaceBase<Scalar> > viewDomainSpace =
359  tpetraVectorSpace<Scalar>(
360  Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
361  tpetraView->getNumVectors(),
362  tpetraView->getMap()->getComm()
363  )
364  );
365 
366  return tpetraMultiVector(
367  tpetraVectorSpace_,
368  viewDomainSpace,
369  tpetraView
370  );
371 }
372 
373 
374 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
377  const ArrayView<const int>& cols_in
378  ) const
379 {
380 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
381  std::cerr << "\nTpetraMultiVector::subView(ArrayView) const called!\n";
382 #endif
383  // Tpetra wants col indices as size_t
384  Array<std::size_t> cols(cols_in.size());
385  for (Array<std::size_t>::size_type i = 0; i < cols.size(); ++i)
386  cols[i] = static_cast<std::size_t>(cols_in[i]);
387 
389  this->getConstTpetraMultiVector()->subView(cols());
390 
391  const RCP<const ScalarProdVectorSpaceBase<Scalar> > viewDomainSpace =
392  tpetraVectorSpace<Scalar>(
393  Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
394  tpetraView->getNumVectors(),
395  tpetraView->getMap()->getComm()
396  )
397  );
398 
399  return constTpetraMultiVector(
400  tpetraVectorSpace_,
401  viewDomainSpace,
402  tpetraView
403  );
404 }
405 
406 
407 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
410  const ArrayView<const int>& cols_in
411  )
412 {
413 #ifdef THYRA_DEFAULT_SPMD_MULTI_VECTOR_VERBOSE_TO_ERROR_OUT
414  std::cerr << "\nTpetraMultiVector::subView(ArrayView) called!\n";
415 #endif
416  // Tpetra wants col indices as size_t
417  Array<std::size_t> cols(cols_in.size());
418  for (Array<std::size_t>::size_type i = 0; i < cols.size(); ++i)
419  cols[i] = static_cast<std::size_t>(cols_in[i]);
420 
422  this->getTpetraMultiVector()->subViewNonConst(cols());
423 
424  const RCP<const ScalarProdVectorSpaceBase<Scalar> > viewDomainSpace =
425  tpetraVectorSpace<Scalar>(
426  Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
427  tpetraView->getNumVectors(),
428  tpetraView->getMap()->getComm()
429  )
430  );
431 
432  return tpetraMultiVector(
433  tpetraVectorSpace_,
434  viewDomainSpace,
435  tpetraView
436  );
437 }
438 
439 
440 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
443  const RTOpPack::RTOpT<Scalar> &primary_op,
444  const ArrayView<const Ptr<const MultiVectorBase<Scalar> > > &multi_vecs,
445  const ArrayView<const Ptr<MultiVectorBase<Scalar> > > &targ_multi_vecs,
446  const ArrayView<const Ptr<RTOpPack::ReductTarget> > &reduct_objs,
447  const Ordinal primary_global_offset
448  ) const
449 {
450 
452  primary_op, multi_vecs, targ_multi_vecs, reduct_objs, primary_global_offset);
453 }
454 
455 
456 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
459  const Range1D &rowRng,
460  const Range1D &colRng,
462  ) const
463 {
465  acquireDetachedMultiVectorViewImpl(rowRng, colRng, sub_mv);
466 }
467 
468 
469 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
472  const Range1D &rowRng,
473  const Range1D &colRng,
475  )
476 {
478  acquireNonconstDetachedMultiVectorViewImpl(rowRng, colRng, sub_mv);
479 }
480 
481 
482 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
486  )
487 {
490 
491 }
492 
493 
494 /* ToDo: Implement these?
495 
496 
497 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
498 RCP<const MultiVectorBase<Scalar> >
499 TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::nonContigSubViewImpl(
500  const ArrayView<const int> &cols
501  ) const
502 {
503  THYRA_DEBUG_ASSERT_MV_COLS("nonContigSubViewImpl(cols)", cols);
504  const int numCols = cols.size();
505  const ArrayRCP<Scalar> localValuesView = createContiguousCopy(cols);
506  return defaultSpmdMultiVector<Scalar>(
507  spmdRangeSpace_,
508  createSmallScalarProdVectorSpaceBase<Scalar>(spmdRangeSpace_, numCols),
509  localValuesView
510  );
511 }
512 
513 
514 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
515 RCP<MultiVectorBase<Scalar> >
516 TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::nonconstNonContigSubViewImpl(
517  const ArrayView<const int> &cols )
518 {
519  THYRA_DEBUG_ASSERT_MV_COLS("nonContigSubViewImpl(cols)", cols);
520  const int numCols = cols.size();
521  const ArrayRCP<Scalar> localValuesView = createContiguousCopy(cols);
522  const Ordinal localSubDim = spmdRangeSpace_->localSubDim();
523  RCP<CopyBackSpmdMultiVectorEntries<Scalar> > copyBackView =
524  copyBackSpmdMultiVectorEntries<Scalar>(cols, localValuesView.getConst(),
525  localSubDim, localValues_.create_weak(), leadingDim_);
526  return Teuchos::rcpWithEmbeddedObjPreDestroy(
527  new TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>(
528  spmdRangeSpace_,
529  createSmallScalarProdVectorSpaceBase<Scalar>(spmdRangeSpace_, numCols),
530  localValuesView),
531  copyBackView
532  );
533 }
534 
535 */
536 
537 
538 // Overridden protected members from SpmdMultiVectorBase
539 
540 
541 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
544 {
545  return tpetraVectorSpace_;
546 }
547 
548 
549 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
551  const Ptr<ArrayRCP<Scalar> > &localValues, const Ptr<Ordinal> &leadingDim
552  )
553 {
554  *localValues = tpetraMultiVector_.getNonconstObj()->get1dViewNonConst();
555  *leadingDim = tpetraMultiVector_->getStride();
556 }
557 
558 
559 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
561  const Ptr<ArrayRCP<const Scalar> > &localValues, const Ptr<Ordinal> &leadingDim
562  ) const
563 {
564  *localValues = tpetraMultiVector_->get1dView();
565  *leadingDim = tpetraMultiVector_->getStride();
566 }
567 
568 
569 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
571  const EOpTransp M_trans,
572  const MultiVectorBase<Scalar> &X,
573  const Ptr<MultiVectorBase<Scalar> > &Y,
574  const Scalar alpha,
575  const Scalar beta
576  ) const
577 {
578  // Try to extract Tpetra objects from X and Y
580  Teuchos::RCP<const TMV> X_tpetra = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(X));
581  Teuchos::RCP<TMV> Y_tpetra = this->getTpetraMultiVector(Teuchos::rcpFromPtr(Y));
582 
583  // If the cast succeeded, call Tpetra directly.
584  // Otherwise, fall back to the default implementation.
585  if (nonnull(X_tpetra) && nonnull(Y_tpetra)) {
587  TEUCHOS_TEST_FOR_EXCEPTION(ST::isComplex && (M_trans == CONJ),
588  std::logic_error,
589  "Error, conjugation without transposition is not allowed for complex scalar types!");
590 
592  switch (M_trans) {
593  case NOTRANS:
594  trans = Teuchos::NO_TRANS;
595  break;
596  case CONJ:
597  trans = Teuchos::NO_TRANS;
598  break;
599  case TRANS:
600  trans = Teuchos::TRANS;
601  break;
602  case CONJTRANS:
603  trans = Teuchos::CONJ_TRANS;
604  break;
605  }
606 
607  Y_tpetra->multiply(trans, Teuchos::NO_TRANS, alpha, *tpetraMultiVector_.getConstObj(), *X_tpetra, beta);
608  Kokkos::fence();
609  } else {
610  SpmdMultiVectorDefaultBase<Scalar>::euclideanApply(M_trans, X, Y, alpha, beta);
611  }
612 
613 }
614 
615 // private
616 
617 
618 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
619 template<class TpetraMultiVector_t>
621  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
622  const RCP<const ScalarProdVectorSpaceBase<Scalar> > &domainSpace,
623  const RCP<TpetraMultiVector_t> &tpetraMultiVector
624  )
625 {
626 #ifdef THYRA_DEBUG
627  TEUCHOS_ASSERT(nonnull(tpetraVectorSpace));
628  TEUCHOS_ASSERT(nonnull(domainSpace));
629  TEUCHOS_ASSERT(nonnull(tpetraMultiVector));
630  // ToDo: Check to make sure that tpetraMultiVector is compatible with
631  // tpetraVectorSpace.
632 #endif
633  tpetraVectorSpace_ = tpetraVectorSpace;
634  domainSpace_ = domainSpace;
635  tpetraMultiVector_.initialize(tpetraMultiVector);
636  this->updateSpmdSpace();
637 }
638 
639 
640 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
641 RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
643 getTpetraMultiVector(const RCP<MultiVectorBase<Scalar> >& mv) const
644 {
645  using Teuchos::rcp_dynamic_cast;
648 
649  RCP<TMV> tmv = rcp_dynamic_cast<TMV>(mv);
650  if (nonnull(tmv)) {
651  return tmv->getTpetraMultiVector();
652  }
653 
654  RCP<TV> tv = rcp_dynamic_cast<TV>(mv);
655  if (nonnull(tv)) {
656  return tv->getTpetraVector();
657  }
658 
659  return Teuchos::null;
660 }
661 
662 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
663 RCP<const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
665 getConstTpetraMultiVector(const RCP<const MultiVectorBase<Scalar> >& mv) const
666 {
667  using Teuchos::rcp_dynamic_cast;
670 
671  RCP<const TMV> tmv = rcp_dynamic_cast<const TMV>(mv);
672  if (nonnull(tmv)) {
673  return tmv->getConstTpetraMultiVector();
674  }
675 
676  RCP<const TV> tv = rcp_dynamic_cast<const TV>(mv);
677  if (nonnull(tv)) {
678  return tv->getConstTpetraVector();
679  }
680 
681  return Teuchos::null;
682 }
683 
684 
685 } // end namespace Thyra
686 
687 
688 #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.