Actual source code: swarm_ex3.c

  1: static char help[] = "Tests DMSwarm with DMShell\n\n";

  3: #include <petscsf.h>
  4: #include <petscdm.h>
  5: #include <petscdmshell.h>
  6: #include <petscdmda.h>
  7: #include <petscdmswarm.h>
  8: #include <petsc/private/dmimpl.h>

 10: PetscErrorCode _DMLocatePoints_DMDARegular_IS(DM dm, Vec pos, IS *iscell)
 11: {
 12:   PetscInt           p, n, bs, npoints, si, sj, milocal, mjlocal, mx, my;
 13:   DM                 dmregular;
 14:   PetscInt          *cellidx;
 15:   const PetscScalar *coor;
 16:   PetscReal          dx, dy;
 17:   PetscMPIInt        rank;

 19:   PetscFunctionBegin;
 20:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
 21:   PetscCall(VecGetLocalSize(pos, &n));
 22:   PetscCall(VecGetBlockSize(pos, &bs));
 23:   npoints = n / bs;

 25:   PetscCall(PetscMalloc1(npoints, &cellidx));
 26:   PetscCall(DMGetApplicationContext(dm, &dmregular));
 27:   PetscCall(DMDAGetCorners(dmregular, &si, &sj, NULL, &milocal, &mjlocal, NULL));
 28:   PetscCall(DMDAGetInfo(dmregular, NULL, &mx, &my, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));

 30:   dx = 2.0 / ((PetscReal)mx);
 31:   dy = 2.0 / ((PetscReal)my);

 33:   PetscCall(VecGetArrayRead(pos, &coor));
 34:   for (p = 0; p < npoints; p++) {
 35:     PetscReal coorx, coory;
 36:     PetscInt  mi, mj;

 38:     coorx = PetscRealPart(coor[2 * p]);
 39:     coory = PetscRealPart(coor[2 * p + 1]);

 41:     mi = (PetscInt)((coorx - (-1.0)) / dx);
 42:     mj = (PetscInt)((coory - (-1.0)) / dy);

 44:     cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;

 46:     if ((mj >= sj) && (mj < sj + mjlocal)) {
 47:       if ((mi >= si) && (mi < si + milocal)) cellidx[p] = (mi - si) + (mj - sj) * milocal;
 48:     }
 49:     if (coorx < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
 50:     if (coorx > 1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
 51:     if (coory < -1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
 52:     if (coory > 1.0) cellidx[p] = DMLOCATEPOINT_POINT_NOT_FOUND;
 53:   }
 54:   PetscCall(VecRestoreArrayRead(pos, &coor));
 55:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, npoints, cellidx, PETSC_OWN_POINTER, iscell));
 56:   PetscFunctionReturn(PETSC_SUCCESS);
 57: }

 59: PetscErrorCode DMLocatePoints_DMDARegular(DM dm, Vec pos, DMPointLocationType ltype, PetscSF cellSF)
 60: {
 61:   IS              iscell;
 62:   PetscSFNode    *cells;
 63:   PetscInt        p, bs, npoints, nfound;
 64:   const PetscInt *boxCells;

 66:   PetscFunctionBegin;
 67:   PetscCall(_DMLocatePoints_DMDARegular_IS(dm, pos, &iscell));
 68:   PetscCall(VecGetLocalSize(pos, &npoints));
 69:   PetscCall(VecGetBlockSize(pos, &bs));
 70:   npoints = npoints / bs;

 72:   PetscCall(PetscMalloc1(npoints, &cells));
 73:   PetscCall(ISGetIndices(iscell, &boxCells));

 75:   for (p = 0; p < npoints; p++) {
 76:     cells[p].rank  = 0;
 77:     cells[p].index = DMLOCATEPOINT_POINT_NOT_FOUND;
 78:     cells[p].index = boxCells[p];
 79:   }
 80:   PetscCall(ISRestoreIndices(iscell, &boxCells));
 81:   PetscCall(ISDestroy(&iscell));
 82:   nfound = npoints;
 83:   PetscCall(PetscSFSetGraph(cellSF, npoints, nfound, NULL, PETSC_OWN_POINTER, cells, PETSC_OWN_POINTER));
 84:   PetscCall(ISDestroy(&iscell));
 85:   PetscFunctionReturn(PETSC_SUCCESS);
 86: }

 88: PetscErrorCode DMGetNeighbors_DMDARegular(DM dm, PetscInt *nneighbors, const PetscMPIInt **neighbors)
 89: {
 90:   DM dmregular;

 92:   PetscFunctionBegin;
 93:   PetscCall(DMGetApplicationContext(dm, &dmregular));
 94:   PetscCall(DMGetNeighbors(dmregular, nneighbors, neighbors));
 95:   PetscFunctionReturn(PETSC_SUCCESS);
 96: }

 98: PetscErrorCode SwarmViewGP(DM dms, const char prefix[])
 99: {
100:   PetscReal  *array;
101:   PetscInt   *iarray;
102:   PetscInt    npoints, p, bs;
103:   FILE       *fp;
104:   char        name[PETSC_MAX_PATH_LEN];
105:   PetscMPIInt rank;

107:   PetscFunctionBegin;
108:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
109:   PetscCall(PetscSNPrintf(name, PETSC_MAX_PATH_LEN - 1, "%s-rank%d.gp", prefix, rank));
110:   fp = fopen(name, "w");
111:   PetscCheck(fp, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "Cannot open file %s", name);
112:   PetscCall(DMSwarmGetLocalSize(dms, &npoints));
113:   PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
114:   PetscCall(DMSwarmGetField(dms, "itag", NULL, NULL, (void **)&iarray));
115:   for (p = 0; p < npoints; p++) fprintf(fp, "%+1.4e %+1.4e %1.4e\n", array[2 * p], array[2 * p + 1], (double)iarray[p]);
116:   PetscCall(DMSwarmRestoreField(dms, "itag", NULL, NULL, (void **)&iarray));
117:   PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
118:   fclose(fp);
119:   PetscFunctionReturn(PETSC_SUCCESS);
120: }

122: /*
123:  Create a DMShell and attach a regularly spaced DMDA for point location
124:  Override methods for point location
125: */
126: PetscErrorCode ex3_1(void)
127: {
128:   DM          dms, dmcell, dmregular;
129:   PetscMPIInt rank;
130:   PetscInt    p, bs, nlocal, overlap, mx, tk;
131:   PetscReal   dx;
132:   PetscReal  *array, dt;
133:   PetscInt   *iarray;
134:   PetscRandom rand;

136:   PetscFunctionBegin;
137:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));

