10 #ifndef IFPACK2_DETAILS_MULTIVECTORLOCALGATHERSCATTER_HPP
11 #define IFPACK2_DETAILS_MULTIVECTORLOCALGATHERSCATTER_HPP
17 #include "Tpetra_MultiVector.hpp"
18 #include "Tpetra_Map.hpp"
51 template <
class MV_in,
class MV_out>
54 typedef typename MV_in::scalar_type InScalar;
55 typedef typename MV_out::scalar_type OutScalar;
56 typedef typename MV_in::local_ordinal_type LO;
57 typedef typename MV_in::global_ordinal_type GO;
58 typedef typename MV_in::node_type NO;
68 const size_t numRows = X_out.getLocalLength();
69 const size_t numVecs = X_in.getNumVectors();
71 for (
size_t j = 0; j < numVecs; ++j) {
72 ArrayRCP<const InScalar> X_in_j = X_in.getData(j);
73 ArrayRCP<OutScalar> X_out_j = X_out.getDataNonConst(j);
74 for (
size_t i = 0; i < numRows; ++i) {
75 const size_t i_perm = perm[i];
76 X_out_j[i] = X_in_j[i_perm];
92 const size_t numBlocks = X_out.getLocalLength() / blockSize;
93 const size_t numVecs = X_in.getNumVectors();
95 for (
size_t j = 0; j < numVecs; ++j) {
96 ArrayRCP<const InScalar> X_in_j = X_in.getData(j);
97 ArrayRCP<OutScalar> X_out_j = X_out.getDataNonConst(j);
98 for (
size_t i = 0; i < numBlocks; ++i) {
99 const size_t i_perm = perm[i];
100 for (LO k = 0; k < blockSize; k++) {
101 X_out_j[i * blockSize + k] = X_in_j[i_perm * blockSize + k];
113 const size_t numRows = X_out.getLocalLength();
114 const size_t numVecs = X_in.getNumVectors();
116 for (
size_t j = 0; j < numVecs; ++j) {
117 ArrayRCP<InScalar> X_in_j = X_in.getDataNonConst(j);
118 ArrayRCP<const OutScalar> X_out_j = X_out.getData(j);
119 for (
size_t i = 0; i < numRows; ++i) {
120 const size_t i_perm = perm[i];
121 X_in_j[i_perm] = X_out_j[i];
132 LO blockSize)
const {
134 const size_t numBlocks = X_out.getLocalLength() / blockSize;
135 const size_t numVecs = X_in.getNumVectors();
137 for (
size_t j = 0; j < numVecs; ++j) {
138 ArrayRCP<const InScalar> X_in_j = X_in.getData(j);
139 ArrayRCP<OutScalar> X_out_j = X_out.getDataNonConst(j);
140 for (
size_t i = 0; i < numBlocks; ++i) {
141 const size_t i_perm = perm[i];
142 for (LO k = 0; k < blockSize; k++) {
143 X_in_j[i_perm * blockSize + k] = X_out_j[i * blockSize + k];
153 template <
typename InView,
typename OutView>
154 void gatherViewToView(OutView X_out,
159 for (
size_t j = 0; j < X_out.extent(1); ++j) {
160 for (
size_t i = 0; i < X_out.extent(0); ++i) {
161 const LO i_perm = perm[i];
162 X_out(i, j) = X_in(i_perm, j);
167 template <
typename InView,
typename OutView>
168 void scatterViewToView(InView X_in,
172 for (
size_t j = 0; j < X_out.extent(1); ++j) {
173 for (
size_t i = 0; i < X_out.extent(0); ++i) {
174 const LO i_perm = perm[i];
175 X_in(i_perm, j) = X_out(i, j);
180 template <
typename InView,
typename OutView>
181 void gatherViewToViewBlock(OutView X_out,
184 LO blockSize)
const {
187 size_t numBlocks = X_out.extent(0) / blockSize;
188 for (
size_t j = 0; j < X_out.extent(1); ++j) {
189 for (
size_t i = 0; i < numBlocks; ++i) {
190 const LO i_perm = perm[i];
191 for (LO k = 0; k < blockSize; k++) {
192 X_out(i * blockSize + k, j) = X_in(i_perm * blockSize + k, j);
198 template <
typename InView,
typename OutView>
199 void scatterViewToViewBlock(InView X_in,
202 LO blockSize)
const {
205 size_t numBlocks = X_out.extent(0) / blockSize;
206 for (
size_t j = 0; j < X_out.extent(1); ++j) {
207 for (
size_t i = 0; i < numBlocks; ++i) {
208 const LO i_perm = perm[i];
209 for (LO k = 0; k < blockSize; k++) {
210 X_in(i_perm * blockSize + k, j) = X_out(i * blockSize + k, j);
219 template <
typename InView>
220 void gatherMVtoView(MV_out X_out,
225 size_t numRows = X_out.getLocalLength();
226 for (
size_t j = 0; j < X_out.getNumVectors(); ++j) {
228 for (
size_t i = 0; i < numRows; ++i) {
229 const LO i_perm = perm[i];
230 X_out_j[i] = X_in(i_perm, j);
235 template <
typename InView>
236 void scatterMVtoView(InView X_in,
239 size_t numRows = X_out.getLocalLength();
241 for (
size_t j = 0; j < X_in.extent(1); ++j) {
243 for (
size_t i = 0; i < numRows; ++i) {
244 const LO i_perm = perm[i];
245 X_in(i_perm, j) = X_out_j[i];
250 template <
typename InView>
251 void gatherMVtoViewBlock(MV_out X_out,
254 LO blockSize)
const {
256 size_t numBlocks = X_out.getLocalLength() / blockSize;
258 for (
size_t j = 0; j < X_out.getNumVectors(); ++j) {
260 for (
size_t i = 0; i < numBlocks; ++i) {
261 const LO i_perm = perm[i];
262 for (LO k = 0; k < blockSize; k++) {
263 X_out_j[i * blockSize + k] = X_in(i_perm * blockSize + k, j);
269 template <
typename InView>
270 void scatterMVtoViewBlock(InView X_in,
273 LO blockSize)
const {
274 size_t numBlocks = X_out.getLocalLength() / blockSize;
276 for (
size_t j = 0; j < X_in.extent(1); ++j) {
278 for (
size_t i = 0; i < numBlocks; ++i) {
279 const LO i_perm = perm[i];
280 for (LO k = 0; k < blockSize; k++) {
281 X_in(i_perm * blockSize + k, j) = X_out_j[i * blockSize + k];
291 #endif // IFPACK2_DETAILS_MULTIVECTORLOCALGATHERSCATTER_HPP
Implementation detail of Ifpack2::Container subclasses.
Definition: Ifpack2_Details_MultiVectorLocalGatherScatter.hpp:52