Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_BlockMultiVector_decl.hpp
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10// clang-format off
11#ifndef TPETRA_BLOCKMULTIVECTOR_DECL_HPP
12#define TPETRA_BLOCKMULTIVECTOR_DECL_HPP
13
16#include "Tpetra_MultiVector.hpp"
17#include <memory>
18#include <utility>
19
20namespace Tpetra {
21
106template<class Scalar,
107 class LO,
108 class GO,
109 class Node>
111 public Tpetra::DistObject<Scalar, LO, GO, Node>
112{
113private:
114 using dist_object_type = Tpetra::DistObject<Scalar, LO, GO, Node>;
115 using STS = Teuchos::ScalarTraits<Scalar>;
116 using packet_type = typename dist_object_type::packet_type;
117
118public:
120
121
126
128 using scalar_type = Scalar;
141
143 using node_type = Node;
144
147
163 typedef Kokkos::View<impl_scalar_type *, device_type> little_vec_type;
164 typedef typename little_vec_type::HostMirror little_host_vec_type;
165
171 typedef Kokkos::View<const impl_scalar_type *, device_type>
173
175 host_device_type;
176 typedef Kokkos::View<const impl_scalar_type*, host_device_type>
177 const_little_host_vec_type;
178
180
182
188
191
194
198
202
205 const Teuchos::DataAccess copyOrView);
206
237 BlockMultiVector (const map_type& meshMap,
238 const LO blockSize,
239 const LO numVecs);
240
245 BlockMultiVector (const map_type& meshMap,
246 const map_type& pointMap,
247 const LO blockSize,
248 const LO numVecs);
249
263 const map_type& meshMap,
264 const LO blockSize);
265
271 const map_type& newMeshMap,
272 const map_type& newPointMap,
273 const size_t offset = 0);
274
280 const map_type& newMeshMap,
281 const size_t offset = 0);
282
284
286
294 static map_type
295 makePointMap (const map_type& meshMap, const LO blockSize);
296
303 static Teuchos::RCP<const map_type>
304 makePointMapRCP (const map_type& meshMap, const LO blockSize);
305
310 const map_type getPointMap () const {
311 return *pointMap_;
312 }
313
315 LO getBlockSize () const {
316 return blockSize_;
317 }
318
320 LO getNumVectors () const {
321 return static_cast<LO> (mv_.getNumVectors ());
322 }
323
329 const mv_type & getMultiVectorView () const;
331
333
335
337 void putScalar (const Scalar& val);
338
340 void scale (const Scalar& val);
341
348 void
349 update (const Scalar& alpha,
351 const Scalar& beta);
352
374 void
375 blockWiseMultiply (const Scalar& alpha,
376 const Kokkos::View<const impl_scalar_type***,
377 device_type, Kokkos::MemoryUnmanaged>& D,
379
410 void
411 blockJacobiUpdate (const Scalar& alpha,
412 const Kokkos::View<const impl_scalar_type***,
413 device_type, Kokkos::MemoryUnmanaged>& D,
416 const Scalar& beta);
417
418
420 template<class TargetMemorySpace>
421 bool need_sync () const {
423 }
424
426 bool need_sync_host() const {
427 return mv_.need_sync_host();
428 }
429
431 bool need_sync_device() const {
432 return mv_.need_sync_device();
433 }
434
436
438
456 bool replaceLocalValues (const LO localRowIndex, const LO colIndex, const Scalar vals[]);
457
468 bool replaceGlobalValues (const GO globalRowIndex, const LO colIndex, const Scalar vals[]);
469
480 bool sumIntoLocalValues (const LO localRowIndex, const LO colIndex, const Scalar vals[]);
481
492 bool sumIntoGlobalValues (const GO globalRowIndex, const LO colIndex, const Scalar vals[]);
493
494
495 const_little_host_vec_type getLocalBlockHost(
496 const LO localRowIndex,
497 const LO colIndex,
498 const Access::ReadOnlyStruct) const;
499
500 little_host_vec_type getLocalBlockHost(
501 const LO localRowIndex,
502 const LO colIndex,
503 const Access::ReadWriteStruct);
504
508 little_host_vec_type getLocalBlockHost(
509 const LO localRowIndex,
510 const LO colIndex,
511 const Access::OverwriteAllStruct);
513
514protected:
520
521
522 // clang-format on
523 virtual bool checkSizes(const Tpetra::SrcDistObject &source) override;
524
525 using dist_object_type::
526 copyAndPermute;
528 // clang-format off
529
530 virtual void
532 (const SrcDistObject& source,
533 const size_t numSameIDs,
534 const Kokkos::DualView<const local_ordinal_type*,
535 buffer_device_type>& permuteToLIDs,
536 const Kokkos::DualView<const local_ordinal_type*,
537 buffer_device_type>& permuteFromLIDs,
538 const CombineMode CM) override;
539
544 // clang-format off
545
546 virtual void
548 (const SrcDistObject& source,
549 const Kokkos::DualView<const local_ordinal_type*,
550 buffer_device_type>& exportLIDs,
551 Kokkos::DualView<packet_type*,
552 buffer_device_type>& exports,
553 Kokkos::DualView<size_t*,
554 buffer_device_type> numPacketsPerLID,
555 size_t& constantNumPackets) override;
556
561 // clang-format off
562
563 virtual void
565 (const Kokkos::DualView<const local_ordinal_type*,
566 buffer_device_type>& importLIDs,
567 Kokkos::DualView<packet_type*,
568 buffer_device_type> imports,
569 Kokkos::DualView<size_t*,
570 buffer_device_type> numPacketsPerLID,
571 const size_t constantNumPackets,
572 const CombineMode combineMode) override;
573 // clang-format off
574
576
577protected:
578
580 size_t getStrideX () const {
581 return static_cast<size_t> (1);
582 }
583
585 size_t getStrideY () const {
586 return mv_.getStride ();
587 }
588
591 bool isValidLocalMeshIndex (const LO meshLocalIndex) const;
592
599
600private:
602 Teuchos::RCP<const map_type> pointMap_;
603
604protected:
607
608private:
609
611 LO blockSize_;
612
614 void
615 replaceLocalValuesImpl (const LO localRowIndex,
616 const LO colIndex,
617 const Scalar vals[]);
619 void
620 sumIntoLocalValuesImpl (const LO localRowIndex,
621 const LO colIndex,
622 const Scalar vals[]);
623
624 static Teuchos::RCP<const mv_type>
625 getMultiVectorFromSrcDistObject (const Tpetra::SrcDistObject&);
626
627 static Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
628 getBlockMultiVectorFromSrcDistObject (const Tpetra::SrcDistObject& src);
629};
630
631} // namespace Tpetra
632
633#endif // TPETRA_BLOCKMULTIVECTOR_DECL_HPP
Forward declaration of Tpetra::BlockCrsMatrix.
Forward declaration of Tpetra::BlockMultiVector.
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.
const mv_type & getMultiVectorView() const
Get a Tpetra::MultiVector that views this BlockMultiVector's data.
BlockMultiVector(BlockMultiVector< Scalar, LO, GO, Node > &&)=default
Move constructor (shallow move).
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 need_sync_host() const
Whether this object needs synchronization to the host.
bool sumIntoLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[])
Sum into all values at the given mesh point, using a local index.
bool need_sync() const
Whether this object needs synchronization to the given memory space.
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.
BlockMultiVector(const map_type &meshMap, const LO blockSize, const LO numVecs)
Constructor that takes a mesh Map, a block size, and a number of vectors (columns).
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 .
BlockMultiVector(const BlockMultiVector< Scalar, LO, GO, Node > &)=default
Copy constructor (shallow copy).
BlockMultiVector< Scalar, LO, GO, Node > & operator=(const BlockMultiVector< Scalar, LO, GO, Node > &)=default
Copy assigment (shallow copy).
little_host_vec_type getLocalBlockHost(const LO localRowIndex, const LO colIndex, const Access::OverwriteAllStruct)
Get a local block on host, with the intent to overwrite all blocks in the BlockMultiVector before acc...
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
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(const BlockMultiVector< Scalar, LO, GO, Node > &X, const map_type &newMeshMap, const size_t offset=0)
View an existing BlockMultiVector using a different mesh Map; compute the new point Map.
BlockMultiVector(const BlockMultiVector< Scalar, LO, GO, Node > &in, const Teuchos::DataAccess copyOrView)
"Copy constructor" with option to deep copy.
const map_type getPointMap() const
Get this BlockMultiVector's (previously computed) point Map.
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
LO getNumVectors() const
Get the number of columns (vectors) in the BlockMultiVector.
Tpetra::MultiVector< Scalar, LO, GO, Node > mv_type
BlockMultiVector(const BlockMultiVector< Scalar, LO, GO, Node > &X, const map_type &newMeshMap, const map_type &newPointMap, const size_t offset=0)
View an existing BlockMultiVector using a different mesh Map, supplying the corresponding point Map.
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
Kokkos::View< impl_scalar_type *, device_type > little_vec_type
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
typename mv_type::device_type device_type
size_t getStrideX() const
Stride between consecutive local entries in the same column.
void scale(const Scalar &val)
Multiply all entries in place by the given value val.
BlockMultiVector(const mv_type &X_mv, const map_type &meshMap, const LO blockSize)
View an existing MultiVector.
Kokkos::View< const impl_scalar_type *, device_type > const_little_vec_type
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 need_sync_device() const
Whether this object needs synchronization to the device.
bool replaceGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[])
Replace all values at the given mesh point, using a global index.
size_t getStrideY() const
Stride between consecutive local entries in the same row.
bool isValidLocalMeshIndex(const LO meshLocalIndex) const
True if and only if meshLocalIndex is a valid local index in the mesh Map.
BlockMultiVector(const map_type &meshMap, const map_type &pointMap, const LO blockSize, const LO numVecs)
Constructor that takes a mesh Map, a point Map, a block size, and a number of vectors (columns).
Base class for distributed Tpetra objects that support data redistribution.
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type
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)
typename ::Kokkos::ArithTraits< Scalar >::val_type packet_type
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)
A parallel distribution of indices over processes.
One or more distributed dense vectors.
typename Kokkos::ArithTraits< Scalar >::val_type impl_scalar_type
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.
CombineMode
Rule for combining data in an Import or Export.