Actual source code: meshBoundary.c

  1: #ifdef PETSC_RCS_HEADER
  2: static char vcid[] = "$Id: meshBoundary.c,v 1.19 2000/10/17 13:48:57 knepley Exp $";
  3: #endif

  5: /*
  6:      Defines the interface to the mesh boundary functions
  7: */

 9:  #include src/mesh/meshimpl.h

 11: /*------------------------------------------- Mesh Boundary Creation ------------------------------------------------*/
 12: /*@
 13:   MeshBoundary2DCreateSimple - This function creates the boundary of a simple square mesh for the mesh generator.

 15:   Input Parameter:
 16: . geomCtx     - The problem geometry information

 18:   Output Parameters in bdCtx:
 19: + numBd       - The number of closed boundaries
 20: . numVertices - The number of vertices
 21: . vertices    - The vertex coordinates
 22: . markers     - The vertex markers
 23: . numSegments - The number of segments
 24: . segments    - The segment endpoints
 25: . segMarkers  - The segment markers
 26: . numHoles    - The number of holes (particles)
 27: - holes       - The hole coordinates

 29:   Level: beginner

 31: .seealso MeshBoundary2DCreate(), MeshBoundary2DDestroy()
 32: .keywords mesh, boundary
 33: @*/
 34: int MeshBoundary2DCreateSimple(MPI_Comm comm, MeshGeometryContext *geomCtx, MeshBoundary2D *bdCtx)
 35: {
 36:   double     length       = geomCtx->size[0];
 37:   double     width        = geomCtx->size[1];
 38:   double     xStart       = geomCtx->start[0];
 39:   double     yStart       = geomCtx->start[1];
 40:   PetscTruth periodicMesh = geomCtx->isPeriodic[0];
 41:   double    *vertices;
 42:   int       *markers;
 43:   int       *segments;
 44:   int       *segMarkers;
 45:   int        rank;
 46:   int        node, seg;
 47:   int        ierr;

 50:   if (periodicMesh == PETSC_TRUE) {
 51:     bdCtx->numBd       = 2;
 52:     bdCtx->numVertices = 12;
 53:     bdCtx->numSegments = 20;
 54:     bdCtx->numHoles    = 0;
 55:   } else {
 56:     bdCtx->numBd       = 1;
 57:     bdCtx->numVertices = 9;
 58:     bdCtx->numSegments = 12;
 59:     bdCtx->numHoles    = 0;
 60:   }
 61:   MPI_Comm_rank(comm, &rank);
 62:   if (rank != 0) return(0);

 64:   /* Setup nodes */
 65:   PetscMalloc((bdCtx->numVertices)*2 * sizeof(double), &vertices);
 66:   PetscMalloc((bdCtx->numVertices)   * sizeof(int),    &markers);
 67:   if (periodicMesh == PETSC_TRUE) {
 68:     vertices[0]  = xStart;                  vertices[1]  = yStart;
 69:     vertices[2]  = xStart +     length/4.0; vertices[3]  = yStart;
 70:     vertices[4]  = xStart +     length/2.0; vertices[5]  = yStart;
 71:     vertices[6]  = xStart + 3.0*length/4.0; vertices[7]  = yStart;
 72:     vertices[8]  = xStart +     length;     vertices[9]  = yStart + width;
 73:     vertices[10] = xStart + 3.0*length/4.0; vertices[11] = yStart + width;
 74:     vertices[12] = xStart +     length/2.0; vertices[13] = yStart + width;
 75:     vertices[14] = xStart +     length/4.0; vertices[15] = yStart + width;
 76:     vertices[16] = xStart;                  vertices[17] = yStart + width/2.0;
 77:     vertices[18] = xStart +     length/4.0; vertices[19] = yStart + width/2.0;
 78:     vertices[20] = xStart +     length/2.0; vertices[21] = yStart + width/2.0;
 79:     vertices[22] = xStart + 3.0*length/4.0; vertices[23] = yStart + width/2.0;
 80:     for(node = 0; node < 4; node++)
 81:       markers[node] = BOTTOM_BD;
 82:     for(node = 4; node < 8; node++)
 83:       markers[node] = TOP_BD;
 84:     for(node = 8; node < 12; node++)
 85:       markers[node] = 0;
 86:   } else {
 87:     vertices[0]  = xStart;              vertices[1]  = yStart;
 88:     vertices[2]  = xStart + length/2.0; vertices[3]  = yStart;
 89:     vertices[4]  = xStart + length;     vertices[5]  = yStart;
 90:     vertices[6]  = xStart + length;     vertices[7]  = yStart + width/2.0;
 91:     vertices[8]  = xStart + length;     vertices[9]  = yStart + width;
 92:     vertices[10] = xStart + length/2.0; vertices[11] = yStart + width;
 93:     vertices[12] = xStart;              vertices[13] = yStart + width;
 94:     vertices[14] = xStart;              vertices[15] = yStart + width/2.0;
 95:     vertices[16] = xStart + length/2.0; vertices[17] = yStart + width/2.0;
 96:     for(node = 0; node < 8; node++)
 97:       markers[node] = OUTER_BD;
 98:     markers[8]      = 0;
 99:   }

101:   /* Setup segments */
102:   PetscMalloc((bdCtx->numSegments)*2 * sizeof(int), &segments);
103:   PetscMalloc((bdCtx->numSegments)   * sizeof(int), &segMarkers);
104:   if (periodicMesh == PETSC_TRUE) {
105:     segments[0]  = 0;  segments[1]  = 1;  segments[2]  = 1;  segments[3]  = 2;
106:     segments[4]  = 2;  segments[5]  = 3;  segments[6]  = 3;  segments[7]  = 0;
107:     segments[8]  = 4;  segments[9]  = 5;  segments[10] = 5;  segments[11] = 6;
108:     segments[12] = 6;  segments[13] = 7;  segments[14] = 7;  segments[15] = 4;
109:     segments[16] = 8;  segments[17] = 9;  segments[18] = 9;  segments[19] = 10;
110:     segments[20] = 10; segments[21] = 11; segments[22] = 11; segments[23] = 8;
111:     segments[24] = 0;  segments[25] = 8;  segments[26] = 4;  segments[27] = 8;
112:     segments[28] = 1;  segments[29] = 9;  segments[30] = 7;  segments[31] = 9;
113:     segments[32] = 2;  segments[33] = 10; segments[34] = 6;  segments[35] = 10;
114:     segments[36] = 3;  segments[37] = 11; segments[38] = 5;  segments[39] = 11;
115:     for(seg = 0; seg < 4; seg++)
116:       segMarkers[seg] = BOTTOM_BD;
117:     for(seg = 4; seg < 8; seg++)
118:       segMarkers[seg] = TOP_BD;
119:     for(seg = 8; seg < 20; seg++)
120:       segMarkers[seg] = 0;
121:   } else {
122:     segments[0]  = 0;  segments[1]  = 1; segments[2]  = 1;  segments[3]  = 2;
123:     segments[4]  = 2;  segments[5]  = 3; segments[6]  = 3;  segments[7]  = 4;
124:     segments[8]  = 4;  segments[9]  = 5; segments[10] = 5;  segments[11] = 6;
125:     segments[12] = 6;  segments[13] = 7; segments[14] = 7;  segments[15] = 0;
126:     segments[16] = 1;  segments[17] = 8; segments[18] = 3;  segments[19] = 8;
127:     segments[20] = 5;  segments[21] = 8; segments[22] = 7;  segments[23] = 8;
128:     for(seg = 0; seg < 8; seg++)
129:       segMarkers[seg] = OUTER_BD;
130:     for(seg = 8; seg < 12; seg++)
131:       segMarkers[seg] = 0;
132:   }

134:   bdCtx->vertices   = vertices;
135:   bdCtx->markers    = markers;
136:   bdCtx->segments   = segments;
137:   bdCtx->segMarkers = segMarkers;
138:   bdCtx->holes      = PETSC_NULL;
139:   return(0);
140: }

