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 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 289dc7b89cSJunchao Zhang Logically collective on da 299dc7b89cSJunchao Zhang 309dc7b89cSJunchao Zhang Input Parameters: 319dc7b89cSJunchao Zhang + da - the distributed array 32*87497f52SBarry 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: 38*87497f52SBarry Smith Call `DMDAVecRestoreKokkosOffsetView()` or `DMDAVecRestoreKokkosOffsetViewWrite()` once you have finished accessing the OffsetView. 399dc7b89cSJunchao Zhang 40*87497f52SBarry Smith If the vector is not of type `VECKOKKOS`, an error will be raised. 419dc7b89cSJunchao Zhang 42*87497f52SBarry 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 44*87497f52SBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 459dc7b89cSJunchao Zhang 46*87497f52SBarry 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. 489dc7b89cSJunchao Zhang Currently, only two memory spaces are supported: Kokkos::HostSpace 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 51*87497f52SBarry 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]), 52*87497f52SBarry Smith where i, j, k are loop variables for the x, y, z dimensions respectively specified in `DMDACreate3d()`, for example. 539dc7b89cSJunchao Zhang 54*87497f52SBarry 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 58*87497f52SBarry 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 59*87497f52SBarry 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 }); 111*87497f52SBarry 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 @*/ 120c6583b63SJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>*); 121c6583b63SJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar*,MemorySpace>*); 122c6583b63SJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewWrite (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar*,MemorySpace>*); 1239dc7b89cSJunchao Zhang 1249dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 1259dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 1269dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewWrite (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 1279dc7b89cSJunchao Zhang 1289dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>*); 1299dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar***,Kokkos::LayoutRight,MemorySpace>*); 1309dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewWrite (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar***,Kokkos::LayoutRight,MemorySpace>*); 1319dc7b89cSJunchao Zhang 1329dc7b89cSJunchao Zhang /*@C 133*87497f52SBarry Smith DMDAVecRestoreKokkosOffsetView - Returns the Kokkos OffsetView that was gotten with `DMDAVecGetKokkosOffsetView()` 1349dc7b89cSJunchao Zhang 1359dc7b89cSJunchao Zhang Synopsis: 1369dc7b89cSJunchao Zhang #include <petscdmda_kokkos.hpp> 137bf34f4f8SJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>* kv); 138bf34f4f8SJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv); 139bf34f4f8SJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv); 1409dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 1419dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 1429dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 1439dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 1449dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 1459dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 1469dc7b89cSJunchao Zhang 1479dc7b89cSJunchao Zhang Logically collective on da 1489dc7b89cSJunchao Zhang 1499dc7b89cSJunchao Zhang Input Parameters: 1509dc7b89cSJunchao Zhang + da - the distributed array 151*87497f52SBarry Smith . v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 1529dc7b89cSJunchao Zhang - kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace 1539dc7b89cSJunchao Zhang 154*87497f52SBarry Smith Note: 155*87497f52SBarry Smith If the vector is not of type `VECKOKKOS`, an error will be raised. 1569dc7b89cSJunchao Zhang 1579dc7b89cSJunchao Zhang Level: intermediate 1589dc7b89cSJunchao Zhang 159db781477SPatrick Sanan .seealso: `DMDAVecGetKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()` 160db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 161db781477SPatrick Sanan `DMStagVecGetArray()` 1629dc7b89cSJunchao Zhang @*/ 1639dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>*); 1649dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar*,MemorySpace>*); 165c6583b63SJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar*,MemorySpace>*); 166c6583b63SJunchao Zhang 1679dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 1689dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 1699dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 170c6583b63SJunchao Zhang 1719dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>*); 1729dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetView (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar***,Kokkos::LayoutRight,MemorySpace>*); 1739dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewWrite(DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar***,Kokkos::LayoutRight,MemorySpace>*); 174c6583b63SJunchao Zhang 1759dc7b89cSJunchao Zhang /*@C 1769dc7b89cSJunchao 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 1779dc7b89cSJunchao Zhang 1789dc7b89cSJunchao Zhang Synopsis: 1799dc7b89cSJunchao Zhang #include <petscdmda_kokkos.hpp> 180bf34f4f8SJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 1819dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 1829dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 1839dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 1849dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 1859dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 1869dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 1879dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 1889dc7b89cSJunchao Zhang PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 1899dc7b89cSJunchao Zhang 1909dc7b89cSJunchao Zhang Logically collective on da 1919dc7b89cSJunchao Zhang 1929dc7b89cSJunchao Zhang Input Parameters: 1939dc7b89cSJunchao Zhang + da - the distributed array 194*87497f52SBarry Smith - v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 1959dc7b89cSJunchao Zhang 1969dc7b89cSJunchao Zhang Output Parameter: 1979dc7b89cSJunchao Zhang . kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace 1989dc7b89cSJunchao Zhang 1999dc7b89cSJunchao Zhang Notes: 200*87497f52SBarry Smith Call `DMDAVecRestoreKokkosOffsetViewDOF()` or `DMDAVecRestoreKokkosOffsetViewDOFWrite()` once you have finished accessing the OffsetView. 2019dc7b89cSJunchao Zhang 202*87497f52SBarry Smith If the vector is not a `VECKOKKOS` an error will be raised. 2039dc7b89cSJunchao Zhang 204*87497f52SBarry Smith If the vector is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is 2059dc7b89cSJunchao Zhang a global vector then the ghost points are not accessible. Of course with the local vector you will have to do the 206*87497f52SBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 2079dc7b89cSJunchao Zhang 208*87497f52SBarry Smith These routines are similar to `DMDAVecGetArrayDOF()` and friends. One can read-only, write-only or read/write access the returned 2099dc7b89cSJunchao Zhang Kokkos OffsetView. Note that passing in a constant OffsetView enables read-only access. 2109dc7b89cSJunchao Zhang Currently, only two memory spaces are supported: Kokkos::HostSpace and Kokkos::DefaultExecutionSpace::memory_space. 2119dc7b89cSJunchao Zhang If needed, a memory copy will be internally called to copy the latest vector data to the given memory space. 2129dc7b89cSJunchao Zhang 213*87497f52SBarry 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]), 214*87497f52SBarry 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()`, 2159dc7b89cSJunchao Zhang for example. 2169dc7b89cSJunchao Zhang 217*87497f52SBarry Smith To give users the same experience as `DMDAVecGetArrayDOF()`, we mandate the returned OffsetView always has Kokkos::LayoutRight (that is, rightest 2189dc7b89cSJunchao 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 2199dc7b89cSJunchao Zhang to use Iterate::Right as IterateInner if one uses Kokkos::MDRangePolicy to access the OffsetView. 2209dc7b89cSJunchao Zhang 221*87497f52SBarry 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 222*87497f52SBarry Smith (rightest) corresponds to DMDA's x direction. 2239dc7b89cSJunchao Zhang 2249dc7b89cSJunchao Zhang If the vector is a global vector, we have 2259dc7b89cSJunchao Zhang .vb 2269dc7b89cSJunchao Zhang kv.extent(0) = zm, kv.begin(0) = zs, kv.end(0) = zs+zm 2279dc7b89cSJunchao Zhang kv.extent(1) = ym, kv.begin(1) = ys, kv.end(1) = ys+ym 2289dc7b89cSJunchao Zhang kv.extent(2) = xm, kv.begin(2) = xs, kv.end(2) = xs+xm 2299dc7b89cSJunchao Zhang kv.extent(3) = dof, kv.begin(3) = 0, kv.end(3) = dof 2309dc7b89cSJunchao Zhang .ve 2319dc7b89cSJunchao Zhang If the vector is a local vector, we have 2329dc7b89cSJunchao Zhang .vb 2339dc7b89cSJunchao Zhang kv.extent(0) = gzm, kv.begin(0) = gzs, kv.end(0) = gzs+gzm 2349dc7b89cSJunchao Zhang kv.extent(1) = gym, kv.begin(1) = gys, kv.end(1) = gys+gym 2359dc7b89cSJunchao Zhang kv.extent(2) = gxm, kv.begin(2) = gxs, kv.end(2) = gxs+gxm 2369dc7b89cSJunchao Zhang kv.extent(3) = dof, kv.begin(3) = 0, kv.end(3) = dof 2379dc7b89cSJunchao Zhang .ve 2389dc7b89cSJunchao Zhang 2399dc7b89cSJunchao Zhang The starts and widths above are obtained by 2409dc7b89cSJunchao Zhang .vb 2419dc7b89cSJunchao Zhang DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm); 2429dc7b89cSJunchao Zhang DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm); 2439dc7b89cSJunchao Zhang .ve 2449dc7b89cSJunchao Zhang 2459dc7b89cSJunchao Zhang For example, to initialize a grid, 2469dc7b89cSJunchao Zhang .vb 2479dc7b89cSJunchao Zhang typedef Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace> PetscScalarKokkosOffsetView4D; 2489dc7b89cSJunchao Zhang 2499dc7b89cSJunchao Zhang PetscScalarKokkosOffsetView4D kv; 2509dc7b89cSJunchao Zhang DMDAVecGetKokkosOffsetViewDOFWrite(da,v,&kv); // v is a global vector 2519dc7b89cSJunchao Zhang DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm); 2529dc7b89cSJunchao Zhang 2539dc7b89cSJunchao Zhang parallel_for ("Label",MDRangePolicy <Rank<4, Iterate::Right, Iterate::Right>>( 2549dc7b89cSJunchao Zhang {zs,ys,xs,0},{zs+zm,ys+ym,xs+xm,dof}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i,PetscInt c) { 2559dc7b89cSJunchao Zhang kv(k,j,i,c) = ...; 2569dc7b89cSJunchao Zhang }); 2579dc7b89cSJunchao Zhang DMDAVecRestoreKokkosOffsetViewDOFWrite(da,v,&kv); 2589dc7b89cSJunchao Zhang .ve 2599dc7b89cSJunchao Zhang 2609dc7b89cSJunchao Zhang Level: intermediate 2619dc7b89cSJunchao Zhang 262db781477SPatrick Sanan .seealso: `DMDAVecRestoreKokkosOffsetViewDOF()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()` 263db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 264db781477SPatrick Sanan `DMStagVecGetArray()` 2659dc7b89cSJunchao Zhang @*/ 2669dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 2679dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 2689dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 2699dc7b89cSJunchao Zhang 2709dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar***, Kokkos::LayoutRight,MemorySpace>*); 2719dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar***, Kokkos::LayoutRight,MemorySpace>*); 2729dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar***, Kokkos::LayoutRight,MemorySpace>*); 2739dc7b89cSJunchao Zhang 2749dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar****, Kokkos::LayoutRight,MemorySpace>*); 2759dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar****, Kokkos::LayoutRight,MemorySpace>*); 2769dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar****, Kokkos::LayoutRight,MemorySpace>*); 2779dc7b89cSJunchao Zhang 2789dc7b89cSJunchao Zhang /*@C 279*87497f52SBarry Smith DMDAVecRestoreKokkosOffsetViewDOF - Returns the Kokkos OffsetView that was gotten from `DMDAVecGetKokkosOffsetViewDOF()` 2809dc7b89cSJunchao Zhang 2819dc7b89cSJunchao Zhang Synopsis: 2829dc7b89cSJunchao Zhang #include <petscdmda_kokkos.hpp> 2839dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 2849dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 2859dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv); 2869dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 2879dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 2889dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv); 2899dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 2909dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 2919dc7b89cSJunchao Zhang PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv); 2929dc7b89cSJunchao Zhang 2939dc7b89cSJunchao Zhang Logically collective on da 2949dc7b89cSJunchao Zhang 2959dc7b89cSJunchao Zhang Input Parameters: 2969dc7b89cSJunchao Zhang + da - the distributed array 297*87497f52SBarry Smith . v - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 2989dc7b89cSJunchao Zhang - kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace 2999dc7b89cSJunchao Zhang 300*87497f52SBarry Smith Note: 301*87497f52SBarry Smith If the vector is not of type `VECKOKKOS`, an error will be raised. 3029dc7b89cSJunchao Zhang 3039dc7b89cSJunchao Zhang Level: intermediate 3049dc7b89cSJunchao Zhang 305db781477SPatrick Sanan .seealso: `DMDAVecGetKokkosOffsetViewDOF()`, `DMDAVecGetKokkosOffsetView()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()` 306db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 307db781477SPatrick Sanan `DMStagVecGetArray()` 3089dc7b89cSJunchao Zhang @*/ 3099dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 3109dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 3119dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar**, Kokkos::LayoutRight,MemorySpace>*); 3129dc7b89cSJunchao Zhang 3139dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar***, Kokkos::LayoutRight,MemorySpace>*); 3149dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar***, Kokkos::LayoutRight,MemorySpace>*); 3159dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar***, Kokkos::LayoutRight,MemorySpace>*); 3169dc7b89cSJunchao Zhang 3179dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView<const PetscScalar****, Kokkos::LayoutRight,MemorySpace>*); 3189dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOF (DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar****, Kokkos::LayoutRight,MemorySpace>*); 3199dc7b89cSJunchao Zhang template<class MemorySpace> PetscErrorCode DMDAVecRestoreKokkosOffsetViewDOFWrite(DM,Vec,Kokkos::Experimental::OffsetView< PetscScalar****, Kokkos::LayoutRight,MemorySpace>*); 320c6583b63SJunchao Zhang #endif 321c6583b63SJunchao Zhang 322c6583b63SJunchao Zhang #endif 323