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
Synopsis#
template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar **, Kokkos::LayoutRight, MemorySpace> *)
Synopsis#
#include <petscdmda_kokkos.hpp>
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOF(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewDOFWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar****,Kokkos::LayoutRight,MemorySpace>* kv);
Logically Collective
Input Parameters#
da - the distributed array
v - the vector, either a vector the same size as one obtained with
DMCreateGlobalVector()
orDMCreateLocalVector()
Output Parameter#
kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace
Notes#
Call DMDAVecRestoreKokkosOffsetViewDOF()
or DMDAVecRestoreKokkosOffsetViewDOFWrite()
once you have finished accessing the OffsetView.
If the vector is not a VECKOKKOS
an error will be raised.
If the vector is a local vector (obtained with DMCreateLocalVector()
etc) then the ghost point locations are accessible. If it is
a global vector then the ghost points are not accessible. Of course with the local vector you will have to do the
appropriate DMGlobalToLocalBegin()
and DMGlobalToLocalEnd()
to have correct values in the ghost locations.
These routines are similar to DMDAVecGetArrayDOF()
and friends. One can read-only, write-only or read/write access the returned
Kokkos OffsetView. Note that passing in a constant OffsetView enables read-only access.
Currently, only two memory spaces are supported: Kokkos::HostSpace and Kokkos::DefaultExecutionSpace::memory_space.
If needed, a memory copy will be internally called to copy the latest vector data to the given memory space.
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]),
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()
,
for example.
To give users the same experience as DMDAVecGetArrayDOF()
, we mandate the returned OffsetView always has Kokkos::LayoutRight (that is, rightest
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
to use Iterate::Right as IterateInner if one uses Kokkos::MDRangePolicy to access the OffsetView.
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
(rightest) corresponds to DMDA’s x direction.
If the vector is a global vector, we have
kv.extent(0) = zm, kv.begin(0) = zs, kv.end(0) = zs+zm
kv.extent(1) = ym, kv.begin(1) = ys, kv.end(1) = ys+ym
kv.extent(2) = xm, kv.begin(2) = xs, kv.end(2) = xs+xm
kv.extent(3) = dof, kv.begin(3) = 0, kv.end(3) = dof
If the vector is a local vector, we have
kv.extent(0) = gzm, kv.begin(0) = gzs, kv.end(0) = gzs+gzm
kv.extent(1) = gym, kv.begin(1) = gys, kv.end(1) = gys+gym
kv.extent(2) = gxm, kv.begin(2) = gxs, kv.end(2) = gxs+gxm
kv.extent(3) = dof, kv.begin(3) = 0, kv.end(3) = dof
The starts and widths above are obtained by
DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
For example, to initialize a grid,
typedef Kokkos::Experimental::OffsetView<const PetscScalar****,Kokkos::LayoutRight,MemorySpace> PetscScalarKokkosOffsetView4D;
PetscScalarKokkosOffsetView4D kv;
DMDAVecGetKokkosOffsetViewDOFWrite(da,v,&kv); // v is a global vector
DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
parallel_for ("Label",MDRangePolicy <Rank<4, Iterate::Right, Iterate::Right>>(
{zs,ys,xs,0},{zs+zm,ys+ym,xs+xm,dof}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i,PetscInt c) {
kv(k,j,i,c) = ...;
});
DMDAVecRestoreKokkosOffsetViewDOFWrite(da,v,&kv);
See Also#
DMDAVecRestoreKokkosOffsetViewDOF()
, DMDAGetGhostCorners()
, DMDAGetCorners()
, VecGetArray()
, VecRestoreArray()
, DMDAVecRestoreArray()
, DMDAVecRestoreArrayDOF()
DMDAVecGetArrayDOF()
, DMDAVecGetArrayWrite()
, DMDAVecRestoreArrayWrite()
, DMDAVecGetArrayRead()
, DMDAVecRestoreArrayRead()
,
DMStagVecGetArray()
Level#
intermediate
Location#
Index of all DMDA routines
Table of Contents for all manual pages
Index of all manual pages