Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tpetra_Experimental_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 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 // @HEADER
41 
42 #ifndef TPETRA_EXPERIMENTAL_BLOCKMULTIVECTOR_DECL_HPP
43 #define TPETRA_EXPERIMENTAL_BLOCKMULTIVECTOR_DECL_HPP
44 
47 #include "Tpetra_MultiVector.hpp"
48 #include <memory>
49 #include <utility>
50 
51 namespace Tpetra {
52 namespace Experimental {
53 
142 template<class Scalar,
143  class LO,
144  class GO,
145  class Node>
147  public Tpetra::DistObject<Scalar, LO, GO, Node>
148 {
149 private:
151  using STS = Teuchos::ScalarTraits<Scalar>;
152  using packet_type = typename dist_object_type::packet_type;
153 
154 public:
156 
157 
162 
164  using scalar_type = Scalar;
172  using local_ordinal_type = LO;
177 
179  using node_type = Node;
180 
183 
199  typedef Kokkos::View<impl_scalar_type*,
200  Kokkos::LayoutRight,
201  device_type,
202  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
204 
210  typedef Kokkos::View<const impl_scalar_type*,
211  Kokkos::LayoutRight,
212  device_type,
213  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
215 
217 
219 
250  BlockMultiVector (const map_type& meshMap,
251  const LO blockSize,
252  const LO numVecs);
253 
258  BlockMultiVector (const map_type& meshMap,
259  const map_type& pointMap,
260  const LO blockSize,
261  const LO numVecs);
262 
275  BlockMultiVector (const mv_type& X_mv,
276  const map_type& meshMap,
277  const LO blockSize);
278 
284  const map_type& newMeshMap,
285  const map_type& newPointMap,
286  const size_t offset = 0);
287 
293  const map_type& newMeshMap,
294  const size_t offset = 0);
295 
300  BlockMultiVector ();
301 
303 
305 
312  static map_type
313  makePointMap (const map_type& meshMap, const LO blockSize);
314 
320  return pointMap_;
321  }
322 
324  LO getBlockSize () const {
325  return blockSize_;
326  }
327 
329  LO getNumVectors () const {
330  return static_cast<LO> (mv_.getNumVectors ());
331  }
332 
338  mv_type getMultiVectorView () const;
339 
341 
343 
345  void putScalar (const Scalar& val);
346 
348  void scale (const Scalar& val);
349 
356  void
357  update (const Scalar& alpha,
359  const Scalar& beta);
360 
382  void
383  blockWiseMultiply (const Scalar& alpha,
384  const Kokkos::View<const impl_scalar_type***,
385  device_type, Kokkos::MemoryUnmanaged>& D,
387 
418  void
419  blockJacobiUpdate (const Scalar& alpha,
420  const Kokkos::View<const impl_scalar_type***,
421  device_type, Kokkos::MemoryUnmanaged>& D,
424  const Scalar& beta);
425 
427 
429 
448  template<class TargetMemorySpace>
449  void sync () {
450  mv_.template sync<typename TargetMemorySpace::memory_space> ();
451  }
452 
454  void sync_host() {
455  mv_.sync_host();
456  }
457 
459  void sync_device() {
460  mv_.sync_device();
461  }
462 
464  template<class TargetMemorySpace>
465  bool need_sync () const {
466  return mv_.template need_sync<typename TargetMemorySpace::memory_space> ();
467  }
468 
470  bool need_sync_host() const {
471  return mv_.need_sync_host();
472  }
473 
475  bool need_sync_device() const {
476  return mv_.need_sync_device();
477  }
478 
484  template<class TargetMemorySpace>
485  void modify () {
486  mv_.template modify<typename TargetMemorySpace::memory_space> ();
487  }
488 
490  void modify_host() {
491  mv_.modify_host();
492  }
493 
495  void modify_device() {
496  mv_.modify_device();
497  }
499 
501 
519  bool replaceLocalValues (const LO localRowIndex, const LO colIndex, const Scalar vals[]) const;
520 
531  bool replaceGlobalValues (const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const;
532 
543  bool sumIntoLocalValues (const LO localRowIndex, const LO colIndex, const Scalar vals[]) const;
544 
555  bool sumIntoGlobalValues (const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const;
556 
567  bool getLocalRowView (const LO localRowIndex, const LO colIndex, Scalar*& vals) const;
568 
579  bool getGlobalRowView (const GO globalRowIndex, const LO colIndex, Scalar*& vals) const;
580 
592  typename little_vec_type::HostMirror
593  getLocalBlock (const LO localRowIndex, const LO colIndex) const;
595 
596 protected:
602 
603 
604  virtual bool checkSizes (const Tpetra::SrcDistObject& source);
605 
606  virtual void
607 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
608  copyAndPermuteNew
609 #else // TPETRA_ENABLE_DEPRECATED_CODE
610  copyAndPermute
611 #endif // TPETRA_ENABLE_DEPRECATED_CODE
612  (const SrcDistObject& source,
613  const size_t numSameIDs,
614  const Kokkos::DualView<const local_ordinal_type*,
615  buffer_device_type>& permuteToLIDs,
616  const Kokkos::DualView<const local_ordinal_type*,
617  buffer_device_type>& permuteFromLIDs);
618 
619  virtual void
620 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
621  packAndPrepareNew
622 #else // TPETRA_ENABLE_DEPRECATED_CODE
623  packAndPrepare
624 #endif // TPETRA_ENABLE_DEPRECATED_CODE
625  (const SrcDistObject& source,
626  const Kokkos::DualView<const local_ordinal_type*,
627  buffer_device_type>& exportLIDs,
628  Kokkos::DualView<packet_type*,
629  buffer_device_type>& exports,
630  Kokkos::DualView<size_t*,
631  buffer_device_type> numPacketsPerLID,
632  size_t& constantNumPackets,
633  Distributor& distor);
634 
635  virtual void
636 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
637  unpackAndCombineNew
638 #else // TPETRA_ENABLE_DEPRECATED_CODE
639  unpackAndCombine
640 #endif // TPETRA_ENABLE_DEPRECATED_CODE
641  (const Kokkos::DualView<const local_ordinal_type*,
642  buffer_device_type>& importLIDs,
643  Kokkos::DualView<packet_type*,
644  buffer_device_type> imports,
645  Kokkos::DualView<size_t*,
646  buffer_device_type> numPacketsPerLID,
647  const size_t constantNumPackets,
648  Distributor& distor,
649  const CombineMode combineMode);
650 
652 
653 protected:
656  return mvData_;
657  }
658 
660  size_t getStrideX () const {
661  return static_cast<size_t> (1);
662  }
663 
665  size_t getStrideY () const {
666  return mv_.getStride ();
667  }
668 
671  bool isValidLocalMeshIndex (const LO meshLocalIndex) const;
672 
673 private:
679  map_type meshMap_;
680 
682  map_type pointMap_;
683 
684 protected:
687 
688 private:
697  impl_scalar_type* mvData_;
698 
700  LO blockSize_;
701 
703  void
704  replaceLocalValuesImpl (const LO localRowIndex,
705  const LO colIndex,
706  const Scalar vals[]) const;
708  void
709  sumIntoLocalValuesImpl (const LO localRowIndex,
710  const LO colIndex,
711  const Scalar vals[]) const;
712 
713  static Teuchos::RCP<const mv_type>
714  getMultiVectorFromSrcDistObject (const Tpetra::SrcDistObject&);
715 
716  static Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
717  getBlockMultiVectorFromSrcDistObject (const Tpetra::SrcDistObject& src);
718 
735  std::pair<int, std::unique_ptr<std::string>>
736  copyAndPermuteImpl
738  const size_t numSameIDs,
739  const Kokkos::DualView<
740  const local_ordinal_type*,
742  >& permuteToLIDs,
743  const Kokkos::DualView<
744  const local_ordinal_type*,
746  >& permuteFromLIDs);
747 
748  std::pair<int, std::unique_ptr<std::string>>
749  packAndPrepareImpl
750  (const Kokkos::DualView<
751  const local_ordinal_type*,
753  >& exportLIDs,
754  Kokkos::DualView<
757  >& exports,
758  Kokkos::DualView<
759  size_t*,
761  > /* numPacketsPerLID */,
762  size_t& constantNumPackets,
763  Distributor& /* distor */ ) const;
764 
765  std::pair<int, std::unique_ptr<std::string>>
766  unpackAndCombineImpl
767  (const Kokkos::DualView<
768  const local_ordinal_type*,
770  >& importLIDs,
771  Kokkos::DualView<
774  > imports,
775  Kokkos::DualView<
776  size_t*,
778  > numPacketsPerLID,
779  const size_t constantNumPackets,
780  Distributor& distor,
781  const CombineMode combineMode);
782 };
783 
784 } // namespace Experimental
785 } // namespace Tpetra
786 
787 #endif // TPETRA_EXPERIMENTAL_BLOCKMULTIVECTOR_DECL_HPP
void sync()
Update data to the given target memory space, only if data in the &quot;other&quot; space have been marked as m...
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.
Tpetra::Map< LO, GO, Node > map_type
The specialization of Tpetra::Map that this class uses.
Forward declaration of Tpetra::Experimental::BlockCrsMatrix.
bool need_sync_host() const
Whether this object needs synchronization to the host.
mv_type mv_
The Tpetra::MultiVector used to represent the data.
void modify_device()
Mark data as modified on the device.
bool need_sync_device() const
Whether this object needs synchronization to the device.
void putScalar(const Scalar &val)
Fill all entries with the given value val.
size_t getNumVectors() const
Number of columns in the multivector.
bool replaceGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using a global index.
bool sumIntoGlobalValues(const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a global index.
MultiVector for multiple degrees of freedom per mesh point.
Forward declaration of Tpetra::Experimental::BlockMultiVector.
bool need_sync() const
Whether this object needs synchronization to the given memory space.
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.
Kokkos::View< impl_scalar_type *, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_vec_type
&quot;Block view&quot; of all degrees of freedom at a mesh point, for a single column of the MultiVector...
little_vec_type::HostMirror getLocalBlock(const LO localRowIndex, const LO colIndex) const
Get a host view of the degrees of freedom at the given mesh point.
void scale(const Scalar &val)
Multiply all entries in place by the given value val.
void update(const Scalar &alpha, const BlockMultiVector< Scalar, LO, GO, Node > &X, const Scalar &beta)
Update: this = beta*this + alpha*X.
bool need_sync_device() const
Whether this MultiVector needs synchronization to the device.
mv_type getMultiVectorView() const
Get a Tpetra::MultiVector that views this BlockMultiVector&#39;s data.
LO getBlockSize() const
Get the number of degrees of freedom per mesh point.
typename::Kokkos::Details::ArithTraits< Scalar >::val_type packet_type
The type of each datum being sent or received in an Import or Export.
impl_scalar_type * getRawPtr() const
Raw pointer to the MultiVector&#39;s data.
void sync_host()
Synchronize to Host.
LO getNumVectors() const
Get the number of columns (vectors) in the BlockMultiVector.
Sets up and executes a communication plan for a Tpetra DistObject.
CombineMode
Rule for combining data in an Import or Export.
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 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 getLocalRowView(const LO localRowIndex, const LO colIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a local index. ...
bool need_sync_host() const
Whether this MultiVector needs synchronization to the host.
void modify_device()
Mark data as modified on the device side.
typename dist_object_type::buffer_device_type buffer_device_type
Kokkos::Device specialization used for communication buffers.
Abstract base class for objects that can be the source of an Import or Export operation.
bool isValidLocalMeshIndex(const LO meshLocalIndex) const
True if and only if meshLocalIndex is a valid local index in the mesh Map.
void modify_host()
Mark data as modified on the host side.
typename map_type::device_type device_type
This class&#39; preferred Kokkos device type.
void modify()
Mark data as modified on the given memory space.
map_type getPointMap() const
Get this BlockMultiVector&#39;s (previously computed) point Map.
typename Kokkos::Details::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
virtual bool checkSizes(const Tpetra::SrcDistObject &source)
Compare the source and target (this) objects for compatibility.
Kokkos::View< const impl_scalar_type *, Kokkos::LayoutRight, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > const_little_vec_type
&quot;Const block view&quot; of all degrees of freedom at a mesh point, for a single column of the MultiVector...
size_t getStrideY() const
Stride between consecutive local entries in the same row.
size_t getStride() const
Stride between columns in the multivector.
typename mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the object.
bool replaceLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using local row and column indices.
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type
size_t getStrideX() const
Stride between consecutive local entries in the same column.
Base class for distributed Tpetra objects that support data redistribution.
bool sumIntoLocalValues(const LO localRowIndex, const LO colIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a local index.
typename mv_type::device_type device_type
The Kokkos Device type.
void sync_device()
Synchronize to Device.
bool getGlobalRowView(const GO globalRowIndex, const LO colIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a global index.