Thyra Version of the Day
Loading...
Searching...
No Matches
Thyra_TpetraVector_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_VECTOR_HPP
11#define THYRA_TPETRA_VECTOR_HPP
12
13
14#include "Thyra_TpetraVector_decl.hpp"
15#include "Thyra_TpetraMultiVector.hpp"
16
17namespace Thyra {
18
19
20// Constructors/initializers/accessors
21
22
23template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
26
27
28template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31 const RCP<Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVector
32 )
33{
34 initializeImpl(tpetraVectorSpace, tpetraVector);
35}
36
37
38template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
41 const RCP<const Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> > &tpetraVector
42 )
43{
44 initializeImpl(tpetraVectorSpace, tpetraVector);
45}
46
47
48template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
51{
52 return tpetraVector_.getNonconstObj();
53}
54
55
56template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
62
63
64// Overridden from VectorDefaultBase
65template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
68{
69 if (domainSpace_.is_null()) {
70 domainSpace_ = tpetraVectorSpace<Scalar>(
71 Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(
72 1,
73 tpetraVector_.getConstObj()->getMap()->getComm()
74 )
75 );
76 }
77 return domainSpace_;
78}
79
80
81// Overridden from SpmdMultiVectorBase
82
83
84template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
87{
88 return tpetraVectorSpace_;
89}
90
91
92// Overridden from SpmdVectorBase
93
94
95template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
97 const Ptr<ArrayRCP<Scalar> > &localValues )
98{
99 *localValues = tpetraVector_.getNonconstObj()->get1dViewNonConst();
100}
101
102
103template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
105 const Ptr<ArrayRCP<const Scalar> > &localValues ) const
106{
107 *localValues = tpetraVector_->get1dView();
108}
109
110
111// Overridden from VectorBase
112
113
114template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
116 Scalar l,
117 Scalar u
118 )
119{
120 // Tpetra randomizes with different seed for each proc, so need a global
121 // reduction to get locally-replicated random vector same on each proc.
122 if (!tpetraVector_.getNonconstObj()->isDistributed()) {
123 auto comm = tpetraVector_.getNonconstObj()->getMap()->getComm();
124 if (tpetraVector_.getConstObj()->getMap()->getComm()->getRank() == 0)
125 tpetraVector_.getNonconstObj()->randomize(l, u);
126 else
127 tpetraVector_.getNonconstObj()->putScalar(Teuchos::ScalarTraits<Scalar>::zero());
128 tpetraVector_.getNonconstObj()->reduce();
129 } else {
130 tpetraVector_.getNonconstObj()->randomize(l, u);
131 }
132}
133
134
135template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
137 const VectorBase<Scalar>& x
138 )
139{
140 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
141
142 // If the cast succeeded, call Tpetra directly.
143 // Otherwise, fall back to the RTOp implementation.
144 if (nonnull(tx)) {
145 tpetraVector_.getNonconstObj()->abs(*tx);
146 } else {
148 }
149}
150
151
152template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
154 const VectorBase<Scalar>& x
155 )
156{
157 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
158
159 // If the cast succeeded, call Tpetra directly.
160 // Otherwise, fall back to the RTOp implementation.
161 if (nonnull(tx)) {
162 tpetraVector_.getNonconstObj()->reciprocal(*tx);
163 } else {
165 }
166}
167
168
169template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
171 const VectorBase<Scalar>& x
172 )
173{
174 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
175
176 // If the cast succeeded, call Tpetra directly.
177 // Otherwise, fall back to the RTOp implementation.
178 if (nonnull(tx)) {
180 tpetraVector_.getNonconstObj()->elementWiseMultiply(
181 ST::one(), *tx, *tpetraVector_.getConstObj(), ST::zero());
182 } else {
184 }
185}
186
187
188template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
191 const VectorBase<Scalar>& x
192 ) const
193{
194 auto tx = this->getConstTpetraVector(Teuchos::rcpFromRef(x));
195
196 // If the cast succeeded, call Tpetra directly.
197 // Otherwise, fall back to the RTOp implementation.
198 if (nonnull(tx)) {
199 // Weighted 2-norm function for Tpetra vector seems to be deprecated...
202 = Tpetra::createVector<Scalar>(tx->getMap());
203 temp->elementWiseMultiply(
204 ST::one(), *tx, *tpetraVector_.getConstObj(), ST::zero());
205 return ST::magnitude(ST::squareroot(tpetraVector_.getConstObj()->dot(*temp)));
206 } else {
208 }
209}
210
211
212template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
214 const RTOpPack::RTOpT<Scalar> &op,
215 const ArrayView<const Ptr<const VectorBase<Scalar> > > &vecs,
216 const ArrayView<const Ptr<VectorBase<Scalar> > > &targ_vecs,
217 const Ptr<RTOpPack::ReductTarget> &reduct_obj,
218 const Ordinal global_offset
219 ) const
220{
221 SpmdVectorDefaultBase<Scalar>::applyOpImpl(op, vecs, targ_vecs, reduct_obj, global_offset);
222}
223
224
225template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
234
235
236template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246
247
248template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
256
257
258// Overridden protected functions from MultiVectorBase
259
260
261template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
263{
264 tpetraVector_.getNonconstObj()->putScalar(alpha);
265}
266
267
268template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
271{
272 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
273
274 // If cast succeeded, call Tpetra directly.
275 // Otherwise, fall back to the RTOp implementation.
276 if (nonnull(tmv)) {
277 tpetraVector_.getNonconstObj()->assign(*tmv);
278 } else {
280 }
281}
282
283
284template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
286{
287 tpetraVector_.getNonconstObj()->scale(alpha);
288}
289
290
291template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
293 Scalar alpha,
295 )
296{
297 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
298
299 // If cast succeeded, call Tpetra directly.
300 // Otherwise, fall back to the RTOp implementation.
302 if (nonnull(tmv)) {
303 tpetraVector_.getNonconstObj()->update(alpha, *tmv, ST::one());
304 } else {
306 }
307}
308
309
310template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
312 const ArrayView<const Scalar>& alpha,
313 const ArrayView<const Ptr<const MultiVectorBase<Scalar> > >& mv,
314 const Scalar& beta
315 )
316{
317#ifdef TEUCHOS_DEBUG
318 TEUCHOS_ASSERT_EQUALITY(alpha.size(), mv.size());
319#endif
320
321 // Try to cast mv to Tpetra objects
324 bool allCastsSuccessful = true;
325 {
326 auto mvIter = mv.begin();
327 auto tmvIter = tmvs.begin();
328 for (; mvIter != mv.end(); ++mvIter, ++tmvIter) {
329 tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromPtr(*mvIter));
330 if (nonnull(tmv)) {
331 *tmvIter = tmv;
332 } else {
333 allCastsSuccessful = false;
334 break;
335 }
336 }
337 }
338
339 // If casts succeeded, or input arrays are size 0, call Tpetra directly.
340 // Otherwise, fall back to the RTOp implementation.
341 auto len = mv.size();
342 if (len == 0) {
343 tpetraVector_.getNonconstObj()->scale(beta);
344 } else if (len == 1 && allCastsSuccessful) {
345 tpetraVector_.getNonconstObj()->update(alpha[0], *tmvs[0], beta);
346 } else if (len == 2 && allCastsSuccessful) {
347 tpetraVector_.getNonconstObj()->update(alpha[0], *tmvs[0], alpha[1], *tmvs[1], beta);
348 } else if (allCastsSuccessful) {
350 auto tmvIter = tmvs.begin();
351 auto alphaIter = alpha.begin();
352
353 // Check if any entry of tmvs aliases this object's wrapped vector.
354 // If so, replace that entry in the array with a copy.
355 tmv = Teuchos::null;
356 for (; tmvIter != tmvs.end(); ++tmvIter) {
357 if (tmvIter->getRawPtr() == tpetraVector_.getConstObj().getRawPtr()) {
358 if (tmv.is_null()) {
359 tmv = Teuchos::rcp(new TpetraMultiVector_t(
360 *tpetraVector_.getConstObj(), Teuchos::Copy));
361 }
362 *tmvIter = tmv;
363 }
364 }
365 tmvIter = tmvs.begin();
366
367 // We add two MVs at a time, so only scale if even num MVs,
368 // and additionally do the first addition if odd num MVs.
369 if ((tmvs.size() % 2) == 0) {
370 tpetraVector_.getNonconstObj()->scale(beta);
371 } else {
372 tpetraVector_.getNonconstObj()->update(*alphaIter, *(*tmvIter), beta);
373 ++tmvIter;
374 ++alphaIter;
375 }
376 for (; tmvIter != tmvs.end(); tmvIter+=2, alphaIter+=2) {
377 tpetraVector_.getNonconstObj()->update(
378 *alphaIter, *(*tmvIter), *(alphaIter+1), *(*(tmvIter+1)), ST::one());
379 }
380 } else {
382 }
383}
384
385
386template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
388 const MultiVectorBase<Scalar>& mv,
389 const ArrayView<Scalar>& prods
390 ) const
391{
392 auto tmv = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(mv));
393
394 // If the cast succeeded, call Tpetra directly.
395 // Otherwise, fall back to the RTOp implementation.
396 if (nonnull(tmv)) {
397 tpetraVector_.getConstObj()->dot(*tmv, prods);
398 } else {
400 }
401}
402
403
404template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
407 ) const
408{
409 tpetraVector_.getConstObj()->norm1(norms);
410}
411
412
413template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
416 ) const
417{
418 tpetraVector_.getConstObj()->norm2(norms);
419}
420
421
422template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
425 ) const
426{
427 tpetraVector_.getConstObj()->normInf(norms);
428}
429
430
431template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
433 const EOpTransp M_trans,
435 const Ptr<MultiVectorBase<Scalar> > &Y,
436 const Scalar alpha,
437 const Scalar beta
438 ) const
439{
440 // Try to extract Tpetra objects from X and Y
441 typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
442 Teuchos::RCP<const TMV> X_tpetra = this->getConstTpetraMultiVector(Teuchos::rcpFromRef(X));
443 Teuchos::RCP<TMV> Y_tpetra = this->getTpetraMultiVector(Teuchos::rcpFromPtr(Y));
444
445 // If the cast succeeded, call Tpetra directly.
446 // Otherwise, fall back to the default implementation.
447 if (nonnull(X_tpetra) && nonnull(Y_tpetra)) {
449 TEUCHOS_TEST_FOR_EXCEPTION(ST::isComplex && (M_trans == CONJ),
450 std::logic_error,
451 "Error, conjugation without transposition is not allowed for complex scalar types!");
452
454 switch (M_trans) {
455 case NOTRANS:
456 trans = Teuchos::NO_TRANS;
457 break;
458 case CONJ:
459 trans = Teuchos::NO_TRANS;
460 break;
461 case TRANS:
462 trans = Teuchos::TRANS;
463 break;
464 case CONJTRANS:
465 trans = Teuchos::CONJ_TRANS;
466 break;
467 }
468
469 Y_tpetra->multiply(trans, Teuchos::NO_TRANS, alpha, *tpetraVector_.getConstObj(), *X_tpetra, beta);
470 Kokkos::fence();
471 } else {
472 VectorDefaultBase<Scalar>::applyImpl(M_trans, X, Y, alpha, beta);
473 }
474
475}
476
477
478// private
479
480
481template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
482template<class TpetraVector_t>
483void TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::initializeImpl(
486 )
487{
488#ifdef TEUCHOS_DEBUG
491#endif
492 tpetraVectorSpace_ = tpetraVectorSpace;
493 tpetraVector_.initialize(tpetraVector);
494 this->updateSpmdSpace();
495}
496
497
498template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
499RCP<Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> >
500TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::
501getTpetraMultiVector(const RCP<MultiVectorBase<Scalar> >& mv) const
502{
503 using Teuchos::rcp_dynamic_cast;
504 typedef TpetraMultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TMV;
505 typedef TpetraVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> TV;
506
507 RCP<TMV> tmv = rcp_dynamic_cast<TMV>(mv);
508 if (nonnull(tmv)) {
509 return tmv->getTpetraMultiVector();
510 }
511
512 RCP<TV> tv = rcp_dynamic_cast<TV>(mv);
513 if (nonnull(tv)) {
514 return tv->getTpetraVector();
515 }
516
517 return Teuchos::null;
518}
519
520
521template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
525{
526 using Teuchos::rcp_dynamic_cast;
529
530 RCP<const TMV> tmv = rcp_dynamic_cast<const TMV>(mv);
531 if (nonnull(tmv)) {
532 return tmv->getConstTpetraMultiVector();
533 }
534
535 RCP<const TV> tv = rcp_dynamic_cast<const TV>(mv);
536 if (nonnull(tv)) {
537 return tv->getConstTpetraVector();
538 }
539
540 return Teuchos::null;
541}
542
543
544template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
547getTpetraVector(const RCP<VectorBase<Scalar> >& v) const
548{
550 RCP<TV> tv = Teuchos::rcp_dynamic_cast<TV>(v);
551 if (nonnull(tv)) {
552 return tv->getTpetraVector();
553 } else {
554 return Teuchos::null;
555 }
556}
557
558
559template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
562getConstTpetraVector(const RCP<const VectorBase<Scalar> >& v) const
563{
565 RCP<const TV> tv = Teuchos::rcp_dynamic_cast<const TV>(v);
566 if (nonnull(tv)) {
567 return tv->getConstTpetraVector();
568 } else {
569 return Teuchos::null;
570 }
571}
572
573
574} // end namespace Thyra
575
576
577#endif // THYRA_TPETRA_VECTOR_HPP
iterator begin() const
size_type size() const
iterator begin()
size_type size() const
iterator end()
bool is_null() const
Interface for a collection of column vectors called a multi-vector.
void norms(const MultiVectorBase< Scalar > &V, const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms)
Column-wise multi-vector natural norm.
virtual void dotsImpl(const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
Default implementation of dots using RTOps.
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.
virtual void assignMultiVecImpl(const MultiVectorBase< Scalar > &mv)
Default implementation of assign(MV) using RTOps.
virtual void updateImpl(Scalar alpha, const MultiVectorBase< Scalar > &mv)
Default implementation of update using RTOps.
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
Calls applyOpImplWithComm(null,op,...).
void acquireNonconstDetachedVectorViewImpl(const Range1D &rng, RTOpPack::SubVectorView< Scalar > *sub_vec)
Implemented through this->getLocalData().
void acquireDetachedVectorViewImpl(const Range1D &rng, RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const
Implemented through this->getLocalData().
void commitNonconstDetachedVectorViewImpl(RTOpPack::SubVectorView< Scalar > *sub_vec)
Implemented through this->commitLocalData().
Concrete implementation of Thyra::MultiVector in terms of Tpetra::MultiVector.
Concrete implementation of an SPMD vector space for Tpetra.
RCP< TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > tpetraVectorSpace(const RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &tpetraMap)
Nonmember constructor that creats a serial vector space.
Concrete Thyra::SpmdVectorBase using Tpetra::Vector.
void acquireNonconstDetachedVectorViewImpl(const Range1D &rng, RTOpPack::SubVectorView< Scalar > *sub_vec)
void getLocalVectorDataImpl(const Ptr< ArrayRCP< const Scalar > > &localValues) const
virtual void assignImpl(Scalar alpha)
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm2WeightedImpl(const VectorBase< Scalar > &x) const
virtual void randomizeImpl(Scalar l, Scalar u)
virtual void eleWiseScaleImpl(const VectorBase< Scalar > &x)
RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetraVector()
Get the embedded non-const Tpetra::Vector.
RCP< const VectorSpaceBase< Scalar > > domain() const
virtual void linearCombinationImpl(const ArrayView< const Scalar > &alpha, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &mv, const Scalar &beta)
RCP< const SpmdVectorSpaceBase< Scalar > > spmdSpaceImpl() const
virtual void dotsImpl(const MultiVectorBase< Scalar > &mv, const ArrayView< Scalar > &prods) const
void getNonconstLocalVectorDataImpl(const Ptr< ArrayRCP< Scalar > > &localValues)
virtual void scaleImpl(Scalar alpha)
TpetraVector()
Construct to uninitialized.
virtual void absImpl(const VectorBase< Scalar > &x)
void constInitialize(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector)
Initialize.
virtual void assignMultiVecImpl(const MultiVectorBase< Scalar > &mv)
virtual void updateImpl(Scalar alpha, const MultiVectorBase< Scalar > &mv)
virtual void normsInfImpl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
void commitNonconstDetachedVectorViewImpl(RTOpPack::SubVectorView< Scalar > *sub_vec)
virtual void norms1Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
virtual void reciprocalImpl(const VectorBase< Scalar > &x)
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
RCP< TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > tpetraVector(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector)
Nonmember constructor for TpetraVector.
void acquireDetachedVectorViewImpl(const Range1D &rng, RTOpPack::ConstSubVectorView< Scalar > *sub_vec) const
void applyImpl(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
virtual void norms2Impl(const ArrayView< typename ScalarTraits< Scalar >::magnitudeType > &norms) const
RCP< const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraVector() const
Get the embedded non-const Tpetra::Vector.
void initialize(const RCP< const TpetraVectorSpace< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVectorSpace, const RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector)
Initialize.
Abstract interface for finite-dimensional dense vectors.
void applyImpl(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
. Applies vector or its adjoint (transpose) as a linear operator.
virtual Teuchos::ScalarTraits< Scalar >::magnitudeType norm2WeightedImpl(const VectorBase< Scalar > &x) const
Default implementation of norm_2 (weighted) using RTOps.
virtual void reciprocalImpl(const VectorBase< Scalar > &x)
Default implementation of reciprocal using RTOps.
virtual void absImpl(const VectorBase< Scalar > &x)
Default implementation of abs using RTOps.
virtual void eleWiseScaleImpl(const VectorBase< Scalar > &x)
Default implementation of ele_wise_scale using RTOps.
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
bool nonnull(const std::shared_ptr< T > &p)
NOTRANS
Type for the dimension of a vector space. `**/ typedef Teuchos::Ordinal Ordinal;.
Teuchos::Range1D Range1D
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)