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