DMDAVecGetKokkosOffsetView#

Gets a Kokkos OffsetView that contains up-to-date data of a vector in the given memory space.

Synopsis#

template <class MemorySpace>
PetscErrorCode DMDAVecGetKokkosOffsetView(DM, Vec, Kokkos::Experimental::OffsetView<const PetscScalar *, MemorySpace> *)

Synopsis#

#include <petscdmda_kokkos.hpp>
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar*,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar*,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar**,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<const PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetView(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);
PetscErrorCode DMDAVecGetKokkosOffsetViewWrite(DM da,Vec v,Kokkos::Experimental::OffsetView<PetscScalar***,Kokkos::LayoutRight,MemorySpace>* kv);

Logically Collective, No Fortran Support

Input Parameters#

Output Parameter#

  • kv - the Kokkos OffsetView with a user-specified template parameter MemorySpace

Notes#

Call DMDAVecRestoreKokkosOffsetView() or DMDAVecRestoreKokkosOffsetViewWrite() once you have finished accessing the OffsetView.

If the vector is not of type 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 DMDAVecGetArray() 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 specified memory space.

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]), where i, j, k are loop variables for the x, y, z dimensions respectively specified in DMDACreate3d(), for example.

To give users the same experience as DMDAVecGetArray(), 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 the OffsetView kv’s first dimension (i.e., the leftest, dim 0) corresponds to the DMDA’s z direction, and its last dimension (rightest) corresponds to DMDA’s x direction.

If the vector is a global vector, we have

      kv.extent(0) = zm*dof,  kv.begin(0) = zs*dof, kv.end(0) = (zs+zm)*dof
      kv.extent(1) = ym*dof,  kv.begin(1) = ys*dof, kv.end(1) = (ys+ym)*dof
      kv.extent(2) = xm*dof,  kv.begin(2) = xs*dof, kv.end(2) = (xs+xm)*dof

If the vector is a local vector, we have

      kv.extent(0) = gzm*dof,  kv.begin(0) = gzs*dof, kv.end(0) = (gzs+gzm)*dof
      kv.extent(1) = gym*dof,  kv.begin(1) = gys*dof, kv.end(1) = (gys+gym)*dof
      kv.extent(2) = gxm*dof,  kv.begin(2) = gxs*dof, kv.end(2) = (gxs+gxm)*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> PetscScalarKokkosOffsetView3D;

    PetscScalarKokkosOffsetView3D kv;
    DMDAVecGetKokkosOffsetViewWrite(da,v,&kv); // v is a global vector and we assume dof is 1
    DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);

    parallel_for ("Label",MDRangePolicy <Rank<3, Iterate::Right, Iterate::Right>>(
      {zs,ys,xs},{zs+zm,ys+ym,xs+xm}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i) {
      kv(k,j,i) = ...;
    });
    DMDAVecRestoreKokkosOffsetViewWrite(da,v,&kv);

For a multi-component problem, one could cast the returned OffsetView to a user’s type. But one has also to shrink the OffsetView’s extent accordingly. For example,

    typedef struct {
      PetscScalar omega,temperature;
    } Node;

    using NodeKokkosOffsetView3D = Kokkos::Experimental::OffsetView<const Node***,Kokkos::LayoutRight,MemorySpace>;
    DMDAVecGetKokkosOffsetViewWrite(da,v,&tv);
    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});

    parallel_for ("Label",MDRangePolicy<Rank<3, Iterate::Right, Iterate::Right>>(
      {zs,ys,xs},{zs+zm,ys+ym,xs+xm}), KOKKOS_LAMBDA (PetscInt k,PetscInt j,PetscInt i) {
      kv(k,j,i).omega = ...;
    });
    DMDAVecRestoreKokkosOffsetViewWrite(da,v,&tv)`;

See Also#

DMDAVecRestoreKokkosOffsetView(), DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF() DMDAVecGetArrayDOF(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecGetArray()

Level#

intermediate

Location#

include/petscdmda_kokkos.hpp

Examples#

src/snes/tutorials/ex55k.kokkos.cxx
src/snes/tutorials/ex3k.kokkos.cxx


Index of all DMDA routines
Table of Contents for all manual pages
Index of all manual pages