40 #ifndef TPETRA_BLOCKMULTIVECTOR_DEF_HPP 41 #define TPETRA_BLOCKMULTIVECTOR_DEF_HPP 45 #include "Teuchos_OrdinalTraits.hpp" 60 template<
class MultiVectorType>
61 struct RawHostPtrFromMultiVector {
62 typedef typename MultiVectorType::impl_scalar_type impl_scalar_type;
64 static impl_scalar_type* getRawPtr (MultiVectorType& X) {
70 auto X_view_host = X.template getLocalView<typename MultiVectorType::dual_view_type::t_host::device_type> ();
71 impl_scalar_type* X_raw = X_view_host.data ();
88 template<
class S,
class LO,
class GO,
class N>
92 return RawHostPtrFromMultiVector<MV>::getRawPtr (X);
99 template<
class Scalar,
class LO,
class GO,
class Node>
107 template<
class Scalar,
class LO,
class GO,
class Node>
108 Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
113 const BMV* src_bmv =
dynamic_cast<const BMV*
> (&src);
114 TEUCHOS_TEST_FOR_EXCEPTION(
115 src_bmv ==
nullptr, std::invalid_argument,
"Tpetra::" 116 "BlockMultiVector: The source object of an Import or Export to a " 117 "BlockMultiVector, must also be a BlockMultiVector.");
118 return Teuchos::rcp (src_bmv,
false);
121 template<
class Scalar,
class LO,
class GO,
class Node>
124 const Teuchos::DataAccess copyOrView) :
126 meshMap_ (in.meshMap_),
127 pointMap_ (in.pointMap_),
128 mv_ (in.mv_, copyOrView),
129 mvData_ (getRawHostPtrFromMultiVector (mv_)),
130 blockSize_ (in.blockSize_)
133 template<
class Scalar,
class LO,
class GO,
class Node>
140 pointMap_ (makePointMap (meshMap, blockSize)),
141 mv_ (
Teuchos::rcpFromRef (pointMap_), numVecs),
142 mvData_ (getRawHostPtrFromMultiVector (mv_)),
143 blockSize_ (blockSize)
146 template<
class Scalar,
class LO,
class GO,
class Node>
154 pointMap_ (pointMap),
155 mv_ (
Teuchos::rcpFromRef (pointMap_), numVecs),
156 mvData_ (getRawHostPtrFromMultiVector (mv_)),
157 blockSize_ (blockSize)
160 template<
class Scalar,
class LO,
class GO,
class Node>
164 const LO blockSize) :
168 blockSize_ (blockSize)
184 RCP<const mv_type> X_view_const;
187 Teuchos::Array<size_t> cols (0);
188 X_view_const = X_mv.
subView (cols ());
190 X_view_const = X_mv.
subView (Teuchos::Range1D (0, numCols-1));
192 TEUCHOS_TEST_FOR_EXCEPTION(
193 X_view_const.is_null (), std::logic_error,
"Tpetra::" 194 "BlockMultiVector constructor: X_mv.subView(...) returned null. This " 195 "should never happen. Please report this bug to the Tpetra developers.");
200 RCP<mv_type> X_view = Teuchos::rcp_const_cast<
mv_type> (X_view_const);
201 TEUCHOS_TEST_FOR_EXCEPTION(
202 X_view->getCopyOrView () != Teuchos::View, std::logic_error,
"Tpetra::" 203 "BlockMultiVector constructor: We just set a MultiVector " 204 "to have view semantics, but it claims that it doesn't have view " 205 "semantics. This should never happen. " 206 "Please report this bug to the Tpetra developers.");
211 Teuchos::RCP<const map_type> pointMap =
mv_.
getMap ();
212 if (! pointMap.is_null ()) {
213 pointMap_ = *pointMap;
215 mvData_ = getRawHostPtrFromMultiVector (
mv_);
218 template<
class Scalar,
class LO,
class GO,
class Node>
223 const size_t offset) :
225 meshMap_ (newMeshMap),
226 pointMap_ (newPointMap),
227 mv_ (X.mv_, newPointMap, offset * X.getBlockSize ()),
228 mvData_ (getRawHostPtrFromMultiVector (mv_)),
229 blockSize_ (X.getBlockSize ())
232 template<
class Scalar,
class LO,
class GO,
class Node>
236 const size_t offset) :
238 meshMap_ (newMeshMap),
239 pointMap_ (makePointMap (newMeshMap, X.getBlockSize ())),
240 mv_ (X.mv_, pointMap_, offset * X.getBlockSize ()),
241 mvData_ (getRawHostPtrFromMultiVector (mv_)),
242 blockSize_ (X.getBlockSize ())
245 template<
class Scalar,
class LO,
class GO,
class Node>
253 template<
class Scalar,
class LO,
class GO,
class Node>
259 typedef typename Teuchos::ArrayView<const GO>::size_type size_type;
261 const GST gblNumMeshMapInds =
263 const size_t lclNumMeshMapIndices =
265 const GST gblNumPointMapInds =
266 gblNumMeshMapInds *
static_cast<GST
> (blockSize);
267 const size_t lclNumPointMapInds =
268 lclNumMeshMapIndices *
static_cast<size_t> (blockSize);
272 return map_type (gblNumPointMapInds, lclNumPointMapInds, indexBase,
280 const size_type lclNumMeshGblInds = lclMeshGblInds.size ();
281 Teuchos::Array<GO> lclPointGblInds (lclNumPointMapInds);
282 for (size_type g = 0; g < lclNumMeshGblInds; ++g) {
283 const GO meshGid = lclMeshGblInds[g];
284 const GO pointGidStart = indexBase +
285 (meshGid - indexBase) * static_cast<GO> (blockSize);
286 const size_type offset = g *
static_cast<size_type
> (blockSize);
287 for (LO k = 0; k < blockSize; ++k) {
288 const GO pointGid = pointGidStart +
static_cast<GO
> (k);
289 lclPointGblInds[offset +
static_cast<size_type
> (k)] = pointGid;
292 return map_type (gblNumPointMapInds, lclPointGblInds (), indexBase,
298 template<
class Scalar,
class LO,
class GO,
class Node>
303 const Scalar vals[])
const 305 auto X_dst = getLocalBlock (localRowIndex, colIndex);
306 typename const_little_vec_type::HostMirror::const_type X_src (reinterpret_cast<const impl_scalar_type*> (vals),
312 template<
class Scalar,
class LO,
class GO,
class Node>
317 const Scalar vals[])
const 319 if (! meshMap_.isNodeLocalElement (localRowIndex)) {
322 replaceLocalValuesImpl (localRowIndex, colIndex, vals);
327 template<
class Scalar,
class LO,
class GO,
class Node>
332 const Scalar vals[])
const 334 const LO localRowIndex = meshMap_.getLocalElement (globalRowIndex);
335 if (localRowIndex == Teuchos::OrdinalTraits<LO>::invalid ()) {
338 replaceLocalValuesImpl (localRowIndex, colIndex, vals);
343 template<
class Scalar,
class LO,
class GO,
class Node>
348 const Scalar vals[])
const 350 auto X_dst = getLocalBlock (localRowIndex, colIndex);
351 typename const_little_vec_type::HostMirror::const_type X_src (reinterpret_cast<const impl_scalar_type*> (vals),
353 AXPY (static_cast<impl_scalar_type> (STS::one ()), X_src, X_dst);
356 template<
class Scalar,
class LO,
class GO,
class Node>
361 const Scalar vals[])
const 363 if (! meshMap_.isNodeLocalElement (localRowIndex)) {
366 sumIntoLocalValuesImpl (localRowIndex, colIndex, vals);
371 template<
class Scalar,
class LO,
class GO,
class Node>
376 const Scalar vals[])
const 378 const LO localRowIndex = meshMap_.getLocalElement (globalRowIndex);
379 if (localRowIndex == Teuchos::OrdinalTraits<LO>::invalid ()) {
382 sumIntoLocalValuesImpl (localRowIndex, colIndex, vals);
387 template<
class Scalar,
class LO,
class GO,
class Node>
392 if (! meshMap_.isNodeLocalElement (localRowIndex)) {
395 auto X_ij = getLocalBlock (localRowIndex, colIndex);
396 vals =
reinterpret_cast<Scalar*
> (X_ij.data ());
401 template<
class Scalar,
class LO,
class GO,
class Node>
406 const LO localRowIndex = meshMap_.getLocalElement (globalRowIndex);
407 if (localRowIndex == Teuchos::OrdinalTraits<LO>::invalid ()) {
410 auto X_ij = getLocalBlock (localRowIndex, colIndex);
411 vals =
reinterpret_cast<Scalar*
> (X_ij.data ());
416 template<
class Scalar,
class LO,
class GO,
class Node>
420 const LO colIndex)
const 434 if (! isValidLocalMeshIndex (localRowIndex)) {
435 return typename little_vec_type::HostMirror ();
437 const size_t blockSize = getBlockSize ();
438 const size_t offset = colIndex * this->getStrideY () +
439 localRowIndex * blockSize;
441 return typename little_vec_type::HostMirror (blockRaw, blockSize);
445 template<
class Scalar,
class LO,
class GO,
class Node>
446 Teuchos::RCP<const typename BlockMultiVector<Scalar, LO, GO, Node>::mv_type>
451 using Teuchos::rcpFromRef;
458 const this_type* srcBlkVec =
dynamic_cast<const this_type*
> (&src);
459 if (srcBlkVec ==
nullptr) {
460 const mv_type* srcMultiVec =
dynamic_cast<const mv_type*
> (&src);
461 if (srcMultiVec ==
nullptr) {
465 return rcp (
new mv_type ());
467 return rcp (srcMultiVec,
false);
470 return rcpFromRef (srcBlkVec->mv_);
474 template<
class Scalar,
class LO,
class GO,
class Node>
478 return ! getMultiVectorFromSrcDistObject (src).is_null ();
481 template<
class Scalar,
class LO,
class GO,
class Node>
485 const size_t numSameIDs,
486 const Kokkos::DualView<
const local_ordinal_type*,
487 buffer_device_type>& permuteToLIDs,
488 const Kokkos::DualView<
const local_ordinal_type*,
489 buffer_device_type>& permuteFromLIDs)
491 TEUCHOS_TEST_FOR_EXCEPTION
492 (
true, std::logic_error,
493 "Tpetra::BlockMultiVector::copyAndPermute: Do NOT use this " 494 "instead, create a point importer using makePointMap function.");
497 template<
class Scalar,
class LO,
class GO,
class Node>
498 void BlockMultiVector<Scalar, LO, GO, Node>::
500 (
const SrcDistObject& src,
501 const Kokkos::DualView<
const local_ordinal_type*,
502 buffer_device_type>& exportLIDs,
503 Kokkos::DualView<packet_type*,
504 buffer_device_type>& exports,
505 Kokkos::DualView<
size_t*,
506 buffer_device_type> numPacketsPerLID,
507 size_t& constantNumPackets,
510 TEUCHOS_TEST_FOR_EXCEPTION
511 (
true, std::logic_error,
512 "Tpetra::BlockMultiVector::copyAndPermute: Do NOT use this; " 513 "instead, create a point importer using makePointMap function.");
516 template<
class Scalar,
class LO,
class GO,
class Node>
517 void BlockMultiVector<Scalar, LO, GO, Node>::
519 (
const Kokkos::DualView<
const local_ordinal_type*,
520 buffer_device_type>& importLIDs,
521 Kokkos::DualView<packet_type*,
522 buffer_device_type> imports,
523 Kokkos::DualView<
size_t*,
524 buffer_device_type> numPacketsPerLID,
525 const size_t constantNumPackets,
529 TEUCHOS_TEST_FOR_EXCEPTION
530 (
true, std::logic_error,
531 "Tpetra::BlockMultiVector::copyAndPermute: Do NOT use this; " 532 "instead, create a point importer using makePointMap function.");
535 template<
class Scalar,
class LO,
class GO,
class Node>
539 return meshLocalIndex != Teuchos::OrdinalTraits<LO>::invalid () &&
540 meshMap_.isNodeLocalElement (meshLocalIndex);
543 template<
class Scalar,
class LO,
class GO,
class Node>
550 template<
class Scalar,
class LO,
class GO,
class Node>
557 template<
class Scalar,
class LO,
class GO,
class Node>
563 mv_.update (alpha, X.
mv_, beta);
568 template <
typename Scalar,
typename ViewY,
typename ViewD,
typename ViewX>
569 struct BlockWiseMultiply {
570 typedef typename ViewD::size_type Size;
573 typedef typename ViewD::device_type Device;
574 typedef typename ViewD::non_const_value_type ImplScalar;
575 typedef Kokkos::MemoryTraits<Kokkos::Unmanaged> Unmanaged;
577 template <
typename View>
578 using UnmanagedView = Kokkos::View<
typename View::data_type,
typename View::array_layout,
579 typename View::device_type, Unmanaged>;
580 typedef UnmanagedView<ViewY> UnMViewY;
581 typedef UnmanagedView<ViewD> UnMViewD;
582 typedef UnmanagedView<ViewX> UnMViewX;
584 const Size block_size_;
591 BlockWiseMultiply (
const Size block_size,
const Scalar& alpha,
592 const ViewY& Y,
const ViewD& D,
const ViewX& X)
593 : block_size_(block_size), alpha_(alpha), Y_(Y), D_(D), X_(X)
596 KOKKOS_INLINE_FUNCTION
597 void operator() (
const Size k)
const {
598 const auto zero = Kokkos::Details::ArithTraits<Scalar>::zero();
599 auto D_curBlk = Kokkos::subview(D_, k, Kokkos::ALL (), Kokkos::ALL ());
600 const auto num_vecs = X_.extent(1);
601 for (Size i = 0; i < num_vecs; ++i) {
602 Kokkos::pair<Size, Size> kslice(k*block_size_, (k+1)*block_size_);
603 auto X_curBlk = Kokkos::subview(X_, kslice, i);
604 auto Y_curBlk = Kokkos::subview(Y_, kslice, i);
613 template <
typename Scalar,
typename ViewY,
typename ViewD,
typename ViewX>
614 inline BlockWiseMultiply<Scalar, ViewY, ViewD, ViewX>
615 createBlockWiseMultiply (
const int block_size,
const Scalar& alpha,
616 const ViewY& Y,
const ViewD& D,
const ViewX& X) {
617 return BlockWiseMultiply<Scalar, ViewY, ViewD, ViewX>(block_size, alpha, Y, D, X);
620 template <
typename ViewY,
624 typename LO =
typename ViewY::size_type>
625 class BlockJacobiUpdate {
627 typedef typename ViewD::device_type Device;
628 typedef typename ViewD::non_const_value_type ImplScalar;
629 typedef Kokkos::MemoryTraits<Kokkos::Unmanaged> Unmanaged;
631 template <
typename ViewType>
632 using UnmanagedView = Kokkos::View<
typename ViewType::data_type,
633 typename ViewType::array_layout,
634 typename ViewType::device_type,
636 typedef UnmanagedView<ViewY> UnMViewY;
637 typedef UnmanagedView<ViewD> UnMViewD;
638 typedef UnmanagedView<ViewZ> UnMViewZ;
648 BlockJacobiUpdate (
const ViewY& Y,
652 const Scalar& beta) :
653 blockSize_ (D.extent (1)),
661 static_assert (static_cast<int> (ViewY::rank) == 1,
662 "Y must have rank 1.");
663 static_assert (static_cast<int> (ViewD::rank) == 3,
"D must have rank 3.");
664 static_assert (static_cast<int> (ViewZ::rank) == 1,
665 "Z must have rank 1.");
671 KOKKOS_INLINE_FUNCTION
void 672 operator() (
const LO& k)
const 675 using Kokkos::subview;
676 typedef Kokkos::pair<LO, LO> range_type;
677 typedef Kokkos::Details::ArithTraits<Scalar> KAT;
681 auto D_curBlk = subview (D_, k, ALL (), ALL ());
682 const range_type kslice (k*blockSize_, (k+1)*blockSize_);
686 auto Z_curBlk = subview (Z_, kslice);
687 auto Y_curBlk = subview (Y_, kslice);
689 if (beta_ == KAT::zero ()) {
692 else if (beta_ != KAT::one ()) {
699 template<
class ViewY,
703 class LO =
typename ViewD::size_type>
705 blockJacobiUpdate (
const ViewY& Y,
711 static_assert (Kokkos::Impl::is_view<ViewY>::value,
"Y must be a Kokkos::View.");
712 static_assert (Kokkos::Impl::is_view<ViewD>::value,
"D must be a Kokkos::View.");
713 static_assert (Kokkos::Impl::is_view<ViewZ>::value,
"Z must be a Kokkos::View.");
714 static_assert (static_cast<int> (ViewY::rank) == static_cast<int> (ViewZ::rank),
715 "Y and Z must have the same rank.");
716 static_assert (static_cast<int> (ViewD::rank) == 3,
"D must have rank 3.");
718 const auto lclNumMeshRows = D.extent (0);
720 #ifdef HAVE_TPETRA_DEBUG 724 const auto blkSize = D.extent (1);
725 const auto lclNumPtRows = lclNumMeshRows * blkSize;
726 TEUCHOS_TEST_FOR_EXCEPTION
727 (Y.extent (0) != lclNumPtRows, std::invalid_argument,
728 "blockJacobiUpdate: Y.extent(0) = " << Y.extent (0) <<
" != " 729 "D.extent(0)*D.extent(1) = " << lclNumMeshRows <<
" * " << blkSize
730 <<
" = " << lclNumPtRows <<
".");
731 TEUCHOS_TEST_FOR_EXCEPTION
732 (Y.extent (0) != Z.extent (0), std::invalid_argument,
733 "blockJacobiUpdate: Y.extent(0) = " << Y.extent (0) <<
" != " 734 "Z.extent(0) = " << Z.extent (0) <<
".");
735 TEUCHOS_TEST_FOR_EXCEPTION
736 (Y.extent (1) != Z.extent (1), std::invalid_argument,
737 "blockJacobiUpdate: Y.extent(1) = " << Y.extent (1) <<
" != " 738 "Z.extent(1) = " << Z.extent (1) <<
".");
739 #endif // HAVE_TPETRA_DEBUG 741 BlockJacobiUpdate<ViewY, Scalar, ViewD, ViewZ, LO> functor (Y, alpha, D, Z, beta);
742 typedef Kokkos::RangePolicy<typename ViewY::execution_space, LO> range_type;
744 range_type range (0, static_cast<LO> (lclNumMeshRows));
745 Kokkos::parallel_for (range, functor);
750 template<
class Scalar,
class LO,
class GO,
class Node>
759 typedef typename device_type::memory_space memory_space;
760 const LO lclNumMeshRows = meshMap_.getNodeNumElements ();
762 if (alpha == STS::zero ()) {
763 this->putScalar (STS::zero ());
766 const LO blockSize = this->getBlockSize ();
768 auto X_lcl = X.
mv_.template getLocalView<memory_space> ();
769 auto Y_lcl = this->mv_.template getLocalView<memory_space> ();
770 auto bwm = Impl::createBlockWiseMultiply (blockSize, alphaImpl, Y_lcl, D, X_lcl);
777 Kokkos::RangePolicy<execution_space, LO> range (0, lclNumMeshRows);
778 Kokkos::parallel_for (range, bwm);
782 template<
class Scalar,
class LO,
class GO,
class Node>
792 using Kokkos::subview;
793 typedef typename device_type::memory_space memory_space;
796 const IST alphaImpl =
static_cast<IST
> (alpha);
797 const IST betaImpl =
static_cast<IST
> (beta);
798 const LO numVecs = mv_.getNumVectors ();
800 auto X_lcl = X.
mv_.template getLocalView<memory_space> ();
801 auto Y_lcl = this->mv_.template getLocalView<memory_space> ();
802 auto Z_lcl = Z.
mv_.template getLocalView<memory_space> ();
804 if (alpha == STS::zero ()) {
808 Z.
update (STS::one (), X, -STS::one ());
809 for (LO j = 0; j < numVecs; ++j) {
810 auto X_lcl_j = subview (X_lcl, ALL (), j);
811 auto Y_lcl_j = subview (Y_lcl, ALL (), j);
812 auto Z_lcl_j = subview (Z_lcl, ALL (), j);
813 Impl::blockJacobiUpdate (Y_lcl_j, alphaImpl, D, Z_lcl_j, betaImpl);
825 #define TPETRA_BLOCKMULTIVECTOR_INSTANT(S,LO,GO,NODE) \ 826 template class BlockMultiVector< S, LO, GO, NODE >; 828 #endif // TPETRA_BLOCKMULTIVECTOR_DEF_HPP little_vec_type::HostMirror getLocalBlock(const LO localRowIndex, const LO colIndex) const
Get a host view of the degrees of freedom at the given mesh point.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
void update(const Scalar &alpha, const BlockMultiVector< Scalar, LO, GO, Node > &X, const Scalar &beta)
Update: this = beta*this + alpha*X.
bool replaceLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using local row and column indices.
mv_type getMultiVectorView() const
Get a Tpetra::MultiVector that views this BlockMultiVector's data.
void putScalar(const Scalar &val)
Fill all entries with the given value val.
void blockJacobiUpdate(const Scalar &alpha, const Kokkos::View< const impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &D, const BlockMultiVector< Scalar, LO, GO, Node > &X, BlockMultiVector< Scalar, LO, GO, Node > &Z, const Scalar &beta)
Block Jacobi update .
void blockWiseMultiply(const Scalar &alpha, const Kokkos::View< const impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &D, const BlockMultiVector< Scalar, LO, GO, Node > &X)
*this := alpha * D * X, where D is a block diagonal matrix.
One or more distributed dense vectors.
bool sumIntoGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a global index.
virtual bool checkSizes(const Tpetra::SrcDistObject &source)
Compare the source and target (this) objects for compatibility.
bool isContiguous() const
True if this Map is distributed contiguously, else false.
Linear algebra kernels for small dense matrices and vectors.
Teuchos::RCP< const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > subView(const Teuchos::Range1D &colRng) const
Return a const MultiVector with const views of selected columns.
KOKKOS_INLINE_FUNCTION void FILL(const ViewType &x, const InputType &val)
Set every entry of x to val.
global_size_t getGlobalNumElements() const
The number of elements in this Map.
static map_type makePointMap(const map_type &meshMap, const LO blockSize)
Create and return the point Map corresponding to the given mesh Map and block size.
MultiVector for multiple degrees of freedom per mesh point.
bool getGlobalRowView(const GO globalRowIndex, const LO colIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a global index.
size_t getNodeNumElements() const
The number of elements belonging to the calling process.
virtual Teuchos::RCP< const map_type > getMap() const
The Map describing the parallel distribution of this object.
KOKKOS_INLINE_FUNCTION void SCAL(const CoefficientType &alpha, const ViewType &x)
x := alpha*x, where x is either rank 1 (a vector) or rank 2 (a matrix).
size_t global_size_t
Global size_t object.
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
Copy the contents of the MultiVector src into dst.
BlockMultiVector()
Default constructor.
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Accessors for the Teuchos::Comm and Kokkos Node objects.
typename device_type::execution_space execution_space
The Kokkos execution space.
global_ordinal_type getIndexBase() const
The index base for this Map.
typename mv_type::device_type device_type
The Kokkos Device type.
CombineMode
Rule for combining data in an Import or Export.
size_t getNumVectors() const
Number of columns in the multivector.
Tpetra::MultiVector< Scalar, LO, GO, Node > mv_type
The specialization of Tpetra::MultiVector that this class uses.
bool replaceGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using a global index.
mv_type mv_
The Tpetra::MultiVector used to represent the data.
Abstract base class for objects that can be the source of an Import or Export operation.
Teuchos::ArrayView< const global_ordinal_type > getNodeElementList() const
Return a NONOWNING view of the global indices owned by this process.
KOKKOS_INLINE_FUNCTION void GEMV(const CoeffType &alpha, const BlkType &A, const VecType1 &x, const VecType2 &y)
y := y + alpha * A * x (dense matrix-vector multiply)
typename Kokkos::Details::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
bool getLocalRowView(const LO localRowIndex, const LO colIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a local index. ...
Teuchos::DataAccess getCopyOrView() const
Get whether this has copy (copyOrView = Teuchos::Copy) or view (copyOrView = Teuchos::View) semantics...
bool isValidLocalMeshIndex(const LO meshLocalIndex) const
True if and only if meshLocalIndex is a valid local index in the mesh Map.
void scale(const Scalar &val)
Multiply all entries in place by the given value val.
typename mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the object.
bool sumIntoLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a local index.
KOKKOS_INLINE_FUNCTION void AXPY(const CoefficientType &alpha, const ViewType1 &x, const ViewType2 &y)
y := y + alpha * x (dense vector or matrix update)
Declaration of Tpetra::Details::Behavior, a class that describes Tpetra's behavior.