Actual source code: ex14.c


  2: static char help[] = "Tests DMCreateDomainDecomposition.\n\n";

  4: /*
  5: Use the options
  6:      -da_grid_x <nx> - number of grid points in x direction, if M < 0
  7:      -da_grid_y <ny> - number of grid points in y direction, if N < 0
  8:      -da_processors_x <MX> number of processors in x directio
  9:      -da_processors_y <MY> number of processors in x direction
 10: */

 12: #include <petscdm.h>
 13: #include <petscdmda.h>

 15: PetscErrorCode FillLocalSubdomain(DM da, Vec gvec)
 16: {
 17:   DMDALocalInfo info;
 18:   PetscMPIInt   rank;
 19:   PetscInt      i, j, k, l;

 22:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
 23:   DMDAGetLocalInfo(da, &info);

 25:   if (info.dim == 3) {
 26:     PetscScalar ***g;
 27:     DMDAVecGetArray(da, gvec, &g);
 28:     /* loop over ghosts */
 29:     for (k = info.zs; k < info.zs + info.zm; k++) {
 30:       for (j = info.ys; j < info.ys + info.ym; j++) {
 31:         for (i = info.xs; i < info.xs + info.xm; i++) {
 32:           g[k][j][info.dof * i + 0] = i;
 33:           g[k][j][info.dof * i + 1] = j;
 34:           g[k][j][info.dof * i + 2] = k;
 35:         }
 36:       }
 37:     }
 38:     DMDAVecRestoreArray(da, gvec, &g);
 39:   }
 40:   if (info.dim == 2) {
 41:     PetscScalar **g;
 42:     DMDAVecGetArray(da, gvec, &g);
 43:     /* loop over ghosts */
 44:     for (j = info.ys; j < info.ys + info.ym; j++) {
 45:       for (i = info.xs; i < info.xs + info.xm; i++) {
 46:         for (l = 0; l < info.dof; l++) {
 47:           g[j][info.dof * i + 0] = i;
 48:           g[j][info.dof * i + 1] = j;
 49:           g[j][info.dof * i + 2] = rank;
 50:         }
 51:       }
 52:     }
 53:     DMDAVecRestoreArray(da, gvec, &g);
 54:   }
 55:   return 0;
 56: }

 58: int main(int argc, char **argv)
 59: {
 60:   DM            da, *subda;
 61:   PetscInt      i, dim = 3;
 62:   PetscInt      M = 25, N = 25, P = 25;
 63:   PetscMPIInt   size, rank;
 64:   Vec           v;
 65:   Vec           slvec, sgvec;
 66:   IS           *ois, *iis;
 67:   VecScatter    oscata;
 68:   VecScatter   *iscat, *oscat, *gscat;
 69:   DMDALocalInfo info;
 70:   PetscBool     patchis_offproc = PETSC_TRUE;

 73:   PetscInitialize(&argc, &argv, (char *)0, help);
 74:   PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL);

 76:   /* Create distributed array and get vectors */
 77:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 78:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
 79:   if (dim == 2) {
 80:     DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, 3, 1, NULL, NULL, &da);
 81:   } else if (dim == 3) {
 82:     DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, P, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 3, 1, NULL, NULL, NULL, &da);
 83:   }
 84:   DMSetFromOptions(da);
 85:   DMSetUp(da);
 86:   DMDAGetLocalInfo(da, &info);

 88:   DMCreateDomainDecomposition(da, NULL, NULL, &iis, &ois, &subda);
 89:   DMCreateDomainDecompositionScatters(da, 1, subda, &iscat, &oscat, &gscat);

 91:   {
 92:     DMDALocalInfo subinfo;
 93:     MatStencil    lower, upper;
 94:     IS            patchis;
 95:     Vec           smallvec;
 96:     Vec           largevec;
 97:     VecScatter    patchscat;

 99:     DMDAGetLocalInfo(subda[0], &subinfo);

101:     lower.i = info.xs;
102:     lower.j = info.ys;
103:     lower.k = info.zs;
104:     upper.i = info.xs + info.xm;
105:     upper.j = info.ys + info.ym;
106:     upper.k = info.zs + info.zm;

108:     /* test the patch IS as a thing to scatter to/from */
109:     DMDACreatePatchIS(da, &lower, &upper, &patchis, patchis_offproc);
110:     DMGetGlobalVector(da, &largevec);

112:     VecCreate(PETSC_COMM_SELF, &smallvec);
113:     VecSetSizes(smallvec, info.dof * (upper.i - lower.i) * (upper.j - lower.j) * (upper.k - lower.k), PETSC_DECIDE);
114:     VecSetFromOptions(smallvec);
115:     VecScatterCreate(smallvec, NULL, largevec, patchis, &patchscat);

117:     FillLocalSubdomain(subda[0], smallvec);
118:     VecSet(largevec, 0);

120:     VecScatterBegin(patchscat, smallvec, largevec, ADD_VALUES, SCATTER_FORWARD);
121:     VecScatterEnd(patchscat, smallvec, largevec, ADD_VALUES, SCATTER_FORWARD);
122:     ISView(patchis, PETSC_VIEWER_STDOUT_WORLD);
123:     VecScatterView(patchscat, PETSC_VIEWER_STDOUT_WORLD);

125:     for (i = 0; i < size; i++) {
126:       if (i == rank) VecView(smallvec, PETSC_VIEWER_STDOUT_SELF);
127:       MPI_Barrier(PETSC_COMM_WORLD);
128:     }

130:     MPI_Barrier(PETSC_COMM_WORLD);
131:     VecView(largevec, PETSC_VIEWER_STDOUT_WORLD);

133:     VecDestroy(&smallvec);
134:     DMRestoreGlobalVector(da, &largevec);
135:     ISDestroy(&patchis);
136:     VecScatterDestroy(&patchscat);
137:   }

