Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 // @HEADER
41 
42 #ifndef TPETRA_BLOCKMULTIVECTOR_DECL_HPP
43 #define TPETRA_BLOCKMULTIVECTOR_DECL_HPP
44 
47 #include "Tpetra_MultiVector.hpp"
48 #include <memory>
49 #include <utility>
50 
51 namespace Tpetra {
52 
141 template<class Scalar,
142  class LO,
143  class GO,
144  class Node>
146  public Tpetra::DistObject<Scalar, LO, GO, Node>
147 {
148 private:
150  using STS = Teuchos::ScalarTraits<Scalar>;
151  using packet_type = typename dist_object_type::packet_type;
152 
153 public:
155 
156 
161 
163  using scalar_type = Scalar;
171  using local_ordinal_type = LO;
176 
178  using node_type = Node;
179 
182 
198  typedef Kokkos::View<impl_scalar_type*,
199  Kokkos::LayoutRight,
200  device_type,
201  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
203 
209  typedef Kokkos::View<const impl_scalar_type*,
210  Kokkos::LayoutRight,
211  device_type,
212  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
214 
216 
218 
223  BlockMultiVector ();
224 
227 
230 
234 
238 
241  const Teuchos::DataAccess copyOrView);
242 
273  BlockMultiVector (const map_type& meshMap,
274  const LO blockSize,
275  const LO numVecs);
276 
281  BlockMultiVector (const map_type& meshMap,
282  const map_type& pointMap,
283  const LO blockSize,
284  const LO numVecs);
285 
298  BlockMultiVector (const mv_type& X_mv,
299  const map_type& meshMap,
300  const LO blockSize);
301 
307  const map_type& newMeshMap,
308  const map_type& newPointMap,
309  const size_t offset = 0);
310 
316  const map_type& newMeshMap,
317  const size_t offset = 0);
318 
320 
322 
329  static map_type
330  makePointMap (const map_type& meshMap, const LO blockSize);
331 
337  return pointMap_;
338  }
339 
341  LO getBlockSize () const {
342  return blockSize_;
343  }
344 
346  LO getNumVectors () const {
347  return static_cast<LO> (mv_.getNumVectors ());
348  }
349 
355  mv_type getMultiVectorView () const;
356 
358 
360 
362  void putScalar (const Scalar& val);
363 
365  void scale (const Scalar& val);
366 
373  void
374  update (const Scalar& alpha,
376  const Scalar& beta);
377 
399  void
400  blockWiseMultiply (const Scalar& alpha,
401  const Kokkos::View<const impl_scalar_type***,
402  device_type, Kokkos::MemoryUnmanaged>& D,
404 
435  void
436  blockJacobiUpdate (const Scalar& alpha,
437  const Kokkos::View<const impl_scalar_type***,
438  device_type, Kokkos::MemoryUnmanaged>& D,
441  const Scalar& beta);
442 
444 
446 
465  template<class TargetMemorySpace>
466  void sync () {
467  mv_.template sync<typename TargetMemorySpace::memory_space> ();
468  }
469 
471  void sync_host() {
472  mv_.sync_host();
473  }
474 
476  void sync_device() {
477  mv_.sync_device();
478  }
479 
481  template<class TargetMemorySpace>
482  bool need_sync () const {
483  return mv_.template need_sync<typename TargetMemorySpace::memory_space> ();
484  }
485 
487  bool need_sync_host() const {
488  return mv_.need_sync_host();
489  }
490 
492  bool need_sync_device() const {
493  return mv_.need_sync_device();
494  }
495 
501  template<class TargetMemorySpace>
502  void modify () {
503  mv_.template modify<typename TargetMemorySpace::memory_space> ();
504  }
505 
507  void modify_host() {
508  mv_.modify_host();
509  }
510 
512  void modify_device() {
513  mv_.modify_device();
514  }
516 
518 
536  bool replaceLocalValues (const LO localRowIndex, const LO colIndex, const Scalar vals[]) const;
537 
548  bool replaceGlobalValues (const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const;
549 
560  bool sumIntoLocalValues (const LO localRowIndex, const LO colIndex, const Scalar vals[]) const;
561 
572  bool sumIntoGlobalValues (const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const;
573 
584  bool getLocalRowView (const LO localRowIndex, const LO colIndex, Scalar*& vals) const;
585 
596  bool getGlobalRowView (const GO globalRowIndex, const LO colIndex, Scalar*& vals) const;
597 
609  typename little_vec_type::HostMirror
610  getLocalBlock (const LO localRowIndex, const LO colIndex) const;
612 
613 protected:
619 
620 
621  virtual bool checkSizes (const Tpetra::SrcDistObject& source);
622 
623  virtual void
624 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
625  copyAndPermuteNew
626 #else // TPETRA_ENABLE_DEPRECATED_CODE
627  copyAndPermute
628 #endif // TPETRA_ENABLE_DEPRECATED_CODE
629  (const SrcDistObject& source,
630  const size_t numSameIDs,
631  const Kokkos::DualView<const local_ordinal_type*,
632  buffer_device_type>& permuteToLIDs,
633  const Kokkos::DualView<const local_ordinal_type*,
634  buffer_device_type>& permuteFromLIDs);
635 
636  virtual void
637 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
638  packAndPrepareNew
639 #else // TPETRA_ENABLE_DEPRECATED_CODE
640  packAndPrepare
641 #endif // TPETRA_ENABLE_DEPRECATED_CODE
642  (const SrcDistObject& source,
643  const Kokkos::DualView<const local_ordinal_type*,
644  buffer_device_type>& exportLIDs,
645  Kokkos::DualView<packet_type*,
646  buffer_device_type>& exports,
647  Kokkos::DualView<size_t*,
648  buffer_device_type> numPacketsPerLID,
649  size_t& constantNumPackets,
650  Distributor& distor);
651 
652  virtual void
653 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
654  unpackAndCombineNew
655 #else // TPETRA_ENABLE_DEPRECATED_CODE
656  unpackAndCombine
657 #endif // TPETRA_ENABLE_DEPRECATED_CODE
658  (const Kokkos::DualView<const local_ordinal_type*,
659  buffer_device_type>& importLIDs,
660  Kokkos::DualView<packet_type*,
661  buffer_device_type> imports,
662  Kokkos::DualView<size_t*,
663  buffer_device_type> numPacketsPerLID,
664  const size_t constantNumPackets,
665  Distributor& distor,
666  const CombineMode combineMode);
667 
669 
670 protected:
673  return mvData_;
674  }
675 
677  size_t getStrideX () const {
678  return static_cast<size_t> (1);
679  }
680 
682  size_t getStrideY () const {
683  return mv_.getStride ();
684  }
685 
688  bool isValidLocalMeshIndex (const LO meshLocalIndex) const;
689 
690 private:
696  map_type meshMap_;
697 
699  map_type pointMap_;
700 
701 protected:
704 
705 private:
714  impl_scalar_type* mvData_;
715 
717  LO blockSize_;
718 
720  void
721  replaceLocalValuesImpl (const LO localRowIndex,
722  const LO colIndex,
723  const Scalar vals[]) const;
725  void
726  sumIntoLocalValuesImpl (const LO localRowIndex,
727  const LO colIndex,
728  const Scalar vals[]) const;
729 
730  static Teuchos::RCP<const mv_type>
731  getMultiVectorFromSrcDistObject (const Tpetra::SrcDistObject&);
732 
733  static Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
734  getBlockMultiVectorFromSrcDistObject (const Tpetra::SrcDistObject& src);
735 };
736 
737 } // namespace Tpetra
738 
739 #endif // TPETRA_BLOCKMULTIVECTOR_DECL_HPP
Node node_type
The Kokkos Node type; a legacy thing that will go away at some point.
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.
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.
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.
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. ...
void sync()
Update data to the given target memory space, only if data in the &quot;other&quot; space have been marked as m...
void putScalar(const Scalar &val)
Fill all entries with the given value val.
bool need_sync_device() const
Whether this object needs synchronization to the device.
LO local_ordinal_type
The type of local indices.
size_t getNumVectors() const
Number of columns in the multivector.
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.
size_t getStrideX() const
Stride between consecutive local entries in the same column.
Forward declaration of Tpetra::BlockCrsMatrix.
GO global_ordinal_type
The type of global indices.
virtual bool checkSizes(const Tpetra::SrcDistObject &source)
Compare the source and target (this) objects for compatibility.
void modify_device()
Mark data as modified on the device.
size_t getStrideY() const
Stride between consecutive local entries in the same row.
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.
bool need_sync_device() const
Whether this MultiVector needs synchronization to the device.
MultiVector for multiple degrees of freedom per mesh point.
bool need_sync_host() const
Whether this object needs synchronization to the host.
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.
void sync_host()
Synchronize to Host.
bool isValidLocalMeshIndex(const LO meshLocalIndex) const
True if and only if meshLocalIndex is a valid local index in the mesh Map.
void modify()
Mark data as modified on the given memory space.
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.
Sets up and executes a communication plan for a Tpetra DistObject.
typename mv_type::device_type device_type
The Kokkos Device type.
CombineMode
Rule for combining data in an Import or Export.
map_type getPointMap() const
Get this BlockMultiVector&#39;s (previously computed) point Map.
void modify_host()
Mark data as modified on the host.
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.
bool need_sync_host() const
Whether this MultiVector needs synchronization to the host.
void modify_device()
Mark data as modified on the device side.
mv_type mv_
The Tpetra::MultiVector used to represent the data.
Abstract base class for objects that can be the source of an Import or Export operation.
Forward declaration of Tpetra::BlockMultiVector.
void modify_host()
Mark data as modified on the host side.
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...
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.
Scalar scalar_type
The type of entries in the object.
typename map_type::device_type device_type
This class&#39; preferred Kokkos device type.
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 need_sync() const
Whether this object needs synchronization to the given memory space.
typename Kokkos::Details::ArithTraits< Scalar >::val_type impl_scalar_type
The type used internally in place of Scalar.
LO getNumVectors() const
Get the number of columns (vectors) in the BlockMultiVector.
mv_type getMultiVectorView() const
Get a Tpetra::MultiVector that views this BlockMultiVector&#39;s data.
typename dist_object_type::buffer_device_type buffer_device_type
Kokkos::Device specialization used for communication buffers.
size_t getStride() const
Stride between columns in the multivector.
impl_scalar_type * getRawPtr() const
Raw pointer to the MultiVector&#39;s data.
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type
void scale(const Scalar &val)
Multiply all entries in place by the given value val.
Base class for distributed Tpetra objects that support data redistribution.
typename mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the object.
void sync_device()
Update data to the device.
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.
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...
void sync_host()
Update data to the host.
void sync_device()
Synchronize to Device.