142: /*@
143:   MeshBoundary2DDestroy - This function frees the memory associated with the boundary of the initial mesh.

145:   Input Parameters in bdCtx:
146: + vertices   - The vertex coordinates
147: . markers    - The vertex markers
148: . segments   - The segment endpoints
149: . segMarkers - The segment markers
150: - holes      - The hole coordinates

152:   Level: beginner

154: .seealso MeshBoundary2DCreateSimple, MeshBoundary2DCreate
155: .keywords mesh, boundary
156: @*/
157: int MeshBoundary2DDestroy(MPI_Comm comm, MeshBoundary2D *bdCtx)
158: {
159:   int rank;

163:   MPI_Comm_rank(comm, &rank);
164:   if (rank == 0) {
165:     PetscFree(bdCtx->vertices);
166:     PetscFree(bdCtx->markers);
167:     PetscFree(bdCtx->segments);
168:     PetscFree(bdCtx->segMarkers);
169:     if (bdCtx->numHoles > 0) {
170:       PetscFree(bdCtx->holes);
171:     }
172:   }
173:   return(0);
174: }

176: /*---------------------------------------- Mesh Boundary Visualization ----------------------------------------------*/
177: static int MeshBoundary2DView_File(Mesh mesh, MeshBoundary2D *bdCtx, PetscViewer viewer)
178: {
179:   FILE *fd;
180:   int   i;
181:   int   ierr;

184:   PetscViewerASCIIGetPointer(viewer, &fd);
185:   PetscFPrintf(PETSC_COMM_WORLD, fd, "Mesh Boundary Object:n");
186:   if (bdCtx->numBd == 1) {
187:     PetscFPrintf(PETSC_COMM_WORLD, fd, "    %d closed boundaryn", bdCtx->numBd);
188:   } else {
189:     PetscFPrintf(PETSC_COMM_WORLD, fd, "    %d closed boundariesn", bdCtx->numBd);
190:   }
191:   PetscFPrintf(PETSC_COMM_WORLD, fd, "      Boundary vertices: %d segments: %d holes: %dn",
192:                bdCtx->numVertices, bdCtx->numSegments, bdCtx->numHoles);
193:   for(i = 0; i < bdCtx->numVertices; i++) {
194:     PetscFPrintf(PETSC_COMM_WORLD, fd, "    %d %g %g %dn", i, bdCtx->vertices[i*2], bdCtx->vertices[i*2+1], bdCtx->markers[i]);
195:   }
196:   for(i = 0; i < bdCtx->numSegments; i++) {
197:     PetscFPrintf(PETSC_COMM_WORLD, fd, "    %d %d %d %dn", i, bdCtx->segments[i*2], bdCtx->segments[i*2+1], bdCtx->segMarkers[i]);
198:   }
199:   for(i = 0; i < bdCtx->numHoles; i++) {
200:     PetscFPrintf(PETSC_COMM_WORLD, fd, "    %d %g %gn", i, bdCtx->holes[i*2], bdCtx->holes[i*2+1]);
201:   }
202:   return(0);
203: }

