forked from petsc/petsc
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathpetscdmda_kokkos.hpp
355 lines (288 loc) · 22.3 KB
/
petscdmda_kokkos.hpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
#pragma once
#include <petscvec_kokkos.hpp>
#include <petscdmda.h>
/* SUBMANSEC = DMDA */
#if defined(PETSC_HAVE_KOKKOS)
#include <Kokkos_Core.hpp>
#include <Kokkos_OffsetView.hpp>
/*@C
DMDAVecGetKokkosOffsetView - Gets a Kokkos OffsetView that contains up-to-date data of a vector in the given memory space.
Synopsis:
#include <petscdmda_kokkos.hpp>
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
Logically Collective
Input Parameters:
+ da - the distributed array
- v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
Output Parameter:
. kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
Notes:
Call `DMDAVecRestoreKokkosOffsetView()` or `DMDAVecRestoreKokkosOffsetViewWrite()` once you have finished accessing the OffsetView.
If the vector is not of type `VECKOKKOS`, an error will be raised.
If the vector is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
a global vector then the ghost points are not accessible. Of course with the local vector you will have to do the
appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
These routines are similar to `DMDAVecGetArray()` and friends. One can read-only, write-only or read/write access the returned
Kokkos OffsetView. Note that passing in a constant OffsetView enables read-only access.
Currently, only two memory spaces are supported: Kokkos::HostSpace and Kokkos::DefaultExecutionSpace::memory_space.
If needed, a memory copy will be internally called to copy the latest vector data to the specified memory space.
In C, to access the returned array of `DMDAVecGetArray()`, the indexing is "backwards", i.e., array[k][j][i] (instead of array[i][j][k]),
where i, j, k are loop variables for the x, y, z dimensions respectively specified in `DMDACreate3d()`, for example.
To give users the same experience as `DMDAVecGetArray()`, we mandate the returned OffsetView always has Kokkos::LayoutRight (that is, rightest
subscript has a stride 1, as in C multi-dimensional arrays), regardless of whether the memory space is host or device. Thus it is important
to use Iterate::Right as IterateInner if one uses Kokkos::MDRangePolicy to access the OffsetView.
Note that the OffsetView kv's first dimension (i.e., the leftest, dim 0) corresponds to the DMDA's z direction, and its last dimension
(rightest) corresponds to DMDA's x direction.
If the vector is a global vector, we have
.vb
kv.extent(0) = zm*dof, kv.begin(0) = zs*dof, kv.end(0) = (zs+zm)*dof
kv.extent(1) = ym*dof, kv.begin(1) = ys*dof, kv.end(1) = (ys+ym)*dof
kv.extent(2) = xm*dof, kv.begin(2) = xs*dof, kv.end(2) = (xs+xm)*dof
.ve
If the vector is a local vector, we have
.vb
kv.extent(0) = gzm*dof, kv.begin(0) = gzs*dof, kv.end(0) = (gzs+gzm)*dof
kv.extent(1) = gym*dof, kv.begin(1) = gys*dof, kv.end(1) = (gys+gym)*dof
kv.extent(2) = gxm*dof, kv.begin(2) = gxs*dof, kv.end(2) = (gxs+gxm)*dof
.ve
The starts and widths above are obtained by
.vb
DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
.ve
For example, to initialize a grid,
.vb
typedef Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace> PetscScalarKokkosOffsetView3D;
PetscScalarKokkosOffsetView3D kv;
DMDAVecGetKokkosOffsetViewWrite(da,v,&kv); // v is a global vector and we assume dof is 1
DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
parallel_for ("Label",MDRangePolicy <Rank<3, Iterate::Right, Iterate::Right>>(
{zs,ys,xs},{zs+zm,ys+ym,xs+xm}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i) {
kv(k,j,i) = ...;
});
DMDAVecRestoreKokkosOffsetViewWrite(da,v,&kv);
.ve
For a multi-component problem, one could cast the returned OffsetView to a user's type. But one has also to shrink
the OffsetView's extent accordingly. For example,
.vb
typedef struct {
PetscScalar omega,temperature;
} Node;
using NodeKokkosOffsetView3D = Kokkos::Experimental::OffsetView<const Node***,Kokkos::LayoutRight,MemorySpace>;
DMDAVecGetKokkosOffsetViewWrite(da,v,&tv);
NodeKokkosOffsetView3D kv(reinterpret_cast<Node*>(tv.data()),{tv.begin(0)/dof,tv.begin(1)/dof,tv.begin(2)/dof}, {tv.end(0)/dof,tv.end(1)/dof,tv.end(2)/dof});
parallel_for ("Label",MDRangePolicy<Rank<3, Iterate::Right, Iterate::Right>>(
{zs,ys,xs},{zs+zm,ys+ym,xs+xm}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i) {
kv(k,j,i).omega = ...;
});
DMDAVecRestoreKokkosOffsetViewWrite(da,v,&tv)`;
.ve
Level: intermediate
.seealso: `DMDAVecRestoreKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
`DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
`DMStagVecGetArray()`
@*/
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar *, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
/*@C
DMDAVecRestoreKokkosOffsetView - Returns the Kokkos OffsetView that was gotten with `DMDAVecGetKokkosOffsetView()`
Synopsis:
#include <petscdmda_kokkos.hpp>
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
Logically Collective
Input Parameters:
+ da - the distributed array
. v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
- kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
Level: intermediate
Note:
If the vector is not of type `VECKOKKOS`, an error will be raised.
.seealso: `DMDAVecGetKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
`DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
`DMStagVecGetArray()`
@*/
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar *, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
/*@C
DMDAVecGetKokkosOffsetViewDOF - Gets a Kokkos OffsetView that contains up-to-date data of a vector in the given memory space, with DOF as the rightest dimension of the OffsetView
Synopsis:
#include <petscdmda_kokkos.hpp>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
Logically Collective
Input Parameters:
+ da - the distributed array
- v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
Output Parameter:
. kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
Notes:
Call `DMDAVecRestoreKokkosOffsetViewDOF()` or `DMDAVecRestoreKokkosOffsetViewDOFWrite()` once you have finished accessing the OffsetView.
If the vector is not a `VECKOKKOS` an error will be raised.
If the vector is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
a global vector then the ghost points are not accessible. Of course with the local vector you will have to do the
appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
These routines are similar to `DMDAVecGetArrayDOF()` and friends. One can read-only, write-only or read/write access the returned
Kokkos OffsetView. Note that passing in a constant OffsetView enables read-only access.
Currently, only two memory spaces are supported: Kokkos::HostSpace and Kokkos::DefaultExecutionSpace::memory_space.
If needed, a memory copy will be internally called to copy the latest vector data to the given memory space.
In C, to access the returned array of `DMDAVecGetArrayDOF()`, the indexing is "backwards", i.e., array[k][j][i][c] (instead of array[c][i][j][k]),
where i, j, k are loop variables for the x, y, z dimensions respectively, and c is the loop variable for DOFs, as specified in `DMDACreate3d()`,
for example.
To give users the same experience as `DMDAVecGetArrayDOF()`, we mandate the returned OffsetView always has Kokkos::LayoutRight (that is, rightest
subscript has a stride 1, as in C multi-dimensional arrays), regardless of whether the memory space is host or device. Thus it is important
to use Iterate::Right as IterateInner if one uses Kokkos::MDRangePolicy to access the OffsetView.
Note that for a 3D `DMDA`, the OffsetView kv's first dimension (i.e., the leftest, dim 0) corresponds to DMDA's z direction, and its second-to-last dimension
(rightest) corresponds to DMDA's x direction.
If the vector is a global vector, we have
.vb
kv.extent(0) = zm, kv.begin(0) = zs, kv.end(0) = zs+zm
kv.extent(1) = ym, kv.begin(1) = ys, kv.end(1) = ys+ym
kv.extent(2) = xm, kv.begin(2) = xs, kv.end(2) = xs+xm
kv.extent(3) = dof, kv.begin(3) = 0, kv.end(3) = dof
.ve
If the vector is a local vector, we have
.vb
kv.extent(0) = gzm, kv.begin(0) = gzs, kv.end(0) = gzs+gzm
kv.extent(1) = gym, kv.begin(1) = gys, kv.end(1) = gys+gym
kv.extent(2) = gxm, kv.begin(2) = gxs, kv.end(2) = gxs+gxm
kv.extent(3) = dof, kv.begin(3) = 0, kv.end(3) = dof
.ve
The starts and widths above are obtained by
.vb
DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
.ve
For example, to initialize a grid,
.vb
typedef Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace> PetscScalarKokkosOffsetView4D;
PetscScalarKokkosOffsetView4D kv;
DMDAVecGetKokkosOffsetViewDOFWrite(da,v,&kv); // v is a global vector
DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
parallel_for ("Label",MDRangePolicy <Rank<4, Iterate::Right, Iterate::Right>>(
{zs,ys,xs,0},{zs+zm,ys+ym,xs+xm,dof}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i,PetscInt c) {
kv(k,j,i,c) = ...;
});
DMDAVecRestoreKokkosOffsetViewDOFWrite(da,v,&kv);
.ve
Level: intermediate
.seealso: `DMDAVecRestoreKokkosOffsetViewDOF()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
`DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
`DMStagVecGetArray()`
@*/
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
/*@C
DMDAVecRestoreKokkosOffsetViewDOF - Returns the Kokkos OffsetView that was gotten from `DMDAVecGetKokkosOffsetViewDOF()`
Synopsis:
#include <petscdmda_kokkos.hpp>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
Logically Collective
Input Parameters:
+ da - the distributed array
. v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
- kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
Level: intermediate
Note:
If the vector is not of type `VECKOKKOS`, an error will be raised.
.seealso: `DMDAVecGetKokkosOffsetViewDOF()`, `DMDAVecGetKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
`DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
`DMStagVecGetArray()`
@*/
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
template <class MemorySpace>
PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
#endif