Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_BlockVector_def.hpp
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef TPETRA_BLOCKVECTOR_DEF_HPP
11#define TPETRA_BLOCKVECTOR_DEF_HPP
12
13namespace Tpetra {
14
15 template<class Scalar, class LO, class GO, class Node>
20
21 template<class Scalar, class LO, class GO, class Node>
27
28 template<class Scalar, class LO, class GO, class Node>
33
34 template<class Scalar, class LO, class GO, class Node>
41
42 template<class Scalar, class LO, class GO, class Node>
45 const map_type& meshMap,
46 const LO blockSize) :
48 {
50 X_mv.getNumVectors () != 1, std::invalid_argument,
51 "Tpetra::BlockVector: Input MultiVector has "
52 << X_mv.getNumVectors () << " != 1 columns.");
53 }
54
55 template<class Scalar, class LO, class GO, class Node>
62
63 template<class Scalar, class LO, class GO, class Node>
71
72 template<class Scalar, class LO, class GO, class Node>
79
80 template<class Scalar, class LO, class GO, class Node>
83 Teuchos::RCP<vec_type> vPtr = this->mv_.getVectorNonConst (0);
84 return *vPtr; // shallow copy
85 }
86
87 template<class Scalar, class LO, class GO, class Node>
88 bool
90 replaceLocalValues (const LO localRowIndex, const Scalar vals[]) {
91 return ((base_type*) this)->replaceLocalValues (localRowIndex, 0, vals);
92 }
93
94 template<class Scalar, class LO, class GO, class Node>
95 bool
98 return ((base_type*) this)->replaceGlobalValues (globalRowIndex, 0, vals);
99 }
100
101 template<class Scalar, class LO, class GO, class Node>
102 bool
104 sumIntoLocalValues (const LO localRowIndex, const Scalar vals[]) {
105 return ((base_type*) this)->sumIntoLocalValues (localRowIndex, 0, vals);
106 }
107
108 template<class Scalar, class LO, class GO, class Node>
109 bool
112 return ((base_type*) this)->sumIntoLocalValues (globalRowIndex, 0, vals);
113 }
114
115 template<class Scalar, class LO, class GO, class Node>
116 typename BlockVector<Scalar, LO, GO, Node>::const_little_host_vec_type
118 getLocalBlockHost (const LO localRowIndex, Access::ReadOnlyStruct) const
119 {
120 return ((const base_type*) this)->getLocalBlockHost(localRowIndex, 0,
121 Access::ReadOnly);
122 }
123
124 template<class Scalar, class LO, class GO, class Node>
125 typename BlockVector<Scalar, LO, GO, Node>::little_host_vec_type
127 getLocalBlockHost (const LO localRowIndex, Access::ReadWriteStruct)
128 {
129 return ((base_type*) this)->getLocalBlockHost(localRowIndex, 0,
130 Access::ReadWrite);
131 }
132
133 template<class Scalar, class LO, class GO, class Node>
134 typename BlockVector<Scalar, LO, GO, Node>::little_host_vec_type
136 getLocalBlockHost (const LO localRowIndex, Access::OverwriteAllStruct)
137 {
138 return ((base_type*) this)->getLocalBlockHost(localRowIndex, 0,
139 Access::OverwriteAll);
140 }
141
142} // namespace Tpetra
143
144//
145// Explicit instantiation macro
146//
147// Must be expanded from within the Tpetra namespace!
148//
149#define TPETRA_BLOCKVECTOR_INSTANT(S,LO,GO,NODE) \
150 template class BlockVector< S, LO, GO, NODE >;
151
152#endif // TPETRA_BLOCKVECTOR_DEF_HPP
const_little_host_vec_type getLocalBlockHost(const LO localRowIndex, Access::ReadOnlyStruct) const
Get a view of the degrees of freedom at the given mesh point, using a local index.
bool sumIntoGlobalValues(const GO globalRowIndex, const Scalar vals[])
Sum into all values at the given mesh point, using a global index.
bool replaceGlobalValues(const GO globalRowIndex, const Scalar vals[])
Replace all values at the given mesh point, using a global index.
BlockVector()
Default constructor.
bool sumIntoLocalValues(const LO localRowIndex, const Scalar vals[])
Sum into all values at the given mesh point, using a local index.
vec_type getVectorView()
Get a Tpetra::Vector that views this BlockVector's data.
bool replaceLocalValues(const LO localRowIndex, const Scalar vals[])
Replace all values at the given mesh point, using a local index.
Struct that holds views of the contents of a CrsMatrix.
A distributed dense vector.
Namespace Tpetra contains the class and methods constituting the Tpetra library.