147c6ae99SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 347c6ae99SBarry Smith 40f99b6f4SBarry Smith /*MC 50f99b6f4SBarry Smith DMDAVecGetArrayF90 - check Fortran Notes at `DMDAVecGetArray()` 671b416bcSMatthew G. Knepley 771b416bcSMatthew G. Knepley Level: intermediate 80f99b6f4SBarry Smith M*/ 90f99b6f4SBarry Smith 1047c6ae99SBarry Smith /*@C 11aa219208SBarry Smith DMDAVecGetArray - Returns a multiple dimension array that shares data with 1247c6ae99SBarry Smith the underlying vector and is indexed using the global dimensions. 1347c6ae99SBarry Smith 1420f4b53cSBarry Smith Logically Collective 1547c6ae99SBarry Smith 16d8d19677SJose E. Roman Input Parameters: 1747c6ae99SBarry Smith + da - the distributed array 18*0af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 1947c6ae99SBarry Smith 2047c6ae99SBarry Smith Output Parameter: 2147c6ae99SBarry Smith . array - the array 2247c6ae99SBarry Smith 23dce8aebaSBarry Smith Level: intermediate 24dce8aebaSBarry Smith 2547c6ae99SBarry Smith Notes: 26dce8aebaSBarry Smith Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries. 2747c6ae99SBarry Smith 2847c6ae99SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 2947c6ae99SBarry Smith 30*0af9b551SBarry Smith If `vec` is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is 31*0af9b551SBarry Smith a global vector then the ghost points are not accessible. Of course, with a local vector you will have had to do the 32dce8aebaSBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 3347c6ae99SBarry Smith 34*0af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 35*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 36*0af9b551SBarry Smith 3795452b02SPatrick Sanan Fortran Notes: 380f99b6f4SBarry Smith Use `DMDAVecGetArrayF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate 39dce8aebaSBarry Smith dimension. For a `DMDA` created with a dof of 1 use the dimension of the `DMDA`, for a `DMDA` created with a dof greater than 1 use one more than the 400f99b6f4SBarry Smith dimension of the `DMDA`. 4147c6ae99SBarry Smith 42*0af9b551SBarry Smith The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise 43*0af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 44*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 4547c6ae99SBarry Smith 46dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()` 47db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 48db781477SPatrick Sanan `DMStagVecGetArray()` 4947c6ae99SBarry Smith @*/ 50d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArray(DM da, Vec vec, void *array) 51d71ae5a4SJacob Faibussowitsch { 5247c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 5347c6ae99SBarry Smith 5447c6ae99SBarry Smith PetscFunctionBegin; 55a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 566db82c94SMatthew G Knepley PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 576db82c94SMatthew G Knepley PetscValidPointer(array, 3); 589566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 599566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 609566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 6147c6ae99SBarry Smith 6247c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 639566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 6447c6ae99SBarry Smith if (N == xm * ym * zm * dof) { 6547c6ae99SBarry Smith gxm = xm; 6647c6ae99SBarry Smith gym = ym; 6747c6ae99SBarry Smith gzm = zm; 6847c6ae99SBarry Smith gxs = xs; 6947c6ae99SBarry Smith gys = ys; 7047c6ae99SBarry Smith gzs = zs; 7163a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof); 7247c6ae99SBarry Smith 7347c6ae99SBarry Smith if (dim == 1) { 749566063dSJacob Faibussowitsch PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 7547c6ae99SBarry Smith } else if (dim == 2) { 769566063dSJacob Faibussowitsch PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 7747c6ae99SBarry Smith } else if (dim == 3) { 789566063dSJacob Faibussowitsch PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 7963a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8147c6ae99SBarry Smith } 8247c6ae99SBarry Smith 830f99b6f4SBarry Smith /*MC 840f99b6f4SBarry Smith DMDAVecRestoreArrayF90 - check Fortran Notes at `DMDAVecRestoreArray()` 8553c0d4aeSBarry Smith 8653c0d4aeSBarry Smith Level: intermediate 870f99b6f4SBarry Smith M*/ 880f99b6f4SBarry Smith 8947c6ae99SBarry Smith /*@ 90dce8aebaSBarry Smith DMDAVecRestoreArray - Restores a multiple dimension array obtained with `DMDAVecGetArray()` 9147c6ae99SBarry Smith 9220f4b53cSBarry Smith Logically Collective 9347c6ae99SBarry Smith 94d8d19677SJose E. Roman Input Parameters: 9547c6ae99SBarry Smith + da - the distributed array 96*0af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 97*0af9b551SBarry Smith - array - the `array` pointer is zeroed 9847c6ae99SBarry Smith 9947c6ae99SBarry Smith Level: intermediate 10047c6ae99SBarry Smith 101dce8aebaSBarry Smith Fortran Note: 102*0af9b551SBarry Smith Use `DMDAVecRestoreArayF90()` 10347c6ae99SBarry Smith 1040f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArayF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, 105db781477SPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 106db781477SPatrick Sanan `DMDStagVecRestoreArray()` 10747c6ae99SBarry Smith @*/ 108d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array) 109d71ae5a4SJacob Faibussowitsch { 11047c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 11147c6ae99SBarry Smith 11247c6ae99SBarry Smith PetscFunctionBegin; 113a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 1146db82c94SMatthew G Knepley PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 1156db82c94SMatthew G Knepley PetscValidPointer(array, 3); 1169566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 1179566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 1189566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 11947c6ae99SBarry Smith 12047c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 1219566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 12247c6ae99SBarry Smith if (N == xm * ym * zm * dof) { 12347c6ae99SBarry Smith gxm = xm; 12447c6ae99SBarry Smith gym = ym; 12547c6ae99SBarry Smith gzm = zm; 12647c6ae99SBarry Smith gxs = xs; 12747c6ae99SBarry Smith gys = ys; 12847c6ae99SBarry Smith gzs = zs; 12963a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof); 13047c6ae99SBarry Smith 13147c6ae99SBarry Smith if (dim == 1) { 1329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 13347c6ae99SBarry Smith } else if (dim == 2) { 1349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 13547c6ae99SBarry Smith } else if (dim == 3) { 1369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 13763a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 1383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13947c6ae99SBarry Smith } 14047c6ae99SBarry Smith 1410f99b6f4SBarry Smith /*MC 1420f99b6f4SBarry Smith DMDAVecGetArrayWriteF90 - check Fortran Notes at `DMDAVecGetArrayWrite()` 14353c0d4aeSBarry Smith 14453c0d4aeSBarry Smith Level: intermediate 1450f99b6f4SBarry Smith M*/ 1460f99b6f4SBarry Smith 14747c6ae99SBarry Smith /*@C 148fdc842d1SBarry Smith DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with 149fdc842d1SBarry Smith the underlying vector and is indexed using the global dimensions. 150fdc842d1SBarry Smith 15120f4b53cSBarry Smith Logically Collective 152fdc842d1SBarry Smith 153d8d19677SJose E. Roman Input Parameters: 154fdc842d1SBarry Smith + da - the distributed array 155*0af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 156fdc842d1SBarry Smith 157fdc842d1SBarry Smith Output Parameter: 158fdc842d1SBarry Smith . array - the array 159fdc842d1SBarry Smith 160dce8aebaSBarry Smith Level: intermediate 161dce8aebaSBarry Smith 162fdc842d1SBarry Smith Notes: 163dce8aebaSBarry Smith Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries. 164fdc842d1SBarry Smith 165fdc842d1SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 166fdc842d1SBarry Smith 167*0af9b551SBarry Smith if `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is 168fdc842d1SBarry Smith a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the 169dce8aebaSBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 170fdc842d1SBarry Smith 171*0af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 172*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 173*0af9b551SBarry Smith 174fdc842d1SBarry Smith Fortran Notes: 1750f99b6f4SBarry Smith From Fortran use `DMDAVecGetArrayWriteF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 176dce8aebaSBarry Smith dimension. For a `DMDA` created with a dof of 1 use the dimension of the `DMDA`, for a `DMDA` created with a dof greater than 1 use one more than the 1770f99b6f4SBarry Smith dimension of the `DMDA`. 178fdc842d1SBarry Smith 179*0af9b551SBarry Smith The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise 180*0af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 181*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 182fdc842d1SBarry Smith 183dce8aebaSBarry Smith Developer Note: 184dce8aebaSBarry Smith This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()` 185fdc842d1SBarry Smith 186dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecRestoreArrayDOF()` 187db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 188fdc842d1SBarry Smith @*/ 189d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayWrite(DM da, Vec vec, void *array) 190d71ae5a4SJacob Faibussowitsch { 191fdc842d1SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 192fdc842d1SBarry Smith 193fdc842d1SBarry Smith PetscFunctionBegin; 194fdc842d1SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 195fdc842d1SBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 196fdc842d1SBarry Smith PetscValidPointer(array, 3); 1971bb6d2a8SBarry Smith if (da->localSection) { 1989566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array)); 1993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 200fdc842d1SBarry Smith } 2019566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 2029566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 2039566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 204fdc842d1SBarry Smith 205fdc842d1SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 2069566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 207fdc842d1SBarry Smith if (N == xm * ym * zm * dof) { 208fdc842d1SBarry Smith gxm = xm; 209fdc842d1SBarry Smith gym = ym; 210fdc842d1SBarry Smith gzm = zm; 211fdc842d1SBarry Smith gxs = xs; 212fdc842d1SBarry Smith gys = ys; 213fdc842d1SBarry Smith gzs = zs; 21463a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof); 215fdc842d1SBarry Smith 216fdc842d1SBarry Smith if (dim == 1) { 2179566063dSJacob Faibussowitsch PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 218fdc842d1SBarry Smith } else if (dim == 2) { 2199566063dSJacob Faibussowitsch PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 220fdc842d1SBarry Smith } else if (dim == 3) { 2219566063dSJacob Faibussowitsch PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 22263a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 224fdc842d1SBarry Smith } 225fdc842d1SBarry Smith 2260f99b6f4SBarry Smith /*MC 2270f99b6f4SBarry Smith DMDAVecRestoreArrayWriteF90 - check Fortran Notes at `DMDAVecRestoreArrayWrite()` 22853c0d4aeSBarry Smith 22953c0d4aeSBarry Smith Level: intermediate 2300f99b6f4SBarry Smith M*/ 2310f99b6f4SBarry Smith 232fdc842d1SBarry Smith /*@ 233dce8aebaSBarry Smith DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayWrite()` 234fdc842d1SBarry Smith 23520f4b53cSBarry Smith Logically Collective 236fdc842d1SBarry Smith 237d8d19677SJose E. Roman Input Parameters: 238fdc842d1SBarry Smith + da - the distributed array 239*0af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 240*0af9b551SBarry Smith - array - the `array` pointer is zeroed 241fdc842d1SBarry Smith 242fdc842d1SBarry Smith Level: intermediate 243fdc842d1SBarry Smith 244dce8aebaSBarry Smith Fortran Note: 2450f99b6f4SBarry Smith Use `DMDAVecRestoreArayWriteF90()` 246fdc842d1SBarry Smith 247dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`, 248db781477SPatrick Sanan `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 249fdc842d1SBarry Smith @*/ 250d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array) 251d71ae5a4SJacob Faibussowitsch { 252fdc842d1SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 253fdc842d1SBarry Smith 254fdc842d1SBarry Smith PetscFunctionBegin; 255fdc842d1SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 256fdc842d1SBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 257fdc842d1SBarry Smith PetscValidPointer(array, 3); 2581bb6d2a8SBarry Smith if (da->localSection) { 2599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec, (PetscScalar **)array)); 2603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 261fdc842d1SBarry Smith } 2629566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 2639566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 2649566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 265fdc842d1SBarry Smith 266fdc842d1SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 2679566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 268fdc842d1SBarry Smith if (N == xm * ym * zm * dof) { 269fdc842d1SBarry Smith gxm = xm; 270fdc842d1SBarry Smith gym = ym; 271fdc842d1SBarry Smith gzm = zm; 272fdc842d1SBarry Smith gxs = xs; 273fdc842d1SBarry Smith gys = ys; 274fdc842d1SBarry Smith gzs = zs; 27563a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof); 276fdc842d1SBarry Smith 277fdc842d1SBarry Smith if (dim == 1) { 2789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 279fdc842d1SBarry Smith } else if (dim == 2) { 2809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 281fdc842d1SBarry Smith } else if (dim == 3) { 2829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 28363a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 2843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 285fdc842d1SBarry Smith } 286fdc842d1SBarry Smith 287fdc842d1SBarry Smith /*@C 288aa219208SBarry Smith DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with 28947c6ae99SBarry Smith the underlying vector and is indexed using the global dimensions. 29047c6ae99SBarry Smith 29120f4b53cSBarry Smith Logically Collective 29247c6ae99SBarry Smith 293d8d19677SJose E. Roman Input Parameters: 29447c6ae99SBarry Smith + da - the distributed array 295*0af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 29647c6ae99SBarry Smith 29747c6ae99SBarry Smith Output Parameter: 298*0af9b551SBarry Smith . array - the `array` pointer is zeroed 29947c6ae99SBarry Smith 300dce8aebaSBarry Smith Level: intermediate 301dce8aebaSBarry Smith 30247c6ae99SBarry Smith Notes: 303dce8aebaSBarry Smith Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries. 30447c6ae99SBarry Smith 305*0af9b551SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF] 306*0af9b551SBarry Smith 307*0af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1][0:ndof-1]` where the values are obtained from 308*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 30947c6ae99SBarry Smith 3100f99b6f4SBarry Smith Fortran Notes: 3110f99b6f4SBarry Smith Use `DMDAVecGetArrayF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 3120f99b6f4SBarry Smith dimension. For a `DMDA` created with a dof of 1 use the dimension of the `DMDA`, for a `DMDA` created with a dof greater than 1 use one more than the 3130f99b6f4SBarry Smith dimension of the `DMDA`. 3140f99b6f4SBarry Smith 315*0af9b551SBarry Smith The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when ndof is 1) otherwise 316*0af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 317*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 3181b82215eSBarry Smith 319dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`, 3204ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()` 32147c6ae99SBarry Smith @*/ 322d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array) 323d71ae5a4SJacob Faibussowitsch { 32447c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 32547c6ae99SBarry Smith 32647c6ae99SBarry Smith PetscFunctionBegin; 3279566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 3289566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 3299566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 33047c6ae99SBarry Smith 33147c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 3329566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 33347c6ae99SBarry Smith if (N == xm * ym * zm * dof) { 33447c6ae99SBarry Smith gxm = xm; 33547c6ae99SBarry Smith gym = ym; 33647c6ae99SBarry Smith gzm = zm; 33747c6ae99SBarry Smith gxs = xs; 33847c6ae99SBarry Smith gys = ys; 33947c6ae99SBarry Smith gzs = zs; 34063a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof); 34147c6ae99SBarry Smith 34247c6ae99SBarry Smith if (dim == 1) { 3439566063dSJacob Faibussowitsch PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 34447c6ae99SBarry Smith } else if (dim == 2) { 3459566063dSJacob Faibussowitsch PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 34647c6ae99SBarry Smith } else if (dim == 3) { 3479566063dSJacob Faibussowitsch PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 34863a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 3493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35047c6ae99SBarry Smith } 35147c6ae99SBarry Smith 35247c6ae99SBarry Smith /*@ 353dce8aebaSBarry Smith DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()` 35447c6ae99SBarry Smith 35520f4b53cSBarry Smith Logically Collective 35647c6ae99SBarry Smith 357d8d19677SJose E. Roman Input Parameters: 35847c6ae99SBarry Smith + da - the distributed array 359*0af9b551SBarry Smith . vec - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 360*0af9b551SBarry Smith - array - the `array` point is zeroed 36147c6ae99SBarry Smith 36247c6ae99SBarry Smith Level: intermediate 36347c6ae99SBarry Smith 3640f99b6f4SBarry Smith Fortran Note: 3650f99b6f4SBarry Smith Use `DMDAVecRestoreArayF90()` 3660f99b6f4SBarry Smith 367dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`, 3684ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 36947c6ae99SBarry Smith @*/ 370d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array) 371d71ae5a4SJacob Faibussowitsch { 37247c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 37347c6ae99SBarry Smith 37447c6ae99SBarry Smith PetscFunctionBegin; 3759566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 3769566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 3779566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 37847c6ae99SBarry Smith 37947c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 3809566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 38147c6ae99SBarry Smith if (N == xm * ym * zm * dof) { 38247c6ae99SBarry Smith gxm = xm; 38347c6ae99SBarry Smith gym = ym; 38447c6ae99SBarry Smith gzm = zm; 38547c6ae99SBarry Smith gxs = xs; 38647c6ae99SBarry Smith gys = ys; 38747c6ae99SBarry Smith gzs = zs; 38847c6ae99SBarry Smith } 38947c6ae99SBarry Smith 39047c6ae99SBarry Smith if (dim == 1) { 3919566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 39247c6ae99SBarry Smith } else if (dim == 2) { 3939566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 39447c6ae99SBarry Smith } else if (dim == 3) { 3959566063dSJacob Faibussowitsch PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 39663a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 3973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39847c6ae99SBarry Smith } 39947c6ae99SBarry Smith 4000f99b6f4SBarry Smith /*MC 4010f99b6f4SBarry Smith DMDAVecGetArrayReadF90 - check Fortran Notes at `DMDAVecGetArrayRead()` 40253c0d4aeSBarry Smith 40353c0d4aeSBarry Smith Level: intermediate 4040f99b6f4SBarry Smith M*/ 4050f99b6f4SBarry Smith 4065edff71fSBarry Smith /*@C 4075edff71fSBarry Smith DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with 4085edff71fSBarry Smith the underlying vector and is indexed using the global dimensions. 4095edff71fSBarry Smith 41020f4b53cSBarry Smith Not Collective 4115edff71fSBarry Smith 412d8d19677SJose E. Roman Input Parameters: 4135edff71fSBarry Smith + da - the distributed array 414*0af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 4155edff71fSBarry Smith 4165edff71fSBarry Smith Output Parameter: 4175edff71fSBarry Smith . array - the array 4185edff71fSBarry Smith 419dce8aebaSBarry Smith Level: intermediate 420dce8aebaSBarry Smith 4215edff71fSBarry Smith Notes: 422dce8aebaSBarry Smith Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries. 4235edff71fSBarry Smith 4245edff71fSBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 4255edff71fSBarry Smith 426*0af9b551SBarry Smith If `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is 4275edff71fSBarry Smith a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the 428dce8aebaSBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 4295edff71fSBarry Smith 430*0af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 431*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 432*0af9b551SBarry Smith 43395452b02SPatrick Sanan Fortran Notes: 4340f99b6f4SBarry Smith Use `DMDAVecGetArrayReadF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate 435dce8aebaSBarry Smith dimension. For a `DMDA` created with a dof of 1 use the dimension of the `DMDA`, for a `DMDA` created with a dof greater than 1 use one more than the 4360f99b6f4SBarry Smith dimension of the `DMDA`. 4370f99b6f4SBarry Smith 438*0af9b551SBarry Smith The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise 439*0af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 440*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 4415edff71fSBarry Smith 4420f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAVecGetArrayReadF90()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`, `DMDAVecRestoreArrayDOF()` 443db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 444db781477SPatrick Sanan `DMStagVecGetArrayRead()` 4455edff71fSBarry Smith @*/ 446d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayRead(DM da, Vec vec, void *array) 447d71ae5a4SJacob Faibussowitsch { 4485edff71fSBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 4495edff71fSBarry Smith 4505edff71fSBarry Smith PetscFunctionBegin; 451a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 4525edff71fSBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 4535edff71fSBarry Smith PetscValidPointer(array, 3); 4549566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 4559566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 4569566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 4575edff71fSBarry Smith 4585edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */ 4599566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 4605edff71fSBarry Smith if (N == xm * ym * zm * dof) { 4615edff71fSBarry Smith gxm = xm; 4625edff71fSBarry Smith gym = ym; 4635edff71fSBarry Smith gzm = zm; 4645edff71fSBarry Smith gxs = xs; 4655edff71fSBarry Smith gys = ys; 4665edff71fSBarry Smith gzs = zs; 46763a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof); 4685edff71fSBarry Smith 4695edff71fSBarry Smith if (dim == 1) { 4709566063dSJacob Faibussowitsch PetscCall(VecGetArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 4715edff71fSBarry Smith } else if (dim == 2) { 4729566063dSJacob Faibussowitsch PetscCall(VecGetArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 4735edff71fSBarry Smith } else if (dim == 3) { 4749566063dSJacob Faibussowitsch PetscCall(VecGetArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 47563a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 4763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4775edff71fSBarry Smith } 4785edff71fSBarry Smith 4790f99b6f4SBarry Smith /*MC 4800f99b6f4SBarry Smith DMDAVecRestoreArrayReadF90 - check Fortran Notes at `DMDAVecRestoreArrayRead()` 48153c0d4aeSBarry Smith 48253c0d4aeSBarry Smith Level: intermediate 4830f99b6f4SBarry Smith M*/ 4840f99b6f4SBarry Smith 4855edff71fSBarry Smith /*@ 486dce8aebaSBarry Smith DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayRead()` 4875edff71fSBarry Smith 48820f4b53cSBarry Smith Not Collective 4895edff71fSBarry Smith 490d8d19677SJose E. Roman Input Parameters: 4915edff71fSBarry Smith + da - the distributed array 492*0af9b551SBarry Smith . vec - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 493*0af9b551SBarry Smith - array - the `array` pointer is zeroed 4945edff71fSBarry Smith 4955edff71fSBarry Smith Level: intermediate 4965edff71fSBarry Smith 497dce8aebaSBarry Smith Fortran Note: 4980f99b6f4SBarry Smith Use `DMDAVecRestoreArrayReadF90()` 4995edff71fSBarry Smith 5000f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAVecRestoreArrayReadF90()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayRead()`, 501db781477SPatrick Sanan `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, 502db781477SPatrick Sanan `DMStagVecRestoreArrayRead()` 5035edff71fSBarry Smith @*/ 504d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayRead(DM da, Vec vec, void *array) 505d71ae5a4SJacob Faibussowitsch { 5065edff71fSBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 5075edff71fSBarry Smith 5085edff71fSBarry Smith PetscFunctionBegin; 509a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 5105edff71fSBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 5115edff71fSBarry Smith PetscValidPointer(array, 3); 5129566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 5139566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 5149566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 5155edff71fSBarry Smith 5165edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */ 5179566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 5185edff71fSBarry Smith if (N == xm * ym * zm * dof) { 5195edff71fSBarry Smith gxm = xm; 5205edff71fSBarry Smith gym = ym; 5215edff71fSBarry Smith gzm = zm; 5225edff71fSBarry Smith gxs = xs; 5235edff71fSBarry Smith gys = ys; 5245edff71fSBarry Smith gzs = zs; 52563a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof); 5265edff71fSBarry Smith 5275edff71fSBarry Smith if (dim == 1) { 5289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 5295edff71fSBarry Smith } else if (dim == 2) { 5309566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 5315edff71fSBarry Smith } else if (dim == 3) { 5329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 53363a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 5343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5355edff71fSBarry Smith } 5365edff71fSBarry Smith 5375edff71fSBarry Smith /*@C 5385edff71fSBarry Smith DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with 5395edff71fSBarry Smith the underlying vector and is indexed using the global dimensions. 5405edff71fSBarry Smith 5415edff71fSBarry Smith Not Collective 5425edff71fSBarry Smith 543d8d19677SJose E. Roman Input Parameters: 5445edff71fSBarry Smith + da - the distributed array 545*0af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 5465edff71fSBarry Smith 5475edff71fSBarry Smith Output Parameter: 5485edff71fSBarry Smith . array - the array 5495edff71fSBarry Smith 550dce8aebaSBarry Smith Level: intermediate 551dce8aebaSBarry Smith 5525edff71fSBarry Smith Notes: 553dce8aebaSBarry Smith Call `DMDAVecRestoreArrayDOFRead()` once you have finished accessing the vector entries. 5545edff71fSBarry Smith 5555edff71fSBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]! 5565edff71fSBarry Smith 557*0af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 558*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 559*0af9b551SBarry Smith 560dce8aebaSBarry Smith Fortran Note: 5610f99b6f4SBarry Smith Use `DMDAVecGetArrayReadF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 5620f99b6f4SBarry Smith dimension. For a `DMDA` created with a dof of 1 use the dimension of the `DMDA`, for a `DMDA` created with a dof greater than 1 use one more than the 5630f99b6f4SBarry Smith dimension of the `DMDA`. 5640f99b6f4SBarry Smith 565*0af9b551SBarry Smith The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise 566*0af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 567*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 5685edff71fSBarry Smith 569dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, 5704ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 5715edff71fSBarry Smith @*/ 572d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array) 573d71ae5a4SJacob Faibussowitsch { 5745edff71fSBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 5755edff71fSBarry Smith 5765edff71fSBarry Smith PetscFunctionBegin; 5779566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 5789566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 5799566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 5805edff71fSBarry Smith 5815edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */ 5829566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 5835edff71fSBarry Smith if (N == xm * ym * zm * dof) { 5845edff71fSBarry Smith gxm = xm; 5855edff71fSBarry Smith gym = ym; 5865edff71fSBarry Smith gzm = zm; 5875edff71fSBarry Smith gxs = xs; 5885edff71fSBarry Smith gys = ys; 5895edff71fSBarry Smith gzs = zs; 59063a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof); 5915edff71fSBarry Smith 5925edff71fSBarry Smith if (dim == 1) { 5939566063dSJacob Faibussowitsch PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 5945edff71fSBarry Smith } else if (dim == 2) { 5959566063dSJacob Faibussowitsch PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 5965edff71fSBarry Smith } else if (dim == 3) { 5979566063dSJacob Faibussowitsch PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 59863a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 5993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6005edff71fSBarry Smith } 6015edff71fSBarry Smith 6025edff71fSBarry Smith /*@ 603dce8aebaSBarry Smith DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFRead()` 6045edff71fSBarry Smith 6055edff71fSBarry Smith Not Collective 6065edff71fSBarry Smith 607d8d19677SJose E. Roman Input Parameters: 6085edff71fSBarry Smith + da - the distributed array 609*0af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 610*0af9b551SBarry Smith - array - the `array` pointer is zeroed 6115edff71fSBarry Smith 6125edff71fSBarry Smith Level: intermediate 6135edff71fSBarry Smith 6140f99b6f4SBarry Smith Fortran Note: 6150f99b6f4SBarry Smith Use `DMDAVecRestoreArrayReadF90()` 6160f99b6f4SBarry Smith 617dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`, 6184ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 6195edff71fSBarry Smith @*/ 620d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array) 621d71ae5a4SJacob Faibussowitsch { 6225edff71fSBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 6235edff71fSBarry Smith 6245edff71fSBarry Smith PetscFunctionBegin; 6259566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 6269566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 6279566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 6285edff71fSBarry Smith 6295edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */ 6309566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 6315edff71fSBarry Smith if (N == xm * ym * zm * dof) { 6325edff71fSBarry Smith gxm = xm; 6335edff71fSBarry Smith gym = ym; 6345edff71fSBarry Smith gzm = zm; 6355edff71fSBarry Smith gxs = xs; 6365edff71fSBarry Smith gys = ys; 6375edff71fSBarry Smith gzs = zs; 6385edff71fSBarry Smith } 6395edff71fSBarry Smith 6405edff71fSBarry Smith if (dim == 1) { 6419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 6425edff71fSBarry Smith } else if (dim == 2) { 6439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 6445edff71fSBarry Smith } else if (dim == 3) { 6459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 64663a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 6473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6485edff71fSBarry Smith } 6495edff71fSBarry Smith 6501e5d2365SBarry Smith /*@C 6511e5d2365SBarry Smith DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with 6521e5d2365SBarry Smith the underlying vector and is indexed using the global dimensions. 6531e5d2365SBarry Smith 6541e5d2365SBarry Smith Not Collective 6551e5d2365SBarry Smith 656d8d19677SJose E. Roman Input Parameters: 6571e5d2365SBarry Smith + da - the distributed array 658*0af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 6591e5d2365SBarry Smith 6601e5d2365SBarry Smith Output Parameter: 6611e5d2365SBarry Smith . array - the array 6621e5d2365SBarry Smith 6630f99b6f4SBarry Smith Level: intermediate 6640f99b6f4SBarry Smith 6651e5d2365SBarry Smith Notes: 666dce8aebaSBarry Smith Call `DMDAVecRestoreArrayDOFWrite()` once you have finished accessing the vector entries. 6671e5d2365SBarry Smith 6681e5d2365SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]! 6691e5d2365SBarry Smith 670*0af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1][0:dof-1]` where the values are obtained from 671*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 672*0af9b551SBarry Smith 6730f99b6f4SBarry Smith Fortran Note: 6740f99b6f4SBarry Smith Use `DMDAVecGetArrayWriteF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 6750f99b6f4SBarry Smith dimension. For a `DMDA` created with a dof of 1 use the dimension of the `DMDA`, for a `DMDA` created with a dof greater than 1 use one more than the 6760f99b6f4SBarry Smith dimension of the `DMDA`. 6771e5d2365SBarry Smith 678*0af9b551SBarry Smith The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise 679*0af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 680*0af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 6811e5d2365SBarry Smith 6820f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArrayWriteF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, 6834ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()` 6841e5d2365SBarry Smith @*/ 685d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array) 686d71ae5a4SJacob Faibussowitsch { 6871e5d2365SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 6881e5d2365SBarry Smith 6891e5d2365SBarry Smith PetscFunctionBegin; 6909566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 6919566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 6929566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 6931e5d2365SBarry Smith 6941e5d2365SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 6959566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 6961e5d2365SBarry Smith if (N == xm * ym * zm * dof) { 6971e5d2365SBarry Smith gxm = xm; 6981e5d2365SBarry Smith gym = ym; 6991e5d2365SBarry Smith gzm = zm; 7001e5d2365SBarry Smith gxs = xs; 7011e5d2365SBarry Smith gys = ys; 7021e5d2365SBarry Smith gzs = zs; 70363a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof); 7041e5d2365SBarry Smith 7051e5d2365SBarry Smith if (dim == 1) { 7069566063dSJacob Faibussowitsch PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 7071e5d2365SBarry Smith } else if (dim == 2) { 7089566063dSJacob Faibussowitsch PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 7091e5d2365SBarry Smith } else if (dim == 3) { 7109566063dSJacob Faibussowitsch PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 71163a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 7123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7131e5d2365SBarry Smith } 7141e5d2365SBarry Smith 7151e5d2365SBarry Smith /*@ 716dce8aebaSBarry Smith DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFWrite()` 7171e5d2365SBarry Smith 7181e5d2365SBarry Smith Not Collective 7191e5d2365SBarry Smith 720d8d19677SJose E. Roman Input Parameters: 7211e5d2365SBarry Smith + da - the distributed array 722*0af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 723*0af9b551SBarry Smith - array - the `array` pointer is zeroed 7241e5d2365SBarry Smith 7251e5d2365SBarry Smith Level: intermediate 7261e5d2365SBarry Smith 7270f99b6f4SBarry Smith Fortran Note: 7280f99b6f4SBarry Smith Use `DMDAVecRestoreArrayWriteF90()` 7290f99b6f4SBarry Smith 730dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`, 7314ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()` 7321e5d2365SBarry Smith @*/ 733d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array) 734d71ae5a4SJacob Faibussowitsch { 7351e5d2365SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 7361e5d2365SBarry Smith 7371e5d2365SBarry Smith PetscFunctionBegin; 7389566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 7399566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 7409566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 7411e5d2365SBarry Smith 7421e5d2365SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 7439566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 7441e5d2365SBarry Smith if (N == xm * ym * zm * dof) { 7451e5d2365SBarry Smith gxm = xm; 7461e5d2365SBarry Smith gym = ym; 7471e5d2365SBarry Smith gzm = zm; 7481e5d2365SBarry Smith gxs = xs; 7491e5d2365SBarry Smith gys = ys; 7501e5d2365SBarry Smith gzs = zs; 7511e5d2365SBarry Smith } 7521e5d2365SBarry Smith 7531e5d2365SBarry Smith if (dim == 1) { 7549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 7551e5d2365SBarry Smith } else if (dim == 2) { 7569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 7571e5d2365SBarry Smith } else if (dim == 3) { 7589566063dSJacob Faibussowitsch PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 75963a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 7603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7611e5d2365SBarry Smith } 762