Tpetra parallel linear algebra  Version of the Day
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 #include "Tpetra_Experimental_BlockVector_decl.hpp"
46 
47 namespace Tpetra {
48 namespace Experimental {
49 
50  template<class Scalar, class LO, class GO, class Node>
52  BlockVector (const map_type& meshMap, const LO blockSize) :
53  base_type (meshMap, blockSize, 1)
54  {}
55 
56  template<class Scalar, class LO, class GO, class Node>
58  BlockVector (const map_type& meshMap,
59  const map_type& pointMap,
60  const LO blockSize) :
61  base_type (meshMap, pointMap, blockSize, 1)
62  {}
63 
64  template<class Scalar, class LO, class GO, class Node>
66  BlockVector (const mv_type& X_mv,
67  const map_type& meshMap,
68  const LO blockSize) :
69  base_type (X_mv, meshMap, blockSize)
70  {
71  TEUCHOS_TEST_FOR_EXCEPTION(
72  X_mv.getNumVectors () != 1, std::invalid_argument,
73  "Tpetra::Experimental::BlockVector: Input MultiVector has "
74  << X_mv.getNumVectors () << " != 1 columns.");
75  }
76 
77  template<class Scalar, class LO, class GO, class Node>
80 
81  template<class Scalar, class LO, class GO, class Node>
84  Teuchos::RCP<vec_type> vPtr = this->mv_.getVectorNonConst (0);
85  vPtr->setCopyOrView (Teuchos::View);
86  return *vPtr; // shallow copy
87  }
88 
89  template<class Scalar, class LO, class GO, class Node>
90  bool
92  replaceLocalValues (const LO localRowIndex, const Scalar vals[]) const {
93  return ((const base_type*) this)->replaceLocalValues (localRowIndex, 0, vals);
94  }
95 
96  template<class Scalar, class LO, class GO, class Node>
97  bool
99  replaceGlobalValues (const GO globalRowIndex, const Scalar vals[]) const {
100  return ((const base_type*) this)->replaceGlobalValues (globalRowIndex, 0, vals);
101  }
102 
103 
104  template<class Scalar, class LO, class GO, class Node>
105  bool
107  sumIntoLocalValues (const LO localRowIndex, const Scalar vals[]) const {
108  return ((const base_type*) this)->sumIntoLocalValues (localRowIndex, 0, vals);
109  }
110 
111  template<class Scalar, class LO, class GO, class Node>
112  bool
114  sumIntoGlobalValues (const GO globalRowIndex, const Scalar vals[]) const {
115  return ((const base_type*) this)->sumIntoLocalValues (globalRowIndex, 0, vals);
116  }
117 
118  template<class Scalar, class LO, class GO, class Node>
119  bool
121  getLocalRowView (const LO localRowIndex, Scalar*& vals) const {
122  return ((const base_type*) this)->getLocalRowView (localRowIndex, 0, vals);
123  }
124 
125  template<class Scalar, class LO, class GO, class Node>
126  bool
128  getGlobalRowView (const GO globalRowIndex, Scalar*& vals) const {
129  return ((const base_type*) this)->getGlobalRowView (globalRowIndex, 0, vals);
130  }
131 
143  template<class Scalar, class LO, class GO, class Node>
146  getLocalBlock (const LO localRowIndex) const
147  {
148  if (! this->isValidLocalMeshIndex (localRowIndex)) {
149  return little_vec_type (NULL, 0, 0);
150  } else {
151  const LO strideX = this->getStrideX ();
152  const size_t blockSize = this->getBlockSize ();
153  const size_t offset = localRowIndex * blockSize * strideX;
154  return little_vec_type (this->getRawPtr () + offset, blockSize, strideX);
155  }
156  }
157 
158 } // namespace Experimental
159 } // namespace Tpetra
160 
161 //
162 // Explicit instantiation macro
163 //
164 // Must be expanded from within the Tpetra namespace!
165 //
166 #define TPETRA_EXPERIMENTAL_BLOCKVECTOR_INSTANT(S,LO,GO,NODE) \
167  template class Experimental::BlockVector< S, LO, GO, NODE >;
168 
169 #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. ...
Namespace Tpetra contains the class and methods constituting the Tpetra library.
mv_type mv_
The Tpetra::MultiVector used to represent the data.
MultiVector for multiple degrees of freedom per mesh point.
LO getBlockSize() const
Get the number of degrees of freedom per mesh point.
bool sumIntoGlobalValues(const GO globalRowIndex, const Scalar vals[]) const
Sum into all values at the given mesh point, using a global index.
impl_scalar_type * getRawPtr() const
Raw pointer to the MultiVector&#39;s data.
LittleVector< impl_scalar_type, LO > little_vec_type
"Block view" of all degrees of freedom at a mesh point.
Nonowning view of a set of degrees of freedom corresponding to a mesh point in a block vector or mult...
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.
bool isValidLocalMeshIndex(const LO meshLocalIndex) const
True if and only if meshLocalIndex is a valid local index in the mesh Map.
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.
size_t getNumVectors() const
Number of columns in the multivector.
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node, classic > > getVectorNonConst(const size_t j)
Return a Vector which is a nonconst view of column j.
size_t getStrideX() const
Stride between consecutive local entries in the same column.
bool replaceLocalValues(const LO localRowIndex, const Scalar vals[]) const
Replace all values at the given mesh point, using a local index.