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#ifndef TPETRA_BLOCKMULTIVECTOR_DECL_HPP
11#define TPETRA_BLOCKMULTIVECTOR_DECL_HPP
12
15#include "Tpetra_MultiVector.hpp"
16#include <memory>
17#include <utility>
18
19namespace Tpetra {
20
105template <class Scalar,
106 class LO,
107 class GO,
108 class Node>
109class BlockMultiVector : public Tpetra::DistObject<Scalar, LO, GO, Node> {
110 private:
112 using STS = Teuchos::ScalarTraits<Scalar>;
113 using packet_type = typename dist_object_type::packet_type;
114
115 public:
117
118
123
138
140 using node_type = Node;
141
144
160 typedef Kokkos::View<impl_scalar_type*, device_type> little_vec_type;
161 typedef typename little_vec_type::host_mirror_type little_host_vec_type;
162
168 typedef Kokkos::View<const impl_scalar_type*, device_type>
170
172 host_device_type;
173 typedef Kokkos::View<const impl_scalar_type*, host_device_type>
174 const_little_host_vec_type;
175
177
179
185
188
191
195
199
202 const Teuchos::DataAccess copyOrView);
203
235 const LO blockSize,
236 const LO numVecs);
237
243 const map_type& pointMap,
244 const LO blockSize,
245 const LO numVecs);
246
260 const map_type& meshMap,
261 const LO blockSize);
262
268 const map_type& newMeshMap,
269 const map_type& newPointMap,
270 const size_t offset = 0);
271
277 const map_type& newMeshMap,
278 const size_t offset = 0);
279
281
283
291 static map_type
292 makePointMap(const map_type& meshMap, const LO blockSize);
293
300 static Teuchos::RCP<const map_type>
301 makePointMapRCP(const map_type& meshMap, const LO blockSize);
302
307 const map_type getPointMap() const {
308 return *pointMap_;
309 }
310
312 LO getBlockSize() const {
313 return blockSize_;
314 }
315
317 LO getNumVectors() const {
318 return static_cast<LO>(mv_.getNumVectors());
319 }
320
326 const mv_type& getMultiVectorView() const;
328
330
332
334 void putScalar(const Scalar& val);
335
337 void scale(const Scalar& val);
338
345 void
346 update(const Scalar& alpha,
348 const Scalar& beta);
349
371 void
373 const Kokkos::View<const impl_scalar_type***,
374 device_type, Kokkos::MemoryUnmanaged>& D,
376
407 void
409 const Kokkos::View<const impl_scalar_type***,
410 device_type, Kokkos::MemoryUnmanaged>& D,
413 const Scalar& beta);
414
416 template <class TargetMemorySpace>
417 bool need_sync() const {
419 }
420
422 bool need_sync_host() const {
423 return mv_.need_sync_host();
424 }
425
427 bool need_sync_device() const {
428 return mv_.need_sync_device();
429 }
430
432
434
452 bool replaceLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[]);
453
464 bool replaceGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[]);
465
476 bool sumIntoLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[]);
477
488 bool sumIntoGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[]);
489
490 const_little_host_vec_type getLocalBlockHost(
491 const LO localRowIndex,
492 const LO colIndex,
493 const Access::ReadOnlyStruct) const;
494
495 little_host_vec_type getLocalBlockHost(
496 const LO localRowIndex,
497 const LO colIndex,
498 const Access::ReadWriteStruct);
499
503 little_host_vec_type getLocalBlockHost(
504 const LO localRowIndex,
505 const LO colIndex,
506 const Access::OverwriteAllStruct);
508
509 protected:
515
516
517 virtual bool checkSizes(const Tpetra::SrcDistObject& source) override;
518
519 using dist_object_type::
522
523 virtual void
525 const size_t numSameIDs,
526 const Kokkos::DualView<const local_ordinal_type*,
528 const Kokkos::DualView<const local_ordinal_type*,
530 const CombineMode CM) override;
531
536
537 virtual void
539 const Kokkos::DualView<const local_ordinal_type*,
541 Kokkos::DualView<packet_type*,
542 buffer_device_type>& exports,
543 Kokkos::DualView<size_t*,
546 size_t& constantNumPackets) override;
547
552
553 virtual void
554 unpackAndCombine(const Kokkos::DualView<const local_ordinal_type*,
556 Kokkos::DualView<packet_type*,
558 imports,
559 Kokkos::DualView<size_t*,
562 const size_t constantNumPackets,
563 const CombineMode combineMode) override;
564
566
567 protected:
569 size_t getStrideX() const {
570 return static_cast<size_t>(1);
571 }
572
574 size_t getStrideY() const {
575 return mv_.getStride();
576 }
577
580 bool isValidLocalMeshIndex(const LO meshLocalIndex) const;
581
588
589 private:
591 Teuchos::RCP<const map_type> pointMap_;
592
593 protected:
596
597 private:
599 LO blockSize_;
600
602 void
603 replaceLocalValuesImpl(const LO localRowIndex,
604 const LO colIndex,
605 const Scalar vals[]);
607 void
608 sumIntoLocalValuesImpl(const LO localRowIndex,
609 const LO colIndex,
610 const Scalar vals[]);
611
612 static Teuchos::RCP<const mv_type>
613 getMultiVectorFromSrcDistObject(const Tpetra::SrcDistObject&);
614
615 static Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
616 getBlockMultiVectorFromSrcDistObject(const Tpetra::SrcDistObject& src);
617};
618
619} // namespace Tpetra
620
621#endif // TPETRA_BLOCKMULTIVECTOR_DECL_HPP
Forward declaration of Tpetra::BlockCrsMatrix.
Forward declaration of Tpetra::BlockMultiVector.
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.
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.
BlockMultiVector< Scalar, LO, GO, Node > & operator=(BlockMultiVector< Scalar, LO, GO, Node > &&)=default
Move assigment (shallow move).
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.
GO global_ordinal_type
The type of global indices.
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.
Kokkos::View< impl_scalar_type *, device_type > little_vec_type
"Block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector.
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 .
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).
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.
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
Kokkos::Device specialization used for communication buffers.
LO getNumVectors() const
Get the number of columns (vectors) in the BlockMultiVector.
Node node_type
The Kokkos Node type; a legacy thing that will go away at some point.
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
typename mv_type::device_type device_type
The Kokkos 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.
LO local_ordinal_type
The type of local indices.
Kokkos::View< const impl_scalar_type *, device_type > const_little_vec_type
"Const block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector.
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.
mv_type mv_
The Tpetra::MultiVector used to represent the data.
Struct that holds views of the contents of a CrsMatrix.
Base class for distributed Tpetra objects that support data redistribution.
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)
Pack data and metadata for communication (sends).
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type
Kokkos::Device specialization for communication buffers.
typename ::Kokkos::ArithTraits< Scalar >::val_type packet_type
The type of each datum being sent or received in an Import or Export.
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)
Perform any unpacking and combining after communication.
size_t getStride() const
Stride between columns in the multivector.
bool need_sync_device() const
Whether this MultiVector needs synchronization to the device.
bool need_sync_host() const
Whether this MultiVector needs synchronization to the host.
typename map_type::device_type device_type
This class' preferred Kokkos device type.
size_t getNumVectors() const
Number of columns in the multivector.
typename Kokkos::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
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.