xref: /petsc/src/dm/tutorials/ex3.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Tests DMCreateInterpolation() for nonuniform DMDA coordinates.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdm.h>
4c4762a1bSJed Brown #include <petscdmda.h>
5c4762a1bSJed Brown 
SetCoordinates1d(DM da)6d71ae5a4SJacob Faibussowitsch PetscErrorCode SetCoordinates1d(DM da)
7d71ae5a4SJacob Faibussowitsch {
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++) {
22ad540459SPierre 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));
283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29c4762a1bSJed Brown }
30c4762a1bSJed Brown 
SetCoordinates2d(DM da)31d71ae5a4SJacob Faibussowitsch PetscErrorCode SetCoordinates2d(DM da)
32d71ae5a4SJacob Faibussowitsch {
33c4762a1bSJed Brown   PetscInt     i, j, mstart, m, nstart, n;
34c4762a1bSJed Brown   Vec          local, global;
35c4762a1bSJed Brown   DMDACoor2d **coors, **coorslocal;
36c4762a1bSJed Brown   DM           cda;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   PetscFunctionBeginUser;
399566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
409566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
419566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &global));
429566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(da, &local));
439566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, global, &coors));
449566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, local, &coorslocal));
459566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(cda, &mstart, &nstart, 0, &m, &n, 0));
46c4762a1bSJed Brown   for (i = mstart; i < mstart + m; i++) {
47c4762a1bSJed Brown     for (j = nstart; j < nstart + n; j++) {
48ad540459SPierre Jolivet       if (i % 2) coors[j][i].x = coorslocal[j][i - 1].x + .1 * (coorslocal[j][i + 1].x - coorslocal[j][i - 1].x);
49ad540459SPierre Jolivet       if (j % 2) coors[j][i].y = coorslocal[j - 1][i].y + .3 * (coorslocal[j + 1][i].y - coorslocal[j - 1][i].y);
50c4762a1bSJed Brown     }
51c4762a1bSJed Brown   }
529566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, global, &coors));
539566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, local, &coorslocal));
54c4762a1bSJed Brown 
559566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(cda, global, INSERT_VALUES, local));
569566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(cda, global, INSERT_VALUES, local));
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
SetCoordinates3d(DM da)60d71ae5a4SJacob Faibussowitsch PetscErrorCode SetCoordinates3d(DM da)
61d71ae5a4SJacob Faibussowitsch {
62c4762a1bSJed Brown   PetscInt      i, j, mstart, m, nstart, n, pstart, p, k;
63c4762a1bSJed Brown   Vec           local, global;
64c4762a1bSJed Brown   DMDACoor3d ***coors, ***coorslocal;
65c4762a1bSJed Brown   DM            cda;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   PetscFunctionBeginUser;
689566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
699566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
709566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &global));
719566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(da, &local));
729566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, global, &coors));
739566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, local, &coorslocal));
749566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(cda, &mstart, &nstart, &pstart, &m, &n, &p));
75c4762a1bSJed Brown   for (i = mstart; i < mstart + m; i++) {
76c4762a1bSJed Brown     for (j = nstart; j < nstart + n; j++) {
77c4762a1bSJed Brown       for (k = pstart; k < pstart + p; k++) {
78ad540459SPierre 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);
79ad540459SPierre 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);
80ad540459SPierre 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);
81c4762a1bSJed Brown       }
82c4762a1bSJed Brown     }
83c4762a1bSJed Brown   }
849566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, global, &coors));
859566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, local, &coorslocal));
869566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(cda, global, INSERT_VALUES, local));
879566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(cda, global, INSERT_VALUES, local));
883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
main(int argc,char ** argv)91d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
92d71ae5a4SJacob Faibussowitsch {
93c4762a1bSJed Brown   PetscInt        M = 5, N = 4, P = 3, m = PETSC_DECIDE, n = PETSC_DECIDE, p = PETSC_DECIDE, dim = 1;
94c4762a1bSJed Brown   DM              dac, daf;
95c4762a1bSJed Brown   DMBoundaryType  bx = DM_BOUNDARY_NONE, by = DM_BOUNDARY_NONE, bz = DM_BOUNDARY_NONE;
96c4762a1bSJed Brown   DMDAStencilType stype = DMDA_STENCIL_BOX;
97c4762a1bSJed Brown   Mat             A;
98c4762a1bSJed Brown 
99327415f7SBarry Smith   PetscFunctionBeginUser;
100*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
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, "-m", &m, NULL));
1059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
1069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-p", &p, NULL));
1079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* Create distributed array and get vectors */
110c4762a1bSJed Brown   if (dim == 1) {
1119566063dSJacob Faibussowitsch     PetscCall(DMDACreate1d(PETSC_COMM_WORLD, bx, M, 1, 1, NULL, &dac));
112c4762a1bSJed Brown   } else if (dim == 2) {
1139566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(PETSC_COMM_WORLD, bx, by, stype, M, N, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &dac));
114c4762a1bSJed Brown   } else if (dim == 3) {
1159566063dSJacob 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));
116c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "dim must be 1,2, or 3");
1179566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dac));
1189566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dac));
119c4762a1bSJed Brown 
1209566063dSJacob Faibussowitsch   PetscCall(DMRefine(dac, PETSC_COMM_WORLD, &daf));
121c4762a1bSJed Brown 
1229566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dac, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
123c4762a1bSJed Brown   if (dim == 1) {
1249566063dSJacob Faibussowitsch     PetscCall(SetCoordinates1d(daf));
125c4762a1bSJed Brown   } else if (dim == 2) {
1269566063dSJacob Faibussowitsch     PetscCall(SetCoordinates2d(daf));
127c4762a1bSJed Brown   } else if (dim == 3) {
1289566063dSJacob Faibussowitsch     PetscCall(SetCoordinates3d(daf));
129c4762a1bSJed Brown   }
1309566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(dac, daf, &A, 0));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* Free memory */
1339566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dac));
1349566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&daf));
1359566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1369566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
137b122ec5aSJacob Faibussowitsch   return 0;
138c4762a1bSJed Brown }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown /*TEST
141c4762a1bSJed Brown 
142c4762a1bSJed Brown    test:
143c4762a1bSJed Brown       nsize: 3
144c4762a1bSJed Brown       args: -mat_view
145c4762a1bSJed Brown 
146c4762a1bSJed Brown    test:
147c4762a1bSJed Brown       suffix: 2
148c4762a1bSJed Brown       nsize: 3
149c4762a1bSJed Brown       args: -mat_view -dim 2
150c4762a1bSJed Brown 
151c4762a1bSJed Brown    test:
152c4762a1bSJed Brown       suffix: 3
153c4762a1bSJed Brown       nsize: 3
154c4762a1bSJed Brown       args: -mat_view -dim 3
155c4762a1bSJed Brown 
156c4762a1bSJed Brown TEST*/
157