1a4963045SJacob Faibussowitsch #pragma once 2c6583b63SJunchao Zhang 311d22bbfSJunchao Zhang #include <petscvec_kokkos.hpp> 4c6583b63SJunchao Zhang #include <petscdmda.h> 5c6583b63SJunchao Zhang 6*ce78bad3SBarry Smith /* MANSEC = DM */ 7ac09b921SBarry Smith /* SUBMANSEC = DMDA */ 8ac09b921SBarry Smith 9c6583b63SJunchao Zhang #if defined(PETSC_HAVE_KOKKOS) 10c6583b63SJunchao Zhang #include <Kokkos_Core.hpp> 11c6583b63SJunchao Zhang #include <Kokkos_OffsetView.hpp> 12c6583b63SJunchao Zhang 139dc7b89cSJunchao Zhang /*@C 149dc7b89cSJunchao Zhang DMDAVecGetKokkosOffsetView - Gets a Kokkos OffsetView that contains up-to-date data of a vector in the given memory space. 159dc7b89cSJunchao Zhang 169dc7b89cSJunchao Zhang Synopsis: 179dc7b89cSJunchao Zhang #include <petscdmda_kokkos.hpp> 18bf34f4f8SJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>* kv); 19bf34f4f8SJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv); 20bf34f4f8SJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv); 219dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 229dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 239dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 249dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 259dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 269dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 279dc7b89cSJunchao Zhang 28cc4c1da9SBarry Smith Logically Collective, No Fortran Support 299dc7b89cSJunchao Zhang 309dc7b89cSJunchao Zhang Input Parameters: 319dc7b89cSJunchao Zhang + da - the distributed array 3287497f52SBarry Smith - v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 339dc7b89cSJunchao Zhang 349dc7b89cSJunchao Zhang Output Parameter: 359dc7b89cSJunchao Zhang . kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace 369dc7b89cSJunchao Zhang 379dc7b89cSJunchao Zhang Notes: 3887497f52SBarry Smith Call `DMDAVecRestoreKokkosOffsetView()` or `DMDAVecRestoreKokkosOffsetViewWrite()` once you have finished accessing the OffsetView. 399dc7b89cSJunchao Zhang 4087497f52SBarry Smith If the vector is not of type `VECKOKKOS`, an error will be raised. 419dc7b89cSJunchao Zhang 4287497f52SBarry Smith If the vector is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is 439dc7b89cSJunchao Zhang a global vector then the ghost points are not accessible. Of course with the local vector you will have to do the 4487497f52SBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 459dc7b89cSJunchao Zhang 4687497f52SBarry Smith These routines are similar to `DMDAVecGetArray()` and friends. One can read-only, write-only or read/write access the returned 479dc7b89cSJunchao Zhang Kokkos OffsetView. Note that passing in a constant OffsetView enables read-only access. 4845402d8aSJunchao Zhang Currently, only two memory spaces are supported: HostMirrorMemorySpace and Kokkos::DefaultExecutionSpace::memory_space. 499dc7b89cSJunchao Zhang If needed, a memory copy will be internally called to copy the latest vector data to the specified memory space. 509dc7b89cSJunchao Zhang 5187497f52SBarry 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]), 5287497f52SBarry Smith where i, j, k are loop variables for the x, y, z dimensions respectively specified in `DMDACreate3d()`, for example. 539dc7b89cSJunchao Zhang 5487497f52SBarry Smith To give users the same experience as `DMDAVecGetArray()`, we mandate the returned OffsetView always has Kokkos::LayoutRight (that is, rightest 559dc7b89cSJunchao 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 569dc7b89cSJunchao Zhang to use Iterate::Right as IterateInner if one uses Kokkos::MDRangePolicy to access the OffsetView. 579dc7b89cSJunchao Zhang 5887497f52SBarry 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 5987497f52SBarry Smith (rightest) corresponds to DMDA's x direction. 609dc7b89cSJunchao Zhang 619dc7b89cSJunchao Zhang If the vector is a global vector, we have 629dc7b89cSJunchao Zhang .vb 639dc7b89cSJunchao Zhang kv.extent(0) = zm*dof, kv.begin(0) = zs*dof, kv.end(0) = (zs+zm)*dof 649dc7b89cSJunchao Zhang kv.extent(1) = ym*dof, kv.begin(1) = ys*dof, kv.end(1) = (ys+ym)*dof 659dc7b89cSJunchao Zhang kv.extent(2) = xm*dof, kv.begin(2) = xs*dof, kv.end(2) = (xs+xm)*dof 669dc7b89cSJunchao Zhang .ve 679dc7b89cSJunchao Zhang If the vector is a local vector, we have 689dc7b89cSJunchao Zhang .vb 699dc7b89cSJunchao Zhang kv.extent(0) = gzm*dof, kv.begin(0) = gzs*dof, kv.end(0) = (gzs+gzm)*dof 709dc7b89cSJunchao Zhang kv.extent(1) = gym*dof, kv.begin(1) = gys*dof, kv.end(1) = (gys+gym)*dof 719dc7b89cSJunchao Zhang kv.extent(2) = gxm*dof, kv.begin(2) = gxs*dof, kv.end(2) = (gxs+gxm)*dof 729dc7b89cSJunchao Zhang .ve 739dc7b89cSJunchao Zhang 749dc7b89cSJunchao Zhang The starts and widths above are obtained by 759dc7b89cSJunchao Zhang .vb 769dc7b89cSJunchao Zhang DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm); 779dc7b89cSJunchao Zhang DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm); 789dc7b89cSJunchao Zhang .ve 799dc7b89cSJunchao Zhang 809dc7b89cSJunchao Zhang For example, to initialize a grid, 819dc7b89cSJunchao Zhang 829dc7b89cSJunchao Zhang .vb 839dc7b89cSJunchao Zhang typedef Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace> PetscScalarKokkosOffsetView3D; 849dc7b89cSJunchao Zhang 859dc7b89cSJunchao Zhang PetscScalarKokkosOffsetView3D kv; 869dc7b89cSJunchao Zhang DMDAVecGetKokkosOffsetViewWrite(da,v,&kv); // v is a global vector and we assume dof is 1 879dc7b89cSJunchao Zhang DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm); 889dc7b89cSJunchao Zhang 899dc7b89cSJunchao Zhang parallel_for ("Label",MDRangePolicy <Rank<3, Iterate::Right, Iterate::Right>>( 909dc7b89cSJunchao Zhang {zs,ys,xs},{zs+zm,ys+ym,xs+xm}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i) { 919dc7b89cSJunchao Zhang kv(k,j,i) = ...; 929dc7b89cSJunchao Zhang }); 939dc7b89cSJunchao Zhang DMDAVecRestoreKokkosOffsetViewWrite(da,v,&kv); 949dc7b89cSJunchao Zhang .ve 959dc7b89cSJunchao Zhang 969dc7b89cSJunchao Zhang For a multi-component problem, one could cast the returned OffsetView to a user's type. But one has also to shrink 979dc7b89cSJunchao Zhang the OffsetView's extent accordingly. For example, 989dc7b89cSJunchao Zhang .vb 999dc7b89cSJunchao Zhang typedef struct { 1009dc7b89cSJunchao Zhang PetscScalar omega,temperature; 1019dc7b89cSJunchao Zhang } Node; 1029dc7b89cSJunchao Zhang 1039dc7b89cSJunchao Zhang using NodeKokkosOffsetView3D = Kokkos::Experimental::OffsetView<const Node***,Kokkos::LayoutRight,MemorySpace>; 1049dc7b89cSJunchao Zhang DMDAVecGetKokkosOffsetViewWrite(da,v,&tv); 1059dc7b89cSJunchao 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}); 1069dc7b89cSJunchao Zhang 1079dc7b89cSJunchao Zhang parallel_for ("Label",MDRangePolicy<Rank<3, Iterate::Right, Iterate::Right>>( 1089dc7b89cSJunchao Zhang {zs,ys,xs},{zs+zm,ys+ym,xs+xm}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i) { 1099dc7b89cSJunchao Zhang kv(k,j,i).omega = ...; 1109dc7b89cSJunchao Zhang }); 11187497f52SBarry Smith DMDAVecRestoreKokkosOffsetViewWrite(da,v,&tv)`; 1129dc7b89cSJunchao Zhang .ve 1139dc7b89cSJunchao Zhang 1149dc7b89cSJunchao Zhang Level: intermediate 1159dc7b89cSJunchao Zhang 116db781477SPatrick Sanan .seealso: `DMDAVecRestoreKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()` 117db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 118db781477SPatrick Sanan `DMStagVecGetArray()` 1199dc7b89cSJunchao Zhang @*/ 1209371c9d4SSatish Balay template <class MemorySpace> 1219371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar *, MemorySpace> *); 1229371c9d4SSatish Balay template <class MemorySpace> 1239371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *); 1249371c9d4SSatish Balay template <class MemorySpace> 1259371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *); 1269dc7b89cSJunchao Zhang 1279371c9d4SSatish Balay template <class MemorySpace> 1289371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 1299371c9d4SSatish Balay template <class MemorySpace> 1309371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 1319371c9d4SSatish Balay template <class MemorySpace> 1329371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 1339dc7b89cSJunchao Zhang 1349371c9d4SSatish Balay template <class MemorySpace> 1359371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 1369371c9d4SSatish Balay template <class MemorySpace> 1379371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 1389371c9d4SSatish Balay template <class MemorySpace> 1399371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 1409dc7b89cSJunchao Zhang 1419dc7b89cSJunchao Zhang /*@C 14287497f52SBarry Smith DMDAVecRestoreKokkosOffsetView - Returns the Kokkos OffsetView that was gotten with `DMDAVecGetKokkosOffsetView()` 1439dc7b89cSJunchao Zhang 1449dc7b89cSJunchao Zhang Synopsis: 1459dc7b89cSJunchao Zhang #include <petscdmda_kokkos.hpp> 146bf34f4f8SJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>* kv); 147bf34f4f8SJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv); 148bf34f4f8SJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv); 1499dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 1509dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 1519dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 1529dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 1539dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 1549dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 1559dc7b89cSJunchao Zhang 156cc4c1da9SBarry Smith Logically Collective, No Fortran Support 1579dc7b89cSJunchao Zhang 1589dc7b89cSJunchao Zhang Input Parameters: 1599dc7b89cSJunchao Zhang + da - the distributed array 16087497f52SBarry Smith . v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 1619dc7b89cSJunchao Zhang - kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace 1629dc7b89cSJunchao Zhang 16320f4b53cSBarry Smith Level: intermediate 16420f4b53cSBarry Smith 16587497f52SBarry Smith Note: 16687497f52SBarry Smith If the vector is not of type `VECKOKKOS`, an error will be raised. 1679dc7b89cSJunchao Zhang 168db781477SPatrick Sanan .seealso: `DMDAVecGetKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()` 169db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 170db781477SPatrick Sanan `DMStagVecGetArray()` 1719dc7b89cSJunchao Zhang @*/ 1729371c9d4SSatish Balay template <class MemorySpace> 1739371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar *, MemorySpace> *); 1749371c9d4SSatish Balay template <class MemorySpace> 1759371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *); 1769371c9d4SSatish Balay template <class MemorySpace> 1779371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar *, MemorySpace> *); 178c6583b63SJunchao Zhang 1799371c9d4SSatish Balay template <class MemorySpace> 1809371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 1819371c9d4SSatish Balay template <class MemorySpace> 1829371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 1839371c9d4SSatish Balay template <class MemorySpace> 1849371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 185c6583b63SJunchao Zhang 1869371c9d4SSatish Balay template <class MemorySpace> 1879371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 1889371c9d4SSatish Balay template <class MemorySpace> 1899371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 1909371c9d4SSatish Balay template <class MemorySpace> 1919371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 192c6583b63SJunchao Zhang 1939dc7b89cSJunchao Zhang /*@C 1949dc7b89cSJunchao 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 1959dc7b89cSJunchao Zhang 1969dc7b89cSJunchao Zhang Synopsis: 1979dc7b89cSJunchao Zhang #include <petscdmda_kokkos.hpp> 198bf34f4f8SJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 1999dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 2009dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 2019dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 2029dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 2039dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 2049dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 2059dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 2069dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 2079dc7b89cSJunchao Zhang 208cc4c1da9SBarry Smith Logically Collective, No Fortran Support 2099dc7b89cSJunchao Zhang 2109dc7b89cSJunchao Zhang Input Parameters: 2119dc7b89cSJunchao Zhang + da - the distributed array 21287497f52SBarry Smith - v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 2139dc7b89cSJunchao Zhang 2149dc7b89cSJunchao Zhang Output Parameter: 2159dc7b89cSJunchao Zhang . kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace 2169dc7b89cSJunchao Zhang 2179dc7b89cSJunchao Zhang Notes: 21887497f52SBarry Smith Call `DMDAVecRestoreKokkosOffsetViewDOF()` or `DMDAVecRestoreKokkosOffsetViewDOFWrite()` once you have finished accessing the OffsetView. 2199dc7b89cSJunchao Zhang 22087497f52SBarry Smith If the vector is not a `VECKOKKOS` an error will be raised. 2219dc7b89cSJunchao Zhang 22287497f52SBarry Smith If the vector is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is 2239dc7b89cSJunchao Zhang a global vector then the ghost points are not accessible. Of course with the local vector you will have to do the 22487497f52SBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 2259dc7b89cSJunchao Zhang 22687497f52SBarry Smith These routines are similar to `DMDAVecGetArrayDOF()` and friends. One can read-only, write-only or read/write access the returned 2279dc7b89cSJunchao Zhang Kokkos OffsetView. Note that passing in a constant OffsetView enables read-only access. 22845402d8aSJunchao Zhang Currently, only two memory spaces are supported: HostMirrorMemorySpace and Kokkos::DefaultExecutionSpace::memory_space. 2299dc7b89cSJunchao Zhang If needed, a memory copy will be internally called to copy the latest vector data to the given memory space. 2309dc7b89cSJunchao Zhang 23187497f52SBarry 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]), 23287497f52SBarry 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()`, 2339dc7b89cSJunchao Zhang for example. 2349dc7b89cSJunchao Zhang 23587497f52SBarry Smith To give users the same experience as `DMDAVecGetArrayDOF()`, we mandate the returned OffsetView always has Kokkos::LayoutRight (that is, rightest 2369dc7b89cSJunchao 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 2379dc7b89cSJunchao Zhang to use Iterate::Right as IterateInner if one uses Kokkos::MDRangePolicy to access the OffsetView. 2389dc7b89cSJunchao Zhang 23987497f52SBarry 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 24087497f52SBarry Smith (rightest) corresponds to DMDA's x direction. 2419dc7b89cSJunchao Zhang 2429dc7b89cSJunchao Zhang If the vector is a global vector, we have 2439dc7b89cSJunchao Zhang .vb 2449dc7b89cSJunchao Zhang kv.extent(0) = zm, kv.begin(0) = zs, kv.end(0) = zs+zm 2459dc7b89cSJunchao Zhang kv.extent(1) = ym, kv.begin(1) = ys, kv.end(1) = ys+ym 2469dc7b89cSJunchao Zhang kv.extent(2) = xm, kv.begin(2) = xs, kv.end(2) = xs+xm 2479dc7b89cSJunchao Zhang kv.extent(3) = dof, kv.begin(3) = 0, kv.end(3) = dof 2489dc7b89cSJunchao Zhang .ve 2499dc7b89cSJunchao Zhang If the vector is a local vector, we have 2509dc7b89cSJunchao Zhang .vb 2519dc7b89cSJunchao Zhang kv.extent(0) = gzm, kv.begin(0) = gzs, kv.end(0) = gzs+gzm 2529dc7b89cSJunchao Zhang kv.extent(1) = gym, kv.begin(1) = gys, kv.end(1) = gys+gym 2539dc7b89cSJunchao Zhang kv.extent(2) = gxm, kv.begin(2) = gxs, kv.end(2) = gxs+gxm 2549dc7b89cSJunchao Zhang kv.extent(3) = dof, kv.begin(3) = 0, kv.end(3) = dof 2559dc7b89cSJunchao Zhang .ve 2569dc7b89cSJunchao Zhang 2579dc7b89cSJunchao Zhang The starts and widths above are obtained by 2589dc7b89cSJunchao Zhang .vb 2599dc7b89cSJunchao Zhang DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm); 2609dc7b89cSJunchao Zhang DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm); 2619dc7b89cSJunchao Zhang .ve 2629dc7b89cSJunchao Zhang 2639dc7b89cSJunchao Zhang For example, to initialize a grid, 2649dc7b89cSJunchao Zhang .vb 2659dc7b89cSJunchao Zhang typedef Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace> PetscScalarKokkosOffsetView4D; 2669dc7b89cSJunchao Zhang 2679dc7b89cSJunchao Zhang PetscScalarKokkosOffsetView4D kv; 2689dc7b89cSJunchao Zhang DMDAVecGetKokkosOffsetViewDOFWrite(da,v,&kv); // v is a global vector 2699dc7b89cSJunchao Zhang DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm); 2709dc7b89cSJunchao Zhang 2719dc7b89cSJunchao Zhang parallel_for ("Label",MDRangePolicy <Rank<4, Iterate::Right, Iterate::Right>>( 2729dc7b89cSJunchao Zhang {zs,ys,xs,0},{zs+zm,ys+ym,xs+xm,dof}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i,PetscInt c) { 2739dc7b89cSJunchao Zhang kv(k,j,i,c) = ...; 2749dc7b89cSJunchao Zhang }); 2759dc7b89cSJunchao Zhang DMDAVecRestoreKokkosOffsetViewDOFWrite(da,v,&kv); 2769dc7b89cSJunchao Zhang .ve 2779dc7b89cSJunchao Zhang 2789dc7b89cSJunchao Zhang Level: intermediate 2799dc7b89cSJunchao Zhang 280db781477SPatrick Sanan .seealso: `DMDAVecRestoreKokkosOffsetViewDOF()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()` 281db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 282db781477SPatrick Sanan `DMStagVecGetArray()` 2839dc7b89cSJunchao Zhang @*/ 2849371c9d4SSatish Balay template <class MemorySpace> 2859371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 2869371c9d4SSatish Balay template <class MemorySpace> 2879371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 2889371c9d4SSatish Balay template <class MemorySpace> 2899371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 2909dc7b89cSJunchao Zhang 2919371c9d4SSatish Balay template <class MemorySpace> 2929371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 2939371c9d4SSatish Balay template <class MemorySpace> 2949371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 2959371c9d4SSatish Balay template <class MemorySpace> 2969371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 2979dc7b89cSJunchao Zhang 2989371c9d4SSatish Balay template <class MemorySpace> 2999371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *); 3009371c9d4SSatish Balay template <class MemorySpace> 3019371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *); 3029371c9d4SSatish Balay template <class MemorySpace> 3039371c9d4SSatish Balay PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *); 3049dc7b89cSJunchao Zhang 3059dc7b89cSJunchao Zhang /*@C 30687497f52SBarry Smith DMDAVecRestoreKokkosOffsetViewDOF - Returns the Kokkos OffsetView that was gotten from `DMDAVecGetKokkosOffsetViewDOF()` 3079dc7b89cSJunchao Zhang 3089dc7b89cSJunchao Zhang Synopsis: 3099dc7b89cSJunchao Zhang #include <petscdmda_kokkos.hpp> 3109dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 3119dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 3129dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 3139dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 3149dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 3159dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 3169dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 3179dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 3189dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 3199dc7b89cSJunchao Zhang 320cc4c1da9SBarry Smith Logically Collective, No Fortran Support 3219dc7b89cSJunchao Zhang 3229dc7b89cSJunchao Zhang Input Parameters: 3239dc7b89cSJunchao Zhang + da - the distributed array 32487497f52SBarry Smith . v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 3259dc7b89cSJunchao Zhang - kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace 3269dc7b89cSJunchao Zhang 32720f4b53cSBarry Smith Level: intermediate 32820f4b53cSBarry Smith 32987497f52SBarry Smith Note: 33087497f52SBarry Smith If the vector is not of type `VECKOKKOS`, an error will be raised. 3319dc7b89cSJunchao Zhang 332db781477SPatrick Sanan .seealso: `DMDAVecGetKokkosOffsetViewDOF()`, `DMDAVecGetKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()` 333db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 334db781477SPatrick Sanan `DMStagVecGetArray()` 3359dc7b89cSJunchao Zhang @*/ 3369371c9d4SSatish Balay template <class MemorySpace> 3379371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 3389371c9d4SSatish Balay template <class MemorySpace> 3399371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 3409371c9d4SSatish Balay template <class MemorySpace> 3419371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar **, Kokkos::LayoutRight, MemorySpace> *); 3429dc7b89cSJunchao Zhang 3439371c9d4SSatish Balay template <class MemorySpace> 3449371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 3459371c9d4SSatish Balay template <class MemorySpace> 3469371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 3479371c9d4SSatish Balay template <class MemorySpace> 3489371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ***, Kokkos::LayoutRight, MemorySpace> *); 3499dc7b89cSJunchao Zhang 3509371c9d4SSatish Balay template <class MemorySpace> 3519371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *); 3529371c9d4SSatish Balay template <class MemorySpace> 3539371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *); 3549371c9d4SSatish Balay template <class MemorySpace> 3559371c9d4SSatish Balay PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM, Vec, Kokkos::Experimental::OffsetView<PetscScalar ****, Kokkos::LayoutRight, MemorySpace> *); 356c6583b63SJunchao Zhang #endif 357