Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tpetra_Experimental_BlockVector_def.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_BLOCKVECTOR_DEF_HPP
43 #define TPETRA_EXPERIMENTAL_BLOCKVECTOR_DEF_HPP
44 
45 namespace Tpetra {
46 namespace Experimental {
47 
48  template<class Scalar, class LO, class GO, class Node>
50  BlockVector (const map_type& meshMap, const LO blockSize) :
51  base_type (meshMap, blockSize, 1)
52  {}
53 
54  template<class Scalar, class LO, class GO, class Node>
56  BlockVector (const map_type& meshMap,
57  const map_type& pointMap,
58  const LO blockSize) :
59  base_type (meshMap, pointMap, blockSize, 1)
60  {}
61 
62  template<class Scalar, class LO, class GO, class Node>
64  BlockVector (const mv_type& X_mv,
65  const map_type& meshMap,
66  const LO blockSize) :
67  base_type (X_mv, meshMap, blockSize)
68  {
69  TEUCHOS_TEST_FOR_EXCEPTION(
70  X_mv.getNumVectors () != 1, std::invalid_argument,
71  "Tpetra::Experimental::BlockVector: Input MultiVector has "
72  << X_mv.getNumVectors () << " != 1 columns.");
73  }
74 
75  template<class Scalar, class LO, class GO, class Node>
77  BlockVector (const vec_type& X_vec,
78  const map_type& meshMap,
79  const LO blockSize) :
80  base_type (X_vec, meshMap, blockSize)
81  {}
82 
83  template<class Scalar, class LO, class GO, class Node>
86  const map_type& newMeshMap,
87  const map_type& newPointMap,
88  const size_t offset) :
89  base_type (X, newMeshMap, newPointMap, offset)
90  {}
91 
92  template<class Scalar, class LO, class GO, class Node>
95  const map_type& newMeshMap,
96  const size_t offset) :
97  base_type (X, newMeshMap, offset)
98  {}
99 
100  template<class Scalar, class LO, class GO, class Node>
103 
104  template<class Scalar, class LO, class GO, class Node>
107  Teuchos::RCP<vec_type> vPtr = this->mv_.getVectorNonConst (0);
108  vPtr->setCopyOrView (Teuchos::View);
109  return *vPtr; // shallow copy
110  }
111 
112  template<class Scalar, class LO, class GO, class Node>
113  bool
115  replaceLocalValues (const LO localRowIndex, const Scalar vals[]) const {
116  return ((const base_type*) this)->replaceLocalValues (localRowIndex, 0, vals);
117  }
118 
119  template<class Scalar, class LO, class GO, class Node>
120  bool
122  replaceGlobalValues (const GO globalRowIndex, const Scalar vals[]) const {
123  return ((const base_type*) this)->replaceGlobalValues (globalRowIndex, 0, vals);
124  }
125 
126 
127  template<class Scalar, class LO, class GO, class Node>
128  bool
130  sumIntoLocalValues (const LO localRowIndex, const Scalar vals[]) const {
131  return ((const base_type*) this)->sumIntoLocalValues (localRowIndex, 0, vals);
132  }
133 
134  template<class Scalar, class LO, class GO, class Node>
135  bool
137  sumIntoGlobalValues (const GO globalRowIndex, const Scalar vals[]) const {
138  return ((const base_type*) this)->sumIntoLocalValues (globalRowIndex, 0, vals);
139  }
140 
141  template<class Scalar, class LO, class GO, class Node>
142  bool
144  getLocalRowView (const LO localRowIndex, Scalar*& vals) const {
145  return ((const base_type*) this)->getLocalRowView (localRowIndex, 0, vals);
146  }
147 
148  template<class Scalar, class LO, class GO, class Node>
149  bool
151  getGlobalRowView (const GO globalRowIndex, Scalar*& vals) const {
152  return ((const base_type*) this)->getGlobalRowView (globalRowIndex, 0, vals);
153  }
154 
166  template<class Scalar, class LO, class GO, class Node>
169  getLocalBlock (const LO localRowIndex) const
170  {
171  if (! this->isValidLocalMeshIndex (localRowIndex)) {
172  return little_vec_type ();
173  } else {
174  const size_t blockSize = this->getBlockSize ();
175  const size_t offset = localRowIndex * blockSize;
176  return little_vec_type (this->getRawPtr () + offset, blockSize);
177  }
178  }
179 
180 } // namespace Experimental
181 } // namespace Tpetra
182 
183 //
184 // Explicit instantiation macro
185 //
186 // Must be expanded from within the Tpetra namespace!
187 //
188 #define TPETRA_EXPERIMENTAL_BLOCKVECTOR_INSTANT(S,LO,GO,NODE) \
189  namespace Experimental { \
190  template class BlockVector< S, LO, GO, NODE >; \
191  }
192 
193 #endif // TPETRA_EXPERIMENTAL_BLOCKVECTOR_DEF_HPP
bool getLocalRowView(const LO localRowIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a local index. ...
size_t getNumVectors() const
Number of columns in the multivector.
MultiVector for multiple degrees of freedom per mesh point.
Vector for multiple degrees of freedom per mesh point.
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...
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.
bool sumIntoGlobalValues(const GO globalRowIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a global index.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getVectorNonConst(const size_t j)
Return a Vector which is a nonconst view of column j.
vec_type getVectorView()
Get a Tpetra::Vector that views this BlockVector&#39;s data.
bool getGlobalRowView(const GO globalRowIndex, Scalar *&vals) const
Get a writeable view of the entries at the given mesh point, using a global index.
little_vec_type getLocalBlock(const LO localRowIndex) const
Get a view of the degrees of freedom at the given mesh point, using a local index.
bool sumIntoLocalValues(const LO localRowIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a local index.
A distributed dense vector.
bool replaceGlobalValues(const GO globalRowIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using a global index.
bool replaceLocalValues(const LO localRowIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using a local index.