xref: /petsc/include/petscdmda_kokkos.hpp (revision db7814771ca77b190574494e87b584e981451db0)
1c6583b63SJunchao Zhang #if !defined(PETSCDMDA_KOKKOS_HPP)
2c6583b63SJunchao Zhang #define PETSCDMDA_KOKKOS_HPP
3c6583b63SJunchao Zhang 
411d22bbfSJunchao Zhang #include <petscvec_kokkos.hpp>
5c6583b63SJunchao Zhang #include <petscdmda.h>
6c6583b63SJunchao Zhang 
7c6583b63SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS)
8c6583b63SJunchao Zhang #include <Kokkos_Core.hpp>
9c6583b63SJunchao Zhang #include <Kokkos_OffsetView.hpp>
10c6583b63SJunchao Zhang 
119dc7b89cSJunchao Zhang /*@C
129dc7b89cSJunchao Zhang    DMDAVecGetKokkosOffsetView - Gets a Kokkos OffsetView that contains up-to-date data of a vector in the given memory space.
139dc7b89cSJunchao Zhang 
149dc7b89cSJunchao Zhang    Synopsis:
159dc7b89cSJunchao Zhang    #include <petscdmda_kokkos.hpp>
169dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar*,Kokkos::LayoutLeft,MemorySpace>* kv);
179dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,Kokkos::LayoutRight,MemorySpace>* kv);
189dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,Kokkos::LayoutRight,MemorySpace>* kv);
199dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
209dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
219dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
229dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
239dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
249dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
259dc7b89cSJunchao Zhang 
269dc7b89cSJunchao Zhang    Logically collective on da
279dc7b89cSJunchao Zhang 
289dc7b89cSJunchao Zhang    Input Parameters:
299dc7b89cSJunchao Zhang +  da - the distributed array
309dc7b89cSJunchao Zhang -  v - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
319dc7b89cSJunchao Zhang 
329dc7b89cSJunchao Zhang    Output Parameter:
339dc7b89cSJunchao Zhang .  kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
349dc7b89cSJunchao Zhang 
359dc7b89cSJunchao Zhang    Notes:
369dc7b89cSJunchao Zhang     Call DMDAVecRestoreKokkosOffsetView() or DMDAVecRestoreKokkosOffsetViewWrite() once you have finished accessing the OffsetView.
379dc7b89cSJunchao Zhang 
389dc7b89cSJunchao Zhang     If the vector is not a Kokkos vector, an error will be raised.
399dc7b89cSJunchao Zhang 
409dc7b89cSJunchao Zhang     If the vector is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
419dc7b89cSJunchao Zhang     a global vector then the ghost points are not accessible. Of course with the local vector you will have to do the
429dc7b89cSJunchao Zhang     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
439dc7b89cSJunchao Zhang 
449dc7b89cSJunchao Zhang     These routines are similar to DMDAVecGetArray() and friends. One can read-only, write-only or read/write access the returned
459dc7b89cSJunchao Zhang     Kokkos OffsetView.  Note that passing in a constant OffsetView enables read-only access.
469dc7b89cSJunchao Zhang     Currently, only two memory spaces are supported: Kokkos::HostSpace and Kokkos::DefaultExecutionSpace::memory_space.
479dc7b89cSJunchao Zhang     If needed, a memory copy will be internally called to copy the latest vector data to the specified memory space.
489dc7b89cSJunchao Zhang 
499dc7b89cSJunchao Zhang     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]),
509dc7b89cSJunchao Zhang     where i, j, k are loop variables for the x, y, z dimensions respectively specified in DMDACreate3d(), for example.
519dc7b89cSJunchao Zhang 
529dc7b89cSJunchao Zhang     To give users the same experience as DMDAVecGetArray(), we mandate the returned OffsetView always has Kokkos::LayoutRight (that is, rightest
539dc7b89cSJunchao 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
549dc7b89cSJunchao Zhang     to use Iterate::Right as IterateInner if one uses Kokkos::MDRangePolicy to access the OffsetView.
559dc7b89cSJunchao Zhang 
569dc7b89cSJunchao Zhang     Note that the OffsetView kv's first dimension (i.e., the leftest, dim 0) corresponds to DA's z direction, and its last dimension
579dc7b89cSJunchao Zhang     (rightest) corresponds to DA's x direction.
589dc7b89cSJunchao Zhang 
599dc7b89cSJunchao Zhang     If the vector is a global vector, we have
609dc7b89cSJunchao Zhang .vb
619dc7b89cSJunchao Zhang       kv.extent(0) = zm*dof,  kv.begin(0) = zs*dof, kv.end(0) = (zs+zm)*dof
629dc7b89cSJunchao Zhang       kv.extent(1) = ym*dof,  kv.begin(1) = ys*dof, kv.end(1) = (ys+ym)*dof
639dc7b89cSJunchao Zhang       kv.extent(2) = xm*dof,  kv.begin(2) = xs*dof, kv.end(2) = (xs+xm)*dof
649dc7b89cSJunchao Zhang .ve
659dc7b89cSJunchao Zhang     If the vector is a local vector, we have
669dc7b89cSJunchao Zhang .vb
679dc7b89cSJunchao Zhang       kv.extent(0) = gzm*dof,  kv.begin(0) = gzs*dof, kv.end(0) = (gzs+gzm)*dof
689dc7b89cSJunchao Zhang       kv.extent(1) = gym*dof,  kv.begin(1) = gys*dof, kv.end(1) = (gys+gym)*dof
699dc7b89cSJunchao Zhang       kv.extent(2) = gxm*dof,  kv.begin(2) = gxs*dof, kv.end(2) = (gxs+gxm)*dof
709dc7b89cSJunchao Zhang .ve
719dc7b89cSJunchao Zhang 
729dc7b89cSJunchao Zhang     The starts and widths above are obtained by
739dc7b89cSJunchao Zhang .vb
749dc7b89cSJunchao Zhang      DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
759dc7b89cSJunchao Zhang      DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
769dc7b89cSJunchao Zhang .ve
779dc7b89cSJunchao Zhang 
789dc7b89cSJunchao Zhang     For example, to initialize a grid,
799dc7b89cSJunchao Zhang 
809dc7b89cSJunchao Zhang .vb
819dc7b89cSJunchao Zhang     typedef Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace> PetscScalarKokkosOffsetView3D;
829dc7b89cSJunchao Zhang 
839dc7b89cSJunchao Zhang     PetscScalarKokkosOffsetView3D kv;
849dc7b89cSJunchao Zhang     DMDAVecGetKokkosOffsetViewWrite(da,v,&kv); // v is a global vector and we assume dof is 1
859dc7b89cSJunchao Zhang     DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
869dc7b89cSJunchao Zhang 
879dc7b89cSJunchao Zhang     parallel_for ("Label",MDRangePolicy <Rank<3, Iterate::Right, Iterate::Right>>(
889dc7b89cSJunchao Zhang       {zs,ys,xs},{zs+zm,ys+ym,xs+xm}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i) {
899dc7b89cSJunchao Zhang       kv(k,j,i) = ...;
909dc7b89cSJunchao Zhang     });
919dc7b89cSJunchao Zhang     DMDAVecRestoreKokkosOffsetViewWrite(da,v,&kv);
929dc7b89cSJunchao Zhang .ve
939dc7b89cSJunchao Zhang 
949dc7b89cSJunchao Zhang     For a multi-component problem, one could cast the returned OffsetView to a user's type. But one has also to shrink
959dc7b89cSJunchao Zhang     the OffsetView's extent accordingly. For example,
969dc7b89cSJunchao Zhang .vb
979dc7b89cSJunchao Zhang     typedef struct {
989dc7b89cSJunchao Zhang       PetscScalar omega,temperature;
999dc7b89cSJunchao Zhang     } Node;
1009dc7b89cSJunchao Zhang 
1019dc7b89cSJunchao Zhang     using NodeKokkosOffsetView3D = Kokkos::Experimental::OffsetView<const Node***,Kokkos::LayoutRight,MemorySpace>;
1029dc7b89cSJunchao Zhang     DMDAVecGetKokkosOffsetViewWrite(da,v,&tv);
1039dc7b89cSJunchao 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});
1049dc7b89cSJunchao Zhang 
1059dc7b89cSJunchao Zhang     parallel_for ("Label",MDRangePolicy<Rank<3, Iterate::Right, Iterate::Right>>(
1069dc7b89cSJunchao Zhang       {zs,ys,xs},{zs+zm,ys+ym,xs+xm}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i) {
1079dc7b89cSJunchao Zhang       kv(k,j,i).omega = ...;
1089dc7b89cSJunchao Zhang     });
1099dc7b89cSJunchao Zhang     DMDAVecRestoreKokkosOffsetViewWrite(da,v,&tv);
1109dc7b89cSJunchao Zhang .ve
1119dc7b89cSJunchao Zhang 
1129dc7b89cSJunchao Zhang   Level: intermediate
1139dc7b89cSJunchao Zhang 
114*db781477SPatrick Sanan .seealso: `DMDAVecRestoreKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
115*db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
116*db781477SPatrick Sanan           `DMStagVecGetArray()`
1179dc7b89cSJunchao Zhang @*/
118c6583b63SJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView         (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>*);
119c6583b63SJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView         (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar*,MemorySpace>*);
120c6583b63SJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewWrite    (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar*,MemorySpace>*);
1219dc7b89cSJunchao Zhang 
1229dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView         (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
1239dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView         (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
1249dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewWrite    (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
1259dc7b89cSJunchao Zhang 
1269dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView         (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>*);
1279dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView         (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar***,Kokkos::LayoutRight,MemorySpace>*);
1289dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewWrite    (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar***,Kokkos::LayoutRight,MemorySpace>*);
1299dc7b89cSJunchao Zhang 
1309dc7b89cSJunchao Zhang /*@C
1319dc7b89cSJunchao Zhang    DMDAVecRestoreKokkosOffsetView - Returns the Kokkos OffsetView that gotten from DMDAVecGetKokkosOffsetView()
1329dc7b89cSJunchao Zhang 
1339dc7b89cSJunchao Zhang    Synopsis:
1349dc7b89cSJunchao Zhang    #include <petscdmda_kokkos.hpp>
1359dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar*,Kokkos::LayoutLeft,MemorySpace>* kv);
1369dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,Kokkos::LayoutRight,MemorySpace>* kv);
1379dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,Kokkos::LayoutRight,MemorySpace>* kv);
1389dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
1399dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
1409dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
1419dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
1429dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
1439dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
1449dc7b89cSJunchao Zhang 
1459dc7b89cSJunchao Zhang    Logically collective on da
1469dc7b89cSJunchao Zhang 
1479dc7b89cSJunchao Zhang    Input Parameters:
1489dc7b89cSJunchao Zhang +  da - the distributed array
1499dc7b89cSJunchao Zhang .  v - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
1509dc7b89cSJunchao Zhang -  kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
1519dc7b89cSJunchao Zhang 
1529dc7b89cSJunchao Zhang    Notes:
1539dc7b89cSJunchao Zhang     If the vector is not of type VECKOKKOS, an error will be raised.
1549dc7b89cSJunchao Zhang 
1559dc7b89cSJunchao Zhang   Level: intermediate
1569dc7b89cSJunchao Zhang 
157*db781477SPatrick Sanan .seealso: `DMDAVecGetKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
158*db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
159*db781477SPatrick Sanan           `DMStagVecGetArray()`
1609dc7b89cSJunchao Zhang @*/
1619dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView     (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>*);
1629dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView     (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar*,MemorySpace>*);
163c6583b63SJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar*,MemorySpace>*);
164c6583b63SJunchao Zhang 
1659dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView     (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
1669dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView     (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
1679dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
168c6583b63SJunchao Zhang 
1699dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView     (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>*);
1709dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView     (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar***,Kokkos::LayoutRight,MemorySpace>*);
1719dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar***,Kokkos::LayoutRight,MemorySpace>*);
172c6583b63SJunchao Zhang 
1739dc7b89cSJunchao Zhang /*@C
1749dc7b89cSJunchao 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
1759dc7b89cSJunchao Zhang 
1769dc7b89cSJunchao Zhang    Synopsis:
1779dc7b89cSJunchao Zhang    #include <petscdmda_kokkos.hpp>
1789dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutLeft,MemorySpace>* kv);
1799dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
1809dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
1819dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
1829dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
1839dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
1849dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
1859dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
1869dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
1879dc7b89cSJunchao Zhang 
1889dc7b89cSJunchao Zhang    Logically collective on da
1899dc7b89cSJunchao Zhang 
1909dc7b89cSJunchao Zhang    Input Parameters:
1919dc7b89cSJunchao Zhang +  da - the distributed array
1929dc7b89cSJunchao Zhang -  v - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
1939dc7b89cSJunchao Zhang 
1949dc7b89cSJunchao Zhang    Output Parameter:
1959dc7b89cSJunchao Zhang .  kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
1969dc7b89cSJunchao Zhang 
1979dc7b89cSJunchao Zhang    Notes:
1989dc7b89cSJunchao Zhang     Call DMDAVecRestoreKokkosOffsetViewDOF() or DMDAVecRestoreKokkosOffsetViewDOFWrite() once you have finished accessing the OffsetView.
1999dc7b89cSJunchao Zhang 
2009dc7b89cSJunchao Zhang     If the vector is not a Kokkos vector, an error will be raised.
2019dc7b89cSJunchao Zhang 
2029dc7b89cSJunchao Zhang     If the vector is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
2039dc7b89cSJunchao Zhang     a global vector then the ghost points are not accessible. Of course with the local vector you will have to do the
2049dc7b89cSJunchao Zhang     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
2059dc7b89cSJunchao Zhang 
2069dc7b89cSJunchao Zhang     These routines are similar to DMDAVecGetArrayDOF() and friends. One can read-only, write-only or read/write access the returned
2079dc7b89cSJunchao Zhang     Kokkos OffsetView.  Note that passing in a constant OffsetView enables read-only access.
2089dc7b89cSJunchao Zhang     Currently, only two memory spaces are supported: Kokkos::HostSpace and Kokkos::DefaultExecutionSpace::memory_space.
2099dc7b89cSJunchao Zhang     If needed, a memory copy will be internally called to copy the latest vector data to the given memory space.
2109dc7b89cSJunchao Zhang 
2119dc7b89cSJunchao Zhang     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]),
2129dc7b89cSJunchao Zhang     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(),
2139dc7b89cSJunchao Zhang     for example.
2149dc7b89cSJunchao Zhang 
2159dc7b89cSJunchao Zhang     To give users the same experience as DMDAVecGetArrayDOF(), we mandate the returned OffsetView always has Kokkos::LayoutRight (that is, rightest
2169dc7b89cSJunchao 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
2179dc7b89cSJunchao Zhang     to use Iterate::Right as IterateInner if one uses Kokkos::MDRangePolicy to access the OffsetView.
2189dc7b89cSJunchao Zhang 
2199dc7b89cSJunchao Zhang     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
2209dc7b89cSJunchao Zhang     (rightest) corresponds to DA's x direction.
2219dc7b89cSJunchao Zhang 
2229dc7b89cSJunchao Zhang     If the vector is a global vector, we have
2239dc7b89cSJunchao Zhang .vb
2249dc7b89cSJunchao Zhang       kv.extent(0) = zm,  kv.begin(0) = zs, kv.end(0) = zs+zm
2259dc7b89cSJunchao Zhang       kv.extent(1) = ym,  kv.begin(1) = ys, kv.end(1) = ys+ym
2269dc7b89cSJunchao Zhang       kv.extent(2) = xm,  kv.begin(2) = xs, kv.end(2) = xs+xm
2279dc7b89cSJunchao Zhang       kv.extent(3) = dof, kv.begin(3) = 0,  kv.end(3) = dof
2289dc7b89cSJunchao Zhang .ve
2299dc7b89cSJunchao Zhang     If the vector is a local vector, we have
2309dc7b89cSJunchao Zhang .vb
2319dc7b89cSJunchao Zhang       kv.extent(0) = gzm, kv.begin(0) = gzs, kv.end(0) = gzs+gzm
2329dc7b89cSJunchao Zhang       kv.extent(1) = gym, kv.begin(1) = gys, kv.end(1) = gys+gym
2339dc7b89cSJunchao Zhang       kv.extent(2) = gxm, kv.begin(2) = gxs, kv.end(2) = gxs+gxm
2349dc7b89cSJunchao Zhang       kv.extent(3) = dof, kv.begin(3) = 0,   kv.end(3) = dof
2359dc7b89cSJunchao Zhang .ve
2369dc7b89cSJunchao Zhang 
2379dc7b89cSJunchao Zhang     The starts and widths above are obtained by
2389dc7b89cSJunchao Zhang .vb
2399dc7b89cSJunchao Zhang      DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
2409dc7b89cSJunchao Zhang      DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
2419dc7b89cSJunchao Zhang .ve
2429dc7b89cSJunchao Zhang 
2439dc7b89cSJunchao Zhang     For example, to initialize a grid,
2449dc7b89cSJunchao Zhang .vb
2459dc7b89cSJunchao Zhang     typedef Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace> PetscScalarKokkosOffsetView4D;
2469dc7b89cSJunchao Zhang 
2479dc7b89cSJunchao Zhang     PetscScalarKokkosOffsetView4D kv;
2489dc7b89cSJunchao Zhang     DMDAVecGetKokkosOffsetViewDOFWrite(da,v,&kv); // v is a global vector
2499dc7b89cSJunchao Zhang     DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
2509dc7b89cSJunchao Zhang 
2519dc7b89cSJunchao Zhang     parallel_for ("Label",MDRangePolicy <Rank<4, Iterate::Right, Iterate::Right>>(
2529dc7b89cSJunchao Zhang       {zs,ys,xs,0},{zs+zm,ys+ym,xs+xm,dof}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i,PetscInt c) {
2539dc7b89cSJunchao Zhang       kv(k,j,i,c) = ...;
2549dc7b89cSJunchao Zhang     });
2559dc7b89cSJunchao Zhang     DMDAVecRestoreKokkosOffsetViewDOFWrite(da,v,&kv);
2569dc7b89cSJunchao Zhang .ve
2579dc7b89cSJunchao Zhang 
2589dc7b89cSJunchao Zhang   Level: intermediate
2599dc7b89cSJunchao Zhang 
260*db781477SPatrick Sanan .seealso: `DMDAVecRestoreKokkosOffsetViewDOF()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
261*db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
262*db781477SPatrick Sanan           `DMStagVecGetArray()`
2639dc7b89cSJunchao Zhang @*/
2649dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF         (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
2659dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF         (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
2669dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite    (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
2679dc7b89cSJunchao Zhang 
2689dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF         (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar***, Kokkos::LayoutRight,MemorySpace>*);
2699dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF         (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar***, Kokkos::LayoutRight,MemorySpace>*);
2709dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite    (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar***, Kokkos::LayoutRight,MemorySpace>*);
2719dc7b89cSJunchao Zhang 
2729dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF         (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar****, Kokkos::LayoutRight,MemorySpace>*);
2739dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF         (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar****, Kokkos::LayoutRight,MemorySpace>*);
2749dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite    (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar****, Kokkos::LayoutRight,MemorySpace>*);
2759dc7b89cSJunchao Zhang 
2769dc7b89cSJunchao Zhang /*@C
2779dc7b89cSJunchao Zhang    DMDAVecRestoreKokkosOffsetViewDOF - Returns the Kokkos OffsetView that gotten from DMDAVecGetKokkosOffsetViewDOF()
2789dc7b89cSJunchao Zhang 
2799dc7b89cSJunchao Zhang    Synopsis:
2809dc7b89cSJunchao Zhang    #include <petscdmda_kokkos.hpp>
2819dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
2829dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
2839dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
2849dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
2859dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
2869dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
2879dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
2889dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
2899dc7b89cSJunchao Zhang    PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
2909dc7b89cSJunchao Zhang 
2919dc7b89cSJunchao Zhang    Logically collective on da
2929dc7b89cSJunchao Zhang 
2939dc7b89cSJunchao Zhang    Input Parameters:
2949dc7b89cSJunchao Zhang +  da - the distributed array
2959dc7b89cSJunchao Zhang .  v - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
2969dc7b89cSJunchao Zhang -  kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
2979dc7b89cSJunchao Zhang 
2989dc7b89cSJunchao Zhang    Notes:
2999dc7b89cSJunchao Zhang     If the vector is not of type VECKOKKOS, an error will be raised.
3009dc7b89cSJunchao Zhang 
3019dc7b89cSJunchao Zhang   Level: intermediate
3029dc7b89cSJunchao Zhang 
303*db781477SPatrick Sanan .seealso: `DMDAVecGetKokkosOffsetViewDOF()`, `DMDAVecGetKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
304*db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
305*db781477SPatrick Sanan           `DMStagVecGetArray()`
3069dc7b89cSJunchao Zhang @*/
3079dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF     (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
3089dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF     (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
3099dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar**, Kokkos::LayoutRight,MemorySpace>*);
3109dc7b89cSJunchao Zhang 
3119dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF     (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar***, Kokkos::LayoutRight,MemorySpace>*);
3129dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF     (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar***, Kokkos::LayoutRight,MemorySpace>*);
3139dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar***, Kokkos::LayoutRight,MemorySpace>*);
3149dc7b89cSJunchao Zhang 
3159dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF     (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar****, Kokkos::LayoutRight,MemorySpace>*);
3169dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF     (DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar****, Kokkos::LayoutRight,MemorySpace>*);
3179dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM,Vec,Kokkos::Experimental::OffsetView<      PetscScalar****, Kokkos::LayoutRight,MemorySpace>*);
318c6583b63SJunchao Zhang #endif
319c6583b63SJunchao Zhang 
320c6583b63SJunchao Zhang #endif
321