Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_BlockMultiVector_decl.hpp
1// @HEADER
2// ***********************************************************************
3//
4// Tpetra: Templated Linear Algebra Services Package
5// Copyright (2008) Sandia Corporation
6//
7// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8// the U.S. Government retains certain rights in this software.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// ************************************************************************
38// @HEADER
39
40#ifndef TPETRA_BLOCKMULTIVECTOR_DECL_HPP
41#define TPETRA_BLOCKMULTIVECTOR_DECL_HPP
42
45#include "Tpetra_MultiVector.hpp"
46#include <memory>
47#include <utility>
48
49namespace Tpetra {
50
135template<class Scalar,
136 class LO,
137 class GO,
138 class Node>
140 public Tpetra::DistObject<Scalar, LO, GO, Node>
141{
142private:
143 using dist_object_type = Tpetra::DistObject<Scalar, LO, GO, Node>;
144 using STS = Teuchos::ScalarTraits<Scalar>;
145 using packet_type = typename dist_object_type::packet_type;
146
147public:
149
150
155
157 using scalar_type = Scalar;
170
172 using node_type = Node;
173
176
192 typedef Kokkos::View<impl_scalar_type *, device_type> little_vec_type;
193 typedef typename little_vec_type::HostMirror little_host_vec_type;
194
200 typedef Kokkos::View<const impl_scalar_type *, device_type>
202
204 host_device_type;
205 typedef Kokkos::View<const impl_scalar_type*, host_device_type>
206 const_little_host_vec_type;
207
209
211
217
220
223
227
231
234 const Teuchos::DataAccess copyOrView);
235
266 BlockMultiVector (const map_type& meshMap,
267 const LO blockSize,
268 const LO numVecs);
269
274 BlockMultiVector (const map_type& meshMap,
275 const map_type& pointMap,
276 const LO blockSize,
277 const LO numVecs);
278
292 const map_type& meshMap,
293 const LO blockSize);
294
300 const map_type& newMeshMap,
301 const map_type& newPointMap,
302 const size_t offset = 0);
303
309 const map_type& newMeshMap,
310 const size_t offset = 0);
311
313
315
322 static map_type
323 makePointMap (const map_type& meshMap, const LO blockSize);
324
330 return pointMap_;
331 }
332
334 LO getBlockSize () const {
335 return blockSize_;
336 }
337
339 LO getNumVectors () const {
340 return static_cast<LO> (mv_.getNumVectors ());
341 }
342
349
351
353
355 void putScalar (const Scalar& val);
356
358 void scale (const Scalar& val);
359
366 void
367 update (const Scalar& alpha,
369 const Scalar& beta);
370
392 void
393 blockWiseMultiply (const Scalar& alpha,
394 const Kokkos::View<const impl_scalar_type***,
395 device_type, Kokkos::MemoryUnmanaged>& D,
397
428 void
429 blockJacobiUpdate (const Scalar& alpha,
430 const Kokkos::View<const impl_scalar_type***,
431 device_type, Kokkos::MemoryUnmanaged>& D,
434 const Scalar& beta);
435
436
438 template<class TargetMemorySpace>
439 bool need_sync () const {
441 }
442
444 bool need_sync_host() const {
445 return mv_.need_sync_host();
446 }
447
449 bool need_sync_device() const {
450 return mv_.need_sync_device();
451 }
452
454
456
474 bool replaceLocalValues (const LO localRowIndex, const LO colIndex, const Scalar vals[]);
475
486 bool replaceGlobalValues (const GO globalRowIndex, const LO colIndex, const Scalar vals[]);
487
498 bool sumIntoLocalValues (const LO localRowIndex, const LO colIndex, const Scalar vals[]);
499
510 bool sumIntoGlobalValues (const GO globalRowIndex, const LO colIndex, const Scalar vals[]);
511
512
513 const_little_host_vec_type getLocalBlockHost(
514 const LO localRowIndex,
515 const LO colIndex,
516 const Access::ReadOnlyStruct) const;
517
518 little_host_vec_type getLocalBlockHost(
519 const LO localRowIndex,
520 const LO colIndex,
521 const Access::ReadWriteStruct);
522
526 little_host_vec_type getLocalBlockHost(
527 const LO localRowIndex,
528 const LO colIndex,
529 const Access::OverwriteAllStruct);
531
532protected:
538
539
540 virtual bool checkSizes (const Tpetra::SrcDistObject& source);
541
542 virtual void
544 (const SrcDistObject& source,
545 const size_t numSameIDs,
546 const Kokkos::DualView<const local_ordinal_type*,
547 buffer_device_type>& permuteToLIDs,
548 const Kokkos::DualView<const local_ordinal_type*,
549 buffer_device_type>& permuteFromLIDs,
550 const CombineMode CM);
551
552 virtual void
554 (const SrcDistObject& source,
555 const Kokkos::DualView<const local_ordinal_type*,
556 buffer_device_type>& exportLIDs,
557 Kokkos::DualView<packet_type*,
558 buffer_device_type>& exports,
559 Kokkos::DualView<size_t*,
560 buffer_device_type> numPacketsPerLID,
561 size_t& constantNumPackets);
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);
573
575
576protected:
577
579 size_t getStrideX () const {
580 return static_cast<size_t> (1);
581 }
582
584 size_t getStrideY () const {
585 return mv_.getStride ();
586 }
587
590 bool isValidLocalMeshIndex (const LO meshLocalIndex) const;
591
598
599private:
601 map_type pointMap_;
602
603protected:
606
607private:
608
610 LO blockSize_;
611
613 void
614 replaceLocalValuesImpl (const LO localRowIndex,
615 const LO colIndex,
616 const Scalar vals[]);
618 void
619 sumIntoLocalValuesImpl (const LO localRowIndex,
620 const LO colIndex,
621 const Scalar vals[]);
622
623 static Teuchos::RCP<const mv_type>
624 getMultiVectorFromSrcDistObject (const Tpetra::SrcDistObject&);
625
626 static Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
627 getBlockMultiVectorFromSrcDistObject (const Tpetra::SrcDistObject& src);
628};
629
630} // namespace Tpetra
631
632#endif // TPETRA_BLOCKMULTIVECTOR_DECL_HPP
Forward declaration of Tpetra::BlockCrsMatrix.
Forward declaration of Tpetra::BlockMultiVector.
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.
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).
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.
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)
Perform copies and permutations that are local to the calling (MPI) process.
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.
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
LO getNumVectors() const
Get the number of columns (vectors) in the BlockMultiVector.
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.
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.
Kokkos::View< impl_scalar_type *, device_type > little_vec_type
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
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.
map_type getPointMap() const
Get this BlockMultiVector's (previously computed) point Map.
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
typename ::Kokkos::Details::ArithTraits< Scalar >::val_type packet_type
A parallel distribution of indices over processes.
One or more distributed dense vectors.
typename Kokkos::Details::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.