205: static int MeshBoundary2DView_Draw_Mesh(Mesh mesh, MeshBoundary2D *bdCtx, PetscDraw draw)
206: {
207:   double *vertices   = bdCtx->vertices;
208:   int    *segments   = bdCtx->segments;
209:   int    *segMarkers = bdCtx->segMarkers;
210:   int     color;
211:   int     seg;
212:   int     ierr;

215:   for(seg = 0; seg < bdCtx->numSegments; seg++) {
216:     if (segMarkers[seg]) {
217:       color = PETSC_DRAW_RED;
218:     } else {
219:       color = PETSC_DRAW_BLACK;
220:     }

222:     MeshDrawLine(mesh, draw, vertices[segments[seg*2]*2], vertices[segments[seg*2]*2+1],
223:                         vertices[segments[seg*2+1]*2], vertices[segments[seg*2+1]*2+1], color);
224: 
225:   }
226:   PetscDrawSynchronizedFlush(draw);
227:   return(0);
228: }

230: static int MeshBoundary2DView_Draw(Mesh mesh, MeshBoundary2D *bdCtx, PetscViewer v)
231: {
232:   PetscReal       startX, endX;
233:   PetscReal       startY, endY;
234:   PetscReal       sizeX,  sizeY;
235:   double         *vertices;
236:   PetscDraw       draw;
237:   PetscTruth      isnull;
238:   PetscDrawButton button;
239:   char            statusLine[256];
240:   int             pause;
241:   PetscReal       xc, yc, scale;
242:   int             xsize, ysize;
243:   int             vert;
244:   int             ierr;

247:   PetscViewerDrawGetDraw(v, 0, &draw);
248:   PetscDrawIsNull(draw, &isnull);
249:   if (isnull) return(0);
250:   MeshGetBoundingBox(mesh, &startX, &startY, PETSC_NULL, &endX, &endY, PETSC_NULL);
251:   vertices = bdCtx->vertices;
252:   for(vert = 0; vert < bdCtx->numVertices; vert++) {
253:     if (startX > vertices[vert*2])   startX = vertices[vert*2];
254:     if (startY > vertices[vert*2+1]) startY = vertices[vert*2+1];
255:     if (endX   < vertices[vert*2])   endX   = vertices[vert*2];
256:     if (endY   < vertices[vert*2+1]) endY   = vertices[vert*2+1];
257:   }
258:   sizeX = endX - startX;
259:   sizeY = endY - startY;
260:   if (sizeX > sizeY) {
261:     ysize  = 300;
262:     xsize  = ysize * (int) (sizeX/sizeY);
263:   } else {
264:     xsize  = 300;
265:     ysize  = xsize * (int) (sizeY/sizeX);
266:   }
267:   PetscDrawResizeWindow(draw, xsize, ysize);
268:   PetscDrawSynchronizedClear(draw);

270:   ierr  = PetscDrawSetCoordinates(draw, startX-0.02*sizeX, startY-0.02*sizeY, endX+0.02*sizeX, endY+0.02*sizeY);
271:   ierr  = MeshBoundary2DView_Draw_Mesh(mesh, bdCtx, draw);
272:   ierr  = PetscDrawGetPause(draw, &pause);
273:   if (pause >= 0) {
274:     PetscSleep(pause);
275:     return(0);
276:   }

278:   /* Allow the mesh to zoom or shrink */
279:   PetscDrawCheckResizedWindow(draw);
280:   PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
281:   while ((button != BUTTON_RIGHT) && (button != BUTTON_CENTER_SHIFT))
282:   {
283:     PetscDrawSynchronizedClear(draw);
284:     statusLine[0] = 0;
285:     switch (button)
286:     {
287:     case BUTTON_LEFT:
288:       scale = 0.5;
289:       break;
290:     case BUTTON_CENTER:
291:       scale = 2.0;
292:       break;
293:     default:
294:       scale = 1.0;
295:       break;
296:     }
297:     if (scale != 1.0) {
298:       startX = scale*(startX + sizeX - xc) + xc - sizeX*scale;
299:       endX   = scale*(endX   - sizeX - xc) + xc + sizeX*scale;
300:       startY = scale*(startY + sizeY - yc) + yc - sizeY*scale;
301:       endY   = scale*(endY   - sizeY - yc) + yc + sizeY*scale;
302:       sizeX *= scale;
303:       sizeY *= scale;
304:     }

306:     PetscDrawSetCoordinates(draw, startX-0.02*sizeX, startY-0.02*sizeY, endX+0.02*sizeX, endY+0.02*sizeY);
307:     PetscDrawString(draw, startX - 0.02*sizeX, startY - 0.02*sizeY, PETSC_DRAW_BLACK, statusLine);
308:     MeshBoundary2DView_Draw_Mesh(mesh, bdCtx, draw);
309:     PetscDrawCheckResizedWindow(draw);
310:     PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
311:   }

313:   return(0);
314: }

