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: MeshBoundary1DCreateSimple - This function creates the boundary of a simple linear 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 MeshBoundary1DCreate(), MeshBoundary1DDestroy()
32: .keywords mesh, boundary
33: @*/
34: int MeshBoundary1DCreateSimple(MPI_Comm comm, MeshGeometryContext *geomCtx, MeshBoundary2D *bdCtx) {
35: double length = geomCtx->size[0];
36: double xStart = geomCtx->start[0];
37: PetscTruth periodicMesh = geomCtx->isPeriodic[0];
38: double *vertices;
39: int *markers;
40: int rank;
41: int ierr;
44: if (periodicMesh == PETSC_TRUE) {
45: bdCtx->numBd = 0;
46: bdCtx->numVertices = 0;
47: } else {
48: bdCtx->numBd = 2;
49: bdCtx->numVertices = 2;
50: }
51: MPI_Comm_rank(comm, &rank);
52: if (rank != 0) return(0);
54: /* Setup nodes */
55: PetscMalloc((bdCtx->numVertices)*2 * sizeof(double), &vertices);
56: PetscMalloc((bdCtx->numVertices) * sizeof(int), &markers);
57: if (periodicMesh == PETSC_FALSE) {
58: vertices[0] = xStart;
59: vertices[1] = xStart + length;
60: markers[0] = BOTTOM_BD;
61: markers[1] = TOP_BD;
62: }
64: bdCtx->vertices = vertices;
65: bdCtx->markers = markers;
66: return(0);
67: }
69: /*@
70: MeshBoundary1DDestroy - This function frees the memory associated with the boundary of the initial mesh.
72: Input Parameters in bdCtx:
73: + vertices - The vertex coordinates
74: . markers - The vertex markers
75: . segments - The segment endpoints
76: . segMarkers - The segment markers
77: - holes - The hole coordinates
79: Level: beginner
81: .seealso MeshBoundary1DCreateSimple, MeshBoundary1DCreate
82: .keywords mesh, boundary
83: @*/
84: int MeshBoundary1DDestroy(MPI_Comm comm, MeshBoundary2D *bdCtx) {
85: int rank;
89: MPI_Comm_rank(comm, &rank);
90: if (rank == 0) {
91: PetscFree(bdCtx->vertices);
92: PetscFree(bdCtx->markers);
93: }
94: return(0);
95: }
97: /*@
98: MeshBoundary2DCreateSimple - This function creates the boundary of a simple square mesh for the mesh generator.
100: Input Parameter:
101: . geomCtx - The problem geometry information
103: Output Parameters in bdCtx:
104: + numBd - The number of closed boundaries
105: . numVertices - The number of vertices
106: . vertices - The vertex coordinates
107: . markers - The vertex markers
108: . numSegments - The number of segments
109: . segments - The segment endpoints
110: . segMarkers - The segment markers
111: . numHoles - The number of holes (particles)
112: - holes - The hole coordinates
114: Level: beginner
116: .seealso MeshBoundary2DCreate(), MeshBoundary2DDestroy()
117: .keywords mesh, boundary
118: @*/
119: int MeshBoundary2DCreateSimple(MPI_Comm comm, MeshGeometryContext *geomCtx, MeshBoundary2D *bdCtx)
120: {
121: double length = geomCtx->size[0];
122: double width = geomCtx->size[1];
123: double xStart = geomCtx->start[0];
124: double yStart = geomCtx->start[1];
125: PetscTruth periodicMesh = geomCtx->isPeriodic[0];
126: double *vertices;
127: int *markers;
128: int *segments;
129: int *segMarkers;
130: int rank;
131: int node, seg;
132: int ierr;
135: if (periodicMesh == PETSC_TRUE) {
136: bdCtx->numBd = 2;
137: bdCtx->numVertices = 12;
138: bdCtx->numSegments = 20;
139: bdCtx->numHoles = 0;
140: } else {
141: bdCtx->numBd = 1;
142: bdCtx->numVertices = 9;
143: bdCtx->numSegments = 12;
144: bdCtx->numHoles = 0;
145: }
146: MPI_Comm_rank(comm, &rank);
147: if (rank != 0) return(0);
149: /* Setup nodes */
150: PetscMalloc((bdCtx->numVertices)*2 * sizeof(double), &vertices);
151: PetscMalloc((bdCtx->numVertices) * sizeof(int), &markers);
152: if (periodicMesh == PETSC_TRUE) {
153: vertices[0] = xStart; vertices[1] = yStart;
154: vertices[2] = xStart + length/4.0; vertices[3] = yStart;
155: vertices[4] = xStart + length/2.0; vertices[5] = yStart;
156: vertices[6] = xStart + 3.0*length/4.0; vertices[7] = yStart;
157: vertices[8] = xStart + length; vertices[9] = yStart + width;
158: vertices[10] = xStart + 3.0*length/4.0; vertices[11] = yStart + width;
159: vertices[12] = xStart + length/2.0; vertices[13] = yStart + width;
160: vertices[14] = xStart + length/4.0; vertices[15] = yStart + width;
161: vertices[16] = xStart; vertices[17] = yStart + width/2.0;
162: vertices[18] = xStart + length/4.0; vertices[19] = yStart + width/2.0;
163: vertices[20] = xStart + length/2.0; vertices[21] = yStart + width/2.0;
164: vertices[22] = xStart + 3.0*length/4.0; vertices[23] = yStart + width/2.0;
165: for(node = 0; node < 4; node++)
166: markers[node] = BOTTOM_BD;
167: for(node = 4; node < 8; node++)
168: markers[node] = TOP_BD;
169: for(node = 8; node < 12; node++)
170: markers[node] = 0;
171: } else {
172: vertices[0] = xStart; vertices[1] = yStart;
173: vertices[2] = xStart + length/2.0; vertices[3] = yStart;
174: vertices[4] = xStart + length; vertices[5] = yStart;
175: vertices[6] = xStart + length; vertices[7] = yStart + width/2.0;
176: vertices[8] = xStart + length; vertices[9] = yStart + width;
177: vertices[10] = xStart + length/2.0; vertices[11] = yStart + width;
178: vertices[12] = xStart; vertices[13] = yStart + width;
179: vertices[14] = xStart; vertices[15] = yStart + width/2.0;
180: vertices[16] = xStart + length/2.0; vertices[17] = yStart + width/2.0;
181: for(node = 0; node < 8; node++)
182: markers[node] = OUTER_BD;
183: markers[8] = 0;
184: }
186: /* Setup segments */
187: PetscMalloc((bdCtx->numSegments)*2 * sizeof(int), &segments);
188: PetscMalloc((bdCtx->numSegments) * sizeof(int), &segMarkers);
189: if (periodicMesh == PETSC_TRUE) {
190: segments[0] = 0; segments[1] = 1; segments[2] = 1; segments[3] = 2;
191: segments[4] = 2; segments[5] = 3; segments[6] = 3; segments[7] = 0;
192: segments[8] = 4; segments[9] = 5; segments[10] = 5; segments[11] = 6;
193: segments[12] = 6; segments[13] = 7; segments[14] = 7; segments[15] = 4;
194: segments[16] = 8; segments[17] = 9; segments[18] = 9; segments[19] = 10;
195: segments[20] = 10; segments[21] = 11; segments[22] = 11; segments[23] = 8;
196: segments[24] = 0; segments[25] = 8; segments[26] = 4; segments[27] = 8;
197: segments[28] = 1; segments[29] = 9; segments[30] = 7; segments[31] = 9;
198: segments[32] = 2; segments[33] = 10; segments[34] = 6; segments[35] = 10;
199: segments[36] = 3; segments[37] = 11; segments[38] = 5; segments[39] = 11;
200: for(seg = 0; seg < 4; seg++)
201: segMarkers[seg] = BOTTOM_BD;
202: for(seg = 4; seg < 8; seg++)
203: segMarkers[seg] = TOP_BD;
204: for(seg = 8; seg < 20; seg++)
205: segMarkers[seg] = 0;
206: } else {
207: segments[0] = 0; segments[1] = 1; segments[2] = 1; segments[3] = 2;
208: segments[4] = 2; segments[5] = 3; segments[6] = 3; segments[7] = 4;
209: segments[8] = 4; segments[9] = 5; segments[10] = 5; segments[11] = 6;
210: segments[12] = 6; segments[13] = 7; segments[14] = 7; segments[15] = 0;
211: segments[16] = 1; segments[17] = 8; segments[18] = 3; segments[19] = 8;
212: segments[20] = 5; segments[21] = 8; segments[22] = 7; segments[23] = 8;
213: for(seg = 0; seg < 8; seg++)
214: segMarkers[seg] = OUTER_BD;
215: for(seg = 8; seg < 12; seg++)
216: segMarkers[seg] = 0;
217: }
219: bdCtx->vertices = vertices;
220: bdCtx->markers = markers;
221: bdCtx->segments = segments;
222: bdCtx->segMarkers = segMarkers;
223: bdCtx->holes = PETSC_NULL;
224: return(0);
225: }
227: /*@
228: MeshBoundary2DDestroy - This function frees the memory associated with the boundary of the initial mesh.
230: Input Parameters in bdCtx:
231: + vertices - The vertex coordinates
232: . markers - The vertex markers
233: . segments - The segment endpoints
234: . segMarkers - The segment markers
235: - holes - The hole coordinates
237: Level: beginner
239: .seealso MeshBoundary2DCreateSimple, MeshBoundary2DCreate
240: .keywords mesh, boundary
241: @*/
242: int MeshBoundary2DDestroy(MPI_Comm comm, MeshBoundary2D *bdCtx)
243: {
244: int rank;
248: MPI_Comm_rank(comm, &rank);
249: if (rank == 0) {
250: PetscFree(bdCtx->vertices);
251: PetscFree(bdCtx->markers);
252: PetscFree(bdCtx->segments);
253: PetscFree(bdCtx->segMarkers);
254: if (bdCtx->numHoles > 0) {
255: PetscFree(bdCtx->holes);
256: }
257: }
258: return(0);
259: }
261: /*---------------------------------------- Mesh Boundary Visualization ----------------------------------------------*/
262: static int MeshBoundary2DView_File(Mesh mesh, MeshBoundary2D *bdCtx, PetscViewer viewer)
263: {
264: FILE *fd;
265: int i;
266: int ierr;
269: PetscViewerASCIIGetPointer(viewer, &fd);
270: PetscFPrintf(PETSC_COMM_WORLD, fd, "Mesh Boundary Object:n");
271: if (bdCtx->numBd == 1) {
272: PetscFPrintf(PETSC_COMM_WORLD, fd, " %d closed boundaryn", bdCtx->numBd);
273: } else {
274: PetscFPrintf(PETSC_COMM_WORLD, fd, " %d closed boundariesn", bdCtx->numBd);
275: }
276: PetscFPrintf(PETSC_COMM_WORLD, fd, " Boundary vertices: %d segments: %d holes: %dn",
277: bdCtx->numVertices, bdCtx->numSegments, bdCtx->numHoles);
278: for(i = 0; i < bdCtx->numVertices; i++) {
279: PetscFPrintf(PETSC_COMM_WORLD, fd, " %d %g %g %dn", i, bdCtx->vertices[i*2], bdCtx->vertices[i*2+1], bdCtx->markers[i]);
280: }
281: for(i = 0; i < bdCtx->numSegments; i++) {
282: PetscFPrintf(PETSC_COMM_WORLD, fd, " %d %d %d %dn", i, bdCtx->segments[i*2], bdCtx->segments[i*2+1], bdCtx->segMarkers[i]);
283: }
284: for(i = 0; i < bdCtx->numHoles; i++) {
285: PetscFPrintf(PETSC_COMM_WORLD, fd, " %d %g %gn", i, bdCtx->holes[i*2], bdCtx->holes[i*2+1]);
286: }
287: return(0);
288: }
290: static int MeshBoundary2DView_Draw_Mesh(Mesh mesh, MeshBoundary2D *bdCtx, PetscDraw draw)
291: {
292: double *vertices = bdCtx->vertices;
293: int *segments = bdCtx->segments;
294: int *segMarkers = bdCtx->segMarkers;
295: int color;
296: int seg;
297: int ierr;
300: for(seg = 0; seg < bdCtx->numSegments; seg++) {
301: if (segMarkers[seg]) {
302: color = PETSC_DRAW_RED;
303: } else {
304: color = PETSC_DRAW_BLACK;
305: }
307: MeshDrawLine(mesh, draw, vertices[segments[seg*2]*2], vertices[segments[seg*2]*2+1],
308: vertices[segments[seg*2+1]*2], vertices[segments[seg*2+1]*2+1], color);
309:
310: }
311: PetscDrawSynchronizedFlush(draw);
312: return(0);
313: }
315: static int MeshBoundary2DView_Draw(Mesh mesh, MeshBoundary2D *bdCtx, PetscViewer v)
316: {
317: PetscReal startX, endX;
318: PetscReal startY, endY;
319: PetscReal sizeX, sizeY;
320: double *vertices;
321: PetscDraw draw;
322: PetscTruth isnull;
323: PetscDrawButton button;
324: char statusLine[256];
325: int pause;
326: PetscReal xc, yc, scale;
327: int xsize, ysize;
328: int vert;
329: int ierr;
332: PetscViewerDrawGetDraw(v, 0, &draw);
333: PetscDrawIsNull(draw, &isnull);
334: if (isnull) return(0);
335: MeshGetBoundingBox(mesh, &startX, &startY, PETSC_NULL, &endX, &endY, PETSC_NULL);
336: vertices = bdCtx->vertices;
337: for(vert = 0; vert < bdCtx->numVertices; vert++) {
338: if (startX > vertices[vert*2]) startX = vertices[vert*2];
339: if (startY > vertices[vert*2+1]) startY = vertices[vert*2+1];
340: if (endX < vertices[vert*2]) endX = vertices[vert*2];
341: if (endY < vertices[vert*2+1]) endY = vertices[vert*2+1];
342: }
343: sizeX = endX - startX;
344: sizeY = endY - startY;
345: if (sizeX > sizeY) {
346: ysize = 300;
347: xsize = ysize * (int) (sizeX/sizeY);
348: } else {
349: xsize = 300;
350: ysize = xsize * (int) (sizeY/sizeX);
351: }
352: PetscDrawResizeWindow(draw, xsize, ysize);
353: PetscDrawSynchronizedClear(draw);
355: ierr = PetscDrawSetCoordinates(draw, startX-0.02*sizeX, startY-0.02*sizeY, endX+0.02*sizeX, endY+0.02*sizeY);
356: ierr = MeshBoundary2DView_Draw_Mesh(mesh, bdCtx, draw);
357: ierr = PetscDrawGetPause(draw, &pause);
358: if (pause >= 0) {
359: PetscSleep(pause);
360: return(0);
361: }
363: /* Allow the mesh to zoom or shrink */
364: PetscDrawCheckResizedWindow(draw);
365: PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
366: while ((button != BUTTON_RIGHT) && (button != BUTTON_CENTER_SHIFT))
367: {
368: PetscDrawSynchronizedClear(draw);
369: statusLine[0] = 0;
370: switch (button)
371: {
372: case BUTTON_LEFT:
373: scale = 0.5;
374: break;
375: case BUTTON_CENTER:
376: scale = 2.0;
377: break;
378: default:
379: scale = 1.0;
380: break;
381: }
382: if (scale != 1.0) {
383: startX = scale*(startX + sizeX - xc) + xc - sizeX*scale;
384: endX = scale*(endX - sizeX - xc) + xc + sizeX*scale;
385: startY = scale*(startY + sizeY - yc) + yc - sizeY*scale;
386: endY = scale*(endY - sizeY - yc) + yc + sizeY*scale;
387: sizeX *= scale;
388: sizeY *= scale;
389: }
391: PetscDrawSetCoordinates(draw, startX-0.02*sizeX, startY-0.02*sizeY, endX+0.02*sizeX, endY+0.02*sizeY);
392: PetscDrawString(draw, startX - 0.02*sizeX, startY - 0.02*sizeY, PETSC_DRAW_BLACK, statusLine);
393: MeshBoundary2DView_Draw_Mesh(mesh, bdCtx, draw);
394: PetscDrawCheckResizedWindow(draw);
395: PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
396: }
398: return(0);
399: }
401: int MeshBoundary2DView(Mesh mesh, MeshBoundary2D *bdCtx, PetscViewer viewer)
402: {
403: PetscTruth isascii, isdraw;
404: int ierr;
407: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &isascii);
408: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_DRAW, &isdraw);
409: if (isascii == PETSC_TRUE) {
410: MeshBoundary2DView_File(mesh, bdCtx, viewer);
411: } else if (isdraw == PETSC_TRUE) {
412: MeshBoundary2DView_Draw(mesh, bdCtx, viewer);
413: }
415: return(0);
416: }
418: /*------------------------------------------- Mesh Boundary Functions ------------------------------------------------*/
419: int MeshBoundary2DSetBoundingBox(Mesh mesh, MeshBoundary2D *bdCtx)
420: {
421: MPI_Comm comm;
422: double startX, startY;
423: double endX, endY;
424: double *vertices;
425: int rank;
426: int p;
427: int ierr;
430: PetscObjectGetComm((PetscObject) mesh, &comm);
431: MPI_Comm_rank(comm, &rank);
432: if (rank == 0) {
433: vertices = bdCtx->vertices;
434: startX = endX = vertices[0];
435: startY = endY = vertices[1];
436: for(p = 0; p < bdCtx->numVertices; p++) {
437: /* Ignore moving boundaries (to accomodate split particles)? This should be done better */
438: if (bdCtx->markers[p] < 0) continue;
439: if (startX > vertices[p*2]) startX = vertices[p*2];
440: if (startY > vertices[p*2+1]) startY = vertices[p*2+1];
441: if (endX < vertices[p*2]) endX = vertices[p*2];
442: if (endY < vertices[p*2+1]) endY = vertices[p*2+1];
443: }
444: }
445: MPI_Bcast(&startX, 1, MPI_DOUBLE, 0, comm);
446: MPI_Bcast(&endX, 1, MPI_DOUBLE, 0, comm);
447: MPI_Bcast(&startY, 1, MPI_DOUBLE, 0, comm);
448: MPI_Bcast(&endY, 1, MPI_DOUBLE, 0, comm);
449: MeshSetBoundingBox(mesh, startX, startY, 0.0, endX, endY, 0.0);
450: return(0);
451: }