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