1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 247c6ae99SBarry Smith 30f99b6f4SBarry Smith /*MC 412b4a537SBarry Smith DMDAVecGetArrayF90 - see Fortran Notes at `DMDAVecGetArray()` 571b416bcSMatthew G. Knepley 671b416bcSMatthew G. Knepley Level: intermediate 70f99b6f4SBarry Smith M*/ 80f99b6f4SBarry Smith 947c6ae99SBarry Smith /*@C 10aa219208SBarry Smith DMDAVecGetArray - Returns a multiple dimension array that shares data with 1112b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA`. 1247c6ae99SBarry Smith 1320f4b53cSBarry Smith Logically Collective 1447c6ae99SBarry Smith 15d8d19677SJose E. Roman Input Parameters: 1672fd0fbdSBarry Smith + da - the `DMDA` 170af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 1847c6ae99SBarry Smith 1947c6ae99SBarry Smith Output Parameter: 2047c6ae99SBarry Smith . array - the array 2147c6ae99SBarry Smith 22dce8aebaSBarry Smith Level: intermediate 23dce8aebaSBarry Smith 2447c6ae99SBarry Smith Notes: 25dce8aebaSBarry Smith Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries. 2647c6ae99SBarry Smith 2747c6ae99SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 2847c6ae99SBarry Smith 290af9b551SBarry Smith If `vec` is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is 300af9b551SBarry Smith a global vector then the ghost points are not accessible. Of course, with a local vector you will have had to do the 31dce8aebaSBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 3247c6ae99SBarry Smith 330af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 340af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 350af9b551SBarry Smith 3695452b02SPatrick Sanan Fortran Notes: 370f99b6f4SBarry Smith Use `DMDAVecGetArrayF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate 38dce8aebaSBarry 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 390f99b6f4SBarry Smith dimension of the `DMDA`. 4047c6ae99SBarry Smith 410af9b551SBarry 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 420af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 430af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 4447c6ae99SBarry Smith 4512b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()` 46db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 47db781477SPatrick Sanan `DMStagVecGetArray()` 4847c6ae99SBarry Smith @*/ 49d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArray(DM da, Vec vec, void *array) 50d71ae5a4SJacob Faibussowitsch { 5147c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 5247c6ae99SBarry Smith 5347c6ae99SBarry Smith PetscFunctionBegin; 54a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 556db82c94SMatthew G Knepley PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 564f572ea9SToby Isaac PetscAssertPointer(array, 3); 579566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 589566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 599566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 6047c6ae99SBarry Smith 6147c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 629566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 6347c6ae99SBarry Smith if (N == xm * ym * zm * dof) { 6447c6ae99SBarry Smith gxm = xm; 6547c6ae99SBarry Smith gym = ym; 6647c6ae99SBarry Smith gzm = zm; 6747c6ae99SBarry Smith gxs = xs; 6847c6ae99SBarry Smith gys = ys; 6947c6ae99SBarry Smith gzs = zs; 7063a3b9bcSJacob 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); 7147c6ae99SBarry Smith 7247c6ae99SBarry Smith if (dim == 1) { 739566063dSJacob Faibussowitsch PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 7447c6ae99SBarry Smith } else if (dim == 2) { 759566063dSJacob Faibussowitsch PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 7647c6ae99SBarry Smith } else if (dim == 3) { 779566063dSJacob Faibussowitsch PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 7863a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8047c6ae99SBarry Smith } 8147c6ae99SBarry Smith 820f99b6f4SBarry Smith /*MC 8312b4a537SBarry Smith DMDAVecRestoreArrayF90 - see Fortran Notes at `DMDAVecRestoreArray()` 8453c0d4aeSBarry Smith 8553c0d4aeSBarry Smith Level: intermediate 860f99b6f4SBarry Smith M*/ 870f99b6f4SBarry Smith 8847c6ae99SBarry Smith /*@ 89dce8aebaSBarry Smith DMDAVecRestoreArray - Restores a multiple dimension array obtained with `DMDAVecGetArray()` 9047c6ae99SBarry Smith 9120f4b53cSBarry Smith Logically Collective 9247c6ae99SBarry Smith 93d8d19677SJose E. Roman Input Parameters: 9472fd0fbdSBarry Smith + da - the `DMDA` 950af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 9612b4a537SBarry Smith - array - the `array` pointer 9747c6ae99SBarry Smith 9847c6ae99SBarry Smith Level: intermediate 9947c6ae99SBarry Smith 10012b4a537SBarry Smith Fortran Note: 101*f8d70eaaSPierre Jolivet Use `DMDAVecRestoreArrayF90()` 10247c6ae99SBarry Smith 103*f8d70eaaSPierre Jolivet .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArrayF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, 104db781477SPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 105*f8d70eaaSPierre Jolivet `DMStagVecRestoreArray()` 10647c6ae99SBarry Smith @*/ 107d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array) 108d71ae5a4SJacob Faibussowitsch { 10947c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 11047c6ae99SBarry Smith 11147c6ae99SBarry Smith PetscFunctionBegin; 112a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 1136db82c94SMatthew G Knepley PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 1144f572ea9SToby Isaac PetscAssertPointer(array, 3); 1159566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 1169566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 1179566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 11847c6ae99SBarry Smith 11947c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 1209566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 12147c6ae99SBarry Smith if (N == xm * ym * zm * dof) { 12247c6ae99SBarry Smith gxm = xm; 12347c6ae99SBarry Smith gym = ym; 12447c6ae99SBarry Smith gzm = zm; 12547c6ae99SBarry Smith gxs = xs; 12647c6ae99SBarry Smith gys = ys; 12747c6ae99SBarry Smith gzs = zs; 12863a3b9bcSJacob 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); 12947c6ae99SBarry Smith 13047c6ae99SBarry Smith if (dim == 1) { 1319566063dSJacob Faibussowitsch PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 13247c6ae99SBarry Smith } else if (dim == 2) { 1339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 13447c6ae99SBarry Smith } else if (dim == 3) { 1359566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 13663a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 1373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13847c6ae99SBarry Smith } 13947c6ae99SBarry Smith 1400f99b6f4SBarry Smith /*MC 14112b4a537SBarry Smith DMDAVecGetArrayWriteF90 - see Fortran Notes at `DMDAVecGetArrayWrite()` 14253c0d4aeSBarry Smith 14353c0d4aeSBarry Smith Level: intermediate 1440f99b6f4SBarry Smith M*/ 1450f99b6f4SBarry Smith 14647c6ae99SBarry Smith /*@C 147fdc842d1SBarry Smith DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with 14812b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA`. 149fdc842d1SBarry Smith 15020f4b53cSBarry Smith Logically Collective 151fdc842d1SBarry Smith 152d8d19677SJose E. Roman Input Parameters: 15372fd0fbdSBarry Smith + da - the `DMDA` 1540af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 155fdc842d1SBarry Smith 156fdc842d1SBarry Smith Output Parameter: 157fdc842d1SBarry Smith . array - the array 158fdc842d1SBarry Smith 159dce8aebaSBarry Smith Level: intermediate 160dce8aebaSBarry Smith 161fdc842d1SBarry Smith Notes: 162dce8aebaSBarry Smith Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries. 163fdc842d1SBarry Smith 164fdc842d1SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 165fdc842d1SBarry Smith 1660af9b551SBarry Smith if `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is 167fdc842d1SBarry Smith a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the 168dce8aebaSBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 169fdc842d1SBarry Smith 1700af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 1710af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 1720af9b551SBarry Smith 173fdc842d1SBarry Smith Fortran Notes: 17412b4a537SBarry Smith Use `DMDAVecGetArrayWriteF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 175dce8aebaSBarry 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 1760f99b6f4SBarry Smith dimension of the `DMDA`. 177fdc842d1SBarry Smith 1780af9b551SBarry 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 1790af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 1800af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 181fdc842d1SBarry Smith 18212b4a537SBarry Smith Developer Note: 183dce8aebaSBarry Smith This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()` 184fdc842d1SBarry Smith 18512b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecRestoreArrayDOF()` 186db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 187fdc842d1SBarry Smith @*/ 188d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayWrite(DM da, Vec vec, void *array) 189d71ae5a4SJacob Faibussowitsch { 190fdc842d1SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 191fdc842d1SBarry Smith 192fdc842d1SBarry Smith PetscFunctionBegin; 193fdc842d1SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 194fdc842d1SBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 1954f572ea9SToby Isaac PetscAssertPointer(array, 3); 1961bb6d2a8SBarry Smith if (da->localSection) { 1979566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array)); 1983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 199fdc842d1SBarry Smith } 2009566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 2019566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 2029566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 203fdc842d1SBarry Smith 204fdc842d1SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 2059566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 206fdc842d1SBarry Smith if (N == xm * ym * zm * dof) { 207fdc842d1SBarry Smith gxm = xm; 208fdc842d1SBarry Smith gym = ym; 209fdc842d1SBarry Smith gzm = zm; 210fdc842d1SBarry Smith gxs = xs; 211fdc842d1SBarry Smith gys = ys; 212fdc842d1SBarry Smith gzs = zs; 21363a3b9bcSJacob 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); 214fdc842d1SBarry Smith 215fdc842d1SBarry Smith if (dim == 1) { 2169566063dSJacob Faibussowitsch PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 217fdc842d1SBarry Smith } else if (dim == 2) { 2189566063dSJacob Faibussowitsch PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 219fdc842d1SBarry Smith } else if (dim == 3) { 2209566063dSJacob Faibussowitsch PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 22163a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 2223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 223fdc842d1SBarry Smith } 224fdc842d1SBarry Smith 2250f99b6f4SBarry Smith /*MC 22612b4a537SBarry Smith DMDAVecRestoreArrayWriteF90 - see Fortran Notes at `DMDAVecRestoreArrayWrite()` 22753c0d4aeSBarry Smith 22853c0d4aeSBarry Smith Level: intermediate 2290f99b6f4SBarry Smith M*/ 2300f99b6f4SBarry Smith 231fdc842d1SBarry Smith /*@ 232dce8aebaSBarry Smith DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayWrite()` 233fdc842d1SBarry Smith 23420f4b53cSBarry Smith Logically Collective 235fdc842d1SBarry Smith 236d8d19677SJose E. Roman Input Parameters: 23772fd0fbdSBarry Smith + da - the `DMDA` 2380af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 23912b4a537SBarry Smith - array - the `array` pointer 240fdc842d1SBarry Smith 241fdc842d1SBarry Smith Level: intermediate 242fdc842d1SBarry Smith 24312b4a537SBarry Smith Fortran Note: 2440f99b6f4SBarry Smith Use `DMDAVecRestoreArayWriteF90()` 245fdc842d1SBarry Smith 24612b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`, 247db781477SPatrick Sanan `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 248fdc842d1SBarry Smith @*/ 249d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array) 250d71ae5a4SJacob Faibussowitsch { 251fdc842d1SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 252fdc842d1SBarry Smith 253fdc842d1SBarry Smith PetscFunctionBegin; 254fdc842d1SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 255fdc842d1SBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 2564f572ea9SToby Isaac PetscAssertPointer(array, 3); 2571bb6d2a8SBarry Smith if (da->localSection) { 2589566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec, (PetscScalar **)array)); 2593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 260fdc842d1SBarry Smith } 2619566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 2629566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 2639566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 264fdc842d1SBarry Smith 265fdc842d1SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 2669566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 267fdc842d1SBarry Smith if (N == xm * ym * zm * dof) { 268fdc842d1SBarry Smith gxm = xm; 269fdc842d1SBarry Smith gym = ym; 270fdc842d1SBarry Smith gzm = zm; 271fdc842d1SBarry Smith gxs = xs; 272fdc842d1SBarry Smith gys = ys; 273fdc842d1SBarry Smith gzs = zs; 27463a3b9bcSJacob 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); 275fdc842d1SBarry Smith 276fdc842d1SBarry Smith if (dim == 1) { 2779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 278fdc842d1SBarry Smith } else if (dim == 2) { 2799566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 280fdc842d1SBarry Smith } else if (dim == 3) { 2819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 28263a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 284fdc842d1SBarry Smith } 285fdc842d1SBarry Smith 286fdc842d1SBarry Smith /*@C 287aa219208SBarry Smith DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with 28812b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA` 28947c6ae99SBarry Smith 29020f4b53cSBarry Smith Logically Collective 29147c6ae99SBarry Smith 292d8d19677SJose E. Roman Input Parameters: 29372fd0fbdSBarry Smith + da - the `DMDA` 2940af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 29547c6ae99SBarry Smith 29647c6ae99SBarry Smith Output Parameter: 29712b4a537SBarry Smith . array - the `array` pointer 29847c6ae99SBarry Smith 299dce8aebaSBarry Smith Level: intermediate 300dce8aebaSBarry Smith 30147c6ae99SBarry Smith Notes: 302dce8aebaSBarry Smith Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries. 30347c6ae99SBarry Smith 3040af9b551SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF] 3050af9b551SBarry Smith 3060af9b551SBarry 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 3070af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 30847c6ae99SBarry Smith 3090f99b6f4SBarry Smith Fortran Notes: 3100f99b6f4SBarry Smith Use `DMDAVecGetArrayF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 3110f99b6f4SBarry 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 3120f99b6f4SBarry Smith dimension of the `DMDA`. 3130f99b6f4SBarry Smith 3140af9b551SBarry 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 3150af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 3160af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 3171b82215eSBarry Smith 31812b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`, 3194ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()` 32047c6ae99SBarry Smith @*/ 321d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array) 322d71ae5a4SJacob Faibussowitsch { 32347c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 32447c6ae99SBarry Smith 32547c6ae99SBarry Smith PetscFunctionBegin; 3269566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 3279566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 3289566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 32947c6ae99SBarry Smith 33047c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 3319566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 33247c6ae99SBarry Smith if (N == xm * ym * zm * dof) { 33347c6ae99SBarry Smith gxm = xm; 33447c6ae99SBarry Smith gym = ym; 33547c6ae99SBarry Smith gzm = zm; 33647c6ae99SBarry Smith gxs = xs; 33747c6ae99SBarry Smith gys = ys; 33847c6ae99SBarry Smith gzs = zs; 33963a3b9bcSJacob 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); 34047c6ae99SBarry Smith 34147c6ae99SBarry Smith if (dim == 1) { 3429566063dSJacob Faibussowitsch PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 34347c6ae99SBarry Smith } else if (dim == 2) { 3449566063dSJacob Faibussowitsch PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 34547c6ae99SBarry Smith } else if (dim == 3) { 3469566063dSJacob Faibussowitsch PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 34763a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 3483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34947c6ae99SBarry Smith } 35047c6ae99SBarry Smith 35147c6ae99SBarry Smith /*@ 352dce8aebaSBarry Smith DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()` 35347c6ae99SBarry Smith 35420f4b53cSBarry Smith Logically Collective 35547c6ae99SBarry Smith 356d8d19677SJose E. Roman Input Parameters: 35772fd0fbdSBarry Smith + da - the `DMDA` 3580af9b551SBarry Smith . vec - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 35912b4a537SBarry Smith - array - the `array` point 36047c6ae99SBarry Smith 36147c6ae99SBarry Smith Level: intermediate 36247c6ae99SBarry Smith 36312b4a537SBarry Smith Fortran Note: 364*f8d70eaaSPierre Jolivet Use `DMDAVecRestoreArrayF90()` 3650f99b6f4SBarry Smith 36612b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, 3674ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 36847c6ae99SBarry Smith @*/ 369d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array) 370d71ae5a4SJacob Faibussowitsch { 37147c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 37247c6ae99SBarry Smith 37347c6ae99SBarry Smith PetscFunctionBegin; 3749566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 3759566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 3769566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 37747c6ae99SBarry Smith 37847c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 3799566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N)); 38047c6ae99SBarry Smith if (N == xm * ym * zm * dof) { 38147c6ae99SBarry Smith gxm = xm; 38247c6ae99SBarry Smith gym = ym; 38347c6ae99SBarry Smith gzm = zm; 38447c6ae99SBarry Smith gxs = xs; 38547c6ae99SBarry Smith gys = ys; 38647c6ae99SBarry Smith gzs = zs; 38747c6ae99SBarry Smith } 38847c6ae99SBarry Smith 38947c6ae99SBarry Smith if (dim == 1) { 3909566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 39147c6ae99SBarry Smith } else if (dim == 2) { 3929566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 39347c6ae99SBarry Smith } else if (dim == 3) { 3949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 39563a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 3963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39747c6ae99SBarry Smith } 39847c6ae99SBarry Smith 3990f99b6f4SBarry Smith /*MC 40012b4a537SBarry Smith DMDAVecGetArrayReadF90 - see Fortran Notes at `DMDAVecGetArrayRead()` 40153c0d4aeSBarry Smith 40253c0d4aeSBarry Smith Level: intermediate 4030f99b6f4SBarry Smith M*/ 4040f99b6f4SBarry Smith 4055edff71fSBarry Smith /*@C 4065edff71fSBarry Smith DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with 40712b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA`. 4085edff71fSBarry Smith 40920f4b53cSBarry Smith Not Collective 4105edff71fSBarry Smith 411d8d19677SJose E. Roman Input Parameters: 41272fd0fbdSBarry Smith + da - the `DMDA` 4130af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 4145edff71fSBarry Smith 4155edff71fSBarry Smith Output Parameter: 4165edff71fSBarry Smith . array - the array 4175edff71fSBarry Smith 418dce8aebaSBarry Smith Level: intermediate 419dce8aebaSBarry Smith 4205edff71fSBarry Smith Notes: 421dce8aebaSBarry Smith Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries. 4225edff71fSBarry Smith 4235edff71fSBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 4245edff71fSBarry Smith 4250af9b551SBarry Smith If `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is 4265edff71fSBarry Smith a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the 427dce8aebaSBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 4285edff71fSBarry Smith 4290af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 4300af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 4310af9b551SBarry Smith 43295452b02SPatrick Sanan Fortran Notes: 4330f99b6f4SBarry Smith Use `DMDAVecGetArrayReadF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate 434dce8aebaSBarry 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 4350f99b6f4SBarry Smith dimension of the `DMDA`. 4360f99b6f4SBarry Smith 4370af9b551SBarry 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 4380af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 4390af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 4405edff71fSBarry Smith 44112b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAVecGetArrayReadF90()`, `DMDAGetGhostCorners()`, 44242747ad1SJacob Faibussowitsch `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`, 44342747ad1SJacob Faibussowitsch `DMDAVecRestoreArrayDOF()`, `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, 44442747ad1SJacob Faibussowitsch `DMDAVecRestoreArray()`, `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); 4534f572ea9SToby Isaac PetscAssertPointer(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 48012b4a537SBarry Smith DMDAVecRestoreArrayReadF90 - see 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: 49172fd0fbdSBarry Smith + da - the `DMDA` 4920af9b551SBarry Smith . vec - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 49312b4a537SBarry Smith - array - the `array` pointer 4945edff71fSBarry Smith 4955edff71fSBarry Smith Level: intermediate 4965edff71fSBarry Smith 49712b4a537SBarry Smith Fortran Note: 4980f99b6f4SBarry Smith Use `DMDAVecRestoreArrayReadF90()` 4995edff71fSBarry Smith 50012b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `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); 5114f572ea9SToby Isaac PetscAssertPointer(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 53912b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA` 5405edff71fSBarry Smith 5415edff71fSBarry Smith Not Collective 5425edff71fSBarry Smith 543d8d19677SJose E. Roman Input Parameters: 54472fd0fbdSBarry Smith + da - the `DMDA` 5450af9b551SBarry 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 5570af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 5580af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 5590af9b551SBarry Smith 56060225df5SJacob Faibussowitsch Fortran Notes: 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 5650af9b551SBarry 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 5660af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 5670af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 5685edff71fSBarry Smith 56912b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `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: 60872fd0fbdSBarry Smith + da - the `DMDA` 6090af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 61012b4a537SBarry Smith - array - the `array` pointer 6115edff71fSBarry Smith 6125edff71fSBarry Smith Level: intermediate 6135edff71fSBarry Smith 61412b4a537SBarry Smith Fortran Note: 6150f99b6f4SBarry Smith Use `DMDAVecRestoreArrayReadF90()` 6160f99b6f4SBarry Smith 61712b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `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 65212b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA` 6531e5d2365SBarry Smith 6541e5d2365SBarry Smith Not Collective 6551e5d2365SBarry Smith 656d8d19677SJose E. Roman Input Parameters: 65772fd0fbdSBarry Smith + da - the `DMDA` 6580af9b551SBarry 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 6700af9b551SBarry 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 6710af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 6720af9b551SBarry Smith 67360225df5SJacob Faibussowitsch Fortran Notes: 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 6780af9b551SBarry 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 6790af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 6800af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 6811e5d2365SBarry Smith 68212b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArrayWriteF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, 68360225df5SJacob Faibussowitsch `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: 72172fd0fbdSBarry Smith + da - the `DMDA` 7220af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 72312b4a537SBarry Smith - array - the `array` pointer 7241e5d2365SBarry Smith 7251e5d2365SBarry Smith Level: intermediate 7261e5d2365SBarry Smith 72712b4a537SBarry Smith Fortran Note: 7280f99b6f4SBarry Smith Use `DMDAVecRestoreArrayWriteF90()` 7290f99b6f4SBarry Smith 73012b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`, 73160225df5SJacob Faibussowitsch `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