Thyra Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_TpetraVector_def.hpp
Go to the documentation of this file.
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_VECTOR_HPP
43 #define THYRA_TPETRA_VECTOR_HPP
44 
45 
48 
49 namespace Thyra {
50 
51 
52 // Constructors/initializers/accessors
53 
54 
55 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
57 {}
58 
59 
60 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
62  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
63  const RCP<Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVector
64  )
65 {
66  initializeImpl(tpetraVectorSpace, tpetraVector);
67 }
68 
69 
70 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
72  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
73  const RCP<const Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVector
74  )
75 {
76  initializeImpl(tpetraVectorSpace, tpetraVector);
77 }
78 
79 
80 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
83 {
84  return tpetraVector_.getNonconstObj();
85 }
86 
87 
88 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
91 {
92  return tpetraVector_;
93 }
94 
95 
96 // Overridden from VectorDefaultBase
97 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
100 {
101  if (domainSpace_.is_null()) {
102  domainSpace_ = tpetraVectorSpace<Scalar>(
103  Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
104  1,
105  tpetraVector_.getConstObj()->getMap()->getComm()
106  )
107  );
108  }
109  return domainSpace_;
110 }
111 
112 
113 // Overridden from SpmdMultiVectorBase
114 
115 
116 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
119 {
120  return tpetraVectorSpace_;
121 }
122 
123 
124 // Overridden from SpmdVectorBase
125 
126 
127 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
129  const Ptr<ArrayRCP<Scalar> > &localValues )
130 {
131  *localValues = tpetraVector_.getNonconstObj()->get1dViewNonConst();
132 }
133 
134 
135 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
137  const Ptr<ArrayRCP<const Scalar> > &localValues ) const
138 {
139  *localValues = tpetraVector_->get1dView();
140 }
141 
142 
143 // Overridden from VectorBase
144 
145 
146 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
148  Scalar l,
149  Scalar u
150  )
151 {
152  // Tpetra randomizes with different seed for each proc, so need a global
153  // reduction to get locally-replicated random vector same on each proc.
154  if (!tpetraVector_.getNonconstObj()->isDistributed()) {
155  auto comm = tpetraVector_.getNonconstObj()->getMap()->getComm();
156  if (tpetraVector_.getConstObj()->getMap()->getComm()->getRank() == 0)
157  tpetraVector_.getNonconstObj()->randomize(l, u);
158  else
159  tpetraVector_.getNonconstObj()->putScalar(Teuchos::ScalarTraits<Scalar>::zero());
160  tpetraVector_.getNonconstObj()->reduce();
161  } else {
162  tpetraVector_.getNonconstObj()->randomize(l, u);
163  }
164 }
165 
166 
167 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
169  const VectorBase<Scalar>& x
170  )
171 {
172  auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
173 
174  // If the cast succeeded, call Tpetra directly.
175  // Otherwise, fall back to the RTOp implementation.
176  if (nonnull(tx)) {
177  tpetraVector_.getNonconstObj()->abs(*tx);
178  } else {
179  // This version will require/modify the host view of this vector.
180  tpetraVector_.getNonconstObj()->sync_host ();
181  tpetraVector_.getNonconstObj()->modify_host ();
182  VectorDefaultBase<Scalar>::absImpl(x);
183  }
184 }
185 
186 
187 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
189  const VectorBase<Scalar>& x
190  )
191 {
192  auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
193 
194  // If the cast succeeded, call Tpetra directly.
195  // Otherwise, fall back to the RTOp implementation.
196  if (nonnull(tx)) {
197  tpetraVector_.getNonconstObj()->reciprocal(*tx);
198  } else {
199  // This version will require/modify the host view of this vector.
200  tpetraVector_.getNonconstObj()->sync_host ();
201  tpetraVector_.getNonconstObj()->modify_host ();
202  VectorDefaultBase<Scalar>::reciprocalImpl(x);
203  }
204 }
205 
206 
207 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
209  const VectorBase<Scalar>& x
210  )
211 {
212  auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
213 
214  // If the cast succeeded, call Tpetra directly.
215  // Otherwise, fall back to the RTOp implementation.
216  if (nonnull(tx)) {
218  tpetraVector_.getNonconstObj()->elementWiseMultiply(
219  ST::one(), *tx, *tpetraVector_.getConstObj(), ST::zero());
220  } else {
221  // This version will require/modify the host view of this vector.
222  tpetraVector_.getNonconstObj()->sync_host ();
223  tpetraVector_.getNonconstObj()->modify_host ();
224  VectorDefaultBase<Scalar>::eleWiseScaleImpl(x);
225  }
226 }
227 
228 
229 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
232  const VectorBase<Scalar>& x
233  ) const
234 {
235  auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
236 
237  // If the cast succeeded, call Tpetra directly.
238  // Otherwise, fall back to the RTOp implementation.
239  if (nonnull(tx)) {
240  // Weighted 2-norm function for Tpetra vector seems to be deprecated...
243  = Tpetra::createVector<Scalar>(tx->getMap());
244  temp->elementWiseMultiply(
245  ST::one(), *tx, *tpetraVector_.getConstObj(), ST::zero());
246  return ST::magnitude(ST::squareroot(tpetraVector_.getConstObj()->dot(*temp)));
247  } else {
248  // This version will require the host view of this vector.
249  tpetraVector_.getNonconstObj()->sync_host ();
250  return VectorDefaultBase<Scalar>::norm2WeightedImpl(x);
251  }
252 }
253 
254 
255 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
257  const RTOpPack::RTOpT<Scalar> &op,
258  const ArrayView<const Ptr<const VectorBase<Scalar> > > &vecs,
259  const ArrayView<const Ptr<VectorBase<Scalar> > > &targ_vecs,
260  const Ptr<RTOpPack::ReductTarget> &reduct_obj,
261  const Ordinal global_offset
262  ) const
263 {
264  // Sync any non-target Tpetra vecs to host space
265  for (auto itr = vecs.begin(); itr != vecs.end(); ++itr) {
266  auto tv = this->getConstTpetraVector(Teuchos::rcpFromPtr(*itr));
267  if (nonnull(tv)) {
268  typedef Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
269  Teuchos::rcp_const_cast<TV>(tv)->sync_host ();
270  }
271  }
272 
273  // Sync any target Tpetra vecs and mark modified on host
274  for (auto itr = targ_vecs.begin(); itr != targ_vecs.end(); ++itr) {
275  auto tv = this->getTpetraVector(Teuchos::rcpFromPtr(*itr));
276  if (nonnull(tv)) {
277  tv->sync_host ();
278  tv->modify_host ();
279  }
280  }
281 
282  SpmdVectorDefaultBase<Scalar>::applyOpImpl(op, vecs, targ_vecs, reduct_obj, global_offset);
283 }
284 
285 
286 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
289  const Range1D& rng,
291  ) const
292 {
293  // Only viewing data, so just sync dual view to host space
294  typedef typename Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
295  Teuchos::rcp_const_cast<TV>(
296  tpetraVector_.getConstObj())->sync_host ();
297 
298  SpmdVectorDefaultBase<Scalar>::acquireDetachedVectorViewImpl(rng, sub_vec);
299 }
300 
301 
302 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
305  const Range1D& rng,
307  )
308 {
309  // Sync to host and mark as modified
310  tpetraVector_.getNonconstObj()->sync_host ();
311  tpetraVector_.getNonconstObj()->modify_host ();
312 
313  SpmdVectorDefaultBase<Scalar>::acquireNonconstDetachedVectorViewImpl(rng, sub_vec);
314 }
315 
316 
317 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
321  )
322 {
323  SpmdVectorDefaultBase<Scalar>::commitNonconstDetachedVectorViewImpl(sub_vec);
324 
325  // Sync changes from host view to execution space
326  typedef typename Tpetra::Vector<
327  Scalar,LocalOrdinal,GlobalOrdinal,Node>::execution_space execution_space;
328  tpetraVector_.getNonconstObj()->template sync<execution_space>();
329 }
330 
331 
332 // Overridden protected functions from MultiVectorBase
333 
334 
335 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
337 {
338  tpetraVector_.getNonconstObj()->putScalar(alpha);
339 }
340 
341 
342 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
344 assignMultiVecImpl(const MultiVectorBase<Scalar>& mv)
345 {
346  auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
347 
348  // If cast succeeded, call Tpetra directly.
349  // Otherwise, fall back to the RTOp implementation.
350  if (nonnull(tmv)) {
351  tpetraVector_.getNonconstObj()->assign(*tmv);
352  } else {
353  // This version will require/modify the host view of this vector.
354  tpetraVector_.getNonconstObj()->sync_host ();
355  tpetraVector_.getNonconstObj()->modify_host ();
356  MultiVectorDefaultBase<Scalar>::assignMultiVecImpl(mv);
357  }
358 }
359 
360 
361 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
363 {
364  tpetraVector_.getNonconstObj()->scale(alpha);
365 }
366 
367 
368 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
370  Scalar alpha,
371  const MultiVectorBase<Scalar>& mv
372  )
373 {
374  auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
375 
376  // If cast succeeded, call Tpetra directly.
377  // Otherwise, fall back to the RTOp implementation.
379  if (nonnull(tmv)) {
380  tpetraVector_.getNonconstObj()->update(alpha, *tmv, ST::one());
381  } else {
382  // This version will require/modify the host view of this vector.
383  tpetraVector_.getNonconstObj()->sync_host();
384  tpetraVector_.getNonconstObj()->modify_host();
385  MultiVectorDefaultBase<Scalar>::updateImpl(alpha, mv);
386  }
387 }
388 
389 
390 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
392  const ArrayView<const Scalar>& alpha,
393  const ArrayView<const Ptr<const MultiVectorBase<Scalar> > >& mv,
394  const Scalar& beta
395  )
396 {
397 #ifdef TEUCHOS_DEBUG
398  TEUCHOS_ASSERT_EQUALITY(alpha.size(), mv.size());
399 #endif
400 
401  // Try to cast mv to Tpetra objects
404  bool allCastsSuccessful = true;
405  {
406  auto mvIter = mv.begin();
407  auto tmvIter = tmvs.begin();
408  for (; mvIter != mv.end(); ++mvIter, ++tmvIter) {
409  tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromPtr(*mvIter));
410  if (nonnull(tmv)) {
411  *tmvIter = tmv;
412  } else {
413  allCastsSuccessful = false;
414  break;
415  }
416  }
417  }
418 
419  // If casts succeeded, or input arrays are size 0, call Tpetra directly.
420  // Otherwise, fall back to the RTOp implementation.
421  auto len = mv.size();
422  if (len == 0) {
423  tpetraVector_.getNonconstObj()->scale(beta);
424  } else if (len == 1 && allCastsSuccessful) {
425  tpetraVector_.getNonconstObj()->update(alpha[0], *tmvs[0], beta);
426  } else if (len == 2 && allCastsSuccessful) {
427  tpetraVector_.getNonconstObj()->update(alpha[0], *tmvs[0], alpha[1], *tmvs[1], beta);
428  } else if (allCastsSuccessful) {
430  auto tmvIter = tmvs.begin();
431  auto alphaIter = alpha.begin();
432 
433  // Check if any entry of tmvs aliases this object's wrapped vector.
434  // If so, replace that entry in the array with a copy.
435  tmv = Teuchos::null;
436  for (; tmvIter != tmvs.end(); ++tmvIter) {
437  if (tmvIter->getRawPtr() == tpetraVector_.getConstObj().getRawPtr()) {
438  if (tmv.is_null()) {
440  *tpetraVector_.getConstObj(), Teuchos::Copy));
441  }
442  *tmvIter = tmv;
443  }
444  }
445  tmvIter = tmvs.begin();
446 
447  // We add two MVs at a time, so only scale if even num MVs,
448  // and additionally do the first addition if odd num MVs.
449  if ((tmvs.size() % 2) == 0) {
450  tpetraVector_.getNonconstObj()->scale(beta);
451  } else {
452  tpetraVector_.getNonconstObj()->update(*alphaIter, *(*tmvIter), beta);
453  ++tmvIter;
454  ++alphaIter;
455  }
456  for (; tmvIter != tmvs.end(); tmvIter+=2, alphaIter+=2) {
457  tpetraVector_.getNonconstObj()->update(
458  *alphaIter, *(*tmvIter), *(alphaIter+1), *(*(tmvIter+1)), ST::one());
459  }
460  } else {
461  // This version will require/modify the host view of this vector.
462  tpetraVector_.getNonconstObj()->sync_host ();
463  tpetraVector_.getNonconstObj()->modify_host ();
464  MultiVectorDefaultBase<Scalar>::linearCombinationImpl(alpha, mv, beta);
465  }
466 }
467 
468 
469 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
471  const MultiVectorBase<Scalar>& mv,
472  const ArrayView<Scalar>& prods
473  ) const
474 {
475  auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
476 
477  // If the cast succeeded, call Tpetra directly.
478  // Otherwise, fall back to the RTOp implementation.
479  if (nonnull(tmv)) {
480  tpetraVector_.getConstObj()->dot(*tmv, prods);
481  } else {
482  // This version will require/modify the host view of this vector.
483  tpetraVector_.getNonconstObj()->sync_host ();
484  tpetraVector_.getNonconstObj()->modify_host ();
485  MultiVectorDefaultBase<Scalar>::dotsImpl(mv, prods);
486  }
487 }
488 
489 
490 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
492  const ArrayView<typename ScalarTraits<Scalar>::magnitudeType>& norms
493  ) const
494 {
495  tpetraVector_.getConstObj()->norm1(norms);
496 }
497 
498 
499 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
501  const ArrayView<typename ScalarTraits<Scalar>::magnitudeType>& norms
502  ) const
503 {
504  tpetraVector_.getConstObj()->norm2(norms);
505 }
506 
507 
508 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
510  const ArrayView<typename ScalarTraits<Scalar>::magnitudeType>& norms
511  ) const
512 {
513  tpetraVector_.getConstObj()->normInf(norms);
514 }
515 
516 
517 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
519  const EOpTransp M_trans,
520  const MultiVectorBase<Scalar> &X,
521  const Ptr<MultiVectorBase<Scalar> > &Y,
522  const Scalar alpha,
523  const Scalar beta
524  ) const
525 {
526  // Try to extract Tpetra objects from X and Y
527  typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
528  typedef Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
529  Teuchos::RCP<const TMV> X_tpetra = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(X));
530  Teuchos::RCP<TMV> Y_tpetra = this->getTpetraMultiVector(Teuchos::rcpFromPtr(Y));
531 
532  // If the cast succeeded, call Tpetra directly.
533  // Otherwise, fall back to the default implementation.
534  if (nonnull(X_tpetra) && nonnull(Y_tpetra)) {
535  // Sync everything to the execution space
536  typedef typename TMV::execution_space execution_space;
537  Teuchos::rcp_const_cast<TMV>(X_tpetra)->template sync<execution_space>();
538  Y_tpetra->template sync<execution_space>();
539  Teuchos::rcp_const_cast<TV>(tpetraVector_.getConstObj())->template sync<execution_space>();
540 
542  TEUCHOS_TEST_FOR_EXCEPTION(ST::isComplex && (M_trans == CONJ),
543  std::logic_error,
544  "Error, conjugation without transposition is not allowed for complex scalar types!");
545 
547  switch (M_trans) {
548  case NOTRANS:
549  trans = Teuchos::NO_TRANS;
550  break;
551  case CONJ:
552  trans = Teuchos::NO_TRANS;
553  break;
554  case TRANS:
555  trans = Teuchos::TRANS;
556  break;
557  case CONJTRANS:
558  trans = Teuchos::CONJ_TRANS;
559  break;
560  }
561 
562  Y_tpetra->template modify<execution_space>();
563  Y_tpetra->multiply(trans, Teuchos::NO_TRANS, alpha, *tpetraVector_.getConstObj(), *X_tpetra, beta);
564  } else {
565  Teuchos::rcp_const_cast<TV>(tpetraVector_.getConstObj())->sync_host ();
566  VectorDefaultBase<Scalar>::applyImpl(M_trans, X, Y, alpha, beta);
567  }
568 
569 }
570 
571 
572 // private
573 
574 
575 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
576 template<class TpetraVector_t>
578  const RCP<const TpetraVectorSpace<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVectorSpace,
579  const RCP<TpetraVector_t> &tpetraVector
580  )
581 {
582 #ifdef TEUCHOS_DEBUG
583  TEUCHOS_ASSERT(nonnull(tpetraVectorSpace));
584  TEUCHOS_ASSERT(nonnull(tpetraVector));
585 #endif
586  tpetraVectorSpace_ = tpetraVectorSpace;
587  tpetraVector_.initialize(tpetraVector);
588  this->updateSpmdSpace();
589 }
590 
591 
592 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
595 getTpetraMultiVector(const RCP<MultiVectorBase<Scalar> >& mv) const
596 {
597  using Teuchos::rcp_dynamic_cast;
600 
601  RCP<TMV> tmv = rcp_dynamic_cast<TMV>(mv);
602  if (nonnull(tmv)) {
603  return tmv->getTpetraMultiVector();
604  }
605 
606  RCP<TV> tv = rcp_dynamic_cast<TV>(mv);
607  if (nonnull(tv)) {
608  return tv->getTpetraVector();
609  }
610 
611  return Teuchos::null;
612 }
613 
614 
615 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
618 getConstTpetraMultiVector(const RCP<const MultiVectorBase<Scalar> >& mv) const
619 {
620  using Teuchos::rcp_dynamic_cast;
623 
624  RCP<const TMV> tmv = rcp_dynamic_cast<const TMV>(mv);
625  if (nonnull(tmv)) {
626  return tmv->getConstTpetraMultiVector();
627  }
628 
629  RCP<const TV> tv = rcp_dynamic_cast<const TV>(mv);
630  if (nonnull(tv)) {
631  return tv->getConstTpetraVector();
632  }
633 
634  return Teuchos::null;
635 }
636 
637 
638 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
641 getTpetraVector(const RCP<VectorBase<Scalar> >& v) const
642 {
644  RCP<TV> tv = Teuchos::rcp_dynamic_cast<TV>(v);
645  if (nonnull(tv)) {
646  return tv->getTpetraVector();
647  } else {
648  return Teuchos::null;
649  }
650 }
651 
652 
653 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
656 getConstTpetraVector(const RCP<const VectorBase<Scalar> >& v) const
657 {
659  RCP<const TV> tv = Teuchos::rcp_dynamic_cast<const TV>(v);
660  if (nonnull(tv)) {
661  return tv->getConstTpetraVector();
662  } else {
663  return Teuchos::null;
664  }
665 }
666 
667 
668 } // end namespace Thyra
669 
670 
671 #endif // THYRA_TPETRA_VECTOR_HPP
void constInitialize(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector)
Initialize.
Concrete implementation of Thyra::MultiVector in terms of Tpetra::MultiVector.
RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetraMultiVector(const RCP< MultiVectorBase< Scalar > > &mv) const
void applyImpl(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
Concrete implementation of an SPMD vector space for Tpetra.
void commitNonconstDetachedVectorViewImpl(RTOpPack::SubVectorView< Scalar > *sub_vec)
void initialize(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector)
Initialize.
RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraMultiVector(const RCP< const MultiVectorBase< Scalar > > &mv) const
iterator begin() const
TpetraVector()
Construct to uninitialized.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void norms1Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
virtual void applyOpImpl(const RTOpPack::RTOpT< Scalar > &op, const ArrayView< const Ptr< const VectorBase< Scalar > > > &vecs, const ArrayView< const Ptr< VectorBase< Scalar > > > &targ_vecs, const Ptr< RTOpPack::ReductTarget > &reduct_obj, const Ordinal global_offset) const
size_type size() const
virtual void norms2Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
void getLocalVectorDataImpl(const Ptr< ArrayRCP< const Scalar > > &localValues) const
virtual void assignMultiVecImpl(const MultiVectorBase< Scalar > &mv)
void acquireNonconstDetachedVectorViewImpl(const Range1D &rng, RTOpPack::SubVectorView< Scalar > *sub_vec)
RCP< const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraVector() const
Get the embedded non-const Tpetra::Vector.
virtual void reciprocalImpl(const VectorBase< Scalar > &x)
void acquireDetachedVectorViewImpl(const Range1D &rng, RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm2WeightedImpl(const VectorBase< Scalar > &x) const
virtual void linearCombinationImpl(const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta)
virtual void absImpl(const VectorBase< Scalar > &x)
void getNonconstLocalVectorDataImpl(const Ptr< ArrayRCP< Scalar > > &localValues)
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraMultiVector_t
RCP< const SpmdVectorSpaceBase< Scalar > > spmdSpaceImpl() const
void initializeImpl(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< TpetraVector_t > &tpetraVector)
Concrete Thyra::SpmdVectorBase using Tpetra::Vector.
virtual void dotsImpl(const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetraVector()
Get the embedded non-const Tpetra::Vector.
bool nonnull(const boost::shared_ptr< T > &p)
virtual void assignImpl(Scalar alpha)
virtual void eleWiseScaleImpl(const VectorBase< Scalar > &x)
virtual void normsInfImpl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
RCP< const VectorSpaceBase< Scalar > > domain() const
virtual void randomizeImpl(Scalar l, Scalar u)
virtual void updateImpl(Scalar alpha, const MultiVectorBase< Scalar > &mv)
virtual void scaleImpl(Scalar alpha)
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
Teuchos_Ordinal Ordinal