40#ifndef TPETRA_BLOCKMULTIVECTOR_DEF_HPP
41#define TPETRA_BLOCKMULTIVECTOR_DEF_HPP
45#include "Teuchos_OrdinalTraits.hpp"
47#ifdef TPETRA_ENABLE_DEPRECATED_CODE
61 template<
class MultiVectorType>
62 struct RawHostPtrFromMultiVector {
63 typedef typename MultiVectorType::impl_scalar_type impl_scalar_type;
65 static impl_scalar_type* getRawPtr (MultiVectorType& X) {
85 template<
class S,
class LO,
class GO,
class N>
89 return RawHostPtrFromMultiVector<MV>::getRawPtr (X);
97template<
class Scalar,
class LO,
class GO,
class Node>
105template<
class Scalar,
class LO,
class GO,
class Node>
106Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
111 const BMV* src_bmv =
dynamic_cast<const BMV*
> (&src);
112 TEUCHOS_TEST_FOR_EXCEPTION(
113 src_bmv ==
nullptr, std::invalid_argument,
"Tpetra::"
114 "BlockMultiVector: The source object of an Import or Export to a "
115 "BlockMultiVector, must also be a BlockMultiVector.");
116 return Teuchos::rcp (src_bmv,
false);
119template<
class Scalar,
class LO,
class GO,
class Node>
122 const Teuchos::DataAccess copyOrView) :
124 meshMap_ (in.meshMap_),
125 pointMap_ (in.pointMap_),
126 mv_ (in.mv_, copyOrView),
127#ifdef TPETRA_ENABLE_DEPRECATED_CODE
128 mvData_ (getRawHostPtrFromMultiVector (mv_)),
130 blockSize_ (in.blockSize_)
133template<
class Scalar,
class LO,
class GO,
class Node>
140 pointMap_ (makePointMap (meshMap, blockSize)),
141 mv_ (Teuchos::rcpFromRef (pointMap_), numVecs),
142#ifdef TPETRA_ENABLE_DEPRECATED_CODE
143 mvData_ (getRawHostPtrFromMultiVector (mv_)),
145 blockSize_ (blockSize)
148template<
class Scalar,
class LO,
class GO,
class Node>
156 pointMap_ (pointMap),
157 mv_ (Teuchos::rcpFromRef (pointMap_), numVecs),
158#ifdef TPETRA_ENABLE_DEPRECATED_CODE
159 mvData_ (getRawHostPtrFromMultiVector (mv_)),
161 blockSize_ (blockSize)
164template<
class Scalar,
class LO,
class GO,
class Node>
168 const LO blockSize) :
171#ifdef TPETRA_ENABLE_DEPRECATED_CODE
174 blockSize_ (blockSize)
190 RCP<const mv_type> X_view_const;
193 Teuchos::Array<size_t> cols (0);
194 X_view_const = X_mv.
subView (cols ());
196 X_view_const = X_mv.
subView (Teuchos::Range1D (0, numCols-1));
198 TEUCHOS_TEST_FOR_EXCEPTION(
199 X_view_const.is_null (), std::logic_error,
"Tpetra::"
200 "BlockMultiVector constructor: X_mv.subView(...) returned null. This "
201 "should never happen. Please report this bug to the Tpetra developers.");
206 RCP<mv_type> X_view = Teuchos::rcp_const_cast<mv_type> (X_view_const);
207 TEUCHOS_TEST_FOR_EXCEPTION(
208 X_view->getCopyOrView () != Teuchos::View, std::logic_error,
"Tpetra::"
209 "BlockMultiVector constructor: We just set a MultiVector "
210 "to have view semantics, but it claims that it doesn't have view "
211 "semantics. This should never happen. "
212 "Please report this bug to the Tpetra developers.");
217 Teuchos::RCP<const map_type> pointMap =
mv_.
getMap ();
218 if (! pointMap.is_null ()) {
219 pointMap_ = *pointMap;
221#ifdef TPETRA_ENABLE_DEPRECATED_CODE
222 mvData_ = getRawHostPtrFromMultiVector (
mv_);
226template<
class Scalar,
class LO,
class GO,
class Node>
231 const size_t offset) :
233 meshMap_ (newMeshMap),
234 pointMap_ (newPointMap),
235 mv_ (X.mv_, newPointMap, offset * X.getBlockSize ()),
236#ifdef TPETRA_ENABLE_DEPRECATED_CODE
237 mvData_ (getRawHostPtrFromMultiVector (mv_)),
239 blockSize_ (X.getBlockSize ())
242template<
class Scalar,
class LO,
class GO,
class Node>
246 const size_t offset) :
248 meshMap_ (newMeshMap),
249 pointMap_ (makePointMap (newMeshMap, X.getBlockSize ())),
250 mv_ (X.mv_, pointMap_, offset * X.getBlockSize ()),
251#ifdef TPETRA_ENABLE_DEPRECATED_CODE
252 mvData_ (getRawHostPtrFromMultiVector (mv_)),
254 blockSize_ (X.getBlockSize ())
257template<
class Scalar,
class LO,
class GO,
class Node>
261#ifdef TPETRA_ENABLE_DEPRECATED_CODE
267template<
class Scalar,
class LO,
class GO,
class Node>
273 typedef typename Teuchos::ArrayView<const GO>::size_type size_type;
275 const GST gblNumMeshMapInds =
277 const size_t lclNumMeshMapIndices =
279 const GST gblNumPointMapInds =
280 gblNumMeshMapInds *
static_cast<GST
> (blockSize);
281 const size_t lclNumPointMapInds =
282 lclNumMeshMapIndices *
static_cast<size_t> (blockSize);
286 return map_type (gblNumPointMapInds, lclNumPointMapInds, indexBase,
294 const size_type lclNumMeshGblInds = lclMeshGblInds.size ();
295 Teuchos::Array<GO> lclPointGblInds (lclNumPointMapInds);
296 for (size_type g = 0; g < lclNumMeshGblInds; ++g) {
297 const GO meshGid = lclMeshGblInds[g];
298 const GO pointGidStart = indexBase +
299 (meshGid - indexBase) *
static_cast<GO
> (blockSize);
300 const size_type offset = g *
static_cast<size_type
> (blockSize);
301 for (LO k = 0; k < blockSize; ++k) {
302 const GO pointGid = pointGidStart +
static_cast<GO
> (k);
303 lclPointGblInds[offset +
static_cast<size_type
> (k)] = pointGid;
306 return map_type (gblNumPointMapInds, lclPointGblInds (), indexBase,
312template<
class Scalar,
class LO,
class GO,
class Node>
319 auto X_dst = getLocalBlockHost (localRowIndex, colIndex, Access::ReadWrite);
320 typename const_little_vec_type::HostMirror::const_type X_src (
reinterpret_cast<const impl_scalar_type*
> (vals),
326template<
class Scalar,
class LO,
class GO,
class Node>
333 if (! meshMap_.isNodeLocalElement (localRowIndex)) {
336 replaceLocalValuesImpl (localRowIndex, colIndex, vals);
341template<
class Scalar,
class LO,
class GO,
class Node>
348 const LO localRowIndex = meshMap_.getLocalElement (globalRowIndex);
349 if (localRowIndex == Teuchos::OrdinalTraits<LO>::invalid ()) {
352 replaceLocalValuesImpl (localRowIndex, colIndex, vals);
357template<
class Scalar,
class LO,
class GO,
class Node>
364 auto X_dst = getLocalBlockHost (localRowIndex, colIndex, Access::ReadWrite);
365 typename const_little_vec_type::HostMirror::const_type X_src (
reinterpret_cast<const impl_scalar_type*
> (vals),
367 AXPY (
static_cast<impl_scalar_type
> (STS::one ()), X_src, X_dst);
370template<
class Scalar,
class LO,
class GO,
class Node>
377 if (! meshMap_.isNodeLocalElement (localRowIndex)) {
380 sumIntoLocalValuesImpl (localRowIndex, colIndex, vals);
385template<
class Scalar,
class LO,
class GO,
class Node>
392 const LO localRowIndex = meshMap_.getLocalElement (globalRowIndex);
393 if (localRowIndex == Teuchos::OrdinalTraits<LO>::invalid ()) {
396 sumIntoLocalValuesImpl (localRowIndex, colIndex, vals);
401#ifdef TPETRA_ENABLE_DEPRECATED_CODE
403template<
class Scalar,
class LO,
class GO,
class Node>
407getLocalRowView (
const LO localRowIndex,
const LO colIndex, Scalar*& vals)
409 if (! meshMap_.isNodeLocalElement (localRowIndex)) {
412 auto X_ij = getLocalBlockHost (localRowIndex, colIndex, Access::ReadWrite);
413 vals =
reinterpret_cast<Scalar*
> (X_ij.data ());
418template<
class Scalar,
class LO,
class GO,
class Node>
421BlockMultiVector<Scalar, LO, GO, Node>::
422getGlobalRowView (
const GO globalRowIndex,
const LO colIndex, Scalar*& vals)
424 const LO localRowIndex = meshMap_.getLocalElement (globalRowIndex);
425 if (localRowIndex == Teuchos::OrdinalTraits<LO>::invalid ()) {
428 auto X_ij = getLocalBlockHost (localRowIndex, colIndex, Access::ReadWrite);
429 vals =
reinterpret_cast<Scalar*
> (X_ij.data ());
434template<
class Scalar,
class LO,
class GO,
class Node>
435typename BlockMultiVector<Scalar, LO, GO, Node>::little_host_vec_type
437BlockMultiVector<Scalar, LO, GO, Node>::
438getLocalBlock (
const LO localRowIndex,
441 if (! isValidLocalMeshIndex (localRowIndex)) {
442 return little_host_vec_type ();
444 const size_t blockSize = getBlockSize ();
445 const size_t offset = colIndex * this->getStrideY () +
446 localRowIndex * blockSize;
447 impl_scalar_type* blockRaw = this->getRawPtr () + offset;
448 return little_host_vec_type (blockRaw, blockSize);
454template<
class Scalar,
class LO,
class GO,
class Node>
455typename BlockMultiVector<Scalar, LO, GO, Node>::const_little_host_vec_type
456BlockMultiVector<Scalar, LO, GO, Node>::
457getLocalBlockHost (
const LO localRowIndex,
459 const Access::ReadOnlyStruct)
const
461 if (!isValidLocalMeshIndex(localRowIndex)) {
462 return const_little_host_vec_type();
464 const size_t blockSize = getBlockSize();
465 auto hostView = mv_.getLocalViewHost(Access::ReadOnly);
466 LO startRow = localRowIndex*blockSize;
467 LO endRow = startRow + blockSize;
468 return Kokkos::subview(hostView, Kokkos::make_pair(startRow, endRow),
473template<
class Scalar,
class LO,
class GO,
class Node>
474typename BlockMultiVector<Scalar, LO, GO, Node>::little_host_vec_type
478 const Access::OverwriteAllStruct)
480 if (!isValidLocalMeshIndex(localRowIndex)) {
481 return little_host_vec_type();
483 const size_t blockSize = getBlockSize();
484 auto hostView = mv_.getLocalViewHost(Access::OverwriteAll);
485 LO startRow = localRowIndex*blockSize;
486 LO endRow = startRow + blockSize;
487 return Kokkos::subview(hostView, Kokkos::make_pair(startRow, endRow),
492template<
class Scalar,
class LO,
class GO,
class Node>
493typename BlockMultiVector<Scalar, LO, GO, Node>::little_host_vec_type
497 const Access::ReadWriteStruct)
499 if (!isValidLocalMeshIndex(localRowIndex)) {
500 return little_host_vec_type();
502 const size_t blockSize = getBlockSize();
503 auto hostView = mv_.getLocalViewHost(Access::ReadWrite);
504 LO startRow = localRowIndex*blockSize;
505 LO endRow = startRow + blockSize;
506 return Kokkos::subview(hostView, Kokkos::make_pair(startRow, endRow),
511template<
class Scalar,
class LO,
class GO,
class Node>
512Teuchos::RCP<const typename BlockMultiVector<Scalar, LO, GO, Node>::mv_type>
513BlockMultiVector<Scalar, LO, GO, Node>::
517 using Teuchos::rcpFromRef;
523 typedef BlockMultiVector<Scalar, LO, GO, Node> this_type;
524 const this_type* srcBlkVec =
dynamic_cast<const this_type*
> (&src);
525 if (srcBlkVec ==
nullptr) {
526 const mv_type* srcMultiVec =
dynamic_cast<const mv_type*
> (&src);
527 if (srcMultiVec ==
nullptr) {
531 return rcp (
new mv_type ());
533 return rcp (srcMultiVec,
false);
536 return rcpFromRef (srcBlkVec->mv_);
540template<
class Scalar,
class LO,
class GO,
class Node>
544 return ! getMultiVectorFromSrcDistObject (src).is_null ();
547template<
class Scalar,
class LO,
class GO,
class Node>
551 const size_t numSameIDs,
552 const Kokkos::DualView<
const local_ordinal_type*,
553 buffer_device_type>& permuteToLIDs,
554 const Kokkos::DualView<
const local_ordinal_type*,
555 buffer_device_type>& permuteFromLIDs,
558 TEUCHOS_TEST_FOR_EXCEPTION
559 (
true, std::logic_error,
560 "Tpetra::BlockMultiVector::copyAndPermute: Do NOT use this "
561 "instead, create a point importer using makePointMap function.");
564template<
class Scalar,
class LO,
class GO,
class Node>
565void BlockMultiVector<Scalar, LO, GO, Node>::
567(
const SrcDistObject& src,
568 const Kokkos::DualView<
const local_ordinal_type*,
569 buffer_device_type>& exportLIDs,
570 Kokkos::DualView<packet_type*,
571 buffer_device_type>& exports,
572 Kokkos::DualView<
size_t*,
573 buffer_device_type> numPacketsPerLID,
574 size_t& constantNumPackets)
576 TEUCHOS_TEST_FOR_EXCEPTION
577 (
true, std::logic_error,
578 "Tpetra::BlockMultiVector::copyAndPermute: Do NOT use this; "
579 "instead, create a point importer using makePointMap function.");
582template<
class Scalar,
class LO,
class GO,
class Node>
583void BlockMultiVector<Scalar, LO, GO, Node>::
585(
const Kokkos::DualView<
const local_ordinal_type*,
586 buffer_device_type>& importLIDs,
587 Kokkos::DualView<packet_type*,
588 buffer_device_type> imports,
589 Kokkos::DualView<
size_t*,
590 buffer_device_type> numPacketsPerLID,
591 const size_t constantNumPackets,
594 TEUCHOS_TEST_FOR_EXCEPTION
595 (
true, std::logic_error,
596 "Tpetra::BlockMultiVector::copyAndPermute: Do NOT use this; "
597 "instead, create a point importer using makePointMap function.");
600template<
class Scalar,
class LO,
class GO,
class Node>
604 return meshLocalIndex != Teuchos::OrdinalTraits<LO>::invalid () &&
605 meshMap_.isNodeLocalElement (meshLocalIndex);
608template<
class Scalar,
class LO,
class GO,
class Node>
615template<
class Scalar,
class LO,
class GO,
class Node>
617scale (
const Scalar& val)
622template<
class Scalar,
class LO,
class GO,
class Node>
624update (
const Scalar& alpha,
628 mv_.update (alpha, X.
mv_, beta);
633template <
typename Scalar,
typename ViewY,
typename ViewD,
typename ViewX>
634struct BlockWiseMultiply {
635 typedef typename ViewD::size_type Size;
638 typedef typename ViewD::device_type Device;
639 typedef typename ViewD::non_const_value_type ImplScalar;
640 typedef Kokkos::MemoryTraits<Kokkos::Unmanaged> Unmanaged;
642 template <
typename View>
643 using UnmanagedView = Kokkos::View<
typename View::data_type,
typename View::array_layout,
644 typename View::device_type, Unmanaged>;
645 typedef UnmanagedView<ViewY> UnMViewY;
646 typedef UnmanagedView<ViewD> UnMViewD;
647 typedef UnmanagedView<ViewX> UnMViewX;
649 const Size block_size_;
656 BlockWiseMultiply (
const Size block_size,
const Scalar& alpha,
657 const ViewY& Y,
const ViewD& D,
const ViewX& X)
658 : block_size_(block_size), alpha_(alpha), Y_(Y), D_(D), X_(X)
661 KOKKOS_INLINE_FUNCTION
662 void operator() (
const Size k)
const {
663 const auto zero = Kokkos::Details::ArithTraits<Scalar>::zero();
664 auto D_curBlk = Kokkos::subview(D_, k, Kokkos::ALL (), Kokkos::ALL ());
665 const auto num_vecs = X_.extent(1);
666 for (Size i = 0; i < num_vecs; ++i) {
667 Kokkos::pair<Size, Size> kslice(k*block_size_, (k+1)*block_size_);
668 auto X_curBlk = Kokkos::subview(X_, kslice, i);
669 auto Y_curBlk = Kokkos::subview(Y_, kslice, i);
678template <
typename Scalar,
typename ViewY,
typename ViewD,
typename ViewX>
679inline BlockWiseMultiply<Scalar, ViewY, ViewD, ViewX>
680createBlockWiseMultiply (
const int block_size,
const Scalar& alpha,
681 const ViewY& Y,
const ViewD& D,
const ViewX& X) {
682 return BlockWiseMultiply<Scalar, ViewY, ViewD, ViewX>(block_size, alpha, Y, D, X);
685template <
typename ViewY,
689 typename LO =
typename ViewY::size_type>
690class BlockJacobiUpdate {
692 typedef typename ViewD::device_type Device;
693 typedef typename ViewD::non_const_value_type ImplScalar;
694 typedef Kokkos::MemoryTraits<Kokkos::Unmanaged> Unmanaged;
696 template <
typename ViewType>
697 using UnmanagedView = Kokkos::View<
typename ViewType::data_type,
698 typename ViewType::array_layout,
699 typename ViewType::device_type,
701 typedef UnmanagedView<ViewY> UnMViewY;
702 typedef UnmanagedView<ViewD> UnMViewD;
703 typedef UnmanagedView<ViewZ> UnMViewZ;
713 BlockJacobiUpdate (
const ViewY& Y,
717 const Scalar& beta) :
718 blockSize_ (D.extent (1)),
726 static_assert (
static_cast<int> (ViewY::rank) == 1,
727 "Y must have rank 1.");
728 static_assert (
static_cast<int> (ViewD::rank) == 3,
"D must have rank 3.");
729 static_assert (
static_cast<int> (ViewZ::rank) == 1,
730 "Z must have rank 1.");
736 KOKKOS_INLINE_FUNCTION
void
737 operator() (
const LO& k)
const
740 using Kokkos::subview;
741 typedef Kokkos::pair<LO, LO> range_type;
742 typedef Kokkos::Details::ArithTraits<Scalar> KAT;
746 auto D_curBlk = subview (D_, k, ALL (), ALL ());
747 const range_type kslice (k*blockSize_, (k+1)*blockSize_);
751 auto Z_curBlk = subview (Z_, kslice);
752 auto Y_curBlk = subview (Y_, kslice);
754 if (beta_ == KAT::zero ()) {
757 else if (beta_ != KAT::one ()) {
768 class LO =
typename ViewD::size_type>
770blockJacobiUpdate (
const ViewY& Y,
776 static_assert (Kokkos::Impl::is_view<ViewY>::value,
"Y must be a Kokkos::View.");
777 static_assert (Kokkos::Impl::is_view<ViewD>::value,
"D must be a Kokkos::View.");
778 static_assert (Kokkos::Impl::is_view<ViewZ>::value,
"Z must be a Kokkos::View.");
779 static_assert (
static_cast<int> (ViewY::rank) ==
static_cast<int> (ViewZ::rank),
780 "Y and Z must have the same rank.");
781 static_assert (
static_cast<int> (ViewD::rank) == 3,
"D must have rank 3.");
783 const auto lclNumMeshRows = D.extent (0);
785#ifdef HAVE_TPETRA_DEBUG
789 const auto blkSize = D.extent (1);
790 const auto lclNumPtRows = lclNumMeshRows * blkSize;
791 TEUCHOS_TEST_FOR_EXCEPTION
792 (Y.extent (0) != lclNumPtRows, std::invalid_argument,
793 "blockJacobiUpdate: Y.extent(0) = " << Y.extent (0) <<
" != "
794 "D.extent(0)*D.extent(1) = " << lclNumMeshRows <<
" * " << blkSize
795 <<
" = " << lclNumPtRows <<
".");
796 TEUCHOS_TEST_FOR_EXCEPTION
797 (Y.extent (0) != Z.extent (0), std::invalid_argument,
798 "blockJacobiUpdate: Y.extent(0) = " << Y.extent (0) <<
" != "
799 "Z.extent(0) = " << Z.extent (0) <<
".");
800 TEUCHOS_TEST_FOR_EXCEPTION
801 (Y.extent (1) != Z.extent (1), std::invalid_argument,
802 "blockJacobiUpdate: Y.extent(1) = " << Y.extent (1) <<
" != "
803 "Z.extent(1) = " << Z.extent (1) <<
".");
806 BlockJacobiUpdate<ViewY, Scalar, ViewD, ViewZ, LO> functor (Y, alpha, D, Z, beta);
807 typedef Kokkos::RangePolicy<typename ViewY::execution_space, LO> range_type;
809 range_type range (0,
static_cast<LO
> (lclNumMeshRows));
810 Kokkos::parallel_for (range, functor);
815template<
class Scalar,
class LO,
class GO,
class Node>
824 const LO lclNumMeshRows = meshMap_.getNodeNumElements ();
826 if (alpha == STS::zero ()) {
827 this->putScalar (STS::zero ());
830 const LO blockSize = this->getBlockSize ();
833 auto Y_lcl = this->mv_.getLocalViewDevice (Access::ReadWrite);
834 auto bwm = Impl::createBlockWiseMultiply (blockSize, alphaImpl, Y_lcl, D, X_lcl);
841 Kokkos::RangePolicy<execution_space, LO> range (0, lclNumMeshRows);
842 Kokkos::parallel_for (range, bwm);
846template<
class Scalar,
class LO,
class GO,
class Node>
856 using Kokkos::subview;
859 const IST alphaImpl =
static_cast<IST
> (alpha);
860 const IST betaImpl =
static_cast<IST
> (beta);
861 const LO numVecs = mv_.getNumVectors ();
863 if (alpha == STS::zero ()) {
867 Z.
update (STS::one (), X, -STS::one ());
868 for (LO j = 0; j < numVecs; ++j) {
869 auto Y_lcl = this->mv_.getLocalViewDevice (Access::ReadWrite);
871 auto Y_lcl_j = subview (Y_lcl, ALL (), j);
872 auto Z_lcl_j = subview (Z_lcl, ALL (), j);
873 Impl::blockJacobiUpdate (Y_lcl_j, alphaImpl, D, Z_lcl_j, betaImpl);
885#define TPETRA_BLOCKMULTIVECTOR_INSTANT(S,LO,GO,NODE) \
886 template class BlockMultiVector< S, LO, GO, NODE >;
Linear algebra kernels for small dense matrices and vectors.
Declaration of Tpetra::Details::Behavior, a class that describes Tpetra's behavior.
MultiVector for multiple degrees of freedom per mesh point.
virtual bool checkSizes(const Tpetra::SrcDistObject &source)
Compare the source and target (this) objects for compatibility.
void putScalar(const Scalar &val)
Fill all entries with the given value val.
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.
bool sumIntoLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[])
Sum into all values at the given mesh point, using a local index.
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 update(const Scalar &alpha, const BlockMultiVector< Scalar, LO, GO, Node > &X, const Scalar &beta)
Update: this = beta*this + alpha*X.
typename mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the object.
bool sumIntoGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[])
Sum into all values at the given mesh point, using a global index.
BlockMultiVector()
Default constructor.
mv_type getMultiVectorView() const
Get a Tpetra::MultiVector that views this BlockMultiVector's data.
bool replaceLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[])
Replace all values at the given mesh point, using local row and column indices.
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.
Tpetra::MultiVector< Scalar, LO, GO, Node > mv_type
The specialization of Tpetra::MultiVector that this class uses.
typename mv_type::device_type device_type
The Kokkos Device type.
void scale(const Scalar &val)
Multiply all entries in place by the given value val.
bool replaceGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[])
Replace all values at the given mesh point, using a global index.
bool isValidLocalMeshIndex(const LO meshLocalIndex) const
True if and only if meshLocalIndex is a valid local index in the mesh Map.
mv_type mv_
The Tpetra::MultiVector used to represent the data.
typename device_type::execution_space execution_space
The Kokkos execution space.
virtual Teuchos::RCP< const map_type > getMap() const
The Map describing the parallel distribution of this object.
Teuchos::ArrayView< const global_ordinal_type > getNodeElementList() const
Return a NONOWNING view of the global indices owned by this process.
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Accessors for the Teuchos::Comm and Kokkos Node objects.
global_ordinal_type getIndexBase() const
The index base for this Map.
global_size_t getGlobalNumElements() const
The number of elements in this Map.
bool isContiguous() const
True if this Map is distributed contiguously, else false.
size_t getNodeNumElements() const
The number of elements belonging to the calling process.
One or more distributed dense 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.
typename Kokkos::Details::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
size_t getNumVectors() const
Number of columns in the multivector.
Teuchos::DataAccess getCopyOrView() const
Get whether this has copy (copyOrView = Teuchos::Copy) or view (copyOrView = Teuchos::View) semantics...
dual_view_type::t_dev::const_type getLocalViewDevice(Access::ReadOnlyStruct) const
Return a read-only, up-to-date view of this MultiVector's local data on device. This requires that th...
Abstract base class for objects that can be the source of an Import or Export operation.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
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.
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)
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).
KOKKOS_INLINE_FUNCTION void FILL(const ViewType &x, const InputType &val)
Set every entry of x to val.
size_t global_size_t
Global size_t object.
KOKKOS_INLINE_FUNCTION void AXPY(const CoefficientType &alpha, const ViewType1 &x, const ViewType2 &y)
y := y + alpha * x (dense vector or matrix update)
CombineMode
Rule for combining data in an Import or Export.