139:   /* Create a regularly spaced DMDA */
140:   mx      = 40;
141:   overlap = 0;
142:   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, mx, mx, PETSC_DECIDE, PETSC_DECIDE, 1, overlap, NULL, NULL, &dmregular));
143:   PetscCall(DMSetFromOptions(dmregular));
144:   PetscCall(DMSetUp(dmregular));

146:   dx = 2.0 / ((PetscReal)mx);
147:   PetscCall(DMDASetUniformCoordinates(dmregular, -1.0 + 0.5 * dx, 1.0 - 0.5 * dx, -1.0 + 0.5 * dx, 1.0 - 0.5 * dx, -1.0, 1.0));

149:   /* Create a DMShell for point location purposes */
150:   PetscCall(DMShellCreate(PETSC_COMM_WORLD, &dmcell));
151:   PetscCall(DMSetApplicationContext(dmcell, dmregular));
152:   dmcell->ops->locatepoints = DMLocatePoints_DMDARegular;
153:   dmcell->ops->getneighbors = DMGetNeighbors_DMDARegular;

155:   /* Create the swarm */
156:   PetscCall(DMCreate(PETSC_COMM_WORLD, &dms));
157:   PetscCall(DMSetType(dms, DMSWARM));
158:   PetscCall(DMSetDimension(dms, 2));
159:   PetscCall(PetscObjectSetName((PetscObject)dms, "Particles"));

161:   PetscCall(DMSwarmSetType(dms, DMSWARM_PIC));
162:   PetscCall(DMSwarmSetCellDM(dms, dmcell));

