Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends 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 // ************************************************************************
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 
49 namespace Tpetra {
50 
139 template<class Scalar,
140  class LO,
141  class GO,
142  class Node>
144  public Tpetra::DistObject<Scalar, LO, GO, Node>
145 {
146 private:
148  using STS = Teuchos::ScalarTraits<Scalar>;
149  using packet_type = typename dist_object_type::packet_type;
150 
151 public:
153 
154 
159 
161  using scalar_type = Scalar;
169  using local_ordinal_type = LO;
174 
176  using node_type = Node;
177 
180 
196  typedef Kokkos::View<impl_scalar_type*,
197  Kokkos::LayoutRight,
198  device_type,
199  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
201 
207  typedef Kokkos::View<const impl_scalar_type*,
208  Kokkos::LayoutRight,
209  device_type,
210  Kokkos::MemoryTraits<Kokkos::Unmanaged> >
212 
214 
216 
221  BlockMultiVector ();
222 
225 
228 
232 
236 
239  const Teuchos::DataAccess copyOrView);
240 
271  BlockMultiVector (const map_type& meshMap,
272  const LO blockSize,
273  const LO numVecs);
274 
279  BlockMultiVector (const map_type& meshMap,
280  const map_type& pointMap,
281  const LO blockSize,
282  const LO numVecs);
283 
296  BlockMultiVector (const mv_type& X_mv,
297  const map_type& meshMap,
298  const LO blockSize);
299 
305  const map_type& newMeshMap,
306  const map_type& newPointMap,
307  const size_t offset = 0);
308 
314  const map_type& newMeshMap,
315  const size_t offset = 0);
316 
318 
320 
327  static map_type
328  makePointMap (const map_type& meshMap, const LO blockSize);
329 
335  return pointMap_;
336  }
337 
339  LO getBlockSize () const {
340  return blockSize_;
341  }
342 
344  LO getNumVectors () const {
345  return static_cast<LO> (mv_.getNumVectors ());
346  }
347 
353  mv_type getMultiVectorView () const;
354 
356 
358 
360  void putScalar (const Scalar& val);
361 
363  void scale (const Scalar& val);
364 
371  void
372  update (const Scalar& alpha,
374  const Scalar& beta);
375 
397  void
398  blockWiseMultiply (const Scalar& alpha,
399  const Kokkos::View<const impl_scalar_type***,
400  device_type, Kokkos::MemoryUnmanaged>& D,
402 
433  void
434  blockJacobiUpdate (const Scalar& alpha,
435  const Kokkos::View<const impl_scalar_type***,
436  device_type, Kokkos::MemoryUnmanaged>& D,
439  const Scalar& beta);
440 
442 
444 
463  template<class TargetMemorySpace>
464  void sync () {
465  mv_.template sync<typename TargetMemorySpace::memory_space> ();
466  }
467 
469  void sync_host() {
470  mv_.sync_host();
471  }
472 
474  void sync_device() {
475  mv_.sync_device();
476  }
477 
479  template<class TargetMemorySpace>
480  bool need_sync () const {
481  return mv_.template need_sync<typename TargetMemorySpace::memory_space> ();
482  }
483 
485  bool need_sync_host() const {
486  return mv_.need_sync_host();
487  }
488 
490  bool need_sync_device() const {
491  return mv_.need_sync_device();
492  }
493 
499  template<class TargetMemorySpace>
500  void modify () {
501  mv_.template modify<typename TargetMemorySpace::memory_space> ();
502  }
503 
505  void modify_host() {
506  mv_.modify_host();
507  }
508 
510  void modify_device() {
511  mv_.modify_device();
512  }
514 
516 
534  bool replaceLocalValues (const LO localRowIndex, const LO colIndex, const Scalar vals[]) const;
535 
546  bool replaceGlobalValues (const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const;
547 
558  bool sumIntoLocalValues (const LO localRowIndex, const LO colIndex, const Scalar vals[]) const;
559 
570  bool sumIntoGlobalValues (const GO globalRowIndex, const LO colIndex, const Scalar vals[]) const;
571 
582  bool getLocalRowView (const LO localRowIndex, const LO colIndex, Scalar*& vals) const;
583 
594  bool getGlobalRowView (const GO globalRowIndex, const LO colIndex, Scalar*& vals) const;
595 
607  typename little_vec_type::HostMirror
608  getLocalBlock (const LO localRowIndex, const LO colIndex) const;
610 
611 protected:
617 
618 
619  virtual bool checkSizes (const Tpetra::SrcDistObject& source);
620 
621  virtual void
622  copyAndPermute
623  (const SrcDistObject& source,
624  const size_t numSameIDs,
625  const Kokkos::DualView<const local_ordinal_type*,
626  buffer_device_type>& permuteToLIDs,
627  const Kokkos::DualView<const local_ordinal_type*,
628  buffer_device_type>& permuteFromLIDs);
629 
630  virtual void
631  packAndPrepare
632  (const SrcDistObject& source,
633  const Kokkos::DualView<const local_ordinal_type*,
634  buffer_device_type>& exportLIDs,
635  Kokkos::DualView<packet_type*,
636  buffer_device_type>& exports,
637  Kokkos::DualView<size_t*,
638  buffer_device_type> numPacketsPerLID,
639  size_t& constantNumPackets,
640  Distributor& distor);
641 
642  virtual void
643  unpackAndCombine
644  (const Kokkos::DualView<const local_ordinal_type*,
645  buffer_device_type>& importLIDs,
646  Kokkos::DualView<packet_type*,
647  buffer_device_type> imports,
648  Kokkos::DualView<size_t*,
649  buffer_device_type> numPacketsPerLID,
650  const size_t constantNumPackets,
651  Distributor& distor,
652  const CombineMode combineMode);
653 
655 
656 protected:
659  return mvData_;
660  }
661 
663  size_t getStrideX () const {
664  return static_cast<size_t> (1);
665  }
666 
668  size_t getStrideY () const {
669  return mv_.getStride ();
670  }
671 
674  bool isValidLocalMeshIndex (const LO meshLocalIndex) const;
675 
676 private:
682  map_type meshMap_;
683 
685  map_type pointMap_;
686 
687 protected:
690 
691 private:
700  impl_scalar_type* mvData_;
701 
703  LO blockSize_;
704 
706  void
707  replaceLocalValuesImpl (const LO localRowIndex,
708  const LO colIndex,
709  const Scalar vals[]) const;
711  void
712  sumIntoLocalValuesImpl (const LO localRowIndex,
713  const LO colIndex,
714  const Scalar vals[]) const;
715 
716  static Teuchos::RCP<const mv_type>
717  getMultiVectorFromSrcDistObject (const Tpetra::SrcDistObject&);
718 
719  static Teuchos::RCP<const BlockMultiVector<Scalar, LO, GO, Node> >
720  getBlockMultiVectorFromSrcDistObject (const Tpetra::SrcDistObject& src);
721 };
722 
723 } // namespace Tpetra
724 
725 #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.