DMStagRestoreProductCoordinateArrays#
restore local array access
Synopsis#
PetscErrorCode DMStagRestoreProductCoordinateArrays(DM dm, void *arrX, void *arrY, void *arrZ)
Logically Collective
Input Parameter#
dm - the
DMSTAG
object
Output Parameters#
arrX - local 1D coordinate arrays for x direction
arrY - local 1D coordinate arrays for y direction
arrZ - local 1D coordinate arrays for z direction
Notes#
This function does not automatically perform a local->global scatter to populate global coordinates from the local coordinates.
Thus, it may be required to explicitly perform these operations in some situations, as in the following partial example#
PetscCall(DMGetCoordinateDM(dm,&cdm));
for (d=0; d<3; ++d) {
DM subdm;
Vec coor,coor_local;
PetscCall(DMProductGetDM(cdm,d,&subdm));
PetscCall(DMGetCoordinates(subdm,&coor));
PetscCall(DMGetCoordinatesLocal(subdm,&coor_local));
PetscCall(DMLocalToGlobal(subdm,coor_local,INSERT_VALUES,coor));
PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Coordinates dim %" PetscInt_FMT ":\n",d));
PetscCall(VecView(coor,PETSC_VIEWER_STDOUT_WORLD));
}
See Also#
DMSTAG: Staggered, Structured Grid, DMSTAG
, DMStagGetProductCoordinateArrays()
, DMStagGetProductCoordinateArraysRead()
Level#
intermediate
Location#
Examples#
src/dm/impls/stag/tutorials/ex6.c
Index of all DMStag routines
Table of Contents for all manual pages
Index of all manual pages