316: int MeshBoundary2DView(Mesh mesh, MeshBoundary2D *bdCtx, PetscViewer viewer)
317: {
318:   PetscTruth isascii, isdraw;
319:   int        ierr;

322:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &isascii);
323:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_DRAW,  &isdraw);
324:   if (isascii == PETSC_TRUE) {
325:     MeshBoundary2DView_File(mesh, bdCtx, viewer);
326:   } else if (isdraw == PETSC_TRUE) {
327:     MeshBoundary2DView_Draw(mesh, bdCtx, viewer);
328:   }

330:   return(0);
331: }

333: /*------------------------------------------- Mesh Boundary Functions ------------------------------------------------*/
334: int MeshBoundary2DSetBoundingBox(Mesh mesh, MeshBoundary2D *bdCtx)
335: {
336:   MPI_Comm comm;
337:   double   startX, startY;
338:   double   endX,   endY;
339:   double  *vertices;
340:   int      rank;
341:   int      p;
342:   int      ierr;

345:   PetscObjectGetComm((PetscObject) mesh, &comm);
346:   MPI_Comm_rank(comm, &rank);
347:   if (rank == 0) {
348:     vertices        = bdCtx->vertices;
349:     startX = endX = vertices[0];
350:     startY = endY = vertices[1];
351:     for(p = 0; p < bdCtx->numVertices; p++) {
352:       /* Ignore moving boundaries (to accomodate split particles)? This should be done better */
353:       if (bdCtx->markers[p] < 0) continue;
354:       if (startX > vertices[p*2])   startX = vertices[p*2];
355:       if (startY > vertices[p*2+1]) startY = vertices[p*2+1];
356:       if (endX   < vertices[p*2])   endX   = vertices[p*2];
357:       if (endY   < vertices[p*2+1]) endY   = vertices[p*2+1];
358:     }
359:   }
360:   MPI_Bcast(&startX, 1, MPI_DOUBLE, 0, comm);
361:   MPI_Bcast(&endX,   1, MPI_DOUBLE, 0, comm);
362:   MPI_Bcast(&startY, 1, MPI_DOUBLE, 0, comm);
363:   MPI_Bcast(&endY,   1, MPI_DOUBLE, 0, comm);
364:   MeshSetBoundingBox(mesh, startX, startY, 0.0, endX, endY, 0.0);
365:   return(0);
366: }