xref: /petsc/src/dm/tutorials/ex3.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 
7d71ae5a4SJacob Faibussowitsch PetscErrorCode SetCoordinates1d(DM da)
8d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown   PetscInt     i, start, m;
10c4762a1bSJed Brown   Vec          local, global;
11c4762a1bSJed Brown   PetscScalar *coors, *coorslocal;
12c4762a1bSJed Brown   DM           cda;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
169566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
179566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &global));
189566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(da, &local));
199566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, global, &coors));
209566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, local, &coorslocal));
219566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(cda, &start, 0, 0, &m, 0, 0));
22c4762a1bSJed Brown   for (i = start; i < start + m; i++) {
23ad540459SPierre Jolivet     if (i % 2) coors[i] = coorslocal[i - 1] + .1 * (coorslocal[i + 1] - coorslocal[i - 1]);
24c4762a1bSJed Brown   }
259566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, global, &coors));
269566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, local, &coorslocal));
279566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(cda, global, INSERT_VALUES, local));
289566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(cda, global, INSERT_VALUES, local));
29*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30c4762a1bSJed Brown }
31c4762a1bSJed Brown 
32d71ae5a4SJacob Faibussowitsch PetscErrorCode SetCoordinates2d(DM da)
33d71ae5a4SJacob Faibussowitsch {
34c4762a1bSJed Brown   PetscInt     i, j, mstart, m, nstart, n;
35c4762a1bSJed Brown   Vec          local, global;
36c4762a1bSJed Brown   DMDACoor2d **coors, **coorslocal;
37c4762a1bSJed Brown   DM           cda;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   PetscFunctionBeginUser;
409566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
419566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
429566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &global));
439566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(da, &local));
449566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, global, &coors));
459566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, local, &coorslocal));
469566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(cda, &mstart, &nstart, 0, &m, &n, 0));
47c4762a1bSJed Brown   for (i = mstart; i < mstart + m; i++) {
48c4762a1bSJed Brown     for (j = nstart; j < nstart + n; j++) {
49ad540459SPierre Jolivet       if (i % 2) coors[j][i].x = coorslocal[j][i - 1].x + .1 * (coorslocal[j][i + 1].x - coorslocal[j][i - 1].x);
50ad540459SPierre Jolivet       if (j % 2) coors[j][i].y = coorslocal[j - 1][i].y + .3 * (coorslocal[j + 1][i].y - coorslocal[j - 1][i].y);
51c4762a1bSJed Brown     }
52c4762a1bSJed Brown   }
539566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, global, &coors));
549566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, local, &coorslocal));
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(cda, global, INSERT_VALUES, local));
579566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(cda, global, INSERT_VALUES, local));
58*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
61d71ae5a4SJacob Faibussowitsch PetscErrorCode SetCoordinates3d(DM da)
62d71ae5a4SJacob Faibussowitsch {
63c4762a1bSJed Brown   PetscInt      i, j, mstart, m, nstart, n, pstart, p, k;
64c4762a1bSJed Brown   Vec           local, global;
65c4762a1bSJed Brown   DMDACoor3d ***coors, ***coorslocal;
66c4762a1bSJed Brown   DM            cda;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   PetscFunctionBeginUser;
699566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
709566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(da, &cda));
719566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &global));
729566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(da, &local));
739566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, global, &coors));
749566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, local, &coorslocal));
759566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(cda, &mstart, &nstart, &pstart, &m, &n, &p));
76c4762a1bSJed Brown   for (i = mstart; i < mstart + m; i++) {
77c4762a1bSJed Brown     for (j = nstart; j < nstart + n; j++) {
78c4762a1bSJed Brown       for (k = pstart; k < pstart + p; k++) {
79ad540459SPierre 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);
80ad540459SPierre 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);
81ad540459SPierre 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);
82c4762a1bSJed Brown       }
83c4762a1bSJed Brown     }
84c4762a1bSJed Brown   }
859566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, global, &coors));
869566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, local, &coorslocal));
879566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(cda, global, INSERT_VALUES, local));
889566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(cda, global, INSERT_VALUES, local));
89*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
90c4762a1bSJed Brown }
91c4762a1bSJed Brown 
92d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
93d71ae5a4SJacob Faibussowitsch {
94c4762a1bSJed Brown   PetscInt        M = 5, N = 4, P = 3, m = PETSC_DECIDE, n = PETSC_DECIDE, p = PETSC_DECIDE, dim = 1;
95c4762a1bSJed Brown   DM              dac, daf;
96c4762a1bSJed Brown   DMBoundaryType  bx = DM_BOUNDARY_NONE, by = DM_BOUNDARY_NONE, bz = DM_BOUNDARY_NONE;
97c4762a1bSJed Brown   DMDAStencilType stype = DMDA_STENCIL_BOX;
98c4762a1bSJed Brown   Mat             A;
99c4762a1bSJed Brown 
100327415f7SBarry Smith   PetscFunctionBeginUser;
1019566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
1039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-N", &N, NULL));
1049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-P", &P, NULL));
1059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
1069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
1079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-p", &p, NULL));
1089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* Create distributed array and get vectors */
111c4762a1bSJed Brown   if (dim == 1) {
1129566063dSJacob Faibussowitsch     PetscCall(DMDACreate1d(PETSC_COMM_WORLD, bx, M, 1, 1, NULL, &dac));
113c4762a1bSJed Brown   } else if (dim == 2) {
1149566063dSJacob Faibussowitsch     PetscCall(DMDACreate2d(PETSC_COMM_WORLD, bx, by, stype, M, N, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, NULL, &dac));
115c4762a1bSJed Brown   } else if (dim == 3) {
1169566063dSJacob 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));
117c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "dim must be 1,2, or 3");
1189566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dac));
1199566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dac));
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch   PetscCall(DMRefine(dac, PETSC_COMM_WORLD, &daf));
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(dac, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
124c4762a1bSJed Brown   if (dim == 1) {
1259566063dSJacob Faibussowitsch     PetscCall(SetCoordinates1d(daf));
126c4762a1bSJed Brown   } else if (dim == 2) {
1279566063dSJacob Faibussowitsch     PetscCall(SetCoordinates2d(daf));
128c4762a1bSJed Brown   } else if (dim == 3) {
1299566063dSJacob Faibussowitsch     PetscCall(SetCoordinates3d(daf));
130c4762a1bSJed Brown   }
1319566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(dac, daf, &A, 0));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* Free memory */
1349566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dac));
1359566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&daf));
1369566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1379566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
138b122ec5aSJacob Faibussowitsch   return 0;
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown /*TEST
142c4762a1bSJed Brown 
143c4762a1bSJed Brown    test:
144c4762a1bSJed Brown       nsize: 3
145c4762a1bSJed Brown       args: -mat_view
146c4762a1bSJed Brown 
147c4762a1bSJed Brown    test:
148c4762a1bSJed Brown       suffix: 2
149c4762a1bSJed Brown       nsize: 3
150c4762a1bSJed Brown       args: -mat_view -dim 2
151c4762a1bSJed Brown 
152c4762a1bSJed Brown    test:
153c4762a1bSJed Brown       suffix: 3
154c4762a1bSJed Brown       nsize: 3
155c4762a1bSJed Brown       args: -mat_view -dim 3
156c4762a1bSJed Brown 
157c4762a1bSJed Brown TEST*/
158