40#ifndef TPETRA_BLOCKMULTIVECTOR_DEF_HPP
41#define TPETRA_BLOCKMULTIVECTOR_DEF_HPP
45#include "Teuchos_OrdinalTraits.hpp"
50template<
class Scalar,
class LO,
class GO,
class Node>
58template<
class Scalar,
class LO,
class GO,
class Node>
59Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
64 const BMV* src_bmv =
dynamic_cast<const BMV*
> (&src);
65 TEUCHOS_TEST_FOR_EXCEPTION(
66 src_bmv ==
nullptr, std::invalid_argument,
"Tpetra::"
67 "BlockMultiVector: The source object of an Import or Export to a "
68 "BlockMultiVector, must also be a BlockMultiVector.");
69 return Teuchos::rcp (src_bmv,
false);
72template<
class Scalar,
class LO,
class GO,
class Node>
75 const Teuchos::DataAccess copyOrView) :
76 dist_object_type (in),
78 pointMap_ (in.pointMap_),
80 blockSize_ (in.blockSize_)
83template<
class Scalar,
class LO,
class GO,
class Node>
88 dist_object_type (Teuchos::rcp (new
map_type (meshMap))),
91 mv_ (pointMap_, numVecs),
92 blockSize_ (blockSize)
95template<
class Scalar,
class LO,
class GO,
class Node>
101 dist_object_type (Teuchos::rcp (new
map_type (meshMap))),
104 mv_ (pointMap_, numVecs),
105 blockSize_ (blockSize)
108template<
class Scalar,
class LO,
class GO,
class Node>
112 const LO blockSize) :
113 dist_object_type (Teuchos::rcp (new
map_type (meshMap))),
115 blockSize_ (blockSize)
131 RCP<const mv_type> X_view_const;
132 const size_t numCols = X_mv.getNumVectors ();
134 Teuchos::Array<size_t> cols (0);
135 X_view_const = X_mv.subView (cols ());
137 X_view_const = X_mv.subView (Teuchos::Range1D (0, numCols-1));
139 TEUCHOS_TEST_FOR_EXCEPTION(
140 X_view_const.is_null (), std::logic_error,
"Tpetra::"
141 "BlockMultiVector constructor: X_mv.subView(...) returned null. This "
142 "should never happen. Please report this bug to the Tpetra developers.");
147 RCP<mv_type> X_view = Teuchos::rcp_const_cast<mv_type> (X_view_const);
148 TEUCHOS_TEST_FOR_EXCEPTION(
149 X_view->getCopyOrView () != Teuchos::View, std::logic_error,
"Tpetra::"
150 "BlockMultiVector constructor: We just set a MultiVector "
151 "to have view semantics, but it claims that it doesn't have view "
152 "semantics. This should never happen. "
153 "Please report this bug to the Tpetra developers.");
158 Teuchos::RCP<const map_type> pointMap = mv_.getMap ();
159 pointMap_ = pointMap;
163template<
class Scalar,
class LO,
class GO,
class Node>
168 const size_t offset) :
169 dist_object_type (Teuchos::rcp (new
map_type (newMeshMap))),
171 pointMap_ (new
map_type(newPointMap)),
176template<
class Scalar,
class LO,
class GO,
class Node>
180 const size_t offset) :
181 dist_object_type (Teuchos::rcp (new
map_type (newMeshMap))),
188template<
class Scalar,
class LO,
class GO,
class Node>
191 dist_object_type (Teuchos::null),
195template<
class Scalar,
class LO,
class GO,
class Node>
201 typedef typename Teuchos::ArrayView<const GO>::size_type size_type;
203 const GST gblNumMeshMapInds =
205 const size_t lclNumMeshMapIndices =
207 const GST gblNumPointMapInds =
208 gblNumMeshMapInds *
static_cast<GST
> (blockSize);
209 const size_t lclNumPointMapInds =
210 lclNumMeshMapIndices *
static_cast<size_t> (blockSize);
214 return map_type (gblNumPointMapInds, lclNumPointMapInds, indexBase,
222 const size_type lclNumMeshGblInds = lclMeshGblInds.size ();
223 Teuchos::Array<GO> lclPointGblInds (lclNumPointMapInds);
224 for (size_type g = 0; g < lclNumMeshGblInds; ++g) {
225 const GO meshGid = lclMeshGblInds[g];
226 const GO pointGidStart = indexBase +
227 (meshGid - indexBase) *
static_cast<GO
> (blockSize);
228 const size_type offset = g *
static_cast<size_type
> (blockSize);
229 for (LO k = 0; k < blockSize; ++k) {
230 const GO pointGid = pointGidStart +
static_cast<GO
> (k);
231 lclPointGblInds[offset +
static_cast<size_type
> (k)] = pointGid;
234 return map_type (gblNumPointMapInds, lclPointGblInds (), indexBase,
241template<
class Scalar,
class LO,
class GO,
class Node>
242Teuchos::RCP<const typename BlockMultiVector<Scalar, LO, GO, Node>::map_type>
247 typedef typename Teuchos::ArrayView<const GO>::size_type size_type;
249 const GST gblNumMeshMapInds =
251 const size_t lclNumMeshMapIndices =
253 const GST gblNumPointMapInds =
254 gblNumMeshMapInds *
static_cast<GST
> (blockSize);
255 const size_t lclNumPointMapInds =
256 lclNumMeshMapIndices *
static_cast<size_t> (blockSize);
260 return Teuchos::rcp(
new map_type (gblNumPointMapInds, lclNumPointMapInds, indexBase,
268 const size_type lclNumMeshGblInds = lclMeshGblInds.size ();
269 Teuchos::Array<GO> lclPointGblInds (lclNumPointMapInds);
270 for (size_type g = 0; g < lclNumMeshGblInds; ++g) {
271 const GO meshGid = lclMeshGblInds[g];
272 const GO pointGidStart = indexBase +
273 (meshGid - indexBase) *
static_cast<GO
> (blockSize);
274 const size_type offset = g *
static_cast<size_type
> (blockSize);
275 for (LO k = 0; k < blockSize; ++k) {
276 const GO pointGid = pointGidStart +
static_cast<GO
> (k);
277 lclPointGblInds[offset +
static_cast<size_type
> (k)] = pointGid;
280 return Teuchos::rcp(
new map_type (gblNumPointMapInds, lclPointGblInds (), indexBase,
286template<
class Scalar,
class LO,
class GO,
class Node>
293 auto X_dst = getLocalBlockHost (localRowIndex, colIndex, Access::ReadWrite);
294 typename const_little_vec_type::HostMirror::const_type X_src (
reinterpret_cast<const impl_scalar_type*
> (vals),
297 using exec_space =
typename device_type::execution_space;
298 Kokkos::deep_copy (exec_space(), X_dst, X_src);
302template<
class Scalar,
class LO,
class GO,
class Node>
309 if (!
meshMap_.isNodeLocalElement (localRowIndex)) {
312 replaceLocalValuesImpl (localRowIndex, colIndex, vals);
317template<
class Scalar,
class LO,
class GO,
class Node>
324 const LO localRowIndex =
meshMap_.getLocalElement (globalRowIndex);
325 if (localRowIndex == Teuchos::OrdinalTraits<LO>::invalid ()) {
328 replaceLocalValuesImpl (localRowIndex, colIndex, vals);
333template<
class Scalar,
class LO,
class GO,
class Node>
340 auto X_dst = getLocalBlockHost (localRowIndex, colIndex, Access::ReadWrite);
341 typename const_little_vec_type::HostMirror::const_type X_src (
reinterpret_cast<const impl_scalar_type*
> (vals),
346template<
class Scalar,
class LO,
class GO,
class Node>
353 if (!
meshMap_.isNodeLocalElement (localRowIndex)) {
356 sumIntoLocalValuesImpl (localRowIndex, colIndex, vals);
361template<
class Scalar,
class LO,
class GO,
class Node>
368 const LO localRowIndex =
meshMap_.getLocalElement (globalRowIndex);
369 if (localRowIndex == Teuchos::OrdinalTraits<LO>::invalid ()) {
372 sumIntoLocalValuesImpl (localRowIndex, colIndex, vals);
378template<
class Scalar,
class LO,
class GO,
class Node>
379typename BlockMultiVector<Scalar, LO, GO, Node>::const_little_host_vec_type
383 const Access::ReadOnlyStruct)
const
385 if (!isValidLocalMeshIndex(localRowIndex)) {
386 return const_little_host_vec_type();
388 const size_t blockSize = getBlockSize();
389 auto hostView = mv_.getLocalViewHost(Access::ReadOnly);
390 LO startRow = localRowIndex*blockSize;
391 LO endRow = startRow + blockSize;
392 return Kokkos::subview(hostView, Kokkos::make_pair(startRow, endRow),
397template<
class Scalar,
class LO,
class GO,
class Node>
398typename BlockMultiVector<Scalar, LO, GO, Node>::little_host_vec_type
402 const Access::OverwriteAllStruct)
405 return little_host_vec_type();
408 auto hostView =
mv_.getLocalViewHost(Access::OverwriteAll);
409 LO startRow = localRowIndex*blockSize;
410 LO endRow = startRow + blockSize;
411 return Kokkos::subview(hostView, Kokkos::make_pair(startRow, endRow),
416template<
class Scalar,
class LO,
class GO,
class Node>
417typename BlockMultiVector<Scalar, LO, GO, Node>::little_host_vec_type
421 const Access::ReadWriteStruct)
423 if (!isValidLocalMeshIndex(localRowIndex)) {
424 return little_host_vec_type();
426 const size_t blockSize = getBlockSize();
427 auto hostView = mv_.getLocalViewHost(Access::ReadWrite);
428 LO startRow = localRowIndex*blockSize;
429 LO endRow = startRow + blockSize;
430 return Kokkos::subview(hostView, Kokkos::make_pair(startRow, endRow),
435template<
class Scalar,
class LO,
class GO,
class Node>
436Teuchos::RCP<const typename BlockMultiVector<Scalar, LO, GO, Node>::mv_type>
441 using Teuchos::rcpFromRef;
448 const this_BMV_type* srcBlkVec =
dynamic_cast<const this_BMV_type*
> (&src);
449 if (srcBlkVec ==
nullptr) {
450 const mv_type* srcMultiVec =
dynamic_cast<const mv_type*
> (&src);
451 if (srcMultiVec ==
nullptr) {
455 return rcp (
new mv_type ());
457 return rcp (srcMultiVec,
false);
460 return rcpFromRef (srcBlkVec->mv_);
464template<
class Scalar,
class LO,
class GO,
class Node>
468 return ! getMultiVectorFromSrcDistObject (src).is_null ();
471template<
class Scalar,
class LO,
class GO,
class Node>
475 const size_t numSameIDs,
482 TEUCHOS_TEST_FOR_EXCEPTION
483 (
true, std::logic_error,
484 "Tpetra::BlockMultiVector::copyAndPermute: Do NOT use this "
485 "instead, create a point importer using makePointMap function.");
488template<
class Scalar,
class LO,
class GO,
class Node>
494 Kokkos::DualView<packet_type*,
496 Kokkos::DualView<
size_t*,
498 size_t& constantNumPackets)
500 TEUCHOS_TEST_FOR_EXCEPTION
501 (
true, std::logic_error,
502 "Tpetra::BlockMultiVector::copyAndPermute: Do NOT use this; "
503 "instead, create a point importer using makePointMap function.");
506template<
class Scalar,
class LO,
class GO,
class Node>
511 Kokkos::DualView<packet_type*,
513 Kokkos::DualView<
size_t*,
515 const size_t constantNumPackets,
518 TEUCHOS_TEST_FOR_EXCEPTION
519 (
true, std::logic_error,
520 "Tpetra::BlockMultiVector::copyAndPermute: Do NOT use this; "
521 "instead, create a point importer using makePointMap function.");
524template<
class Scalar,
class LO,
class GO,
class Node>
528 return meshLocalIndex != Teuchos::OrdinalTraits<LO>::invalid () &&
529 meshMap_.isNodeLocalElement (meshLocalIndex);
532template<
class Scalar,
class LO,
class GO,
class Node>
539template<
class Scalar,
class LO,
class GO,
class Node>
541scale (
const Scalar& val)
546template<
class Scalar,
class LO,
class GO,
class Node>
548update (
const Scalar& alpha,
552 mv_.update (alpha, X.
mv_, beta);
557template <
typename Scalar,
typename ViewY,
typename ViewD,
typename ViewX>
558struct BlockWiseMultiply {
559 typedef typename ViewD::size_type Size;
562 typedef typename ViewD::device_type Device;
563 typedef typename ViewD::non_const_value_type ImplScalar;
564 typedef Kokkos::MemoryTraits<Kokkos::Unmanaged> Unmanaged;
566 template <
typename View>
567 using UnmanagedView = Kokkos::View<
typename View::data_type,
typename View::array_layout,
568 typename View::device_type, Unmanaged>;
569 typedef UnmanagedView<ViewY> UnMViewY;
570 typedef UnmanagedView<ViewD> UnMViewD;
571 typedef UnmanagedView<ViewX> UnMViewX;
573 const Size block_size_;
580 BlockWiseMultiply (
const Size block_size,
const Scalar& alpha,
581 const ViewY& Y,
const ViewD& D,
const ViewX& X)
582 : block_size_(block_size), alpha_(alpha), Y_(Y), D_(D), X_(X)
585 KOKKOS_INLINE_FUNCTION
586 void operator() (
const Size k)
const {
587 const auto zero = Kokkos::ArithTraits<Scalar>::zero();
588 auto D_curBlk = Kokkos::subview(D_, k, Kokkos::ALL (), Kokkos::ALL ());
589 const auto num_vecs = X_.extent(1);
590 for (Size i = 0; i < num_vecs; ++i) {
591 Kokkos::pair<Size, Size> kslice(k*block_size_, (k+1)*block_size_);
592 auto X_curBlk = Kokkos::subview(X_, kslice, i);
593 auto Y_curBlk = Kokkos::subview(Y_, kslice, i);
602template <
typename Scalar,
typename ViewY,
typename ViewD,
typename ViewX>
603inline BlockWiseMultiply<Scalar, ViewY, ViewD, ViewX>
604createBlockWiseMultiply (
const int block_size,
const Scalar& alpha,
605 const ViewY& Y,
const ViewD& D,
const ViewX& X) {
606 return BlockWiseMultiply<Scalar, ViewY, ViewD, ViewX>(block_size, alpha, Y, D, X);
609template <
typename ViewY,
613 typename LO =
typename ViewY::size_type>
614class BlockJacobiUpdate {
616 typedef typename ViewD::device_type Device;
617 typedef typename ViewD::non_const_value_type ImplScalar;
618 typedef Kokkos::MemoryTraits<Kokkos::Unmanaged> Unmanaged;
620 template <
typename ViewType>
621 using UnmanagedView = Kokkos::View<
typename ViewType::data_type,
622 typename ViewType::array_layout,
623 typename ViewType::device_type,
625 typedef UnmanagedView<ViewY> UnMViewY;
626 typedef UnmanagedView<ViewD> UnMViewD;
627 typedef UnmanagedView<ViewZ> UnMViewZ;
637 BlockJacobiUpdate (
const ViewY& Y,
641 const Scalar& beta) :
642 blockSize_ (D.extent (1)),
650 static_assert (
static_cast<int> (ViewY::rank) == 1,
651 "Y must have rank 1.");
652 static_assert (
static_cast<int> (ViewD::rank) == 3,
"D must have rank 3.");
653 static_assert (
static_cast<int> (ViewZ::rank) == 1,
654 "Z must have rank 1.");
660 KOKKOS_INLINE_FUNCTION
void
661 operator() (
const LO& k)
const
664 using Kokkos::subview;
665 typedef Kokkos::pair<LO, LO> range_type;
666 typedef Kokkos::ArithTraits<Scalar> KAT;
670 auto D_curBlk = subview (D_, k, ALL (), ALL ());
671 const range_type kslice (k*blockSize_, (k+1)*blockSize_);
675 auto Z_curBlk = subview (Z_, kslice);
676 auto Y_curBlk = subview (Y_, kslice);
678 if (beta_ == KAT::zero ()) {
681 else if (beta_ != KAT::one ()) {
692 class LO =
typename ViewD::size_type>
694blockJacobiUpdate (
const ViewY& Y,
700 static_assert (Kokkos::is_view<ViewY>::value,
"Y must be a Kokkos::View.");
701 static_assert (Kokkos::is_view<ViewD>::value,
"D must be a Kokkos::View.");
702 static_assert (Kokkos::is_view<ViewZ>::value,
"Z must be a Kokkos::View.");
703 static_assert (
static_cast<int> (ViewY::rank) ==
static_cast<int> (ViewZ::rank),
704 "Y and Z must have the same rank.");
705 static_assert (
static_cast<int> (ViewD::rank) == 3,
"D must have rank 3.");
707 const auto lclNumMeshRows = D.extent (0);
709#ifdef HAVE_TPETRA_DEBUG
713 const auto blkSize = D.extent (1);
714 const auto lclNumPtRows = lclNumMeshRows * blkSize;
715 TEUCHOS_TEST_FOR_EXCEPTION
716 (Y.extent (0) != lclNumPtRows, std::invalid_argument,
717 "blockJacobiUpdate: Y.extent(0) = " << Y.extent (0) <<
" != "
718 "D.extent(0)*D.extent(1) = " << lclNumMeshRows <<
" * " << blkSize
719 <<
" = " << lclNumPtRows <<
".");
720 TEUCHOS_TEST_FOR_EXCEPTION
721 (Y.extent (0) != Z.extent (0), std::invalid_argument,
722 "blockJacobiUpdate: Y.extent(0) = " << Y.extent (0) <<
" != "
723 "Z.extent(0) = " << Z.extent (0) <<
".");
724 TEUCHOS_TEST_FOR_EXCEPTION
725 (Y.extent (1) != Z.extent (1), std::invalid_argument,
726 "blockJacobiUpdate: Y.extent(1) = " << Y.extent (1) <<
" != "
727 "Z.extent(1) = " << Z.extent (1) <<
".");
730 BlockJacobiUpdate<ViewY, Scalar, ViewD, ViewZ, LO> functor (Y, alpha, D, Z, beta);
731 typedef Kokkos::RangePolicy<typename ViewY::execution_space, LO> range_type;
733 range_type range (0,
static_cast<LO
> (lclNumMeshRows));
734 Kokkos::parallel_for (range, functor);
739template<
class Scalar,
class LO,
class GO,
class Node>
747 typedef typename device_type::execution_space exec_space;
748 const LO lclNumMeshRows =
meshMap_.getLocalNumElements ();
750 if (alpha == STS::zero ()) {
757 auto Y_lcl = this->
mv_.getLocalViewDevice (Access::ReadWrite);
758 auto bwm = Impl::createBlockWiseMultiply (blockSize, alphaImpl, Y_lcl, D, X_lcl);
765 Kokkos::RangePolicy<exec_space, LO> range (0, lclNumMeshRows);
766 Kokkos::parallel_for (range, bwm);
770template<
class Scalar,
class LO,
class GO,
class Node>
780 using Kokkos::subview;
783 const IST alphaImpl =
static_cast<IST
> (alpha);
784 const IST betaImpl =
static_cast<IST
> (beta);
785 const LO numVecs =
mv_.getNumVectors ();
787 if (alpha == STS::zero ()) {
791 Z.
update (STS::one (), X, -STS::one ());
792 for (LO j = 0; j < numVecs; ++j) {
793 auto Y_lcl = this->
mv_.getLocalViewDevice (Access::ReadWrite);
795 auto Y_lcl_j = subview (Y_lcl, ALL (), j);
796 auto Z_lcl_j = subview (Z_lcl, ALL (), j);
797 Impl::blockJacobiUpdate (Y_lcl_j, alphaImpl, D, Z_lcl_j, betaImpl);
809#define TPETRA_BLOCKMULTIVECTOR_INSTANT(S,LO,GO,NODE) \
810 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) override
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.
static Teuchos::RCP< const map_type > makePointMapRCP(const map_type &meshMap, const LO blockSize)
Create and return an owning RCP to the point Map corresponding to the given mesh Map and block size.
map_type meshMap_
Mesh Map given to constructor.
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.
typename dist_object_type::buffer_device_type buffer_device_type
Kokkos::Device specialization used for communication buffers.
Tpetra::MultiVector< Scalar, LO, GO, Node > mv_type
The specialization of Tpetra::MultiVector that this class uses.
virtual void copyAndPermute(const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM) override
virtual void packAndPrepare(const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets) override
Tpetra::Map< LO, GO, Node > map_type
The specialization of Tpetra::Map 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.
LO local_ordinal_type
The type of local indices.
virtual void unpackAndCombine(const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode) override
LO getBlockSize() const
Get the number of degrees of freedom per mesh point.
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.
Teuchos::ArrayView< const global_ordinal_type > getLocalElementList() 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 getLocalNumElements() const
The number of elements belonging to the calling process.
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 for new Tpetra features that are not ready for public release, but are ready for evaluation...
Namespace Tpetra contains the class and methods constituting the Tpetra library.
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 AXPY(const CoefficientType &alpha, const ViewType1 &x, const ViewType2 &y)
y := y + alpha * x (dense vector or matrix update)
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 SCAL(const CoefficientType &alpha, const ViewType &x)
x := alpha*x, where x is either rank 1 (a vector) or rank 2 (a matrix).
CombineMode
Rule for combining data in an Import or Export.