xref: /petsc/include/petscdmda_kokkos.hpp (revision 45402d8af0650df8cd128ec935cf42f37f0f84dd)
1a4963045SJacob Faibussowitsch #pragma once
2c6583b63SJunchao Zhang 
311d22bbfSJunchao Zhang #include <petscvec_kokkos.hpp>
4c6583b63SJunchao Zhang #include <petscdmda.h>
5c6583b63SJunchao Zhang 
6ac09b921SBarry Smith /* SUBMANSEC = DMDA */
7ac09b921SBarry Smith 
8c6583b63SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
9c6583b63SJunchao Zhang   #include <Kokkos_Core.hpp>
10c6583b63SJunchao Zhang   #include <Kokkos_OffsetView.hpp>
11c6583b63SJunchao Zhang 
129dc7b89cSJunchao Zhang /*@C
139dc7b89cSJunchao Zhang    DMDAVecGetKokkosOffsetView - Gets a Kokkos OffsetView that contains up-to-date data of a vector in the given memory space.
149dc7b89cSJunchao Zhang 
159dc7b89cSJunchao Zhang    Synopsis:
169dc7b89cSJunchao Zhang    #include <petscdmda_kokkos.hpp>
17bf34f4f8SJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>* kv);
18bf34f4f8SJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv);
19bf34f4f8SJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv);
209dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
219dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
229dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
239dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
249dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
259dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
269dc7b89cSJunchao Zhang 
27cc4c1da9SBarry Smith    Logically Collective, No Fortran Support
289dc7b89cSJunchao Zhang 
299dc7b89cSJunchao Zhang    Input Parameters:
309dc7b89cSJunchao Zhang +  da - the distributed array
3187497f52SBarry Smith -  v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
329dc7b89cSJunchao Zhang 
339dc7b89cSJunchao Zhang    Output Parameter:
349dc7b89cSJunchao Zhang .  kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
359dc7b89cSJunchao Zhang 
369dc7b89cSJunchao Zhang    Notes:
3787497f52SBarry Smith     Call `DMDAVecRestoreKokkosOffsetView()` or `DMDAVecRestoreKokkosOffsetViewWrite()` once you have finished accessing the OffsetView.
389dc7b89cSJunchao Zhang 
3987497f52SBarry Smith     If the vector is not of type `VECKOKKOS`, an error will be raised.
409dc7b89cSJunchao Zhang 
4187497f52SBarry Smith     If the vector is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
429dc7b89cSJunchao Zhang     a global vector then the ghost points are not accessible. Of course with the local vector you will have to do the
4387497f52SBarry Smith     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
449dc7b89cSJunchao Zhang 
4587497f52SBarry Smith     These routines are similar to `DMDAVecGetArray()` and friends. One can read-only, write-only or read/write access the returned
469dc7b89cSJunchao Zhang     Kokkos OffsetView.  Note that passing in a constant OffsetView enables read-only access.
47*45402d8aSJunchao Zhang     Currently, only two memory spaces are supported: HostMirrorMemorySpace and Kokkos::DefaultExecutionSpace::memory_space.
489dc7b89cSJunchao Zhang     If needed, a memory copy will be internally called to copy the latest vector data to the specified memory space.
499dc7b89cSJunchao Zhang 
5087497f52SBarry Smith     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]),
5187497f52SBarry Smith     where i, j, k are loop variables for the x, y, z dimensions respectively specified in `DMDACreate3d()`, for example.
529dc7b89cSJunchao Zhang 
5387497f52SBarry Smith     To give users the same experience as `DMDAVecGetArray()`, we mandate the returned OffsetView always has Kokkos::LayoutRight (that is, rightest
549dc7b89cSJunchao Zhang     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
559dc7b89cSJunchao Zhang     to use Iterate::Right as IterateInner if one uses Kokkos::MDRangePolicy to access the OffsetView.
569dc7b89cSJunchao Zhang 
5787497f52SBarry Smith     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
5887497f52SBarry Smith     (rightest) corresponds to DMDA's x direction.
599dc7b89cSJunchao Zhang 
609dc7b89cSJunchao Zhang     If the vector is a global vector, we have
619dc7b89cSJunchao Zhang .vb
629dc7b89cSJunchao Zhang       kv.extent(0) = zm*dof,  kv.begin(0) = zs*dof, kv.end(0) = (zs+zm)*dof
639dc7b89cSJunchao Zhang       kv.extent(1) = ym*dof,  kv.begin(1) = ys*dof, kv.end(1) = (ys+ym)*dof
649dc7b89cSJunchao Zhang       kv.extent(2) = xm*dof,  kv.begin(2) = xs*dof, kv.end(2) = (xs+xm)*dof
659dc7b89cSJunchao Zhang .ve
669dc7b89cSJunchao Zhang     If the vector is a local vector, we have
679dc7b89cSJunchao Zhang .vb
689dc7b89cSJunchao Zhang       kv.extent(0) = gzm*dof,  kv.begin(0) = gzs*dof, kv.end(0) = (gzs+gzm)*dof
699dc7b89cSJunchao Zhang       kv.extent(1) = gym*dof,  kv.begin(1) = gys*dof, kv.end(1) = (gys+gym)*dof
709dc7b89cSJunchao Zhang       kv.extent(2) = gxm*dof,  kv.begin(2) = gxs*dof, kv.end(2) = (gxs+gxm)*dof
719dc7b89cSJunchao Zhang .ve
729dc7b89cSJunchao Zhang 
739dc7b89cSJunchao Zhang     The starts and widths above are obtained by
749dc7b89cSJunchao Zhang .vb
759dc7b89cSJunchao Zhang      DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
769dc7b89cSJunchao Zhang      DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
779dc7b89cSJunchao Zhang .ve
789dc7b89cSJunchao Zhang 
799dc7b89cSJunchao Zhang     For example, to initialize a grid,
809dc7b89cSJunchao Zhang 
819dc7b89cSJunchao Zhang .vb
829dc7b89cSJunchao Zhang     typedef Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace> PetscScalarKokkosOffsetView3D;
839dc7b89cSJunchao Zhang 
849dc7b89cSJunchao Zhang     PetscScalarKokkosOffsetView3D kv;
859dc7b89cSJunchao Zhang     DMDAVecGetKokkosOffsetViewWrite(da,v,&kv); // v is a global vector and we assume dof is 1
869dc7b89cSJunchao Zhang     DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
879dc7b89cSJunchao Zhang 
889dc7b89cSJunchao Zhang     parallel_for ("Label",MDRangePolicy <Rank<3, Iterate::Right, Iterate::Right>>(
899dc7b89cSJunchao Zhang       {zs,ys,xs},{zs+zm,ys+ym,xs+xm}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i) {
909dc7b89cSJunchao Zhang       kv(k,j,i) = ...;
919dc7b89cSJunchao Zhang     });
929dc7b89cSJunchao Zhang     DMDAVecRestoreKokkosOffsetViewWrite(da,v,&kv);
939dc7b89cSJunchao Zhang .ve
949dc7b89cSJunchao Zhang 
959dc7b89cSJunchao Zhang     For a multi-component problem, one could cast the returned OffsetView to a user's type. But one has also to shrink
969dc7b89cSJunchao Zhang     the OffsetView's extent accordingly. For example,
979dc7b89cSJunchao Zhang .vb
989dc7b89cSJunchao Zhang     typedef struct {
999dc7b89cSJunchao Zhang       PetscScalar omega,temperature;
1009dc7b89cSJunchao Zhang     } Node;
1019dc7b89cSJunchao Zhang 
1029dc7b89cSJunchao Zhang     using NodeKokkosOffsetView3D = Kokkos::Experimental::OffsetView<const Node***,Kokkos::LayoutRight,MemorySpace>;
1039dc7b89cSJunchao Zhang     DMDAVecGetKokkosOffsetViewWrite(da,v,&tv);
1049dc7b89cSJunchao Zhang     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});
1059dc7b89cSJunchao Zhang 
1069dc7b89cSJunchao Zhang     parallel_for ("Label",MDRangePolicy<Rank<3, Iterate::Right, Iterate::Right>>(
1079dc7b89cSJunchao Zhang       {zs,ys,xs},{zs+zm,ys+ym,xs+xm}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i) {
1089dc7b89cSJunchao Zhang       kv(k,j,i).omega = ...;
1099dc7b89cSJunchao Zhang     });
11087497f52SBarry Smith     DMDAVecRestoreKokkosOffsetViewWrite(da,v,&tv)`;
1119dc7b89cSJunchao Zhang .ve
1129dc7b89cSJunchao Zhang 
1139dc7b89cSJunchao Zhang   Level: intermediate
1149dc7b89cSJunchao Zhang 
115db781477SPatrick Sanan .seealso: `DMDAVecRestoreKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
116db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
117db781477SPatrick Sanan           `DMStagVecGetArray()`
1189dc7b89cSJunchao Zhang @*/
1199371c9d4SSatish Balay template <class MemorySpace>
1209371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar *, MemorySpace> *);
1219371c9d4SSatish Balay template <class MemorySpace>
1229371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *);
1239371c9d4SSatish Balay template <class MemorySpace>
1249371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *);
1259dc7b89cSJunchao Zhang 
1269371c9d4SSatish Balay template <class MemorySpace>
1279371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
1289371c9d4SSatish Balay template <class MemorySpace>
1299371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
1309371c9d4SSatish Balay template <class MemorySpace>
1319371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
1329dc7b89cSJunchao Zhang 
1339371c9d4SSatish Balay template <class MemorySpace>
1349371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
1359371c9d4SSatish Balay template <class MemorySpace>
1369371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
1379371c9d4SSatish Balay template <class MemorySpace>
1389371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
1399dc7b89cSJunchao Zhang 
1409dc7b89cSJunchao Zhang /*@C
14187497f52SBarry Smith    DMDAVecRestoreKokkosOffsetView - Returns the Kokkos OffsetView that was gotten with `DMDAVecGetKokkosOffsetView()`
1429dc7b89cSJunchao Zhang 
1439dc7b89cSJunchao Zhang    Synopsis:
1449dc7b89cSJunchao Zhang    #include <petscdmda_kokkos.hpp>
145bf34f4f8SJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>* kv);
146bf34f4f8SJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv);
147bf34f4f8SJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv);
1489dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
1499dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
1509dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
1519dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
1529dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
1539dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
1549dc7b89cSJunchao Zhang 
155cc4c1da9SBarry Smith    Logically Collective, No Fortran Support
1569dc7b89cSJunchao Zhang 
1579dc7b89cSJunchao Zhang    Input Parameters:
1589dc7b89cSJunchao Zhang +  da - the distributed array
15987497f52SBarry Smith .  v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
1609dc7b89cSJunchao Zhang -  kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
1619dc7b89cSJunchao Zhang 
16220f4b53cSBarry Smith   Level: intermediate
16320f4b53cSBarry Smith 
16487497f52SBarry Smith    Note:
16587497f52SBarry Smith     If the vector is not of type `VECKOKKOS`, an error will be raised.
1669dc7b89cSJunchao Zhang 
167db781477SPatrick Sanan .seealso: `DMDAVecGetKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
168db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
169db781477SPatrick Sanan           `DMStagVecGetArray()`
1709dc7b89cSJunchao Zhang @*/
1719371c9d4SSatish Balay template <class MemorySpace>
1729371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar *, MemorySpace> *);
1739371c9d4SSatish Balay template <class MemorySpace>
1749371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *);
1759371c9d4SSatish Balay template <class MemorySpace>
1769371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *);
177c6583b63SJunchao Zhang 
1789371c9d4SSatish Balay template <class MemorySpace>
1799371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
1809371c9d4SSatish Balay template <class MemorySpace>
1819371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
1829371c9d4SSatish Balay template <class MemorySpace>
1839371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
184c6583b63SJunchao Zhang 
1859371c9d4SSatish Balay template <class MemorySpace>
1869371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
1879371c9d4SSatish Balay template <class MemorySpace>
1889371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
1899371c9d4SSatish Balay template <class MemorySpace>
1909371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
191c6583b63SJunchao Zhang 
1929dc7b89cSJunchao Zhang /*@C
1939dc7b89cSJunchao Zhang    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
1949dc7b89cSJunchao Zhang 
1959dc7b89cSJunchao Zhang    Synopsis:
1969dc7b89cSJunchao Zhang    #include <petscdmda_kokkos.hpp>
197bf34f4f8SJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
1989dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
1999dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
2009dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
2019dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
2029dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
2039dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
2049dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
2059dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
2069dc7b89cSJunchao Zhang 
207cc4c1da9SBarry Smith    Logically Collective, No Fortran Support
2089dc7b89cSJunchao Zhang 
2099dc7b89cSJunchao Zhang    Input Parameters:
2109dc7b89cSJunchao Zhang +  da - the distributed array
21187497f52SBarry Smith -  v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
2129dc7b89cSJunchao Zhang 
2139dc7b89cSJunchao Zhang    Output Parameter:
2149dc7b89cSJunchao Zhang .  kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
2159dc7b89cSJunchao Zhang 
2169dc7b89cSJunchao Zhang    Notes:
21787497f52SBarry Smith     Call `DMDAVecRestoreKokkosOffsetViewDOF()` or `DMDAVecRestoreKokkosOffsetViewDOFWrite()` once you have finished accessing the OffsetView.
2189dc7b89cSJunchao Zhang 
21987497f52SBarry Smith     If the vector is not a `VECKOKKOS` an error will be raised.
2209dc7b89cSJunchao Zhang 
22187497f52SBarry Smith     If the vector is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
2229dc7b89cSJunchao Zhang     a global vector then the ghost points are not accessible. Of course with the local vector you will have to do the
22387497f52SBarry Smith     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
2249dc7b89cSJunchao Zhang 
22587497f52SBarry Smith     These routines are similar to `DMDAVecGetArrayDOF()` and friends. One can read-only, write-only or read/write access the returned
2269dc7b89cSJunchao Zhang     Kokkos OffsetView.  Note that passing in a constant OffsetView enables read-only access.
227*45402d8aSJunchao Zhang     Currently, only two memory spaces are supported: HostMirrorMemorySpace and Kokkos::DefaultExecutionSpace::memory_space.
2289dc7b89cSJunchao Zhang     If needed, a memory copy will be internally called to copy the latest vector data to the given memory space.
2299dc7b89cSJunchao Zhang 
23087497f52SBarry Smith     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]),
23187497f52SBarry Smith     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()`,
2329dc7b89cSJunchao Zhang     for example.
2339dc7b89cSJunchao Zhang 
23487497f52SBarry Smith     To give users the same experience as `DMDAVecGetArrayDOF()`, we mandate the returned OffsetView always has Kokkos::LayoutRight (that is, rightest
2359dc7b89cSJunchao Zhang     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
2369dc7b89cSJunchao Zhang     to use Iterate::Right as IterateInner if one uses Kokkos::MDRangePolicy to access the OffsetView.
2379dc7b89cSJunchao Zhang 
23887497f52SBarry Smith     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
23987497f52SBarry Smith     (rightest) corresponds to DMDA's x direction.
2409dc7b89cSJunchao Zhang 
2419dc7b89cSJunchao Zhang     If the vector is a global vector, we have
2429dc7b89cSJunchao Zhang .vb
2439dc7b89cSJunchao Zhang       kv.extent(0) = zm,  kv.begin(0) = zs, kv.end(0) = zs+zm
2449dc7b89cSJunchao Zhang       kv.extent(1) = ym,  kv.begin(1) = ys, kv.end(1) = ys+ym
2459dc7b89cSJunchao Zhang       kv.extent(2) = xm,  kv.begin(2) = xs, kv.end(2) = xs+xm
2469dc7b89cSJunchao Zhang       kv.extent(3) = dof, kv.begin(3) = 0,  kv.end(3) = dof
2479dc7b89cSJunchao Zhang .ve
2489dc7b89cSJunchao Zhang     If the vector is a local vector, we have
2499dc7b89cSJunchao Zhang .vb
2509dc7b89cSJunchao Zhang       kv.extent(0) = gzm, kv.begin(0) = gzs, kv.end(0) = gzs+gzm
2519dc7b89cSJunchao Zhang       kv.extent(1) = gym, kv.begin(1) = gys, kv.end(1) = gys+gym
2529dc7b89cSJunchao Zhang       kv.extent(2) = gxm, kv.begin(2) = gxs, kv.end(2) = gxs+gxm
2539dc7b89cSJunchao Zhang       kv.extent(3) = dof, kv.begin(3) = 0,   kv.end(3) = dof
2549dc7b89cSJunchao Zhang .ve
2559dc7b89cSJunchao Zhang 
2569dc7b89cSJunchao Zhang     The starts and widths above are obtained by
2579dc7b89cSJunchao Zhang .vb
2589dc7b89cSJunchao Zhang      DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
2599dc7b89cSJunchao Zhang      DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
2609dc7b89cSJunchao Zhang .ve
2619dc7b89cSJunchao Zhang 
2629dc7b89cSJunchao Zhang     For example, to initialize a grid,
2639dc7b89cSJunchao Zhang .vb
2649dc7b89cSJunchao Zhang     typedef Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace> PetscScalarKokkosOffsetView4D;
2659dc7b89cSJunchao Zhang 
2669dc7b89cSJunchao Zhang     PetscScalarKokkosOffsetView4D kv;
2679dc7b89cSJunchao Zhang     DMDAVecGetKokkosOffsetViewDOFWrite(da,v,&kv); // v is a global vector
2689dc7b89cSJunchao Zhang     DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
2699dc7b89cSJunchao Zhang 
2709dc7b89cSJunchao Zhang     parallel_for ("Label",MDRangePolicy <Rank<4, Iterate::Right, Iterate::Right>>(
2719dc7b89cSJunchao Zhang       {zs,ys,xs,0},{zs+zm,ys+ym,xs+xm,dof}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i,PetscInt c) {
2729dc7b89cSJunchao Zhang       kv(k,j,i,c) = ...;
2739dc7b89cSJunchao Zhang     });
2749dc7b89cSJunchao Zhang     DMDAVecRestoreKokkosOffsetViewDOFWrite(da,v,&kv);
2759dc7b89cSJunchao Zhang .ve
2769dc7b89cSJunchao Zhang 
2779dc7b89cSJunchao Zhang   Level: intermediate
2789dc7b89cSJunchao Zhang 
279db781477SPatrick Sanan .seealso: `DMDAVecRestoreKokkosOffsetViewDOF()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
280db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
281db781477SPatrick Sanan           `DMStagVecGetArray()`
2829dc7b89cSJunchao Zhang @*/
2839371c9d4SSatish Balay template <class MemorySpace>
2849371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
2859371c9d4SSatish Balay template <class MemorySpace>
2869371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
2879371c9d4SSatish Balay template <class MemorySpace>
2889371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
2899dc7b89cSJunchao Zhang 
2909371c9d4SSatish Balay template <class MemorySpace>
2919371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
2929371c9d4SSatish Balay template <class MemorySpace>
2939371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
2949371c9d4SSatish Balay template <class MemorySpace>
2959371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
2969dc7b89cSJunchao Zhang 
2979371c9d4SSatish Balay template <class MemorySpace>
2989371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
2999371c9d4SSatish Balay template <class MemorySpace>
3009371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
3019371c9d4SSatish Balay template <class MemorySpace>
3029371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
3039dc7b89cSJunchao Zhang 
3049dc7b89cSJunchao Zhang /*@C
30587497f52SBarry Smith    DMDAVecRestoreKokkosOffsetViewDOF - Returns the Kokkos OffsetView that was gotten from `DMDAVecGetKokkosOffsetViewDOF()`
3069dc7b89cSJunchao Zhang 
3079dc7b89cSJunchao Zhang    Synopsis:
3089dc7b89cSJunchao Zhang    #include <petscdmda_kokkos.hpp>
3099dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
3109dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
3119dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
3129dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
3139dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
3149dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
3159dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
3169dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
3179dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
3189dc7b89cSJunchao Zhang 
319cc4c1da9SBarry Smith    Logically Collective, No Fortran Support
3209dc7b89cSJunchao Zhang 
3219dc7b89cSJunchao Zhang    Input Parameters:
3229dc7b89cSJunchao Zhang +  da - the distributed array
32387497f52SBarry Smith .  v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
3249dc7b89cSJunchao Zhang -  kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
3259dc7b89cSJunchao Zhang 
32620f4b53cSBarry Smith   Level: intermediate
32720f4b53cSBarry Smith 
32887497f52SBarry Smith    Note:
32987497f52SBarry Smith     If the vector is not of type `VECKOKKOS`, an error will be raised.
3309dc7b89cSJunchao Zhang 
331db781477SPatrick Sanan .seealso: `DMDAVecGetKokkosOffsetViewDOF()`, `DMDAVecGetKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
332db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
333db781477SPatrick Sanan           `DMStagVecGetArray()`
3349dc7b89cSJunchao Zhang @*/
3359371c9d4SSatish Balay template <class MemorySpace>
3369371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
3379371c9d4SSatish Balay template <class MemorySpace>
3389371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
3399371c9d4SSatish Balay template <class MemorySpace>
3409371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *);
3419dc7b89cSJunchao Zhang 
3429371c9d4SSatish Balay template <class MemorySpace>
3439371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
3449371c9d4SSatish Balay template <class MemorySpace>
3459371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
3469371c9d4SSatish Balay template <class MemorySpace>
3479371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *);
3489dc7b89cSJunchao Zhang 
3499371c9d4SSatish Balay template <class MemorySpace>
3509371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
3519371c9d4SSatish Balay template <class MemorySpace>
3529371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
3539371c9d4SSatish Balay template <class MemorySpace>
3549371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *);
355c6583b63SJunchao Zhang #endif
356