139:   /* view the various parts */
140:   {
141:     for (i = 0; i < size; i++) {
142:       if (i == rank) {
143:         PetscPrintf(PETSC_COMM_SELF, "Processor %d: \n", i);
144:         DMView(subda[0], PETSC_VIEWER_STDOUT_SELF);
145:       }
146:       MPI_Barrier(PETSC_COMM_WORLD);
147:     }

149:     DMGetLocalVector(subda[0], &slvec);
150:     DMGetGlobalVector(subda[0], &sgvec);
151:     DMGetGlobalVector(da, &v);

153:     /* test filling outer between the big DM and the small ones with the IS scatter*/
154:     VecScatterCreate(v, ois[0], sgvec, NULL, &oscata);

156:     FillLocalSubdomain(subda[0], sgvec);

158:     VecScatterBegin(oscata, sgvec, v, ADD_VALUES, SCATTER_REVERSE);
159:     VecScatterEnd(oscata, sgvec, v, ADD_VALUES, SCATTER_REVERSE);

161:     /* test the local-to-local scatter */

163:     /* fill up the local subdomain and then add them together */
164:     FillLocalSubdomain(da, v);

166:     VecScatterBegin(gscat[0], v, slvec, ADD_VALUES, SCATTER_FORWARD);
167:     VecScatterEnd(gscat[0], v, slvec, ADD_VALUES, SCATTER_FORWARD);

169:     VecView(v, PETSC_VIEWER_STDOUT_WORLD);

171:     /* test ghost scattering backwards */

173:     VecSet(v, 0);

175:     VecScatterBegin(gscat[0], slvec, v, ADD_VALUES, SCATTER_REVERSE);
176:     VecScatterEnd(gscat[0], slvec, v, ADD_VALUES, SCATTER_REVERSE);

178:     VecView(v, PETSC_VIEWER_STDOUT_WORLD);

180:     /* test overlap scattering backwards */

182:     DMLocalToGlobalBegin(subda[0], slvec, ADD_VALUES, sgvec);
183:     DMLocalToGlobalEnd(subda[0], slvec, ADD_VALUES, sgvec);

185:     VecSet(v, 0);

187:     VecScatterBegin(oscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE);
188:     VecScatterEnd(oscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE);

190:     VecView(v, PETSC_VIEWER_STDOUT_WORLD);

192:     /* test interior scattering backwards */

194:     VecSet(v, 0);

196:     VecScatterBegin(iscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE);
197:     VecScatterEnd(iscat[0], sgvec, v, ADD_VALUES, SCATTER_REVERSE);

199:     VecView(v, PETSC_VIEWER_STDOUT_WORLD);

201:     /* test matrix allocation */
202:     for (i = 0; i < size; i++) {
203:       if (i == rank) {
204:         Mat m;
205:         PetscPrintf(PETSC_COMM_SELF, "Processor %d: \n", i);
206:         DMSetMatType(subda[0], MATAIJ);
207:         DMCreateMatrix(subda[0], &m);
208:         MatView(m, PETSC_VIEWER_STDOUT_SELF);
209:         MatDestroy(&m);
210:       }
211:       MPI_Barrier(PETSC_COMM_WORLD);
212:     }
213:     DMRestoreLocalVector(subda[0], &slvec);
214:     DMRestoreGlobalVector(subda[0], &sgvec);
215:     DMRestoreGlobalVector(da, &v);
216:   }

218:   DMDestroy(&subda[0]);
219:   ISDestroy(&ois[0]);
220:   ISDestroy(&iis[0]);

222:   VecScatterDestroy(&iscat[0]);
223:   VecScatterDestroy(&oscat[0]);
224:   VecScatterDestroy(&gscat[0]);
225:   VecScatterDestroy(&oscata);

227:   PetscFree(iscat);
228:   PetscFree(oscat);
229:   PetscFree(gscat);
230:   PetscFree(oscata);

232:   PetscFree(subda);
233:   PetscFree(ois);
234:   PetscFree(iis);

236:   DMDestroy(&da);
237:   PetscFinalize();
238:   return 0;
239: }