1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests DMCreateInterpolation() for nonuniform DMDA coordinates.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown 79371c9d4SSatish Balay PetscErrorCode SetCoordinates1d(DM da) { 8c4762a1bSJed Brown PetscInt i, start, m; 9c4762a1bSJed Brown Vec local, global; 10c4762a1bSJed Brown PetscScalar *coors, *coorslocal; 11c4762a1bSJed Brown DM cda; 12c4762a1bSJed Brown 13c4762a1bSJed Brown PetscFunctionBeginUser; 149566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 159566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda)); 169566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &global)); 179566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(da, &local)); 189566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda, global, &coors)); 199566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, local, &coorslocal)); 209566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &start, 0, 0, &m, 0, 0)); 21c4762a1bSJed Brown for (i = start; i < start + m; i++) { 22*ad540459SPierre Jolivet if (i % 2) coors[i] = coorslocal[i - 1] + .1 * (coorslocal[i + 1] - coorslocal[i - 1]); 23c4762a1bSJed Brown } 249566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda, global, &coors)); 259566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, local, &coorslocal)); 269566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(cda, global, INSERT_VALUES, local)); 279566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(cda, global, INSERT_VALUES, local)); 28c4762a1bSJed Brown PetscFunctionReturn(0); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown 319371c9d4SSatish Balay PetscErrorCode SetCoordinates2d(DM da) { 32c4762a1bSJed Brown PetscInt i, j, mstart, m, nstart, n; 33c4762a1bSJed Brown Vec local, global; 34c4762a1bSJed Brown DMDACoor2d **coors, **coorslocal; 35c4762a1bSJed Brown DM cda; 36c4762a1bSJed Brown 37c4762a1bSJed Brown PetscFunctionBeginUser; 389566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda)); 409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &global)); 419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(da, &local)); 429566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda, global, &coors)); 439566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, local, &coorslocal)); 449566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &mstart, &nstart, 0, &m, &n, 0)); 45c4762a1bSJed Brown for (i = mstart; i < mstart + m; i++) { 46c4762a1bSJed Brown for (j = nstart; j < nstart + n; j++) { 47*ad540459SPierre Jolivet if (i % 2) coors[j][i].x = coorslocal[j][i - 1].x + .1 * (coorslocal[j][i + 1].x - coorslocal[j][i - 1].x); 48*ad540459SPierre Jolivet if (j % 2) coors[j][i].y = coorslocal[j - 1][i].y + .3 * (coorslocal[j + 1][i].y - coorslocal[j - 1][i].y); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown } 519566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda, global, &coors)); 529566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, local, &coorslocal)); 53c4762a1bSJed Brown 549566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(cda, global, INSERT_VALUES, local)); 559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(cda, global, INSERT_VALUES, local)); 56c4762a1bSJed Brown PetscFunctionReturn(0); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 599371c9d4SSatish Balay PetscErrorCode SetCoordinates3d(DM da) { 60c4762a1bSJed Brown PetscInt i, j, mstart, m, nstart, n, pstart, p, k; 61c4762a1bSJed Brown Vec local, global; 62c4762a1bSJed Brown DMDACoor3d ***coors, ***coorslocal; 63c4762a1bSJed Brown DM cda; 64c4762a1bSJed Brown 65c4762a1bSJed Brown PetscFunctionBeginUser; 669566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 679566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(da, &cda)); 689566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da, &global)); 699566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(da, &local)); 709566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(cda, global, &coors)); 719566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(cda, local, &coorslocal)); 729566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(cda, &mstart, &nstart, &pstart, &m, &n, &p)); 73c4762a1bSJed Brown for (i = mstart; i < mstart + m; i++) { 74c4762a1bSJed Brown for (j = nstart; j < nstart + n; j++) { 75c4762a1bSJed Brown for (k = pstart; k < pstart + p; k++) { 76*ad540459SPierre Jolivet if (i % 2) coors[k][j][i].x = coorslocal[k][j][i - 1].x + .1 * (coorslocal[k][j][i + 1].x - coorslocal[k][j][i - 1].x); 77*ad540459SPierre Jolivet if (j % 2) coors[k][j][i].y = coorslocal[k][j - 1][i].y + .3 * (coorslocal[k][j + 1][i].y - coorslocal[k][j - 1][i].y); 78*ad540459SPierre Jolivet if (k % 2) coors[k][j][i].z = coorslocal[k - 1][j][i].z + .4 * (coorslocal[k + 1][j][i].z - coorslocal[k - 1][j][i].z); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown } 81c4762a1bSJed Brown } 829566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(cda, global, &coors)); 839566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(cda, local, &coorslocal)); 849566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(cda, global, INSERT_VALUES, local)); 859566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(cda, global, INSERT_VALUES, local)); 86c4762a1bSJed Brown PetscFunctionReturn(0); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 899371c9d4SSatish Balay int main(int argc, char **argv) { 90c4762a1bSJed Brown PetscInt M = 5, N = 4, P = 3, m = PETSC_DECIDE, n = PETSC_DECIDE, p = PETSC_DECIDE, dim = 1; 91c4762a1bSJed Brown DM dac, daf; 92c4762a1bSJed Brown DMBoundaryType bx = DM_BOUNDARY_NONE, by = DM_BOUNDARY_NONE, bz = DM_BOUNDARY_NONE; 93c4762a1bSJed Brown DMDAStencilType stype = DMDA_STENCIL_BOX; 94c4762a1bSJed Brown Mat A; 95c4762a1bSJed Brown 96327415f7SBarry Smith PetscFunctionBeginUser; 979566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL)); 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL)); 1009566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-P", &P, NULL)); 1019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL)); 1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-p", &p, NULL)); 1049566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL)); 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* Create distributed array and get vectors */ 107c4762a1bSJed Brown if (dim == 1) { 1089566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, bx, M, 1, 1, NULL, &dac)); 109c4762a1bSJed Brown } else if (dim == 2) { 1109566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, bx, by, stype, M, N, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &dac)); 111c4762a1bSJed Brown } else if (dim == 3) { 1129566063dSJacob Faibussowitsch PetscCall(DMDACreate3d(PETSC_COMM_WORLD, bx, by, bz, stype, M, N, P, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, NULL, &dac)); 113c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "dim must be 1,2, or 3"); 1149566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dac)); 1159566063dSJacob Faibussowitsch PetscCall(DMSetUp(dac)); 116c4762a1bSJed Brown 1179566063dSJacob Faibussowitsch PetscCall(DMRefine(dac, PETSC_COMM_WORLD, &daf)); 118c4762a1bSJed Brown 1199566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(dac, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0)); 120c4762a1bSJed Brown if (dim == 1) { 1219566063dSJacob Faibussowitsch PetscCall(SetCoordinates1d(daf)); 122c4762a1bSJed Brown } else if (dim == 2) { 1239566063dSJacob Faibussowitsch PetscCall(SetCoordinates2d(daf)); 124c4762a1bSJed Brown } else if (dim == 3) { 1259566063dSJacob Faibussowitsch PetscCall(SetCoordinates3d(daf)); 126c4762a1bSJed Brown } 1279566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dac, daf, &A, 0)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* Free memory */ 1309566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dac)); 1319566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daf)); 1329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1339566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 134b122ec5aSJacob Faibussowitsch return 0; 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /*TEST 138c4762a1bSJed Brown 139c4762a1bSJed Brown test: 140c4762a1bSJed Brown nsize: 3 141c4762a1bSJed Brown args: -mat_view 142c4762a1bSJed Brown 143c4762a1bSJed Brown test: 144c4762a1bSJed Brown suffix: 2 145c4762a1bSJed Brown nsize: 3 146c4762a1bSJed Brown args: -mat_view -dim 2 147c4762a1bSJed Brown 148c4762a1bSJed Brown test: 149c4762a1bSJed Brown suffix: 3 150c4762a1bSJed Brown nsize: 3 151c4762a1bSJed Brown args: -mat_view -dim 3 152c4762a1bSJed Brown 153c4762a1bSJed Brown TEST*/ 154