164:   PetscCall(DMSwarmRegisterPetscDatatypeField(dms, "itag", 1, PETSC_INT));
165:   PetscCall(DMSwarmFinalizeFieldRegister(dms));
166:   {
167:     PetscInt           si, sj, milocal, mjlocal;
168:     const PetscScalar *LA_coors;
169:     Vec                coors;
170:     PetscInt           cnt;

172:     PetscCall(DMDAGetCorners(dmregular, &si, &sj, NULL, &milocal, &mjlocal, NULL));
173:     PetscCall(DMGetCoordinates(dmregular, &coors));
174:     PetscCall(VecGetArrayRead(coors, &LA_coors));
175:     PetscCall(DMSwarmSetLocalSizes(dms, milocal * mjlocal, 4));
176:     PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
177:     PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
178:     cnt = 0;
179:     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
180:     PetscCall(PetscRandomSetInterval(rand, -1.0, 1.0));
181:     for (p = 0; p < nlocal; p++) {
182:       PetscReal px, py, rx, ry, r2;

184:       PetscCall(PetscRandomGetValueReal(rand, &rx));
185:       PetscCall(PetscRandomGetValueReal(rand, &ry));

187:       px = PetscRealPart(LA_coors[2 * p + 0]) + 0.1 * rx * dx;
188:       py = PetscRealPart(LA_coors[2 * p + 1]) + 0.1 * ry * dx;

190:       r2 = px * px + py * py;
191:       if (r2 < 0.75 * 0.75) {
192:         array[bs * cnt + 0] = px;
193:         array[bs * cnt + 1] = py;
194:         cnt++;
195:       }
196:     }
197:     PetscCall(PetscRandomDestroy(&rand));
198:     PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
199:     PetscCall(VecRestoreArrayRead(coors, &LA_coors));
200:     PetscCall(DMSwarmSetLocalSizes(dms, cnt, 4));

202:     PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
203:     PetscCall(DMSwarmGetField(dms, "itag", &bs, NULL, (void **)&iarray));
204:     for (p = 0; p < nlocal; p++) iarray[p] = (PetscInt)rank;
205:     PetscCall(DMSwarmRestoreField(dms, "itag", &bs, NULL, (void **)&iarray));
206:   }

208:   PetscCall(DMView(dms, PETSC_VIEWER_STDOUT_WORLD));
209:   PetscCall(SwarmViewGP(dms, "step0"));

211:   dt = 0.1;
212:   for (tk = 1; tk < 20; tk++) {
213:     char prefix[PETSC_MAX_PATH_LEN];
214:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Step %" PetscInt_FMT " \n", tk));
215:     /* push points */
216:     PetscCall(DMSwarmGetLocalSize(dms, &nlocal));
217:     PetscCall(DMSwarmGetField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));
218:     for (p = 0; p < nlocal; p++) {
219:       PetscReal cx, cy, vx, vy;

221:       cx = array[2 * p];
222:       cy = array[2 * p + 1];
223:       vx = cy;
224:       vy = -cx;

226:       array[2 * p] += dt * vx;
227:       array[2 * p + 1] += dt * vy;
228:     }
229:     PetscCall(DMSwarmRestoreField(dms, DMSwarmPICField_coor, &bs, NULL, (void **)&array));

231:     /* migrate points */
232:     PetscCall(DMSwarmMigrate(dms, PETSC_TRUE));
233:     /* view points */
234:     PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN - 1, "step%" PetscInt_FMT, tk));
235:     /* should use the regular SwarmView() api, not one for a particular type */
236:     PetscCall(SwarmViewGP(dms, prefix));
237:   }
238:   PetscCall(DMDestroy(&dmregular));
239:   PetscCall(DMDestroy(&dmcell));
240:   PetscCall(DMDestroy(&dms));
241:   PetscFunctionReturn(PETSC_SUCCESS);
242: }

244: int main(int argc, char **argv)
245: {
246:   PetscFunctionBeginUser;
247:   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
248:   PetscCall(ex3_1());
249:   PetscCall(PetscFinalize());
250:   return 0;
251: }

253: /*TEST

255:    build:
256:       requires: double !complex

258:    test:
259:       filter: grep -v atomic
260:       filter_output: grep -v atomic
261: TEST*/