Actual source code: tri2d.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: tri2d.c,v 1.38 2000/10/17 13:48:55 knepley Exp $";
3: #endif
5: /* Implements 2d triangular grids */
6: #include "src/mesh/impls/triangular/2d/2dimpl.h" /*I "mesh.h" I*/
7: #include "tri2dView.h"
8: #include "tri2dCSR.h"
9: #ifdef PETSC_HAVE_TRIANGLE
10: #include "tri2d_Triangle.h"
11: #endif
13: extern int PartitionGhostNodeIndex_Private(Partition, int, int *);
14: EXTERN_C_BEGIN
15: int MeshSerialize_Triangular_2D(MPI_Comm, Mesh *, PetscViewer, PetscTruth);
16: EXTERN_C_END
18: /*--------------------------------------------- Public Functions ----------------------------------------------------*/
19: /*@C
20: MeshTriangular2DCalcAreas - Calculate and store the area of each element
22: Collective on Mesh
24: Input Parameters:
25: + mesh - The mesh
26: - create - The flag which creates the storage if it does not already exist
28: Note:
29: If areas have not been calculated before, and 'create' is PETSC_FALSE,
30: then no action taken.
32: Level: beginner
34: .keywords: mesh, element, area
35: .seealso: MeshTriangular2DCalcAspectRatios()
36: @*/
37: int MeshTriangular2DCalcAreas(Mesh mesh, PetscTruth create)
38: {
39: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
40: int numFaces = mesh->numFaces;
41: int numCorners = mesh->numCorners;
42: int *faces = tri->faces;
43: double *nodes = tri->nodes;
44: double x21, x31, y21, y31;
45: int elem;
46: int ierr;
49: if (tri->areas || create == PETSC_TRUE)
50: {
51: if (tri->areas) {
52: PetscFree(tri->areas);
53: }
54: PetscMalloc(numFaces * sizeof(double), &tri->areas);
55: for(elem = 0; elem < numFaces; elem++)
56: {
57: if (mesh->isPeriodic == PETSC_TRUE) {
58: x21 = MeshPeriodicDiffX(mesh, nodes[faces[elem*numCorners+1]*2] - nodes[faces[elem*numCorners]*2]);
59: y21 = MeshPeriodicDiffY(mesh, nodes[faces[elem*numCorners+1]*2+1] - nodes[faces[elem*numCorners]*2+1]);
60: x31 = MeshPeriodicDiffX(mesh, nodes[faces[elem*numCorners+2]*2] - nodes[faces[elem*numCorners]*2]);
61: y31 = MeshPeriodicDiffY(mesh, nodes[faces[elem*numCorners+2]*2+1] - nodes[faces[elem*numCorners]*2+1]);
62: } else {
63: x21 = nodes[faces[elem*numCorners+1]*2] - nodes[faces[elem*numCorners]*2];
64: y21 = nodes[faces[elem*numCorners+1]*2+1] - nodes[faces[elem*numCorners]*2+1];
65: x31 = nodes[faces[elem*numCorners+2]*2] - nodes[faces[elem*numCorners]*2];
66: y31 = nodes[faces[elem*numCorners+2]*2+1] - nodes[faces[elem*numCorners]*2+1];
67: }
68: tri->areas[elem] = PetscAbsReal(x21*y31 - x31*y21)/2.0;
69: }
70: PetscLogFlops(numFaces*8);
71: }
72: return(0);
73: }
75: /*@C
76: MeshTriangular2DCalcAspectRatios - Calculate and store the aspect ratio of each element
78: Collective on Mesh
80: Input Parameters:
81: + mesh - The mesh
82: - create - The flag which creates the storage if it does not already exist
84: Note:
85: If aspect ratios have not been calculated before, and 'create' is PETSC_FALSE,
86: then no action taken.
88: Level: beginner
90: .keywords: mesh, element, aspect ratio
91: .seealso: MeshTriangular2DCalcAreas()
92: @*/
93: int MeshTriangular2DCalcAspectRatios(Mesh mesh, PetscTruth create)
94: {
95: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
96: int numFaces = mesh->numFaces;
97: int numCorners = mesh->numCorners;
98: int *faces = tri->faces;
99: double *nodes = tri->nodes;
100: double x21, x31, x32, y21, y31, y32;
101: double area;
102: int elem;
103: int ierr;
106: if (tri->aspectRatios || create == PETSC_TRUE)
107: {
108: if (tri->aspectRatios) {
109: PetscFree(tri->aspectRatios);
110: }
111: PetscMalloc(numFaces * sizeof(double), &tri->aspectRatios);
112: for(elem = 0; elem < numFaces; elem++)
113: {
114: if (mesh->isPeriodic == PETSC_TRUE) {
115: x21 = MeshPeriodicDiffX(mesh, nodes[faces[elem*numCorners+1]*2] - nodes[faces[elem*numCorners]*2]);
116: y21 = MeshPeriodicDiffY(mesh, nodes[faces[elem*numCorners+1]*2+1] - nodes[faces[elem*numCorners]*2+1]);
117: x31 = MeshPeriodicDiffX(mesh, nodes[faces[elem*numCorners+2]*2] - nodes[faces[elem*numCorners]*2]);
118: y31 = MeshPeriodicDiffY(mesh, nodes[faces[elem*numCorners+2]*2+1] - nodes[faces[elem*numCorners]*2+1]);
119: x32 = MeshPeriodicDiffX(mesh, nodes[faces[elem*numCorners+2]*2] - nodes[faces[elem*numCorners+1]*2]);
120: y32 = MeshPeriodicDiffY(mesh, nodes[faces[elem*numCorners+2]*2+1] - nodes[faces[elem*numCorners+1]*2+1]);
121: } else {
122: x21 = nodes[faces[elem*numCorners+1]*2] - nodes[faces[elem*numCorners]*2];
123: y21 = nodes[faces[elem*numCorners+1]*2+1] - nodes[faces[elem*numCorners]*2+1];
124: x31 = nodes[faces[elem*numCorners+2]*2] - nodes[faces[elem*numCorners]*2];
125: y31 = nodes[faces[elem*numCorners+2]*2+1] - nodes[faces[elem*numCorners]*2+1];
126: x32 = nodes[faces[elem*numCorners+2]*2] - nodes[faces[elem*numCorners+1]*2];
127: y32 = nodes[faces[elem*numCorners+2]*2+1] - nodes[faces[elem*numCorners+1]*2+1];
128: }
129: area = PetscAbsReal(x21*y31 - x31*y21);
130: tri->aspectRatios[elem] = PetscMax(x32*x32 + y32*y32, PetscMax(x21*x21 + y21*y21, x31*x31 + y31*y31)) / area;
131: }
132: PetscLogFlops(numFaces*19);
133: }
134: return(0);
135: }
137: /*--------------------------------------------- Private Functions ---------------------------------------------------*/
138: static int MeshDestroy_Triangular_2D(Mesh mesh)
139: {
140: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
141: int ierr;
144: if (tri->nodes != PETSC_NULL) {
145: PetscFree(tri->nodes);
146: PetscFree(tri->markers);
147: PetscFree(tri->degrees);
148: }
149: if (tri->nodesOld != PETSC_NULL) {
150: PetscFree(tri->nodesOld);
151: }
152: if (tri->edges != PETSC_NULL) {
153: PetscFree(tri->edges);
154: PetscFree(tri->edgemarkers);
155: }
156: if (tri->faces != PETSC_NULL) {
157: PetscFree(tri->faces);
158: }
159: if (tri->neighbors != PETSC_NULL) {
160: PetscFree(tri->neighbors);
161: }
162: if (tri->bdNodes != PETSC_NULL) {
163: PetscFree(tri->bdNodes);
164: PetscFree(tri->bdEdges);
165: PetscFree(tri->bdMarkers);
166: PetscFree(tri->bdBegin);
167: PetscFree(tri->bdEdgeBegin);
168: }
169: #if 0
170: if (tri->bdCtx != PETSC_NULL) {
171: MeshBoundaryDestroy(tri->bdCtx);
172: }
173: if (tri->bdCtxNew != PETSC_NULL) {
174: MeshBoundaryDestroy(tri->bdCtxNew);
175: }
176: #endif
177: if (tri->areas != PETSC_NULL) {
178: PetscFree(tri->areas);
179: }
180: if (tri->aspectRatios != PETSC_NULL) {
181: PetscFree(tri->aspectRatios);
182: }
183: PetscFree(tri);
184: if (mesh->support != PETSC_NULL) {
185: PetscFree(mesh->support);
186: }
187: PartitionDestroy(mesh->part);
188: if (mesh->nodeOrdering != PETSC_NULL) {
189: AODestroy(mesh->nodeOrdering);
190: }
191: if (mesh->coarseMap != PETSC_NULL) {
192: AODestroy(mesh->coarseMap);
193: }
195: return(0);
196: }
198: int MeshDebug_Triangular_2D(Mesh mesh, PetscTruth distributed)
199: {
200: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
201: Partition p = mesh->part;
202: Partition_Triangular_2D *q = PETSC_NULL;
203: int numBd = mesh->numBd;
204: int numElements = mesh->numFaces;
205: int numCorners = mesh->numCorners;
206: int numNodes = mesh->numNodes;
207: int *elements = tri->faces;
208: int *markers = tri->markers;
209: int *bdMarkers = tri->bdMarkers;
210: int degree;
211: int *support;
212: int *degrees;
213: int *supports;
214: int *sizes;
215: int numProcs, rank;
216: PetscTruth isMidnode;
217: int proc, bd, elem, nElem, sElem, sElem2, corner, node, newNode, size;
218: int ierr;
221: MPI_Comm_size(mesh->comm, &numProcs);
222: MPI_Comm_rank(mesh->comm, &rank);
223: if (p != PETSC_NULL) {
224: q = (Partition_Triangular_2D *) p->data;
225: }
227: /* Check basic sizes */
228: PetscMalloc(numProcs * sizeof(int), &sizes);
229: MPI_Allgather(&mesh->numBd, 1, MPI_INT, sizes, 1, MPI_INT, mesh->comm);
230: for(proc = 0, 0; proc < numProcs-1; proc++) {
231: if (sizes[proc] != sizes[proc+1]) {
232: PetscPrintf(mesh->comm, "The number of boundaries is different on proc %d(%d) and proc %d(%d)n",
233: proc, sizes[proc], proc+1, sizes[proc+1]);
234: 1;
235: }
236: }
237:
238: #if 0
239: /* I made vertices distributed in linear meshes */
240: MPI_Allgather(&tri->numVertices, 1, MPI_INT, sizes, 1, MPI_INT, mesh->comm);
241: for(proc = 0, 0; proc < numProcs-1; proc++) {
242: if (sizes[proc] != sizes[proc+1]) {
243: PetscPrintf(mesh->comm, "The number of vertices is different on proc %d(%d) and proc %d(%d)n",
244: proc, sizes[proc], proc+1, sizes[proc+1]);
245: 1;
246: }
247: }
248:
249: #endif
250: MPI_Allgather(&mesh->numCorners, 1, MPI_INT, sizes, 1, MPI_INT, mesh->comm);
251: for(proc = 0, 0; proc < numProcs-1; proc++) {
252: if (sizes[proc] != sizes[proc+1]) {
253: PetscPrintf(mesh->comm, "The number of face corners is different on proc %d(%d) and proc %d(%d)n",
254: proc, sizes[proc], proc+1, sizes[proc+1]);
255: 1;
256: }
257: }
258:
259: MPI_Allgather(&mesh->numBdEdges, 1, MPI_INT, sizes, 1, MPI_INT, mesh->comm);
260: for(proc = 0, 0; proc < numProcs-1; proc++) {
261: if (sizes[proc] != sizes[proc+1]) {
262: PetscPrintf(mesh->comm, "The number of boundary edges is different on proc %d(%d) and proc %d(%d)n",
263: proc, sizes[proc], proc+1, sizes[proc+1]);
264: 1;
265: }
266: }
267:
268: if (distributed == PETSC_FALSE) {
269: MPI_Allgather(&mesh->numNodes, 1, MPI_INT, sizes, 1, MPI_INT, mesh->comm);
270: for(proc = 0, 0; proc < numProcs-1; proc++) {
271: if (sizes[proc] != sizes[proc+1]) {
272: PetscPrintf(mesh->comm, "The number of nodes is different on proc %d(%d) and proc %d(%d)n",
273: proc, sizes[proc], proc+1, sizes[proc+1]);
274: 1;
275: }
276: }
277:
278: MPI_Allgather(&mesh->numBdNodes, 1, MPI_INT, sizes, 1, MPI_INT, mesh->comm);
279: for(proc = 0, 0; proc < numProcs-1; proc++) {
280: if (sizes[proc] != sizes[proc+1]) {
281: PetscPrintf(mesh->comm, "The number of boundary nodes is different on proc %d(%d) and proc %d(%d)n",
282: proc, sizes[proc], proc+1, sizes[proc+1]);
283: 1;
284: }
285: }
286:
287: MPI_Allgather(&mesh->numEdges, 1, MPI_INT, sizes, 1, MPI_INT, mesh->comm);
288: for(proc = 0, 0; proc < numProcs-1; proc++) {
289: if (sizes[proc] != sizes[proc+1]) {
290: PetscPrintf(mesh->comm, "The number of edges is different on proc %d(%d) and proc %d(%d)n",
291: proc, sizes[proc], proc+1, sizes[proc+1]);
292: 1;
293: }
294: }
295:
296: MPI_Allgather(&mesh->numFaces, 1, MPI_INT, sizes, 1, MPI_INT, mesh->comm);
297: for(proc = 0, 0; proc < numProcs-1; proc++) {
298: if (sizes[proc] != sizes[proc+1]) {
299: PetscPrintf(mesh->comm, "The number of faces is different on proc %d(%d) and proc %d(%d)n",
300: proc, sizes[proc], proc+1, sizes[proc+1]);
301: 1;
302: }
303: }
304:
305: } else {
306: MPI_Allreduce(&mesh->numNodes, &size, 1, MPI_INT, MPI_SUM, mesh->comm);
307: if (size != q->numNodes) {
308: SETERRQ2(PETSC_ERR_PLIB, "The total number of nodes %d should be %dn", size, q->numNodes);
309: }
310: #if 0
311: MPI_Allreduce(&tri->numBdNodes, &size, 1, MPI_INT, MPI_SUM, mesh->comm);
312: #else
313: size = mesh->numBdNodes;
314: #endif
315: if (size != q->numBdNodes) {
316: SETERRQ2(PETSC_ERR_PLIB, "The total number of boundary nodes %d should be %dn", size, q->numBdNodes);
317: }
318: MPI_Allreduce(&mesh->numEdges, &size, 1, MPI_INT, MPI_SUM, mesh->comm);
319: if (size != q->numEdges) {
320: SETERRQ2(PETSC_ERR_PLIB, "The total number of edges %d should be %dn", size, q->numEdges);
321: }
322: MPI_Allreduce(&mesh->numFaces, &size, 1, MPI_INT, MPI_SUM, mesh->comm);
323: if (size != mesh->part->numElements) {
324: SETERRQ2(PETSC_ERR_PLIB, "The total number of faces %d should be %dn", size, mesh->part->numElements);
325: }
326: }
327: PetscFree(sizes);
329: if (distributed == PETSC_FALSE) {
330: /* Check markers */
331: for(node = 0; node < numNodes; node++) {
332: if (!markers[node]) continue;
334: for(bd = 0; bd < numBd; bd++)
335: if(bdMarkers[bd] == markers[node])
336: break;
337: if (bd == numBd) {
338: SETERRQ2(PETSC_ERR_PLIB, "The marker %d for node %d is invalidn", markers[node], node);
339: }
340: }
341: /* Check mesh connectivity */
342: for(node = 0; node < numNodes; node++) {
343: MeshGetNodeSupport(mesh, node, 0, °ree, &support);
345: /* Check node degree */
346: PetscMalloc(numProcs * sizeof(int), °rees);
347: MPI_Allgather(°ree, 1, MPI_INT, degrees, 1, MPI_INT, mesh->comm);
348: for(proc = 0, 0; proc < numProcs-1; proc++)
349: if (degrees[proc] != degrees[proc+1]) {
350: PetscPrintf(mesh->comm, "Degree of node %d is different on proc %d(%d) and proc %d(%d)n",
351: node, proc, degrees[proc], proc+1, degrees[proc+1]);
352: PetscPrintf(mesh->comm, "Node Information:n");
353: PetscSequentialPhaseBegin(mesh->comm, 1);
354: PetscPrintf(PETSC_COMM_SELF, " Processor %dn", rank);
355: if ((rank == proc) || (rank == proc+1))
356: for(sElem = 0; sElem < degree; sElem++)
357: PetscPrintf(PETSC_COMM_SELF, " Support : %dn", support[sElem]);
358: PetscPrintf(PETSC_COMM_SELF, " Marker : %dn", tri->markers[node]);
359: PetscPrintf(PETSC_COMM_SELF, " Location : (%g,%g)n", tri->nodes[node*2], tri->nodes[node*2+1]);
360: PetscPrintf(PETSC_COMM_SELF, " In Elements:");
361: for(elem = 0, isMidnode = PETSC_FALSE; elem < numElements; elem++)
362: for(corner = 0; corner < numCorners; corner++)
363: if (elements[elem*numCorners+corner] == node) {
364: PetscPrintf(PETSC_COMM_SELF, " %d", elem);
365: if (corner >= 3)
366: isMidnode = PETSC_TRUE;
367: }
368: PetscPrintf(PETSC_COMM_SELF, "n");
369: if (isMidnode == PETSC_TRUE)
370: PetscPrintf(PETSC_COMM_SELF, " Midnoden");
371: else
372: PetscPrintf(PETSC_COMM_SELF, " Vertexn");
373: PetscSequentialPhaseEnd(mesh->comm, 1);
374: 1;
375: }
376:
377: PetscFree(degrees);
379: /* Check support elements */
380: PetscMalloc(degree*numProcs * sizeof(int), &supports);
381: MPI_Allgather(support, degree, MPI_INT, supports, degree, MPI_INT, mesh->comm);
382: for(sElem = 0, 0; sElem < degree; sElem++) {
383: nElem = supports[sElem];
384: for(proc = 1; proc < numProcs; proc++) {
385: for(sElem2 = 0; sElem2 < degree; sElem2++)
386: if (supports[proc*degree+sElem2] == nElem)
387: break;
388: if (sElem2 == degree) {
389: PetscPrintf(mesh->comm, "Support element %d of node %d on proc %d(%d) is not present on proc %dn",
390: sElem, node, 0, supports[sElem], proc);
391: PetscPrintf(mesh->comm, " Support of node %d on proc %d:n ", node, proc);
392: for(sElem2 = 0; sElem2 < degree; sElem2++)
393: PetscPrintf(mesh->comm, " %d", supports[proc*degree+sElem2]);
394: PetscPrintf(mesh->comm, "n");
395: 1;
396: }
397: }
398: }
399:
400: PetscFree(supports);
402: /* Check that node only appears inside elements in the support */
403: for(elem = 0, 0; elem < numElements; elem++)
404: for(corner = 0; corner < numCorners; corner++) {
405: newNode = elements[elem*numCorners+corner];
406: if (node == newNode) {
407: for(sElem = 0; sElem < degree; sElem++)
408: if (support[sElem] == elem)
409: break;
410: if (sElem == degree) {
411: PetscPrintf(mesh->comm, "Node %d found in element %d which is not present in the supportn",
412: node, elem);
413: 1;
414: }
415: }
416: }
417:
419: MeshRestoreNodeSupport(mesh, node, 0, °ree, &support);
420: }
421: } else {
422: /* Check markers */
423: for(node = 0; node < q->numOverlapNodes; node++) {
424: if (!markers[node]) continue;
426: for(bd = 0; bd < numBd; bd++)
427: if(bdMarkers[bd] == markers[node])
428: break;
429: if (bd == numBd) {
430: SETERRQ2(PETSC_ERR_PLIB, "The marker %d for node %d is invalidn", markers[node], node);
431: }
432: }
433: /* Check mesh connectivity */
434: for(node = 0; node < q->numLocNodes; node++) {
435: MeshGetNodeSupport(mesh, node, 0, °ree, &support);
437: /* Check that node only appears inside elements in the support */
438: for(elem = 0, 0; elem < p->numOverlapElements; elem++) {
439: for(corner = 0; corner < numCorners; corner++) {
440: newNode = elements[elem*numCorners+corner];
441: if (node == newNode) {
442: for(sElem = 0; sElem < degree; sElem++) {
443: if (support[sElem] == elem) break;
444: }
445: if (sElem == degree) {
446: PetscPrintf(PETSC_COMM_SELF, "[%d]Node %d found in element %d which is not present in the supportn",
447: p->rank, node, elem);
448: 1;
449: }
450: }
451: }
452: }
453:
455: MeshRestoreNodeSupport(mesh, node, 0, °ree, &support);
456: }
457: }
458: return(0);
459: }
461: int MeshSetupSupport_Triangular_2D(Mesh mesh)
462: {
463: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
464: int edge, node;
465: int ierr;
468: /* Calculate maximum degree of vertices */
469: if (mesh->numNodes > 0) {
470: PetscMalloc(mesh->numNodes * sizeof(int), &tri->degrees);
471: }
472: PetscMemzero(tri->degrees, mesh->numNodes * sizeof(int));
473: for(edge = 0; edge < mesh->numEdges*2; edge++) {
474: tri->degrees[tri->edges[edge]]++;
475: }
476: for(node = 0, mesh->maxDegree = 0; node < mesh->numNodes; node++) {
477: mesh->maxDegree = PetscMax(mesh->maxDegree, tri->degrees[node]);
478: }
479: if (mesh->maxDegree > 0) {
480: PetscMalloc(mesh->maxDegree * sizeof(int), &mesh->support);
481: }
482: mesh->supportTaken = PETSC_FALSE;
483: return(0);
484: }
486: int MeshCheckBoundary_Triangular_2D(Mesh mesh) {
487: MeshBoundary2D *bdCtx = mesh->bdCtx;
488: int *markers = bdCtx->markers;
489: int *segments = bdCtx->segments;
490: int *segMarkers = bdCtx->segMarkers;
491: int numBd = 0;
492: int numBdVertices = 0;
493: int numBdSegments = 0;
494: int *bdMarkers;
495: int *bdBegin;
496: int *bdSegmentBegin;
497: int bd, vertex, segment, marker, rank;
498: int ierr;
502: MPI_Comm_rank(mesh->comm, &rank);
503: if (rank != 0) return(0);
504: PetscMalloc(bdCtx->numBd * sizeof(int), &bdMarkers);
505: PetscMalloc((bdCtx->numBd+1) * sizeof(int), &bdBegin);
506: PetscMalloc((bdCtx->numBd+1) * sizeof(int), &bdSegmentBegin);
507: PetscMemzero(bdBegin, (bdCtx->numBd+1) * sizeof(int));
508: PetscMemzero(bdSegmentBegin, (bdCtx->numBd+1) * sizeof(int));
509: for(vertex = 0; vertex < bdCtx->numVertices; vertex++) {
510: if (markers[vertex] == 0) continue;
511: numBdVertices++;
512: /* Check for new marker */
513: for(bd = 0; bd < numBd; bd++) {
514: if (markers[vertex] == bdMarkers[bd]) break;
515: }
516: if (bd == numBd) {
517: if (numBd >= bdCtx->numBd) SETERRQ1(PETSC_ERR_ARG_CORRUPT, "More markers present than declared: %d", bdCtx->numBd);
518: /* Insert new marker */
519: for(bd = 0; bd < numBd; bd++) {
520: if (markers[vertex] < bdMarkers[bd]) break;
521: }
522: if (bd < numBd) {
523: PetscMemmove(&bdMarkers[bd+1], &bdMarkers[bd], (numBd-bd) * sizeof(int));
524: PetscMemmove(&bdBegin[bd+2], &bdBegin[bd+1], (numBd-bd) * sizeof(int));
525: bdBegin[bd+1] = 0;
526: bdSegmentBegin[bd+1] = 0;
527: }
528: bdMarkers[bd] = markers[vertex];
529: numBd++;
530: bdBegin[bd+1]++;
531: } else {
532: bdBegin[bd+1]++;
533: }
534: }
535: for(segment = 0; segment < bdCtx->numSegments; segment++) {
536: if (segMarkers[segment] == 0) continue;
537: marker = segMarkers[segment];
538: for(bd = 0; bd < numBd; bd++) {
539: if (marker == bdMarkers[bd]) break;
540: }
541: if (bd == numBd) SETERRQ2(PETSC_ERR_ARG_CORRUPT, "Invalid segment marker %d on segment %d", marker, segment);
542: numBdSegments++;
543: bdSegmentBegin[bd+1]++;
544: if (markers[segments[segment*2]] != marker) {
545: SETERRQ4(PETSC_ERR_PLIB, "Marker %d on left vertex %d does not match marker %d on its segment %d",
546: markers[segments[segment*2]], segments[segment*2], marker, segment);
547: }
548: if (markers[segments[segment*2+1]] != marker) {
549: SETERRQ4(PETSC_ERR_PLIB, "Marker %d on right vertex %d does not match marker %d on its segment %d",
550: markers[segments[segment*2+1]], segments[segment*2+1], marker, segment);
551: }
552: }
554: /* Do prefix sums to get position offsets */
555: for(bd = 2; bd <= numBd; bd++) {
556: bdBegin[bd] = bdBegin[bd-1] + bdBegin[bd];
557: bdSegmentBegin[bd] = bdSegmentBegin[bd-1] + bdSegmentBegin[bd];
558: }
560: if (numBd != bdCtx->numBd) {
561: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of boundaries %d should be %d", numBd, bdCtx->numBd);
562: }
563: if (bdBegin[numBd] != numBdVertices) {
564: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of boundary vertices %d should be %d", bdBegin[numBd], numBdVertices);
565: }
566: if (bdSegmentBegin[numBd] != numBdSegments) {
567: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of boundary segments %d should be %d", bdSegmentBegin[numBd], numBdSegments);
568: }
569: PetscFree(bdMarkers);
570: PetscFree(bdBegin);
571: PetscFree(bdSegmentBegin);
572: return(0);
573: }
575: int MeshSetupBoundary_Triangular_2D(Mesh mesh)
576: {
577: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
578: int numBd = mesh->numBd;
579: int numNodes = mesh->numNodes;
580: int *markers = tri->markers;
581: int numEdges = mesh->numEdges;
582: int *edges = tri->edges;
583: int *edgemarkers = tri->edgemarkers;
584: int numFaces = mesh->numFaces;
585: int numCorners = mesh->numCorners;
586: int *faces = tri->faces;
587: int newNumBd; /* Current number of different boundary markers */
588: int numBdEdges; /* Current offset into bdEdges[] */
589: int *bdNodes; /* Node numbers for boundary nodes ordered by boundary */
590: int *bdEdges; /* Node numbers for boundary edges ordered by boundary */
591: int *bdBeginOff; /* Current offset into the bdNodes or bdEdges array */
592: int *bdSeen; /* Flags for boundaries already processed */
593: PetscTruth closed; /* Indicates whether a boundary has been closed */
594: int degree;
595: int *support;
596: int rank, bd, marker, node, nextNode, midnode, bdNode, nextBdNode, midBdNode, edge, bdEdge;
597: int elem, cElem, supElem, corner, supCorner, tmp;
598: int ierr;
601: MPI_Comm_rank(mesh->comm, &rank);
602: PetscLogEventBegin(MESH_SetupBoundary, mesh, 0, 0, 0);
603: tri->areas = PETSC_NULL;
604: tri->aspectRatios = PETSC_NULL;
605: mesh->numBdNodes = 0;
606: mesh->numBdEdges = 0;
607: if (numBd == 0) {
608: tri->bdMarkers = PETSC_NULL;
609: tri->bdBegin = PETSC_NULL;
610: tri->bdEdgeBegin = PETSC_NULL;
611: tri->bdNodes = PETSC_NULL;
612: tri->bdEdges = PETSC_NULL;
613: return(0);
614: }
615: PetscMalloc(numBd * sizeof(int), &tri->bdMarkers);
616: PetscMalloc((numBd+1) * sizeof(int), &tri->bdBegin);
617: PetscMalloc((numBd+1) * sizeof(int), &tri->bdEdgeBegin);
618: PetscMalloc((numBd+1) * sizeof(int), &bdBeginOff);
619: PetscLogObjectMemory(mesh, (numBd + (numBd+1)*2) * sizeof(int));
620: PetscMemzero(tri->bdMarkers, numBd * sizeof(int));
621: PetscMemzero(tri->bdBegin, (numBd+1) * sizeof(int));
622: PetscMemzero(tri->bdEdgeBegin, (numBd+1) * sizeof(int));
623: PetscMemzero(bdBeginOff, (numBd+1) * sizeof(int));
625: for(node = 0, newNumBd = 0; node < numNodes; node++) {
626: /* Get number of boundary nodes and markers */
627: if (markers[node]) {
628: mesh->numBdNodes++;
629: /* Check for new marker */
630: for(bd = 0; bd < newNumBd; bd++) {
631: if (markers[node] == tri->bdMarkers[bd]) break;
632: }
633: if (bd == newNumBd) {
634: /* Insert new marker */
635: for(bd = 0; bd < newNumBd; bd++) {
636: if (markers[node] < tri->bdMarkers[bd]) break;
637: }
638: if (bd < newNumBd) {
639: PetscMemmove(&tri->bdMarkers[bd+1], &tri->bdMarkers[bd], (newNumBd-bd) * sizeof(int));
640: PetscMemmove(&tri->bdBegin[bd+2], &tri->bdBegin[bd+1], (newNumBd-bd) * sizeof(int));
641: tri->bdBegin[bd+1] = 0;
642: }
643: tri->bdMarkers[bd] = markers[node];
644: newNumBd++;
645: tri->bdBegin[bd+1]++;
646: } else {
647: tri->bdBegin[bd+1]++;
648: }
649: }
650: }
652: /* Do prefix sums to get position offsets */
653: for(bd = 2; bd <= numBd; bd++) {
654: tri->bdBegin[bd] = tri->bdBegin[bd-1] + tri->bdBegin[bd];
655: }
657: /* This is shameful -- I will write the code to look over edges soon */
658: if (numCorners == 3) {
659: PetscMemcpy(tri->bdEdgeBegin, tri->bdBegin, (numBd+1) * sizeof(int));
660: } else if (numCorners == 6) {
661: /* We know that there are an even number of nodes in every boundary */
662: for(bd = 0; bd <= numBd; bd++) {
663: if (tri->bdBegin[bd]%2) {
664: SETERRQ2(PETSC_ERR_PLIB, "There are %d nodes on boundary %d (not divisible by two)",
665: tri->bdBegin[bd]-tri->bdBegin[bd-1], tri->bdMarkers[bd]);
666: } else {
667: tri->bdEdgeBegin[bd] = tri->bdBegin[bd]/2;
668: }
669: }
670: } else {
671: SETERRQ1(PETSC_ERR_SUP, "Number of local nodes %d not supported", numCorners);
672: }
674: /* Get number of boundary edges */
675: for(edge = 0; edge < numEdges; edge++) {
676: marker = edgemarkers[edge];
677: if (marker) {
678: mesh->numBdEdges++;
679: MeshGetMidnodeFromEdge(mesh, edge, &midnode);
680: if (markers[edges[edge*2]] != marker) {
681: if ((markers[edges[edge*2+1]] == marker) &&
682: ((markers[midnode] == markers[edges[edge*2]]) || (markers[midnode] == markers[edges[edge*2+1]]))) {
683: /* I assume here that the generator mistakenly included an edge between two boundaries */
684: PetscPrintf(PETSC_COMM_SELF, "[%d]Removing edge %d between boundaries %d and %d from boundaryn",
685: rank, edge, markers[edges[edge*2]], markers[edges[edge*2+1]]);
686: edgemarkers[edge] = 0;
687: markers[midnode] = 0;
688: continue;
689: } else {
690: SETERRQ5(PETSC_ERR_PLIB, "Marker %d on left node %d does not match marker %d on its edge %d, right marker is %d",
691: markers[edges[edge*2]], edges[edge*2], marker, edge, markers[edges[edge*2+1]]);
692: }
693: }
694: if (markers[edges[edge*2+1]] != marker) {
695: if ((markers[edges[edge*2]] == marker) &&
696: ((markers[midnode] == markers[edges[edge*2]]) || (markers[midnode] == markers[edges[edge*2+1]]))) {
697: /* I assume here that the generator mistakenly included an edge between two boundaries */
698: PetscPrintf(PETSC_COMM_SELF, "[%d]Removing edge %d between boundaries %d and %d from boundaryn",
699: rank, edge, markers[edges[edge*2]], markers[edges[edge*2+1]]);
700: edgemarkers[edge] = 0;
701: markers[midnode] = 0;
702: continue;
703: } else {
704: SETERRQ5(PETSC_ERR_PLIB, "Marker %d on right node %d does not match marker %d on its edge %d, left marker is %d",
705: markers[edges[edge*2+1]], edges[edge*2+1], marker, edge, markers[edges[edge*2]]);
706: }
707: }
708: if (markers[midnode] != marker) {
709: SETERRQ5(PETSC_ERR_PLIB, "Marker %d on midnode %d does not match marker %d on its edge %d, left marker is %d",
710: markers[midnode], midnode, marker, edge, markers[edges[edge*2]]);
711: }
712: }
713: }
715: /* Check boundary information consistency */
716: if (newNumBd != numBd) {
717: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of boundaries %d should be %d", newNumBd, numBd);
718: }
719: if (tri->bdBegin[numBd] != mesh->numBdNodes) {
720: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of boundary nodes %d should be %d", tri->bdBegin[numBd], mesh->numBdNodes);
721: }
722: if (tri->bdEdgeBegin[numBd] != mesh->numBdEdges) {
723: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of boundary edges %d should be %d", tri->bdEdgeBegin[numBd], mesh->numBdEdges);
724: }
726: PetscMalloc(mesh->numBdNodes * sizeof(int), &bdNodes);
727: PetscMalloc(mesh->numBdEdges * sizeof(int), &bdEdges);
728: PetscLogObjectMemory(mesh, (mesh->numBdNodes + mesh->numBdEdges) * sizeof(int));
730: /* Split nodes by marker */
731: PetscMemcpy(bdBeginOff, tri->bdBegin, (numBd+1) * sizeof(int));
732: for(node = 0; node < numNodes; node++) {
733: for(bd = 0; bd < numBd; bd++) {
734: if (markers[node] == tri->bdMarkers[bd]) {
735: #ifdef MESH_TRACE_BOUNDARY_CREATE
736: PetscPrintf(PETSC_COMM_WORLD, "bd: %d bdNode[%d] = %dn", bd, bdBeginOff[bd], node);
737: #endif
738: bdNodes[bdBeginOff[bd]++] = node;
739: }
740: }
741: }
742: for(bd = 0; bd < numBd; bd++) {
743: if (tri->bdBegin[bd+1] != bdBeginOff[bd])
744: SETERRQ(PETSC_ERR_PLIB, "Invalid boundary node marker information");
745: }
747: /* Split edges by marker */
748: PetscMemcpy(bdBeginOff, tri->bdEdgeBegin, (numBd+1) * sizeof(int));
749: for(edge = 0; edge < numEdges; edge++) {
750: for(bd = 0; bd < numBd; bd++) {
751: if (edgemarkers[edge] == tri->bdMarkers[bd]) {
752: #ifdef MESH_TRACE_BOUNDARY_CREATE
753: PetscPrintf(PETSC_COMM_WORLD, "bd: %d bdEdge[%d] = %dn", bd, bdBeginOff[bd], edge);
754: #endif
755: bdEdges[bdBeginOff[bd]++] = edge;
756: }
757: }
758: }
759: for(bd = 0; bd < numBd; bd++) {
760: if (tri->bdEdgeBegin[bd+1] != bdBeginOff[bd])
761: SETERRQ(PETSC_ERR_PLIB, "Invalid boundary edge marker information");
762: }
763: PetscFree(bdBeginOff);
765: PetscMalloc(numBd * sizeof(int), &bdSeen);
766: PetscMemzero(bdSeen, numBd * sizeof(int));
767: /* Loop over elements */
768: for(elem = 0; elem < numFaces; elem++) {
769: /* Find an element with a node on the boundary */
770: for(corner = 0; corner < 3; corner++) {
771: if (markers[faces[elem*numCorners+corner]]) {
772: /* Find boundary index */
773: cElem = elem;
774: node = faces[elem*numCorners+corner];
775: for(bd = 0; bd < numBd; bd++)
776: if (markers[node] == tri->bdMarkers[bd])
777: break;
778: if (bd == numBd) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary marker");
779: /* Check for processed boundaries */
780: if (bdSeen[bd])
781: continue;
782: numBdEdges = tri->bdEdgeBegin[bd];
784: /* Start building this boundary */
785: #ifdef MESH_TRACE_BOUNDARY_CREATE
786: PetscPrintf(mesh->comm, "Starting boundary %d beginning at %d ending before %dn",
787: bd, tri->bdBegin[bd], tri->bdBegin[bd+1]);
788: #endif
789: /* Find initial node */
790: for(bdNode = tri->bdBegin[bd]; bdNode < tri->bdBegin[bd+1]; bdNode++) {
791: if (bdNodes[bdNode] == node) {
792: /* Move this node to the head of the list */
793: #ifdef MESH_TRACE_BOUNDARY_CREATE
794: PetscPrintf(mesh->comm, " moving node %d from %d to %dn", bdNodes[bdNode], bdNode, tri->bdBegin[bd]);
795: #endif
796: tmp = bdNodes[tri->bdBegin[bd]];
797: bdNodes[tri->bdBegin[bd]] = bdNodes[bdNode];
798: bdNodes[bdNode] = tmp;
799: break;
800: }
801: }
803: /* Order edges counterclockwise around a boundary */
804: /* I do not currently check the orientation of the constructed boundary */
805: for(bdNode = tri->bdBegin[bd], closed = PETSC_FALSE; bdNode < tri->bdBegin[bd+1]; bdNode++) {
806: #ifdef MESH_TRACE_BOUNDARY_CREATE
807: PetscPrintf(mesh->comm, "n At boundary node %d x:%lf y: %lfn",
808: bdNode, tri->nodes[bdNodes[bdNode]*2], tri->nodes[bdNodes[bdNode]*2+1]);
809: #ifdef MESH_TRACE_BOUNDARY_CREATE_DETAIL
810: for(node = tri->bdBegin[bd]; node < tri->bdBegin[bd+1]; node++)
811: PetscPrintf(mesh->comm, " bdNode[%d]: %dn", node, bdNodes[node]);
812: #endif
813: #endif
815: /* Find a neighbor of the point -- Could maybe do better than linear search */
816: for(bdEdge = numBdEdges; bdEdge < tri->bdEdgeBegin[bd+1]; bdEdge++) {
817: edge = bdEdges[bdEdge];
818: if ((edgemarkers[edge] == tri->bdMarkers[bd]) &&
819: (((edges[edge*2] == bdNodes[bdNode]) && (markers[edges[edge*2+1]] == tri->bdMarkers[bd])) ||
820: ((edges[edge*2+1] == bdNodes[bdNode]) && (markers[edges[edge*2]] == tri->bdMarkers[bd]))))
821: {
822: /* Get neighboring node number */
823: if (edges[edge*2] == bdNodes[bdNode]) {
824: node = edges[edge*2+1];
825: } else {
826: node = edges[edge*2];
827: }
829: /* Find neighboring node in bdNodes[] */
830: for(nextBdNode = tri->bdBegin[bd]; nextBdNode < tri->bdBegin[bd+1]; nextBdNode++)
831: if (bdNodes[nextBdNode] == node) break;
832: #ifdef MESH_TRACE_BOUNDARY_CREATE
833: PetscPrintf(mesh->comm, " found connection along edge %d to bd node %d = node %dn",
834: edge, nextBdNode, node);
835: #endif
837: if (nextBdNode > bdNode) {
838: /* Insert midnode */
839: if (numCorners == 6) {
840: /* Move this node next to the connected one */
841: #ifdef MESH_TRACE_BOUNDARY_CREATE
842: PetscPrintf(mesh->comm, " moving node %d from %d to %dn", bdNodes[nextBdNode], nextBdNode, bdNode+2);
843: #endif
844: tmp = bdNodes[bdNode+2];
845: bdNodes[bdNode+2] = bdNodes[nextBdNode];
846: bdNodes[nextBdNode] = tmp;
847: nextBdNode = bdNode+2;
849: /* Walk around the node looking for a boundary edge and reset containing element */
850: node = bdNodes[bdNode];
851: nextNode = bdNodes[nextBdNode];
852: ierr = MeshGetNodeSupport(mesh, node, cElem, °ree, &support);
853: for(supElem = 0, midnode = -1; supElem < degree; supElem++) {
854: for(supCorner = 0; supCorner < 3; supCorner++) {
855: cElem = support[supElem];
856: if ((faces[cElem*numCorners+supCorner] == node) &&
857: (faces[cElem*numCorners+((supCorner+1)%3)] == nextNode))
858: {
859: midnode = faces[cElem*numCorners+((((supCorner*2+1)%3)*2)%3 + 3)];
860: supElem = degree;
861: break;
862: }
863: else if ((faces[cElem*numCorners+supCorner] == node) &&
864: (faces[cElem*numCorners+((supCorner+2)%3)] == nextNode))
865: {
866: midnode = faces[cElem*numCorners+((((supCorner*2+2)%3)*2)%3 + 3)];
867: supElem = degree;
868: break;
869: }
870: }
871: }
872: ierr = MeshRestoreNodeSupport(mesh, node, cElem, °ree, &support);
873: /* Find midnode in bdNodes[] */
874: for(midBdNode = tri->bdBegin[bd]; midBdNode < tri->bdBegin[bd+1]; midBdNode++)
875: if (bdNodes[midBdNode] == midnode)
876: break;
877: if (midBdNode == tri->bdBegin[bd+1]) SETERRQ(PETSC_ERR_PLIB, "Unable to locate midnode on boundary");
878: #ifdef MESH_TRACE_BOUNDARY_CREATE
879: PetscPrintf(mesh->comm, " moving midnode %d in elem %d from %d to %dn",
880: midnode, cElem, midBdNode, bdNode+1);
881: #endif
882: tmp = bdNodes[bdNode+1];
883: bdNodes[bdNode+1] = bdNodes[midBdNode];
884: bdNodes[midBdNode] = tmp;
885: bdNode++;
886: } else {
887: /* numCorners == 3 */
888: /* Move this node next to the connected one */
889: #ifdef MESH_TRACE_BOUNDARY_CREATE
890: PetscPrintf(mesh->comm, " moving node %d from %d to %dn", bdNodes[nextBdNode], nextBdNode, bdNode+1);
891: #endif
892: tmp = bdNodes[bdNode+1];
893: bdNodes[bdNode+1] = bdNodes[nextBdNode];
894: bdNodes[nextBdNode] = tmp;
895: }
897: /* Reorder bdEdges[] */
898: #ifdef MESH_TRACE_BOUNDARY_CREATE
899: PetscPrintf(mesh->comm, " moving edge %d from %d to %dn", edge, bdEdge, numBdEdges);
900: #endif
901: tmp = bdEdges[numBdEdges];
902: bdEdges[numBdEdges] = bdEdges[bdEdge];
903: bdEdges[bdEdge] = tmp;
904: numBdEdges++;
905: break;
906: }
907: /* Check that first and last nodes are connected (closed boundary) */
908: else if ((nextBdNode == tri->bdBegin[bd]) && (bdNode != (numCorners == 6 ? nextBdNode+2 : nextBdNode+1)))
909: {
910: if (bdEdges[numBdEdges++] != edge) SETERRQ(PETSC_ERR_PLIB, "Invalid closing edge");
911: closed = PETSC_TRUE;
912: /* End the loop since only a midnode is left */
913: bdNode = tri->bdBegin[bd+1];
914: #ifdef MESH_TRACE_BOUNDARY_CREATE
915: PetscPrintf(mesh->comm, " adding edge %d total: %dn", edge, numBdEdges);
916: #endif
917: break;
918: }
919: }
920: }
921: if (bdEdge == tri->bdEdgeBegin[bd+1]) SETERRQ(PETSC_ERR_PLIB, "No connection");
922: }
924: /* Check for closed boundary */
925: if (closed == PETSC_FALSE) {
926: PetscPrintf(PETSC_COMM_SELF, "Boundary %d with marker %d is not closedn", bd, tri->bdMarkers[bd]);
927: }
929: /* Check size of boundary */
930: if (tri->bdBegin[bd+1] != numBdEdges*(numCorners == 3 ? 1 : 2)) {
931: PetscPrintf(PETSC_COMM_SELF, "On boundary %d with marker %d, edge count %d x %d should equal %d (was %d)n",
932: bd, tri->bdMarkers[bd], numBdEdges, (numCorners == 3 ? 1 : 2), tri->bdBegin[bd+1], tri->bdBegin[bd]);
933: SETERRQ(PETSC_ERR_PLIB, "Invalid boundary edge marker information");
934: }
935: #ifdef MESH_TRACE_BOUNDARY_CREATE
936: /* Check boundary nodes */
937: for(bdNode = tri->bdBegin[bd]; bdNode < tri->bdBegin[bd+1]; bdNode++)
938: PetscPrintf(mesh->comm, "bd: %4d bdNode: %4d node: %5dn",
939: bd, bdNode, bdNodes[bdNode], bdEdges[bdNode], edges[bdEdges[bdNode]*2], edges[bdEdges[bdNode]*2+1]);
940: /* Check boundary edges */
941: for(bdEdge = tri->bdEdgeBegin[bd]; bdEdge < tri->bdEdgeBegin[bd+1]; bdEdge++)
942: PetscPrintf(mesh->comm, "bd: %4d bdEdge: %4d edge: %5d start: %5d end: %5dn",
943: bd, bdEdge, bdEdges[bdEdge], edges[bdEdges[bdEdge]*2], edges[bdEdges[bdEdge]*2+1]);
944: #endif
945: bdSeen[bd] = 1;
946: }
947: }
948: }
949: PetscFree(bdSeen);
951: /* Set fields */
952: tri->bdNodes = bdNodes;
953: tri->bdEdges = bdEdges;
954: /* Diagnostics */
955: #ifdef MESH_TRACE_BOUNDARY_REORDER
956: int minNode, minBdNode, minEdge, minBdEdge, dir, nextBdEdge;
958: PetscMalloc(mesh->numBdNodes * sizeof(int), &bdNodes);
959: PetscMalloc(mesh->numBdEdges * sizeof(int), &bdEdges);
960: for(bd = 0; bd < numBd; bd++) {
961: /* Find smallest node */
962: for(bdNode = tri->bdBegin[bd], minNode = tri->numNodes, minBdNode = -1; bdNode < tri->bdBegin[bd+1]; bdNode++)
963: if (tri->bdNodes[bdNode] < minNode) {
964: minNode = tri->bdNodes[bdNode];
965: minBdNode = bdNode;
966: }
967: /* Proceed in the direction of the smaller node */
968: if (minBdNode == tri->bdBegin[bd]) {
969: if (tri->bdNodes[minBdNode+1] < tri->bdNodes[tri->bdBegin[bd+1]-1])
970: dir = 1;
971: else
972: dir = 0;
973: } else if (minBdNode == tri->bdBegin[bd+1]-1) {
974: if (tri->bdNodes[tri->bdBegin[bd]] < tri->bdNodes[minBdNode-1])
975: dir = 0;
976: else
977: dir = 1;
978: } else if (tri->bdNodes[minBdNode+1] < tri->bdNodes[minBdNode-1]) {
979: dir = 1;
980: } else {
981: dir = 0;
982: }
983: if (dir)
984: {
985: for(nextBdNode = tri->bdBegin[bd], bdNode = minBdNode; bdNode < tri->bdBegin[bd+1]; nextBdNode++, bdNode++)
986: bdNodes[nextBdNode] = tri->bdNodes[bdNode];
987: for(bdNode = tri->bdBegin[bd]; bdNode < minBdNode; nextBdNode++, bdNode++)
988: bdNodes[nextBdNode] = tri->bdNodes[bdNode];
989: }
990: else
991: {
992: for(nextBdNode = tri->bdBegin[bd], bdNode = minBdNode; bdNode >= tri->bdBegin[bd]; nextBdNode++, bdNode--)
993: bdNodes[nextBdNode] = tri->bdNodes[bdNode];
994: for(bdNode = tri->bdBegin[bd+1]-1; bdNode > minBdNode; nextBdNode++, bdNode--)
995: bdNodes[nextBdNode] = tri->bdNodes[bdNode];
996: }
997: }
998: for(bd = 0; bd < numBd; bd++) {
999: /* Find smallest edge */
1000: for(bdEdge = tri->bdEdgeBegin[bd], minEdge = tri->numEdges, minBdEdge = -1; bdEdge < tri->bdEdgeBegin[bd+1]; bdEdge++)
1001: if (tri->bdEdges[bdEdge] < minEdge) {
1002: minEdge = tri->bdEdges[bdEdge];
1003: minBdEdge = bdEdge;
1004: }
1005: /* Proceed in the direction of the smaller edge */
1006: if (minBdEdge == tri->bdEdgeBegin[bd]) {
1007: if (tri->bdEdges[minBdEdge+1] < tri->bdEdges[tri->bdEdgeBegin[bd+1]-1])
1008: dir = 1;
1009: else
1010: dir = 0;
1011: } else if (minBdEdge == tri->bdEdgeBegin[bd+1]-1) {
1012: if (tri->bdEdges[tri->bdEdgeBegin[bd]] < tri->bdEdges[minBdEdge-1])
1013: dir = 0;
1014: else
1015: dir = 1;
1016: } else if (tri->bdEdges[minBdEdge+1] < tri->bdEdges[minBdEdge-1]) {
1017: dir = 1;
1018: } else {
1019: dir = 0;
1020: }
1021: if (dir)
1022: {
1023: for(nextBdEdge = tri->bdEdgeBegin[bd], bdEdge = minBdEdge; bdEdge < tri->bdEdgeBegin[bd+1]; nextBdEdge++, bdEdge++)
1024: bdEdges[nextBdEdge] = tri->bdEdges[bdEdge];
1025: for(bdEdge = tri->bdEdgeBegin[bd]; bdEdge < minBdEdge; nextBdEdge++, bdEdge++)
1026: bdEdges[nextBdEdge] = tri->bdEdges[bdEdge];
1027: }
1028: else
1029: {
1030: for(nextBdEdge = tri->bdEdgeBegin[bd], bdEdge = minBdEdge; bdEdge >= tri->bdEdgeBegin[bd]; nextBdEdge++, bdEdge--)
1031: bdEdges[nextBdEdge] = tri->bdEdges[bdEdge];
1032: for(bdEdge = tri->bdEdgeBegin[bd+1]-1; bdEdge > minBdEdge; nextBdEdge++, bdEdge--)
1033: bdEdges[nextBdEdge] = tri->bdEdges[bdEdge];
1034: }
1035: }
1036: PetscMemcpy(tri->bdNodes, bdNodes, mesh->numBdNodes * sizeof(int));
1037: PetscMemcpy(tri->bdEdges, bdEdges, mesh->numBdEdges * sizeof(int));
1038: PetscFree(bdNodes);
1039: PetscFree(bdEdges);
1040: #endif
1041: #ifdef MESH_TRACE_BOUNDARY
1042: int minNode, minBdNode, minEdge, minBdEdge, dir;
1044: for(bd = 0; bd < numBd; bd++) {
1045: /* Find smallest node */
1046: for(bdNode = tri->bdBegin[bd], minNode = tri->numNodes, minBdNode = -1; bdNode < tri->bdBegin[bd+1]; bdNode++)
1047: if (tri->bdNodes[bdNode] < minNode) {
1048: minNode = tri->bdNodes[bdNode];
1049: minBdNode = bdNode;
1050: }
1051: /* Proceed in the direction of the smaller node */
1052: if (minBdNode == tri->bdBegin[bd]) {
1053: if (tri->bdNodes[minBdNode+1] < tri->bdNodes[tri->bdBegin[bd+1]-1])
1054: dir = 1;
1055: else
1056: dir = 0;
1057: } else if (minBdNode == tri->bdBegin[bd+1]-1) {
1058: if (tri->bdNodes[tri->bdBegin[bd]] < tri->bdNodes[minBdNode-1])
1059: dir = 0;
1060: else
1061: dir = 1;
1062: } else if (tri->bdNodes[minBdNode+1] < tri->bdNodes[minBdNode-1]) {
1063: dir = 1;
1064: } else {
1065: dir = 0;
1066: }
1067: if (dir)
1068: {
1069: for(bdNode = minBdNode; bdNode < tri->bdBegin[bd+1]; bdNode++)
1070: PetscPrintf(mesh->comm, "bd: %4d node: %5dn", bd, tri->bdNodes[bdNode]);
1071: for(bdNode = tri->bdBegin[bd]; bdNode < minBdNode; bdNode++)
1072: PetscPrintf(mesh->comm, "bd: %4d node: %5dn", bd, tri->bdNodes[bdNode]);
1073: }
1074: else
1075: {
1076: for(bdNode = minBdNode; bdNode >= tri->bdBegin[bd]; bdNode--)
1077: PetscPrintf(mesh->comm, "bd: %4d node: %5dn", bd, tri->bdNodes[bdNode]);
1078: for(bdNode = tri->bdBegin[bd+1]-1; bdNode > minBdNode; bdNode--)
1079: PetscPrintf(mesh->comm, "bd: %4d node: %5dn", bd, tri->bdNodes[bdNode]);
1080: }
1081: }
1082: for(bd = 0; bd < numBd; bd++) {
1083: /* Find smallest edge */
1084: for(bdEdge = tri->bdEdgeBegin[bd], minEdge = tri->numEdges, minBdEdge = -1; bdEdge < tri->bdEdgeBegin[bd+1]; bdEdge++)
1085: if (tri->bdEdges[bdEdge] < minEdge) {
1086: minEdge = tri->bdEdges[bdEdge];
1087: minBdEdge = bdEdge;
1088: }
1089: /* Proceed in the direction of the smaller edge */
1090: if (minBdEdge == tri->bdEdgeBegin[bd]) {
1091: if (tri->bdEdges[minBdEdge+1] < tri->bdEdges[tri->bdEdgeBegin[bd+1]-1])
1092: dir = 1;
1093: else
1094: dir = 0;
1095: } else if (minBdEdge == tri->bdEdgeBegin[bd+1]-1) {
1096: if (tri->bdEdges[tri->bdEdgeBegin[bd]] < tri->bdEdges[minBdEdge-1])
1097: dir = 0;
1098: else
1099: dir = 1;
1100: } else if (tri->bdEdges[minBdEdge+1] < tri->bdEdges[minBdEdge-1]) {
1101: dir = 1;
1102: } else {
1103: dir = 0;
1104: }
1105: if (dir)
1106: {
1107: for(bdEdge = minBdEdge; bdEdge < tri->bdEdgeBegin[bd+1]; bdEdge++)
1108: PetscPrintf(mesh->comm, "bd: %4d edge: %5d start: %5d end: %5dn",
1109: bd, tri->bdEdges[bdEdge], edges[tri->bdEdges[bdEdge]*2], edges[tri->bdEdges[bdEdge]*2+1]);
1110: for(bdEdge = tri->bdEdgeBegin[bd]; bdEdge < minBdEdge; bdEdge++)
1111: PetscPrintf(mesh->comm, "bd: %4d edge: %5d start: %5d end: %5dn",
1112: bd, tri->bdEdges[bdEdge], edges[tri->bdEdges[bdEdge]*2], edges[tri->bdEdges[bdEdge]*2+1]);
1113: }
1114: else
1115: {
1116: for(bdEdge = minBdEdge; bdEdge >= tri->bdEdgeBegin[bd]; bdEdge--)
1117: PetscPrintf(mesh->comm, "bd: %4d edge: %5d start: %5d end: %5dn",
1118: bd, tri->bdEdges[bdEdge], edges[tri->bdEdges[bdEdge]*2], edges[tri->bdEdges[bdEdge]*2+1]);
1119: for(bdEdge = tri->bdEdgeBegin[bd+1]-1; bdEdge > minBdEdge; bdEdge--)
1120: PetscPrintf(mesh->comm, "bd: %4d edge: %5d start: %5d end: %5dn",
1121: bd, tri->bdEdges[bdEdge], edges[tri->bdEdges[bdEdge]*2], edges[tri->bdEdges[bdEdge]*2+1]);
1122: }
1123: }
1124: #endif
1126: PetscLogEventEnd(MESH_SetupBoundary, mesh, 0, 0, 0);
1127: return(0);
1128: }
1130: /*
1131: MeshAssemble_Private - This function assembles the mesh entirely on the first processor.
1133: Collective on Mesh
1135: Input Parameters:
1136: . mesh - The mesh being refined
1138: Output Parameter:
1139: + nodes - The node coordinates
1140: . markers - The node markers
1141: . faces - The nodes for each element
1142: - edges - The nodes for each edge
1144: Level: developer
1146: .keywords mesh, assemble
1147: .seealso MeshInitRefineInput_Triangle()
1148: */
1149: int MeshAssemble_Private(Mesh mesh, double **nodes, int **markers, int **faces, int **edges)
1150: {
1151: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1152: Partition p = mesh->part;
1153: Partition_Triangular_2D *q = (Partition_Triangular_2D *) p->data;
1154: int numNodes = q->numNodes;
1155: int numFaces = p->numElements;
1156: int numEdges = q->numEdges;
1157: int numCorners = mesh->numCorners;
1158: int numProcs = p->numProcs;
1159: int rank = p->rank;
1160: int *locFaces = PETSC_NULL;
1161: int *locEdges = PETSC_NULL;
1162: int *numRecvNodes, *numRecvMarkers, *numRecvFaces, *numRecvEdges;
1163: int *nodeOffsets, *markerOffsets, *faceOffsets, *edgeOffsets;
1164: int proc, elem, edge, size;
1165: int ierr;
1168: /* Allocate global arrays */
1169: if (rank == 0) {
1170: PetscMalloc(numNodes*2 * sizeof(double), nodes);
1171: PetscMalloc(numNodes * sizeof(int), markers);
1172: PetscMalloc(numFaces*numCorners * sizeof(int), faces);
1173: PetscMalloc(numEdges*2 * sizeof(int), edges);
1174: }
1176: if (numProcs > 1) {
1177: if (mesh->numFaces > 0) {
1178: PetscMalloc(mesh->numFaces*numCorners * sizeof(int), &locFaces);
1179: }
1180: if (mesh->numEdges > 0) {
1181: PetscMalloc(mesh->numEdges*2 * sizeof(int), &locEdges);
1182: }
1184: /* Calculate offsets */
1185: PetscMalloc(numProcs * sizeof(int), &numRecvNodes);
1186: PetscMalloc(numProcs * sizeof(int), &numRecvMarkers);
1187: PetscMalloc(numProcs * sizeof(int), &numRecvEdges);
1188: PetscMalloc(numProcs * sizeof(int), &numRecvFaces);
1189: PetscMalloc(numProcs * sizeof(int), &nodeOffsets);
1190: PetscMalloc(numProcs * sizeof(int), &markerOffsets);
1191: PetscMalloc(numProcs * sizeof(int), &edgeOffsets);
1192: PetscMalloc(numProcs * sizeof(int), &faceOffsets);
1193: for(proc = 0; proc < numProcs; proc++) {
1194: numRecvNodes[proc] = (q->firstNode[proc+1] - q->firstNode[proc])*2;
1195: numRecvMarkers[proc] = (q->firstNode[proc+1] - q->firstNode[proc]);
1196: numRecvEdges[proc] = (q->firstEdge[proc+1] - q->firstEdge[proc])*2;
1197: numRecvFaces[proc] = (p->firstElement[proc+1] - p->firstElement[proc])*numCorners;
1198: }
1199: nodeOffsets[0] = 0;
1200: markerOffsets[0] = 0;
1201: edgeOffsets[0] = 0;
1202: faceOffsets[0] = 0;
1203: for(proc = 1; proc < numProcs; proc++) {
1204: nodeOffsets[proc] = numRecvNodes[proc-1] + nodeOffsets[proc-1];
1205: markerOffsets[proc] = numRecvMarkers[proc-1] + markerOffsets[proc-1];
1206: edgeOffsets[proc] = numRecvEdges[proc-1] + edgeOffsets[proc-1];
1207: faceOffsets[proc] = numRecvFaces[proc-1] + faceOffsets[proc-1];
1208: }
1210: /* Local to global node number conversion */
1211: for(elem = 0; elem < mesh->numFaces*numCorners; elem++) {
1212: PartitionLocalToGlobalNodeIndex(p, tri->faces[elem], &locFaces[elem]);
1213: }
1214: for(edge = 0; edge < mesh->numEdges*2; edge++) {
1215: PartitionLocalToGlobalNodeIndex(p, tri->edges[edge], &locEdges[edge]);
1216: }
1218: /* Collect global arrays */
1219: size = mesh->numNodes*2;
1220: MPI_Gatherv(tri->nodes, size, MPI_DOUBLE, *nodes, numRecvNodes, nodeOffsets, MPI_DOUBLE, 0, p->comm);
1221:
1222: size = mesh->numNodes;
1223: MPI_Gatherv(tri->markers, size, MPI_INT, *markers, numRecvMarkers, markerOffsets, MPI_INT, 0, p->comm);
1224:
1225: size = mesh->numEdges*2;
1226: MPI_Gatherv(locEdges, size, MPI_INT, *edges, numRecvEdges, edgeOffsets, MPI_INT, 0, p->comm);
1227:
1228: size = mesh->numFaces*numCorners;
1229: MPI_Gatherv(locFaces, size, MPI_INT, *faces, numRecvFaces, faceOffsets, MPI_INT, 0, p->comm);
1230:
1232: /* Cleanup */
1233: if (locFaces != PETSC_NULL) {
1234: PetscFree(locFaces);
1235: }
1236: if (locEdges != PETSC_NULL) {
1237: PetscFree(locEdges);
1238: }
1239: PetscFree(numRecvNodes);
1240: PetscFree(numRecvMarkers);
1241: PetscFree(numRecvEdges);
1242: PetscFree(numRecvFaces);
1243: PetscFree(nodeOffsets);
1244: PetscFree(markerOffsets);
1245: PetscFree(edgeOffsets);
1246: PetscFree(faceOffsets);
1247: } else {
1248: /* Uniprocessor case */
1249: PetscMemcpy(*nodes, tri->nodes, numNodes*2 * sizeof(double));
1250: PetscMemcpy(*markers, tri->markers, numNodes * sizeof(int));
1251: PetscMemcpy(*faces, tri->faces, numFaces*numCorners * sizeof(int));
1252: PetscMemcpy(*edges, tri->edges, numEdges*2 * sizeof(int));
1253: }
1254: return(0);
1255: }
1257: /*
1258: MeshCoarsen_Triangular_2D - Coarsens a two dimensional unstructured mesh using area constraints
1260: Collective on Mesh
1262: Input Parameters:
1263: + mesh - The mesh begin coarsened
1264: - area - A function which gives an area constraint when evaluated inside an element
1266: Output Parameters:
1267: . newmesh - The coarse mesh
1269: Note:
1270: If PETSC_NULL is given as the 'area' argument, the mesh consisting only of vertices is returned.
1272: Level: developer
1274: .keywords unstructured grid
1275: .seealso GridCoarsen(), MeshRefine(), MeshReform()
1276: */
1277: static int MeshCoarsen_Triangular_2D(Mesh mesh, PointFunction area, Mesh *newmesh)
1278: {
1279: ParameterDict dict;
1280: Mesh m;
1281: Mesh_Triangular *tri;
1282: Partition p;
1283: Partition_Triangular_2D *q;
1284: Mesh_Triangular *triOld = (Mesh_Triangular *) mesh->data;
1285: Partition pOld = mesh->part;
1286: Partition_Triangular_2D *qOld = (Partition_Triangular_2D *) pOld->data;
1287: int (*refine)(Mesh, PointFunction, Mesh);
1288: int *newNodes; /* newNodes[fine node #] = coarse node # or -1 */
1289: int *FineMapping; /* FineMapping[n] = node n in the fine mesh numbering */
1290: int *CoarseMapping; /* CoarseMapping[n] = node n in the coarse mesh numbering */
1291: PetscTruth hasIndex;
1292: int *nodePerm, *temp;
1293: int numGhostElements, numGhostNodes;
1294: int proc, bd, elem, edge, corner, node, newNode, gNode, ghostNode, bdNode, count;
1295: int ierr;
1298: if (area == PETSC_NULL) {
1299: if (mesh->numCorners == 3) {
1300: /* Return the mesh and increment the reference count instead of copying */
1301: PetscObjectReference((PetscObject) mesh);
1302: *newmesh = mesh;
1303: return(0);
1304: } else if (mesh->numCorners != 6) {
1305: SETERRQ1(PETSC_ERR_SUP, "Number of local nodes %d not supported", mesh->numCorners);
1306: }
1308: /* We would like to preserve the partition of vertices, although it may be somewhat
1309: unbalanced, since this makes mapping from the coarse mesh back to the original
1310: much easier.
1311: */
1312: MeshCreate(mesh->comm, &m);
1313: ParameterDictCreate(mesh->comm, &dict);
1314: ParameterDictSetObject(dict, "bdCtx", mesh->bdCtx);
1315: ParameterDictSetInteger(dict, "dim", 2);
1316: ParameterDictSetInteger(dict, "numLocNodes", mesh->numCorners);
1317: PetscObjectSetParameterDict((PetscObject) m, dict);
1318: ParameterDictDestroy(dict);
1319: PetscNew(Mesh_Triangular, &tri);
1320: PetscLogObjectMemory(m, sizeof(Mesh_Triangular));
1321: PetscMemcpy(m->ops, mesh->ops, sizeof(struct _MeshOps));
1322: PetscStrallocpy(mesh->type_name, &m->type_name);
1323: PetscStrallocpy(mesh->serialize_name, &m->serialize_name);
1324: m->data = (void *) tri;
1325: m->dim = 2;
1326: m->startX = mesh->startX;
1327: m->endX = mesh->endX;
1328: m->sizeX = mesh->sizeX;
1329: m->startY = mesh->startY;
1330: m->endY = mesh->endY;
1331: m->sizeY = mesh->sizeY;
1332: m->isPeriodic = mesh->isPeriodic;
1333: m->isPeriodicDim[0] = mesh->isPeriodicDim[0];
1334: m->isPeriodicDim[1] = mesh->isPeriodicDim[1];
1335: m->isPeriodicDim[2] = mesh->isPeriodicDim[2];
1336: m->nodeOrdering = PETSC_NULL;
1337: m->partitioned = 1;
1338: m->highlightElement = mesh->highlightElement;
1339: m->usr = mesh->usr;
1341: /* Copy function list */
1342: PetscObjectQueryFunction((PetscObject) mesh, "MeshTriangular2D_Refine_C", (void (**)(void)) &refine);
1343: if (refine != PETSC_NULL) {
1344: #ifdef PETSC_HAVE_TRIANGLE
1345: PetscObjectComposeFunction((PetscObject) m, "MeshTriangular2D_Refine_C", "MeshRefine_Triangle",
1346: (void (*)(void)) MeshRefine_Triangle);
1347: #else
1348: /* The query needs to return the name also */
1349: PetscObjectComposeFunction((PetscObject) m, "MeshTriangular2D_Refine_C", "", (void (*)(void)) refine);
1350: #endif
1351:
1352: }
1354: /* Create the partition object */
1355: PartitionCreate(m, &p);
1356: PetscNew(Partition_Triangular_2D, &q);
1357: PetscLogObjectMemory(p, sizeof(Partition_Triangular_2D));
1358: PetscMemcpy(p->ops, pOld->ops, sizeof(struct _PartitionOps));
1359: PetscObjectChangeTypeName((PetscObject) p, pOld->type_name);
1360: PetscObjectChangeSerializeName((PetscObject) p, pOld->serialize_name);
1361: p->data = (void *) q;
1362: m->part = p;
1363: PetscLogObjectParent(m, p);
1365: /* Construct element partition */
1366: p->numProcs = pOld->numProcs;
1367: p->rank = pOld->rank;
1368: p->ordering = PETSC_NULL;
1369: p->numLocElements = pOld->numLocElements;
1370: p->numElements = pOld->numElements;
1371: p->numOverlapElements = pOld->numOverlapElements;
1372: PetscMalloc((p->numProcs+1) * sizeof(int), &p->firstElement);
1373: PetscMemcpy(p->firstElement, pOld->firstElement, (p->numProcs+1) * sizeof(int));
1374: p->ghostElements = PETSC_NULL;
1375: p->ghostElementProcs = PETSC_NULL;
1376: numGhostElements = pOld->numOverlapElements - pOld->numLocElements;
1377: if (numGhostElements > 0) {
1378: PetscMalloc(numGhostElements * sizeof(int), &p->ghostElements);
1379: PetscMalloc(numGhostElements * sizeof(int), &p->ghostElementProcs);
1380: PetscLogObjectMemory(p, numGhostElements*2 * sizeof(int));
1381: PetscMemcpy(p->ghostElements, pOld->ghostElements, numGhostElements * sizeof(int));
1382: PetscMemcpy(p->ghostElementProcs, pOld->ghostElementProcs, numGhostElements * sizeof(int));
1383: }
1384: q->nodeOrdering = PETSC_NULL;
1385: q->edgeOrdering = PETSC_NULL;
1386: q->numLocEdges = qOld->numLocEdges;
1387: q->numEdges = qOld->numEdges;
1388: PetscMalloc((p->numProcs+1) * sizeof(int), &q->firstEdge);
1389: PetscMemcpy(q->firstEdge, qOld->firstEdge, (p->numProcs+1) * sizeof(int));
1390: PetscLogObjectMemory(p, (p->numProcs+1)*2 * sizeof(int));
1392: /* We must handle ghost members since we manually construct the partition */
1393: m->numBd = mesh->numBd;
1394: m->numFaces = mesh->numFaces;
1395: m->numCorners = 3;
1396: m->numEdges = mesh->numEdges;
1397: PetscMalloc(p->numOverlapElements*m->numCorners * sizeof(int), &tri->faces);
1398: PetscMalloc(p->numOverlapElements*3 * sizeof(int), &tri->neighbors);
1399: PetscMalloc(m->numEdges*2 * sizeof(int), &tri->edges);
1400: PetscMalloc(m->numEdges * sizeof(int), &tri->edgemarkers);
1401: PetscMalloc(qOld->numOverlapNodes * sizeof(int), &newNodes);
1402: PetscLogObjectMemory(m, (p->numOverlapElements*(m->numCorners + 3) + m->numEdges*3) * sizeof(int));
1403: for(node = 0; node < qOld->numOverlapNodes; node++)
1404: newNodes[node] = -1;
1405: /* Renumber the vertices and construct the face list */
1406: for(elem = 0, newNode = 0, numGhostNodes = 0; elem < p->numOverlapElements; elem++) {
1407: for(corner = 0; corner < 3; corner++) {
1408: node = triOld->faces[elem*mesh->numCorners+corner];
1409: if (newNodes[node] == -1) {
1410: if (node >= mesh->numNodes) {
1411: /* Mark them as ghost nodes */
1412: numGhostNodes++;
1413: newNodes[node] = -2;
1414: } else {
1415: newNodes[node] = newNode++;
1416: }
1417: }
1418: tri->faces[elem*m->numCorners+corner] = node;
1419: }
1420: }
1421: PetscMemcpy(tri->neighbors, triOld->neighbors, p->numOverlapElements*3 * sizeof(int));
1422: PetscMemcpy(tri->edges, triOld->edges, m->numEdges*2 * sizeof(int));
1423: PetscMemcpy(tri->edgemarkers, triOld->edgemarkers, m->numEdges * sizeof(int));
1425: /* We may have extra ghost nodes from the edges */
1426: for(edge = 0; edge < q->numLocEdges; edge++) {
1427: for(corner = 0; corner < 2; corner++) {
1428: node = tri->edges[edge*2+corner];
1429: if (newNodes[node] == -1) {
1430: if (node >= mesh->numNodes) {
1431: /* Mark them as ghost nodes */
1432: numGhostNodes++;
1433: newNodes[node] = -2;
1434: }
1435: }
1436: }
1437: }
1439: /* Compute vertex sizes */
1440: m->numNodes = newNode;
1441: m->numVertices = m->numNodes;
1442: MPI_Allreduce(&m->numNodes, &q->numNodes, 1, MPI_INT, MPI_SUM, m->comm);
1443: q->numLocNodes = newNode;
1444: q->numOverlapNodes = newNode + numGhostNodes;
1446: /* Compute vertex partition */
1447: PetscMalloc((p->numProcs+1) * sizeof(int), &q->firstNode);
1448: PetscLogObjectMemory(p, (p->numProcs+1) * sizeof(int));
1449: MPI_Allgather(&q->numLocNodes, 1, MPI_INT, &q->firstNode[1], 1, MPI_INT, p->comm);
1450: for(proc = 1, q->firstNode[0] = 0; proc <= p->numProcs; proc++) {
1451: q->firstNode[proc] += q->firstNode[proc-1];
1452: }
1453: q->ghostNodes = PETSC_NULL;
1454: q->ghostNodeProcs = PETSC_NULL;
1455: if (numGhostNodes > 0) {
1456: PetscMalloc(numGhostNodes * sizeof(int), &q->ghostNodes);
1457: PetscMalloc(numGhostNodes * sizeof(int), &q->ghostNodeProcs);
1458: PetscLogObjectMemory(p, numGhostNodes*2 * sizeof(int));
1459: }
1461: /* Create the mapping from coarse nodes to vertices in the original mesh */
1462: PetscMalloc(m->numNodes * sizeof(int), &FineMapping);
1463: PetscMalloc(m->numNodes * sizeof(int), &CoarseMapping);
1464: for(node = 0, count = 0; node < mesh->numNodes; node++) {
1465: newNode = newNodes[node];
1466: if (newNode >= 0) {
1467: /* These are all interior nodes */
1468: PartitionLocalToGlobalNodeIndex(pOld, node, &FineMapping[count]);
1469: PartitionLocalToGlobalNodeIndex(p, newNode, &CoarseMapping[count]);
1470: count++;
1471: }
1472: }
1473: if (count != m->numNodes) SETERRQ2(PETSC_ERR_PLIB, "Invalid coarse map size %d should be %d", count, m->numNodes);
1474: AOCreateMapping(m->comm, m->numNodes, FineMapping, CoarseMapping, &m->coarseMap);
1475: PetscFree(FineMapping);
1476: PetscFree(CoarseMapping);
1478: /* Setup ghost nodes */
1479: for(ghostNode = mesh->numNodes, count = 0; ghostNode < qOld->numOverlapNodes; ghostNode++) {
1480: if (newNodes[ghostNode] == -2) {
1481: PartitionLocalToGlobalNodeIndex(pOld, ghostNode, &gNode);
1482: AOApplicationToPetsc(m->coarseMap, 1, &gNode);
1483: q->ghostNodes[count] = gNode;
1484: q->ghostNodeProcs[count] = qOld->ghostNodeProcs[ghostNode-mesh->numNodes];
1485: count++;
1486: }
1487: }
1488: if (count != numGhostNodes) SETERRQ2(PETSC_ERR_PLIB, "Invalid number of ghost nodes %d should be %d", count, numGhostNodes);
1489: /* Cleanup */
1490: PetscFree(newNodes);
1492: /* Resort ghost nodes */
1493: if (numGhostNodes > 0) {
1494: PetscMalloc(numGhostNodes*2 * sizeof(int), &nodePerm);
1495: temp = nodePerm + numGhostNodes;
1496: for(node = 0; node < numGhostNodes; node++)
1497: nodePerm[node] = node;
1498: PetscSortIntWithPermutation(numGhostNodes, q->ghostNodes, nodePerm);
1499: for(node = 0; node < numGhostNodes; node++)
1500: temp[node] = q->ghostNodes[nodePerm[node]];
1501: for(node = 0; node < numGhostNodes; node++)
1502: q->ghostNodes[node] = temp[node];
1503: for(node = 0; node < numGhostNodes; node++)
1504: temp[node] = q->ghostNodeProcs[nodePerm[node]];
1505: for(node = 0; node < numGhostNodes; node++)
1506: q->ghostNodeProcs[node] = temp[node];
1507: PetscFree(nodePerm);
1508: }
1510: /* Copy vertex information */
1511: PetscMalloc(q->numOverlapNodes*2 * sizeof(double), &tri->nodes);
1512: PetscMalloc(q->numOverlapNodes * sizeof(int), &tri->markers);
1513: PetscLogObjectMemory(m, q->numOverlapNodes*2 * sizeof(double));
1514: PetscLogObjectMemory(m, q->numOverlapNodes * sizeof(int));
1515: for(newNode = 0; newNode < q->numOverlapNodes; newNode++) {
1516: PartitionLocalToGlobalNodeIndex(p, newNode, &gNode);
1517: AOPetscToApplication(m->coarseMap, 1, &gNode);
1518: PartitionGlobalToLocalNodeIndex(pOld, gNode, &node);
1519: tri->nodes[newNode*2] = triOld->nodes[node*2];
1520: tri->nodes[newNode*2+1] = triOld->nodes[node*2+1];
1521: tri->markers[newNode] = triOld->markers[node];
1522: }
1524: /* Renumber nodes */
1525: for(elem = 0; elem < p->numOverlapElements*m->numCorners; elem++) {
1526: PartitionLocalToGlobalNodeIndex(pOld, tri->faces[elem], &tri->faces[elem]);
1527: }
1528: for(edge = 0; edge < m->numEdges*2; edge++) {
1529: PartitionLocalToGlobalNodeIndex(pOld, tri->edges[edge], &tri->edges[edge]);
1530: }
1531: AOApplicationToPetsc(m->coarseMap, p->numOverlapElements*m->numCorners, tri->faces);
1532: AOApplicationToPetsc(m->coarseMap, m->numEdges*2, tri->edges);
1533: for(elem = 0; elem < p->numOverlapElements*m->numCorners; elem++) {
1534: PartitionGlobalToLocalNodeIndex(p, tri->faces[elem], &tri->faces[elem]);
1535: }
1536: for(edge = 0; edge < m->numEdges*2; edge++) {
1537: PartitionGlobalToLocalNodeIndex(p, tri->edges[edge], &tri->edges[edge]);
1538: }
1540: /* Construct derived and boundary information */
1541: tri->areas = PETSC_NULL;
1542: tri->aspectRatios = PETSC_NULL;
1543: m->numBdEdges = mesh->numBdEdges;
1544: PetscMalloc(m->numBd * sizeof(int), &tri->bdMarkers);
1545: PetscMalloc((m->numBd+1) * sizeof(int), &tri->bdBegin);
1546: PetscMalloc((m->numBd+1) * sizeof(int), &tri->bdEdgeBegin);
1547: PetscMalloc(m->numBdEdges * sizeof(int), &tri->bdEdges);
1548: PetscLogObjectMemory(m, (m->numBd + (m->numBd+1)*2 + m->numBdEdges) * sizeof(int));
1549: PetscMemcpy(tri->bdMarkers, triOld->bdMarkers, m->numBd * sizeof(int));
1550: PetscMemcpy(tri->bdEdgeBegin, triOld->bdEdgeBegin, (m->numBd+1) * sizeof(int));
1551: PetscMemcpy(tri->bdEdges, triOld->bdEdges, m->numBdEdges * sizeof(int));
1552: tri->bdBegin[0] = 0;
1553: for(bd = 0; bd < m->numBd; bd++) {
1554: tri->bdBegin[bd+1] = tri->bdBegin[bd];
1555: for(bdNode = triOld->bdBegin[bd]; bdNode < triOld->bdBegin[bd+1]; bdNode++) {
1556: node = triOld->bdNodes[bdNode];
1557: AOMappingHasApplicationIndex(m->coarseMap, node, &hasIndex);
1558: if (hasIndex == PETSC_TRUE)
1559: tri->bdBegin[bd+1]++;
1560: }
1561: }
1562: m->numBdNodes = tri->bdBegin[m->numBd];
1563: /* Assume closed boundaries */
1564: if (m->numBdNodes != m->numBdEdges) SETERRQ(PETSC_ERR_PLIB, "Invalid mesh boundary");
1565: #ifdef PETSC_USE_BOPT_g
1566: MPI_Scan(&m->numBdNodes, &node, 1, MPI_INT, MPI_SUM, m->comm);
1567: if (node != m->numBdNodes*(p->rank+1)) SETERRQ(PETSC_ERR_PLIB, "Invalid mesh boundary");
1568: #endif
1569: PetscMalloc(m->numBdNodes * sizeof(int), &tri->bdNodes);
1570: PetscLogObjectMemory(m, m->numBdNodes * sizeof(int));
1571: for(bd = 0, count = 0; bd < m->numBd; bd++) {
1572: for(bdNode = triOld->bdBegin[bd]; bdNode < triOld->bdBegin[bd+1]; bdNode++) {
1573: node = triOld->bdNodes[bdNode];
1574: AOMappingHasApplicationIndex(m->coarseMap, node, &hasIndex);
1575: if (hasIndex == PETSC_TRUE) {
1576: AOApplicationToPetsc(m->coarseMap, 1, &node);
1577: tri->bdNodes[count++] = node;
1578: }
1579: }
1580: }
1581: if (count != m->numBdNodes) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary node partition");
1583: /* Partition boundary nodes */
1584: PetscMalloc((p->numProcs+1) * sizeof(int), &q->firstBdNode);
1585: PetscLogObjectMemory(p, (p->numProcs+1) * sizeof(int));
1586: q->numLocBdNodes = 0;
1587: for(bdNode = 0; bdNode < m->numBdNodes; bdNode++) {
1588: newNode = tri->bdNodes[bdNode];
1589: if ((newNode >= q->firstNode[p->rank]) && (newNode < q->firstNode[p->rank+1]))
1590: q->numLocBdNodes++;
1591: }
1592: MPI_Allgather(&q->numLocBdNodes, 1, MPI_INT, &q->firstBdNode[1], 1, MPI_INT, p->comm);
1593: for(proc = 1, q->firstBdNode[0] = 0; proc <= p->numProcs; proc++) {
1594: q->firstBdNode[proc] += q->firstBdNode[proc-1];
1595: }
1596: q->numBdNodes = q->firstBdNode[p->numProcs];
1598: /* Process ghost boundary nodes */
1599: q->numOverlapBdNodes = q->numLocBdNodes;
1600: for(node = 0; node < numGhostNodes; node++) {
1601: if (tri->markers[m->numNodes+node] != 0)
1602: q->numOverlapBdNodes++;
1603: }
1604: q->ghostBdNodes = PETSC_NULL;
1605: if (q->numOverlapBdNodes > q->numLocBdNodes) {
1606: PetscMalloc((q->numOverlapBdNodes - q->numLocBdNodes) * sizeof(int), &q->ghostBdNodes);
1607: for(node = 0, bdNode = 0; node < numGhostNodes; node++) {
1608: if (tri->markers[m->numNodes+node] != 0)
1609: q->ghostBdNodes[bdNode++] = node;
1610: }
1611: if (bdNode != q->numOverlapBdNodes - q->numLocBdNodes) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary node partition");
1612: }
1614: /* Copy holes from previous mesh */
1615: m->numHoles = mesh->numHoles;
1616: m->holes = PETSC_NULL;
1617: if (m->numHoles > 0) {
1618: PetscMalloc(m->numHoles*2 * sizeof(double), &m->holes);
1619: PetscMemcpy(m->holes, mesh->holes, m->numHoles*2 * sizeof(double));
1620: PetscLogObjectMemory(m, m->numHoles*2 * sizeof(double));
1621: }
1623: /* Calculate maximum degree of vertices */
1624: m->maxDegree = mesh->maxDegree;
1625: PetscMalloc(m->numNodes * sizeof(int), &tri->degrees);
1626: PetscMalloc(m->maxDegree * sizeof(int), &m->support);
1627: for(newNode = 0; newNode < m->numNodes; newNode++) {
1628: PartitionLocalToGlobalNodeIndex(p, newNode, &gNode);
1629: AOPetscToApplication(m->coarseMap, 1, &gNode);
1630: PartitionGlobalToLocalNodeIndex(pOld, gNode, &node);
1631: tri->degrees[newNode] = triOld->degrees[node];
1632: }
1633: m->supportTaken = PETSC_FALSE;
1635: #ifdef PETSC_USE_BOPT_g
1636: /* Check mesh integrity */
1637: MeshDebug_Triangular_2D(m, PETSC_TRUE);
1638: #endif
1640: /* Initialize save space */
1641: tri->nodesOld = PETSC_NULL;
1643: } else {
1644: SETERRQ(PETSC_ERR_SUP, "No coarsening strategies currently supported");
1645: }
1647: PetscObjectSetName((PetscObject) m, "Coarse Mesh");
1648: *newmesh = m;
1649: return(0);
1650: }
1652: /*
1653: MeshRefine_Triangular_2D - Refines a two dimensional unstructured mesh using area constraints
1655: Collective on Mesh
1657: Input Parameters:
1658: + mesh - The mesh begin refined
1659: - area - A function which gives an area constraint when evaluated inside an element
1661: Output Parameter:
1662: . newmesh - The refined mesh
1664: Level: developer
1666: .keywords unstructured grid
1667: .seealso GridRefine(), MeshCoarsen(), MeshReform()
1668: */
1669: static int MeshRefine_Triangular_2D(Mesh mesh, PointFunction area, Mesh *newmesh)
1670: {
1671: ParameterDict dict;
1672: Mesh m;
1673: int (*f)(Mesh, PointFunction, Mesh);
1674: int ierr;
1677: MeshCreate(mesh->comm, &m);
1678: ParameterDictCreate(mesh->comm, &dict);
1679: ParameterDictSetObject(dict, "bdCtx", mesh->bdCtx);
1680: ParameterDictSetInteger(dict, "dim", 2);
1681: ParameterDictSetInteger(dict, "numLocNodes", mesh->numCorners);
1682: PetscObjectSetParameterDict((PetscObject) m, dict);
1683: ParameterDictDestroy(dict);
1685: PetscObjectQueryFunction((PetscObject) mesh, "MeshTriangular2D_Refine_C", (void (**)(void)) &f);
1686: if (f == PETSC_NULL) SETERRQ(PETSC_ERR_SUP, "Mesh has no refinement function");
1687: (*f)(mesh, area, m);
1688: PetscObjectSetName((PetscObject) m, "Refined Mesh");
1689: *newmesh = m;
1690: return(0);
1691: }
1693: static int MeshResetNodes_Triangular_2D(Mesh mesh, PetscTruth resetBd)
1694: {
1695: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1696: int *elements = tri->faces;
1697: int *markers = tri->markers;
1698: double *nodes = tri->nodes;
1699: int numCorners = mesh->numCorners;
1700: int elem, corner, node1, node2, midNode, midCorner;
1703: if (numCorners == 3)
1704: return(0);
1706: if (mesh->isPeriodic == PETSC_TRUE) {
1707: for(elem = 0; elem < mesh->part->numOverlapElements; elem++) {
1708: for(corner = 0; corner < 3; corner++) {
1709: node1 = elements[elem*numCorners+corner];
1710: node2 = elements[elem*numCorners+((corner+1)%3)];
1711: if ((markers[node1] != markers[node2]) || (markers[node1] == 0) || (resetBd == PETSC_TRUE)) {
1712: midCorner = (corner+2)%3 + 3;
1713: midNode = elements[elem*numCorners+midCorner];
1714: nodes[midNode*2] = MeshPeriodicX(mesh, 0.5*(nodes[node1*2] +
1715: MeshPeriodicRelativeX(mesh, nodes[node2*2], nodes[node1*2])));
1716: nodes[midNode*2+1] = MeshPeriodicY(mesh, 0.5*(nodes[node1*2+1] +
1717: MeshPeriodicRelativeY(mesh, nodes[node2*2+1], nodes[node1*2+1])));
1718: }
1719: }
1720: }
1721: } else {
1722: for(elem = 0; elem < mesh->part->numOverlapElements; elem++) {
1723: for(corner = 0; corner < 3; corner++) {
1724: node1 = elements[elem*numCorners+corner];
1725: node2 = elements[elem*numCorners+((corner+1)%3)];
1726: if ((markers[node1] != markers[node2]) || (markers[node1] == 0) || (resetBd == PETSC_TRUE)) {
1727: midCorner = (corner+2)%3 + 3;
1728: midNode = elements[elem*numCorners+midCorner];
1729: nodes[midNode*2] = 0.5*(nodes[node1*2] + nodes[node2*2]);
1730: nodes[midNode*2+1] = 0.5*(nodes[node1*2+1] + nodes[node2*2+1]);
1731: }
1732: }
1733: }
1734: }
1735: return(0);
1736: }
1738: int MeshPartition_Triangular_2D(Mesh mesh)
1739: {
1743: PartitionCreateTriangular2D(mesh, &mesh->part);
1744: mesh->partitioned = 1;
1745: PetscFunctionReturn(ierr);
1746: }
1748: int MeshUpdateBoundingBox_Triangular_2D(Mesh mesh)
1749: {
1750: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1751: double *nodes = tri->nodes;
1752: Partition part;
1753: int dim, numOverlapNodes;
1754: int node;
1755: int ierr;
1757: /* Calculate the local bounding box */
1759: MeshGetPartition(mesh, &part);
1760: MeshGetDimension(mesh, &dim);
1761: PartitionGetNumOverlapNodes(part, &numOverlapNodes);
1762: if (numOverlapNodes > 0) {
1763: mesh->locStartX = mesh->locEndX = nodes[0];
1764: mesh->locStartY = mesh->locEndY = nodes[1];
1765: } else {
1766: mesh->locStartX = mesh->locEndX = 0.0;
1767: mesh->locStartY = mesh->locEndY = 0.0;
1768: }
1769: for(node = 0; node < numOverlapNodes; node++) {
1770: if (nodes[node*dim] < mesh->locStartX) mesh->locStartX = nodes[node*dim];
1771: if (nodes[node*dim] > mesh->locEndX) mesh->locEndX = nodes[node*dim];
1772: if (nodes[node*dim+1] < mesh->locStartY) mesh->locStartY = nodes[node*dim+1];
1773: if (nodes[node*dim+1] > mesh->locEndY) mesh->locEndY = nodes[node*dim+1];
1774: }
1775: mesh->locSizeX = mesh->locEndX - mesh->locStartX;
1776: mesh->locSizeY = mesh->locEndY - mesh->locStartY;
1777: return(0);
1778: }
1780: int MeshGetElemNeighbor_Triangular_2D(Mesh mesh, int elem, int corner, int *neighbor)
1781: {
1782: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1783: #ifdef PETSC_USE_BOPT_g
1784: int numOverlapElements;
1785: int ierr;
1786: #endif
1789: #ifdef PETSC_USE_BOPT_g
1790: PartitionGetNumOverlapElements(mesh->part, &numOverlapElements);
1791: if ((elem < 0) || (elem >= numOverlapElements)) {
1792: SETERRQ2(PETSC_ERR_ARG_WRONG, "Invalid element %d should be in [0,%d)", elem, numOverlapElements);
1793: }
1794: if ((corner < 0) || (corner >= mesh->numCorners)) {
1795: SETERRQ2(PETSC_ERR_ARG_WRONG, "Invalid corner %d should be in [0,%d)", corner, mesh->numCorners);
1796: }
1797: #endif
1798: *neighbor = tri->neighbors[elem*3+corner];
1799: return(0);
1800: }
1802: int MeshGetNodeCoords_Triangular_2D(Mesh mesh, int node, double *x, double *y, double *z)
1803: {
1804: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1805: #ifdef PETSC_USE_BOPT_g
1806: int numOverlapNodes;
1807: int ierr;
1808: #endif
1811: #ifdef PETSC_USE_BOPT_g
1812: ierr = PartitionGetNumOverlapNodes(mesh->part, &numOverlapNodes);
1813: if ((node < 0) || (node >= numOverlapNodes)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Invalid node specified");
1814: #endif
1815: if (x != PETSC_NULL) *x = tri->nodes[node*2];
1816: if (y != PETSC_NULL) *y = tri->nodes[node*2+1];
1817: if (z != PETSC_NULL) *z = 0.0;
1818: return(0);
1819: }
1821: int MeshSetNodeCoords_Triangular_2D(Mesh mesh, int node, double x, double y, double z)
1822: {
1823: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1824: #ifdef PETSC_USE_BOPT_g
1825: int numOverlapNodes;
1826: int ierr;
1827: #endif
1830: #ifdef PETSC_USE_BOPT_g
1831: ierr = PartitionGetNumOverlapNodes(mesh->part, &numOverlapNodes);
1832: if ((node < 0) || (node >= numOverlapNodes)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Invalid node specified");
1833: #endif
1834: tri->nodes[node*2] = MeshPeriodicX(mesh, x);
1835: tri->nodes[node*2+1] = MeshPeriodicY(mesh, y);
1836: return(0);
1837: }
1839: int MeshGetNodeCoordsSaved_Triangular_2D(Mesh mesh, int node, double *x, double *y, double *z)
1840: {
1841: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1842: #ifdef PETSC_USE_BOPT_g
1843: int numOverlapNodes;
1844: int ierr;
1845: #endif
1848: #ifdef PETSC_USE_BOPT_g
1849: ierr = PartitionGetNumOverlapNodes(mesh->part, &numOverlapNodes);
1850: if ((node < 0) || (node >= numOverlapNodes)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Invalid node specified");
1851: #endif
1852: if (tri->nodesOld != PETSC_NULL) {
1853: if (x != PETSC_NULL) *x = tri->nodesOld[node*2];
1854: if (y != PETSC_NULL) *y = tri->nodesOld[node*2+1];
1855: if (z != PETSC_NULL) *z = 0.0;
1856: } else {
1857: if (x != PETSC_NULL) *x = tri->nodes[node*2];
1858: if (y != PETSC_NULL) *y = tri->nodes[node*2+1];
1859: if (z != PETSC_NULL) *z = 0.0;
1860: }
1861: return(0);
1862: }
1864: int MeshGetNodeBoundary_Triangular_2D(Mesh mesh, int node, int *bd)
1865: {
1866: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1867: #ifdef PETSC_USE_BOPT_g
1868: int numOverlapNodes;
1869: int ierr;
1870: #endif
1873: #ifdef PETSC_USE_BOPT_g
1874: ierr = PartitionGetNumOverlapNodes(mesh->part, &numOverlapNodes);
1875: if ((node < 0) || (node >= numOverlapNodes)) {
1876: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid node %d should be in [0,%d)", node, numOverlapNodes);
1877: }
1878: #endif
1879: *bd = tri->markers[node];
1880: return(0);
1881: }
1883: int MeshGetNodeFromElement_Triangular_2D(Mesh mesh, int elem, int corner, int *node)
1884: {
1885: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1886: int numCorners = mesh->numCorners;
1889: *node = tri->faces[elem*numCorners+corner];
1890: return(0);
1891: }
1893: int MeshGetNodeFromEdge_Triangular_2D(Mesh mesh, int edge, int corner, int *node)
1894: {
1895: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1898: *node = tri->edges[edge*2+corner];
1899: return(0);
1900: }
1902: int MeshNodeIsVertex_Triangular_2D(Mesh mesh, int node, PetscTruth *isVertex)
1903: {
1904: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1905: #ifdef PETSC_USE_BOPT_g
1906: int numOverlapNodes;
1907: int ierr;
1908: #endif
1911: #ifdef PETSC_USE_BOPT_g
1912: ierr = PartitionGetNumOverlapNodes(mesh->part, &numOverlapNodes);
1913: if ((node < 0) || (node >= numOverlapNodes)) {
1914: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid node %d should be in [0,%d)", node, numOverlapNodes);
1915: }
1916: #endif
1917: if (tri->degrees[node] == 0)
1918: *isVertex = PETSC_FALSE;
1919: else
1920: *isVertex = PETSC_TRUE;
1921: return(0);
1922: }
1924: int MeshGetBoundarySize_Triangular_2D(Mesh mesh, int boundary, int *size)
1925: {
1926: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1927: int b; /* Canonical boundary number */
1930: /* Find canonical boundary number */
1931: for(b = 0; b < mesh->numBd; b++)
1932: if (tri->bdMarkers[b] == boundary) break;
1933: if (b == mesh->numBd) SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid boundary specified");
1934: *size = tri->bdBegin[b+1] - tri->bdBegin[b];
1935: return(0);
1936: }
1938: int MeshGetBoundaryIndex_Triangular_2D(Mesh mesh, int boundary, int *index)
1939: {
1940: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1941: int b; /* Canonical boundary number */
1944: /* Find canonical boundary number */
1945: for(b = 0; b < mesh->numBd; b++)
1946: if (tri->bdMarkers[b] == boundary) break;
1947: if (b == mesh->numBd) SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid boundary marker %d", boundary);
1948: *index = b;
1949: return(0);
1950: }
1952: int MeshGetBoundaryStart_Triangular_2D(Mesh mesh, int boundary, PetscTruth ghost, int *node)
1953: {
1954: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1955: int b; /* Canonical boundary number */
1956: int ierr;
1959: /* Find canonical boundary number */
1960: for(b = 0; b < mesh->numBd; b++)
1961: if (tri->bdMarkers[b] == boundary) break;
1962: if (b == mesh->numBd) SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid boundary %d specified", boundary);
1963: if (mesh->activeBd != -1) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Already iterating over boundary %d", mesh->activeBd);
1964: /* Find first boundary node */
1965: mesh->activeBd = b;
1966: mesh->activeBdOld = b;
1967: mesh->activeBdNode = tri->bdBegin[b] - 1;
1968: MeshGetBoundaryNext(mesh, boundary, ghost, node);
1969: PetscFunctionReturn(ierr);
1970: }
1972: int MeshGetBoundaryNext_Triangular_2D(Mesh mesh, int boundary, PetscTruth ghost, int *node)
1973: {
1974: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1975: Partition p = mesh->part;
1976: Partition_Triangular_2D *q = (Partition_Triangular_2D *) p->data;
1977: int numBd = mesh->numBd;
1978: int numNodes = mesh->numNodes;
1979: int b; /* Canonical boundary number */
1980: int tempNode; /* Canonical local node number */
1981: int ierr;
1984: /* Find canonical boundary number */
1985: for(b = 0; b < numBd; b++) {
1986: if (tri->bdMarkers[b] == boundary) break;
1987: }
1988: if (b == numBd) SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid boundary %d specified", boundary);
1989: /* Find next boundary node */
1990: if (mesh->activeBd != b) SETERRQ1(PETSC_ERR_ARG_WRONG, "Boundary %d not active", boundary);
1991: do {
1992: mesh->activeBdNode++;
1993: if (mesh->activeBdNode == tri->bdBegin[b+1]) {
1994: mesh->activeBd = mesh->activeBdNode = -1;
1995: *node = -1;
1996: return(0);
1997: }
1998: tempNode = tri->bdNodes[mesh->activeBdNode] - q->firstNode[p->rank];
1999: /* Translate ghost nodes */
2000: if (ghost == PETSC_TRUE) {
2001: if ((tempNode < 0) || (tempNode >= numNodes)) {
2002: PartitionGhostNodeIndex_Private(p, tri->bdNodes[mesh->activeBdNode], &tempNode);
2003: if (ierr == 0) {
2004: tempNode += numNodes;
2005: break;
2006: }
2007: }
2008: }
2009: }
2010: while((tempNode < 0) || (tempNode >= numNodes));
2011: *node = tempNode;
2012: return(0);
2013: }
2015: int MeshGetActiveBoundary_Triangular_2D(Mesh mesh, int *boundary)
2016: {
2017: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
2020: if (mesh->activeBdOld >= 0) {
2021: *boundary = tri->bdMarkers[mesh->activeBdOld];
2022: } else {
2023: *boundary = 0;
2024: }
2025: return(0);
2026: }
2028: int MeshGetNodeSupport_Triangular_2D(Mesh mesh, int node, int startElem, int *degree, int **support)
2029: {
2030: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
2031: int numOverlapElements = mesh->numFaces;
2032: int numCorners = mesh->numCorners;
2033: int *elements = tri->faces;
2034: int *neighbors = tri->neighbors;
2035: int deg = -1;
2036: PetscTruth found;
2037: int firstElem;
2038: int prevElem, nextElem;
2039: int elem, neighbor, corner;
2042: /* First check suggested element */
2043: firstElem = -1;
2044: for(corner = 0; corner < numCorners; corner++) {
2045: if (elements[startElem*numCorners+corner] == node) {
2046: firstElem = startElem;
2047: if (corner >= 3) {
2048: deg = 1;
2049: } else {
2050: deg = 0;
2051: }
2052: break;
2053: }
2054: }
2056: /* Locate an element containing the node */
2057: if (mesh->part != PETSC_NULL) {
2058: numOverlapElements = mesh->part->numOverlapElements;
2059: }
2060: for(elem = 0; (elem < numOverlapElements) && (firstElem < 0); elem++) {
2061: for(corner = 0; corner < numCorners; corner++) {
2062: if (elements[elem*numCorners+corner] == node) {
2063: firstElem = elem;
2064: if (corner >= 3) {
2065: deg = 1;
2066: } else {
2067: deg = 0;
2068: }
2069: break;
2070: }
2071: }
2072: }
2073: if (firstElem < 0) {
2074: SETERRQ1(PETSC_ERR_ARG_WRONG, "Node %d not found in mesh", node);
2075: }
2077: /* Check for midnode */
2078: if (deg != 0) {
2079: mesh->support[0] = firstElem;
2080: mesh->support[1] = -1;
2081: /* Find other element containing node */
2082: for(neighbor = 0; neighbor < 3; neighbor++) {
2083: elem = neighbors[firstElem*3+neighbor];
2084: if (elem < 0)
2085: continue;
2087: for(corner = 3; corner < numCorners; corner++) {
2088: if (elements[elem*numCorners+corner] == node) {
2089: mesh->support[1] = elem;
2090: deg++;
2091: break;
2092: }
2093: }
2094: if (corner < numCorners) break;
2095: }
2096: #ifdef PETSC_USE_BOPT_g
2097: if (mesh->support[1] == -1) {
2098: if (tri->markers[node] == 0) SETERRQ1(PETSC_ERR_PLIB, "Interior node %d with incomplete support", node);
2099: }
2100: #endif
2101: *support = mesh->support;
2102: *degree = deg;
2103: return(0);
2104: }
2106: /* Walk around node */
2107: prevElem = -1;
2108: nextElem = -1;
2109: found = PETSC_TRUE;
2110: while((nextElem != firstElem) && (found == PETSC_TRUE)) {
2111: /* First time through */
2112: if (prevElem == -1)
2113: nextElem = firstElem;
2115: /* Add element to the list */
2116: mesh->support[deg++] = nextElem;
2118: /* Look for a neighbor containing the node */
2119: found = PETSC_FALSE;
2120: for(neighbor = 0; neighbor < 3; neighbor++) {
2121: /* Reject the element we just came from (we could reject the opposite face also) */
2122: elem = neighbors[nextElem*3+neighbor];
2123: if ((elem == prevElem) || (elem < 0)) continue;
2125: for(corner = 0; corner < 3; corner++) {
2126: if (elements[elem*numCorners+corner] == node) {
2127: prevElem = nextElem;
2128: nextElem = elem;
2129: found = PETSC_TRUE;
2130: break;
2131: }
2132: }
2133: if (corner < 3) break;
2134: }
2136: if (found == PETSC_FALSE) {
2137: #ifdef PETSC_USE_BOPT_g
2138: /* Make sure that this is a boundary node */
2139: if (tri->markers[node] == 0) SETERRQ(PETSC_ERR_PLIB, "Interior node with incomplete support");
2140: #endif
2142: /* Also walk the other way for a boundary node --
2143: We could be nice and reorder the elements afterwards to all be adjacent */
2144: if (deg > 1) {
2145: prevElem = mesh->support[1];
2146: nextElem = mesh->support[0];
2147: found = PETSC_TRUE;
2148: while(found == PETSC_TRUE) {
2149: /* Look for a neighbor containing the node */
2150: found = PETSC_FALSE;
2151: for(neighbor = 0; neighbor < 3; neighbor++) {
2152: /* Reject the element we just came from (we could reject the opposite face also) */
2153: elem = neighbors[nextElem*3+neighbor];
2154: if ((elem == prevElem) || (elem < 0)) continue;
2156: for(corner = 0; corner < 3; corner++) {
2157: if (elements[elem*numCorners+corner] == node) {
2158: prevElem = nextElem;
2159: nextElem = elem;
2160: found = PETSC_TRUE;
2162: /* Add element to the list */
2163: mesh->support[deg++] = nextElem;
2164: break;
2165: }
2166: }
2167: if (corner < 3) break;
2168: }
2169: }
2170: }
2171: }
2173: #ifdef PETSC_USE_BOPT_g
2174: /* Check for overflow */
2175: if (deg > mesh->maxDegree) SETERRQ1(PETSC_ERR_MEM, "Maximum degree %d exceeded", mesh->maxDegree);
2176: #endif
2177: }
2179: *support = mesh->support;
2180: *degree = deg;
2181: return(0);
2182: }
2184: #if 0
2185: int MeshLocatePoint_Triangular_2D_Directed(Mesh mesh, double xx, double yy, double zz, int *containingElem)
2186: {
2187: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
2188: double start_x, start_y; /* Coordinates of a vertex */
2189: double v_x, v_y; /* Coordinates of an edge vector */
2190: double p_x, p_y; /* Coordinates of the vector from a vertex (start) to a point */
2191: double cross; /* The magnitude of the cross product of two vectors */
2192: PetscTruth inside[3]; /* Are we in the correct halfplane? */
2193: PetscTruth found; /* We have found the element containing a point */
2194: int nextElem; /* The next element to search */
2195: int prevElem; /* The last element searched */
2196: int numCorners = mesh->numCorners;
2197: int numElements = mesh->numFaces;
2198: double *nodes = tri->nodes;
2199: int *elements = tri->faces;
2200: int *neighbors = tri->neighbors;
2201: double x = xx;
2202: double y = yy;
2203: double coords[6];
2204: int neighbor;
2205: int elem, edge;
2208: *containingElem = -1;
2209: /* Quick bounding rectangle check */
2210: if ((mesh->isPeriodic == PETSC_FALSE) &&
2211: ((mesh->locStartX > x) || (mesh->locEndX < x) || (mesh->locStartY > y) || (mesh->locEndY < y)))
2212: return(0);
2213: /* Simple search */
2214: for(elem = 0, found = PETSC_FALSE; elem < numElements; elem++) {
2215: /* A point lies within a triangle is it is on the inner side of every edge --
2216: We take the cross product of the edge oriented counterclockwise with
2217: the vector from the edge start to the point:
2219: 2
2220: |
2221: |
2222: |
2223: 0---1----> v_0
2224:
2225: theta
2226:
2227: P
2229: If
2231: |v_i cross P| equiv |v_i| |P| sintheta > 0
2233: then the point lies in the correct half-plane. This has the advantage of
2234: only using elementary arithmetic. Robust predicates could be used to catch
2235: points near a boundary.
2236: */
2237: if (mesh->isPeriodic == PETSC_TRUE) {
2238: coords[0] = nodes[elements[elem*numCorners]*2];
2239: coords[1] = nodes[elements[elem*numCorners]*2+1];
2240: coords[2] = MeshPeriodicRelativeX(mesh, nodes[elements[elem*numCorners+1]*2], nodes[elements[elem*numCorners]*2]);
2241: coords[3] = MeshPeriodicRelativeY(mesh, nodes[elements[elem*numCorners+1]*2+1], nodes[elements[elem*numCorners]*2+1]);
2242: coords[4] = MeshPeriodicRelativeX(mesh, nodes[elements[elem*numCorners+2]*2], nodes[elements[elem*numCorners]*2]);
2243: coords[5] = MeshPeriodicRelativeY(mesh, nodes[elements[elem*numCorners+2]*2+1], nodes[elements[elem*numCorners]*2+1]);
2244: x = MeshPeriodicRelativeX(mesh, xx, nodes[elements[elem*numCorners]*2]);
2245: y = MeshPeriodicRelativeY(mesh, yy, nodes[elements[elem*numCorners]*2+1]);
2246: }
2247: for(edge = 0; edge < 3; edge++) {
2248: /* First edge 0--1 */
2249: if (mesh->isPeriodic == PETSC_TRUE) {
2250: start_x = coords[edge*2];
2251: start_y = coords[edge*2+1];
2252: v_x = coords[((edge+1)%3)*2] - start_x;
2253: v_y = coords[((edge+1)%3)*2+1] - start_y;
2254: } else {
2255: start_x = nodes[elements[elem*numCorners+edge]*2];
2256: start_y = nodes[elements[elem*numCorners+edge]*2+1];
2257: v_x = nodes[elements[elem*numCorners+((edge+1)%3)]*2] - start_x;
2258: v_y = nodes[elements[elem*numCorners+((edge+1)%3)]*2+1] - start_y;
2259: }
2260: p_x = x - start_x;
2261: p_y = y - start_y;
2262: cross = v_x*p_y - v_y*p_x;
2264: /* We will use this routine to find integration points which lie on the
2265: boundary so we must give a little slack to the test. I do not think
2266: this will bias our results. */
2267: if ((cross >= 0) || (PetscAbsScalar(cross) < PETSC_MACHINE_EPSILON)) {
2268: inside[edge] = PETSC_TRUE;
2269: } else {
2270: inside[edge] = PETSC_FALSE;
2271: }
2272: }
2274: /* This information could be used for a directed search --
2276: 4 T|
2277: F|
2278: F|
2279: |
2280: |
2281: 2
2282: T | T
2283: 5 T | F 3
2284: F |T T
2285: ----0---1------
2286: F | F F
2287: 6 T | T 1 F 2
2288: F | T T
2290: and clearly all F's means that the area of the triangle is 0.
2291: The regions are ordered so that all the Hamming distances are
2292: 1 (Gray code).
2293: */
2294: if (inside[2] == PETSC_TRUE) {
2295: if (inside[1] == PETSC_TRUE) {
2296: if (inside[0] == PETSC_TRUE) {
2297: /* Region 0 -- TTT */
2298: found = PETSC_TRUE;
2299: } else {
2300: /* Region 1 -- FTT */
2301: nextElem = neighbors[elem*3+2];
2302: }
2303: } else {
2304: if (inside[0] == PETSC_TRUE) {
2305: /* Region 3 -- TFT */
2306: nextElem = neighbors[elem*3];
2307: } else {
2308: /* Region 2 -- FFT */
2309: if (neighbors[elem*3] >= 0) {
2310: nextElem = neighbors[elem*3];
2311: } else {
2312: nextElem = neighbors[elem*3+2];
2313: }
2314: }
2315: }
2316: } else {
2317: if (inside[1] == PETSC_TRUE) {
2318: if (inside[0] == PETSC_TRUE) {
2319: /* Region 5 -- TTF */
2320: nextElem = neighbors[elem*3+1];
2321: } else {
2322: /* Region 6 -- FTF */
2323: if (neighbors[elem*3+1] >= 0) {
2324: nextElem = neighbors[elem*3+1];
2325: } else {
2326: nextElem = neighbors[elem*3+2];
2327: }
2328: }
2329: } else {
2330: if (inside[0] == PETSC_TRUE) {
2331: /* Region 4 -- TFF */
2332: if (neighbors[elem*3] >= 0) {
2333: nextElem = neighbors[elem*3];
2334: } else {
2335: nextElem = neighbors[elem*3+1];
2336: }
2337: } else {
2338: /* Region 7 -- FFF */
2339: PetscPrintf(PETSC_COMM_SELF, "Element %d: (%g, %g)-(%g, %g)-(%g, %g) and (%g, %g)n", elem,
2340: nodes[elements[elem*numCorners]*2], nodes[elements[elem*numCorners]*2+1],
2341: nodes[elements[elem*numCorners+1]*2], nodes[elements[elem*numCorners+1]*2+1],
2342: nodes[elements[elem*numCorners+2]*2], nodes[elements[elem*numCorners+2]*2+1], x, y);
2343: SETERRQ1(PETSC_ERR_MESH_NULL_ELEM, "Element %d has no area", elem);
2344: }
2345: }
2346: }
2348: /* We found the point */
2349: if (found) {
2350: *containingElem = elem;
2351: break;
2352: }
2354: /* Choose next direction */
2355: if (nextElem < 0) {
2356: /* Choose randomly */
2357: for(neighbor = 0; neighbor < 3; neighbor++) {
2358: tempElem = neighbors[elem*3+neighbor];
2359: if ((tempElem >= 0) && (tempElem != prevElem)) {
2360: nextElem = tempElem;
2361: break;
2362: }
2363: }
2364: if (nextElem < 0) {
2365: SETERRQ1(PETSC_ERR_MAX_ITER, "Element search hit dead end in element %d", elem);
2366: }
2367: }
2369: prevElem = elem;
2370: }
2372: return(0);
2373: }
2374: #endif
2376: int MeshLocatePoint_Triangular_2D(Mesh mesh, double xx, double yy, double zz, int *containingElem) {
2377: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
2378: double start_x, start_y; /* Coordinates of a vertex */
2379: double v_x, v_y; /* Coordinates of an edge vector */
2380: double p_x, p_y; /* Coordinates of the vector from a vertex (start) to a point */
2381: double cross; /* The magnitude of the cross product of two vectors */
2382: PetscTruth inside[3]; /* Are we in the correct halfplane? */
2383: PetscTruth found; /* We have found the element containing a point */
2384: int numCorners = mesh->numCorners;
2385: int numElements = mesh->numFaces;
2386: double *nodes = tri->nodes;
2387: int *markers = tri->markers;
2388: int *elements = tri->faces;
2389: int *neighbors = tri->neighbors;
2390: double x = xx;
2391: double y = yy;
2392: double coords[6];
2393: double pNormSq;
2394: int elem, edge;
2397: *containingElem = -1;
2398: /* Quick bounding rectangle check */
2399: if ((mesh->isPeriodic == PETSC_FALSE) &&
2400: ((mesh->locStartX > x) || (mesh->locEndX < x) || (mesh->locStartY > y) || (mesh->locEndY < y)))
2401: return(0);
2402: /* Simple search */
2403: for(elem = 0, found = PETSC_FALSE; elem < numElements; elem++) {
2404: /* A point lies within a triangle is it is on the inner side of every edge --
2405: We take the cross product of the edge oriented counterclockwise with
2406: the vector from the edge start to the point:
2408: 2
2409: |
2410: |
2411: |
2412: 0---1----> v_0
2413:
2414: theta
2415:
2416: P
2418: If
2420: |v_i cross P| equiv |v_i| |P| sintheta > 0
2422: then the point lies in the correct half-plane. This has the advantage of
2423: only using elementary arithmetic. Robust predicates could be used to catch
2424: points near a boundary.
2425: */
2426: if (mesh->isPeriodic == PETSC_TRUE) {
2427: coords[0] = nodes[elements[elem*numCorners]*2];
2428: coords[1] = nodes[elements[elem*numCorners]*2+1];
2429: coords[2] = MeshPeriodicRelativeX(mesh, nodes[elements[elem*numCorners+1]*2], nodes[elements[elem*numCorners]*2]);
2430: coords[3] = MeshPeriodicRelativeY(mesh, nodes[elements[elem*numCorners+1]*2+1], nodes[elements[elem*numCorners]*2+1]);
2431: coords[4] = MeshPeriodicRelativeX(mesh, nodes[elements[elem*numCorners+2]*2], nodes[elements[elem*numCorners]*2]);
2432: coords[5] = MeshPeriodicRelativeY(mesh, nodes[elements[elem*numCorners+2]*2+1], nodes[elements[elem*numCorners]*2+1]);
2433: x = MeshPeriodicRelativeX(mesh, xx, nodes[elements[elem*numCorners]*2]);
2434: y = MeshPeriodicRelativeY(mesh, yy, nodes[elements[elem*numCorners]*2+1]);
2435: }
2436: for(edge = 0; edge < 3; edge++) {
2437: /* First edge 0--1 */
2438: if (mesh->isPeriodic == PETSC_TRUE) {
2439: start_x = coords[edge*2];
2440: start_y = coords[edge*2+1];
2441: v_x = coords[((edge+1)%3)*2] - start_x;
2442: v_y = coords[((edge+1)%3)*2+1] - start_y;
2443: } else {
2444: start_x = nodes[elements[elem*numCorners+edge]*2];
2445: start_y = nodes[elements[elem*numCorners+edge]*2+1];
2446: v_x = nodes[elements[elem*numCorners+((edge+1)%3)]*2] - start_x;
2447: v_y = nodes[elements[elem*numCorners+((edge+1)%3)]*2+1] - start_y;
2448: }
2449: p_x = x - start_x;
2450: p_y = y - start_y;
2451: cross = v_x*p_y - v_y*p_x;
2452: pNormSq = p_x*p_x+p_y*p_y;
2454: /* Check for identical points, being looser at the boundary */
2455: if (pNormSq < PETSC_MACHINE_EPSILON) {
2456: *containingElem = elem;
2457: return(0);
2458: } else if (markers[elements[elem*numCorners+edge]] != 0) {
2459: if (pNormSq < 0.001*PetscMax(v_x*v_x, v_y*v_y)) {
2460: *containingElem = elem;
2461: return(0);
2462: }
2463: }
2465: if (cross >= 0) {
2466: inside[edge] = PETSC_TRUE;
2467: } else if ((neighbors[elem*3+((edge+2)%3)] < 0) && (cross*cross < 0.087155743*(v_x*v_x+v_y*v_y)*pNormSq)) {
2468: /* We are losing precision here (I think) for large problems (since points are turning up inside particles)
2469: and therefore I am putting a 5 degree limit on this check (sin 5 = 0.087155743), which is only active for
2470: boundary edges. */
2471: inside[edge] = PETSC_TRUE;
2472: } else {
2473: inside[edge] = PETSC_FALSE;
2474: }
2475: }
2477: /* This information could be used for a directed search --
2479: 4 T|
2480: F|
2481: F|
2482: |
2483: |
2484: 2
2485: T | T
2486: 5 T | F 3
2487: F |T T
2488: ----0---1------
2489: F | F F
2490: 6 T | T 1 F 2
2491: F | T T
2493: and clearly all F's means that the area of the triangle is 0.
2494: The regions are ordered so that all the Hamming distances are
2495: 1 (Gray code).
2496: */
2497: if ((inside[2] == PETSC_TRUE) && (inside[1] == PETSC_TRUE) && (inside[0] == PETSC_TRUE)) {
2498: double minX = coords[0];
2499: double maxX = coords[0];
2500: double minY = coords[1];
2501: double maxY = coords[1];
2502: double offX, offY;
2504: if (coords[2] < minX) minX = coords[2];
2505: if (coords[2] > maxX) maxX = coords[2];
2506: if (coords[3] < minY) minY = coords[3];
2507: if (coords[3] > maxY) maxY = coords[3];
2508: if (coords[4] < minX) minX = coords[4];
2509: if (coords[4] > maxX) maxX = coords[4];
2510: if (coords[5] < minY) minY = coords[5];
2511: if (coords[5] > maxY) maxY = coords[5];
2512: offX = 0.01*(maxX - minX);
2513: offY = 0.01*(maxY - minY);
2515: if ((x >= minX-offX) && (x <= maxX+offX) && (y >= minY-offY) && (y <= maxY+offY)) {
2516: /* Region 0 -- TTT */
2517: found = PETSC_TRUE;
2518: }
2519: } else if ((inside[2] == PETSC_FALSE) && (inside[1] == PETSC_FALSE) && (inside[0] == PETSC_FALSE)) {
2520: /* Region 7 -- FFF */
2521: PetscPrintf(PETSC_COMM_SELF, "Element %d: (%g, %g)-(%g, %g)-(%g, %g) and (%g, %g)n", elem, coords[0], coords[1], coords[2],
2522: coords[3], coords[4], coords[5], x, y);
2523: SETERRQ1(PETSC_ERR_MESH_NULL_ELEM, "Element %d has no area or is inverted", elem);
2524: }
2526: /* We found the point */
2527: if (found) {
2528: *containingElem = elem;
2529: break;
2530: }
2531: }
2533: return(0);
2534: }
2536: int MeshNearestNode_Triangular_2D(Mesh mesh, double x, double y, double z, PetscTruth outsideDomain, int *node)
2537: {
2538: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
2539: Partition p = mesh->part;
2540: Partition_Triangular_2D *q = (Partition_Triangular_2D *) p->data;
2541: int numCorners = mesh->numCorners;
2542: int *elements = tri->faces;
2543: int numNodes = mesh->numNodes;
2544: double *nodes = tri->nodes;
2545: int *firstNode = q->firstNode;
2546: int rank = p->rank;
2547: int elem;
2548: double minDist, dist, allMinDist;
2549: int minNode, globalMinNode, allMinNode;
2550: int corner, n;
2551: int ierr;
2554: if (outsideDomain == PETSC_FALSE) {
2555: MeshLocatePoint(mesh, x, y, z, &elem);
2556: if (elem >= 0) {
2557: minDist = PetscSqr(MeshPeriodicDiffX(mesh, nodes[elements[elem*numCorners]*2] - x)) +
2558: PetscSqr(MeshPeriodicDiffY(mesh, nodes[elements[elem*numCorners]*2+1] - y));
2559: minNode = elements[elem*numCorners];
2560: for(corner = 1; corner < numCorners; corner++) {
2561: dist = PetscSqr(MeshPeriodicDiffX(mesh, nodes[elements[elem*numCorners+corner]*2] - x)) +
2562: PetscSqr(MeshPeriodicDiffY(mesh, nodes[elements[elem*numCorners+corner]*2+1] - y));
2563: if (dist < minDist) {
2564: minDist = dist;
2565: minNode = elements[elem*numCorners+corner];
2566: }
2567: }
2568: PartitionLocalToGlobalNodeIndex(p, minNode, &globalMinNode);
2569: } else {
2570: minNode = -1;
2571: globalMinNode = -1;
2572: }
2574: /* The minimum node might still be a ghost node */
2575: MPI_Allreduce(&globalMinNode, &allMinNode, 1, MPI_INT, MPI_MAX, p->comm);
2576: if ((allMinNode >= firstNode[rank]) && (allMinNode < firstNode[rank+1])) {
2577: *node = allMinNode - firstNode[rank];
2578: } else {
2579: *node = -1;
2580: }
2581: if (allMinNode < 0)
2582: PetscFunctionReturn(1);
2583: } else {
2584: /* Brute force check */
2585: minNode = 0;
2586: minDist = PetscSqr(MeshPeriodicDiffX(mesh, nodes[0*2] - x)) +
2587: PetscSqr(MeshPeriodicDiffY(mesh, nodes[0*2+1] - y));
2588: for(n = 1; n < numNodes; n++) {
2589: dist = PetscSqr(MeshPeriodicDiffX(mesh, nodes[n*2] - x)) +
2590: PetscSqr(MeshPeriodicDiffY(mesh, nodes[n*2+1] - y));
2591: if (dist < minDist) {
2592: minDist = dist;
2593: minNode = n;
2594: }
2595: }
2597: /* Find the minimum distance */
2598: MPI_Allreduce(&minDist, &allMinDist, 1, MPI_DOUBLE, MPI_MIN, p->comm);
2599: if (minDist == allMinDist) {
2600: PartitionLocalToGlobalNodeIndex(p, minNode, &globalMinNode);
2601: } else {
2602: globalMinNode = -1;
2603: }
2605: /* We might still have ties */
2606: MPI_Allreduce(&globalMinNode, &allMinNode, 1, MPI_INT, MPI_MAX, p->comm);
2607: if ((allMinNode >= firstNode[rank]) && (allMinNode < firstNode[rank+1])) {
2608: *node = allMinNode - firstNode[rank];
2609: } else {
2610: *node = -1;
2611: }
2612: if (allMinNode < 0)
2613: PetscFunctionReturn(1);
2614: }
2616: return(0);
2617: }
2619: int MeshSaveMesh_Triangular_2D(Mesh mesh) {
2620: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
2621: int numOverlapNodes;
2622: int ierr;
2625: ierr = PartitionGetNumOverlapNodes(mesh->part, &numOverlapNodes);
2626: if (tri->nodesOld == PETSC_NULL) {
2627: PetscMalloc(numOverlapNodes*2 * sizeof(double), &tri->nodesOld);
2628: PetscLogObjectMemory(mesh, numOverlapNodes*2 * sizeof(double));
2629: }
2630: PetscMemcpy(tri->nodesOld, tri->nodes, numOverlapNodes*2 * sizeof(double));
2631: return(0);
2632: }
2634: int MeshRestoreMesh_Triangular_2D(Mesh mesh) {
2635: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
2636: int numOverlapNodes;
2637: int ierr;
2640: ierr = PartitionGetNumOverlapNodes(mesh->part, &numOverlapNodes);
2641: if (tri->nodesOld != PETSC_NULL) {
2642: PetscMemcpy(tri->nodes, tri->nodesOld, numOverlapNodes*2 * sizeof(double));
2643: }
2644: return(0);
2645: }
2647: int MeshIsDistorted_Triangular_2D(Mesh mesh, PetscTruth *flag) {
2648: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
2649: double *nodes = tri->nodes;
2650: int *faces = tri->faces;
2651: int numCorners = mesh->numCorners;
2652: int numElements = mesh->numFaces;
2653: double maxAspectRatio = mesh->maxAspectRatio;
2654: double area; /* Twice the area of the triangle */
2655: double aspectRatio; /* (Longest edge)^2/(Twice the area of the triangle) */
2656: double max;
2657: PetscTruth locFlag;
2658: double x21, y21, x31, y31, x32, y32;
2659: int elem;
2660: int locRetVal;
2661: int retVal;
2662: int ierr;
2665: for(elem = 0, locFlag = PETSC_FALSE, locRetVal = 0, max = 0.0; elem < numElements; elem++) {
2666: /* Find the area */
2667: if (mesh->isPeriodic == PETSC_TRUE) {
2668: x21 = MeshPeriodicDiffX(mesh, nodes[faces[elem*numCorners+1]*2] - nodes[faces[elem*numCorners]*2]);
2669: y21 = MeshPeriodicDiffY(mesh, nodes[faces[elem*numCorners+1]*2+1] - nodes[faces[elem*numCorners]*2+1]);
2670: x31 = MeshPeriodicDiffX(mesh, nodes[faces[elem*numCorners+2]*2] - nodes[faces[elem*numCorners]*2]);
2671: y31 = MeshPeriodicDiffY(mesh, nodes[faces[elem*numCorners+2]*2+1] - nodes[faces[elem*numCorners]*2+1]);
2672: x32 = MeshPeriodicDiffX(mesh, nodes[faces[elem*numCorners+2]*2] - nodes[faces[elem*numCorners+1]*2]);
2673: y32 = MeshPeriodicDiffY(mesh, nodes[faces[elem*numCorners+2]*2+1] - nodes[faces[elem*numCorners+1]*2+1]);
2674: } else {
2675: x21 = nodes[faces[elem*numCorners+1]*2] - nodes[faces[elem*numCorners]*2];
2676: y21 = nodes[faces[elem*numCorners+1]*2+1] - nodes[faces[elem*numCorners]*2+1];
2677: x31 = nodes[faces[elem*numCorners+2]*2] - nodes[faces[elem*numCorners]*2];
2678: y31 = nodes[faces[elem*numCorners+2]*2+1] - nodes[faces[elem*numCorners]*2+1];
2679: x32 = nodes[faces[elem*numCorners+2]*2] - nodes[faces[elem*numCorners+1]*2];
2680: y32 = nodes[faces[elem*numCorners+2]*2+1] - nodes[faces[elem*numCorners+1]*2+1];
2681: }
2682: area = x21*y31 - x31*y21;
2683: /* Check for inverted elements */
2684: if (area < 0.0) {
2685: locFlag = PETSC_TRUE;
2686: locRetVal = 1;
2687: break;
2688: } else if (area == 0.0) {
2689: SETERRQ1(PETSC_ERR_ARG_CORRUPT, "Element %d has no area", elem);
2690: }
2692: /* Find the aspect ratio = (longest edge length)^2 / (2 area) */
2693: aspectRatio = PetscMax(x21*x21 + y21*y21, x31*x31 + y31*y31);
2694: aspectRatio = PetscMax(x32*x32 + y32*y32, aspectRatio);
2695: aspectRatio /= area;
2696: /* We cannot bail here since we must check for inverted elements */
2697: max = PetscMax(max, aspectRatio);
2698: }
2699: if (max > maxAspectRatio) {
2700: PetscLogInfo(mesh, "Aspect ratio too large in element %d: %g > %gn", elem, max, maxAspectRatio);
2701: locFlag = PETSC_TRUE;
2702: } else {
2703: PetscLogInfo(mesh, "Maximum aspect ratio in element %d: %gn", elem, max);
2704: }
2705: PetscLogFlops(elem*19);
2706: MPI_Allreduce(&locFlag, flag, 1, MPI_INT, MPI_LOR, mesh->comm);
2707: MPI_Allreduce(&locRetVal, &retVal, 1, MPI_INT, MPI_LOR, mesh->comm);
2708: PetscFunctionReturn(retVal);
2709: }
2711: static struct _MeshOps MeshOps = {/* Generic Operations */
2712: PETSC_NULL/* MeshSetup */,
2713: PETSC_NULL/* MeshSetFromOptions */,
2714: MeshView_Triangular_2D,
2715: PETSC_NULL/* MeshCopy */,
2716: PETSC_NULL/* MeshDuplicate */,
2717: MeshDestroy_Triangular_2D,
2718: /* Mesh-Specific Operations */
2719: MeshPartition_Triangular_2D,
2720: MeshCoarsen_Triangular_2D,
2721: MeshRefine_Triangular_2D,
2722: MeshResetNodes_Triangular_2D,
2723: MeshSaveMesh_Triangular_2D,
2724: MeshRestoreMesh_Triangular_2D,
2725: /* Mesh Query Functions */
2726: MeshUpdateBoundingBox_Triangular_2D,
2727: MeshIsDistorted_Triangular_2D,
2728: /* Mesh Boundary Query Functions */
2729: MeshGetBoundarySize_Triangular_2D,
2730: MeshGetBoundaryIndex_Triangular_2D,
2731: MeshGetBoundaryStart_Triangular_2D,
2732: MeshGetBoundaryNext_Triangular_2D,
2733: MeshGetActiveBoundary_Triangular_2D,
2734: /* Mesh Node Query Functions */
2735: MeshGetNodeBoundary_Triangular_2D,
2736: MeshNodeIsVertex_Triangular_2D,
2737: MeshGetNodeCoords_Triangular_2D,
2738: MeshSetNodeCoords_Triangular_2D,
2739: MeshGetNodeCoordsSaved_Triangular_2D,
2740: MeshNearestNode_Triangular_2D,
2741: MeshGetNodeSupport_Triangular_2D,
2742: /* Mesh Element Query Functions */
2743: MeshGetElemNeighbor_Triangular_2D,
2744: MeshLocatePoint_Triangular_2D,
2745: /* Mesh Embedding Query Functions */
2746: MeshGetNodeFromElement_Triangular_2D,
2747: MeshGetNodeFromEdge_Triangular_2D,
2748: /* CSR Support Functions */
2749: MeshCreateLocalCSR_Triangular_2D,
2750: MeshCreateFullCSR_Triangular_2D,
2751: MeshCreateDualCSR_Triangular_2D};
2753: /*@
2754: MeshCreateTriangular2DCSR - Creates a basic two dimensional unstructured grid from an existing CSR representation.
2756: Collective on MPI_Comm
2758: Input Parameters:
2759: + comm - The communicator that shares the grid
2760: . numNodes - The number of mesh nodes, here identical to vertices
2761: . numFaces - The number of elements in the mesh
2762: . nodes - The coordinates for each node
2763: . offsets - The offset of each node's adjacency list
2764: . adj - The adjacency list for the mesh
2765: - faces - The list of nodes for each element
2767: Output Parameter:
2768: . mesh - The mesh
2770: Level: beginner
2772: .keywords unstructured mesh
2773: .seealso MeshCreateTriangular2D(), MeshCreateTriangular2DPeriodic()
2774: @*/
2775: int MeshCreateTriangular2DCSR(MPI_Comm comm, int numNodes, int numFaces, double *nodes, int *offsets, int *adj, int *faces, Mesh *mesh)
2776: {
2777: ParameterDict dict;
2778: MeshBoundary2D *bdCtx;
2779: int ierr;
2782: /* REWORK: We are really stretching this interface, but I need this yesterday */
2783: PetscMalloc(sizeof(MeshBoundary2D), &bdCtx);
2784: bdCtx->numBd = 0;
2785: bdCtx->numVertices = numNodes;
2786: bdCtx->vertices = nodes;
2787: bdCtx->markers = offsets;
2788: bdCtx->segMarkers = adj;
2789: bdCtx->numSegments = numFaces;
2790: bdCtx->segments = faces;
2791: bdCtx->numHoles = 0;
2792: bdCtx->holes = PETSC_NULL;
2793: MeshCreate(comm, mesh);
2794: ParameterDictCreate(comm, &dict);
2795: ParameterDictSetObject(dict, "bdCtx", bdCtx);
2796: ParameterDictSetInteger(dict, "numLocNodes", 3);
2797: PetscObjectSetParameterDict((PetscObject) *mesh, dict);
2798: ParameterDictDestroy(dict);
2799: MeshSetDimension(*mesh, 2);
2800: /* Setup special mechanisms */
2801: PetscObjectComposeFunction((PetscObject) *mesh, "MeshTriangular2D_Create_C", "MeshCreate_CSR",
2802: (void (*)(void)) MeshCreate_CSR);
2803:
2804: PetscObjectComposeFunction((PetscObject) *mesh, "MeshTriangular2D_Refine_C", "MeshRefine_CSR",
2805: (void (*)(void)) MeshRefine_CSR);
2806:
2807: MeshSetType(*mesh, MESH_TRIANGULAR_2D);
2809: PetscObjectSetName((PetscObject) *mesh, "Constructed Mesh");
2810: PetscFree(bdCtx);
2811: return(0);
2812: }
2814: EXTERN_C_BEGIN
2815: int MeshSerialize_Triangular_2D(MPI_Comm comm, Mesh *mesh, PetscViewer viewer, PetscTruth store)
2816: {
2817: Mesh m;
2818: Mesh_Triangular *tri;
2819: int fd;
2820: int zero = 0;
2821: int one = 1;
2822: int numNodes, numEdges, numFaces, numCorners, numBd, numBdNodes, numBdEdges, old;
2823: int ierr;
2826: PetscViewerBinaryGetDescriptor(viewer, &fd);
2827: if (store) {
2828: m = *mesh;
2829: PetscBinaryWrite(fd, &m->highlightElement, 1, PETSC_INT, 0);
2830: PetscBinaryWrite(fd, &m->maxDegree, 1, PETSC_INT, 0);
2831: PetscBinaryWrite(fd, &m->startX, 1, PETSC_DOUBLE, 0);
2832: PetscBinaryWrite(fd, &m->endX, 1, PETSC_DOUBLE, 0);
2833: PetscBinaryWrite(fd, &m->sizeX, 1, PETSC_DOUBLE, 0);
2834: PetscBinaryWrite(fd, &m->startY, 1, PETSC_DOUBLE, 0);
2835: PetscBinaryWrite(fd, &m->endY, 1, PETSC_DOUBLE, 0);
2836: PetscBinaryWrite(fd, &m->sizeY, 1, PETSC_DOUBLE, 0);
2837: PetscBinaryWrite(fd, &m->startZ, 1, PETSC_DOUBLE, 0);
2838: PetscBinaryWrite(fd, &m->endZ, 1, PETSC_DOUBLE, 0);
2839: PetscBinaryWrite(fd, &m->sizeZ, 1, PETSC_DOUBLE, 0);
2840: PetscBinaryWrite(fd, &m->locStartX, 1, PETSC_DOUBLE, 0);
2841: PetscBinaryWrite(fd, &m->locEndX, 1, PETSC_DOUBLE, 0);
2842: PetscBinaryWrite(fd, &m->locSizeX, 1, PETSC_DOUBLE, 0);
2843: PetscBinaryWrite(fd, &m->locStartY, 1, PETSC_DOUBLE, 0);
2844: PetscBinaryWrite(fd, &m->locEndY, 1, PETSC_DOUBLE, 0);
2845: PetscBinaryWrite(fd, &m->locSizeY, 1, PETSC_DOUBLE, 0);
2846: PetscBinaryWrite(fd, &m->locStartZ, 1, PETSC_DOUBLE, 0);
2847: PetscBinaryWrite(fd, &m->locEndZ, 1, PETSC_DOUBLE, 0);
2848: PetscBinaryWrite(fd, &m->locSizeZ, 1, PETSC_DOUBLE, 0);
2849: PetscBinaryWrite(fd, &m->isPeriodic, 1, PETSC_INT, 0);
2850: PetscBinaryWrite(fd, &m->isPeriodicDim, 3, PETSC_INT, 0);
2852: PetscBinaryWrite(fd, &m->numBd, 1, PETSC_INT, 0);
2853: PetscBinaryWrite(fd, &m->numVertices, 1, PETSC_INT, 0);
2854: PetscBinaryWrite(fd, &m->numNodes, 1, PETSC_INT, 0);
2855: PetscBinaryWrite(fd, &m->numBdNodes, 1, PETSC_INT, 0);
2856: PetscBinaryWrite(fd, &m->numEdges, 1, PETSC_INT, 0);
2857: PetscBinaryWrite(fd, &m->numBdEdges, 1, PETSC_INT, 0);
2858: PetscBinaryWrite(fd, &m->numFaces, 1, PETSC_INT, 0);
2859: PetscBinaryWrite(fd, &m->numBdFaces, 1, PETSC_INT, 0);
2860: PetscBinaryWrite(fd, &m->numCells, 1, PETSC_INT, 0);
2861: PetscBinaryWrite(fd, &m->numCorners, 1, PETSC_INT, 0);
2862: PetscBinaryWrite(fd, &m->numCellCorners, 1, PETSC_INT, 0);
2863: PetscBinaryWrite(fd, &m->numHoles, 1, PETSC_INT, 0);
2865: PartitionSerialize(m, &m->part, viewer, store);
2867: PartitionGetNumOverlapElements(m->part, &numFaces);
2868: PartitionGetNumOverlapNodes(m->part, &numNodes);
2869: numEdges = m->numEdges;
2870: numCorners = m->numCorners;
2871: numBd = m->numBd;
2872: numBdNodes = m->numBdNodes;
2873: numBdEdges = m->numBdEdges;
2875: tri = (Mesh_Triangular *) (*mesh)->data;
2876: PetscBinaryWrite(fd, tri->nodes, numNodes*2, PETSC_DOUBLE, 0);
2877: if (tri->nodesOld == PETSC_NULL) {
2878: PetscBinaryWrite(fd, &zero, 1, PETSC_INT, 0);
2879: } else {
2880: PetscBinaryWrite(fd, &one, 1, PETSC_INT, 0);
2881: PetscBinaryWrite(fd, tri->nodesOld, numNodes*2, PETSC_DOUBLE, 0);
2882: }
2883: PetscBinaryWrite(fd, tri->markers, numNodes, PETSC_INT, 0);
2884: PetscBinaryWrite(fd, tri->degrees, numNodes, PETSC_INT, 0);
2885: PetscBinaryWrite(fd, tri->edges, numEdges*2, PETSC_INT, 0);
2886: PetscBinaryWrite(fd, tri->edgemarkers, numEdges, PETSC_INT, 0);
2887: PetscBinaryWrite(fd, tri->faces, numFaces*numCorners, PETSC_INT, 0);
2888: PetscBinaryWrite(fd, tri->neighbors, numFaces*3, PETSC_INT, 0);
2889: PetscBinaryWrite(fd, tri->bdNodes, numBdNodes, PETSC_INT, 0);
2890: PetscBinaryWrite(fd, tri->bdEdges, numBdEdges, PETSC_INT, 0);
2891: PetscBinaryWrite(fd, tri->bdMarkers, numBd, PETSC_INT, 0);
2892: PetscBinaryWrite(fd, tri->bdBegin, numBd+1, PETSC_INT, 0);
2893: PetscBinaryWrite(fd, tri->bdEdgeBegin, numBd+1, PETSC_INT, 0);
2894: PetscBinaryWrite(fd, m->holes, m->numHoles*2, PETSC_DOUBLE, 0);
2896: PetscBinaryWrite(fd, &m->maxAspectRatio, 1, PETSC_DOUBLE, 0);
2897: } else {
2898: /* Create the mesh context */
2899: MeshCreate(comm, &m);
2900: PetscNew(Mesh_Triangular, &tri);
2901: PetscLogObjectMemory(m, sizeof(Mesh_Triangular));
2902: PetscMemcpy(m->ops, &MeshOps, sizeof(struct _MeshOps));
2903: PetscStrallocpy(MESH_TRIANGULAR_2D, &m->type_name);
2904: m->data = (void *) tri;
2905: m->dim = 2;
2907: PetscBinaryRead(fd, &m->highlightElement, 1, PETSC_INT);
2908: PetscBinaryRead(fd, &m->maxDegree, 1, PETSC_INT);
2909: PetscMalloc(m->maxDegree * sizeof(int), &m->support);
2911: PetscBinaryRead(fd, &m->startX, 1, PETSC_DOUBLE);
2912: PetscBinaryRead(fd, &m->endX, 1, PETSC_DOUBLE);
2913: PetscBinaryRead(fd, &m->sizeX, 1, PETSC_DOUBLE);
2914: PetscBinaryRead(fd, &m->startY, 1, PETSC_DOUBLE);
2915: PetscBinaryRead(fd, &m->endY, 1, PETSC_DOUBLE);
2916: PetscBinaryRead(fd, &m->sizeY, 1, PETSC_DOUBLE);
2917: PetscBinaryRead(fd, &m->startZ, 1, PETSC_DOUBLE);
2918: PetscBinaryRead(fd, &m->endZ, 1, PETSC_DOUBLE);
2919: PetscBinaryRead(fd, &m->sizeZ, 1, PETSC_DOUBLE);
2920: PetscBinaryRead(fd, &m->locStartX, 1, PETSC_DOUBLE);
2921: PetscBinaryRead(fd, &m->locEndX, 1, PETSC_DOUBLE);
2922: PetscBinaryRead(fd, &m->locSizeX, 1, PETSC_DOUBLE);
2923: PetscBinaryRead(fd, &m->locStartY, 1, PETSC_DOUBLE);
2924: PetscBinaryRead(fd, &m->locEndY, 1, PETSC_DOUBLE);
2925: PetscBinaryRead(fd, &m->locSizeY, 1, PETSC_DOUBLE);
2926: PetscBinaryRead(fd, &m->locStartZ, 1, PETSC_DOUBLE);
2927: PetscBinaryRead(fd, &m->locEndZ, 1, PETSC_DOUBLE);
2928: PetscBinaryRead(fd, &m->locSizeZ, 1, PETSC_DOUBLE);
2929: PetscBinaryRead(fd, &m->isPeriodic, 1, PETSC_INT);
2930: PetscBinaryRead(fd, &m->isPeriodicDim, 3, PETSC_INT);
2932: m->activeBd = -1;
2933: m->activeBdOld = -1;
2934: m->activeBdNode = -1;
2935: tri->areas = PETSC_NULL;
2936: tri->aspectRatios = PETSC_NULL;
2938: PetscBinaryRead(fd, &m->numBd, 1, PETSC_INT);
2939: PetscBinaryRead(fd, &m->numVertices, 1, PETSC_INT);
2940: PetscBinaryRead(fd, &m->numNodes, 1, PETSC_INT);
2941: PetscBinaryRead(fd, &m->numBdNodes, 1, PETSC_INT);
2942: PetscBinaryRead(fd, &m->numEdges, 1, PETSC_INT);
2943: PetscBinaryRead(fd, &m->numBdEdges, 1, PETSC_INT);
2944: PetscBinaryRead(fd, &m->numFaces, 1, PETSC_INT);
2945: PetscBinaryRead(fd, &m->numBdFaces, 1, PETSC_INT);
2946: PetscBinaryRead(fd, &m->numCells, 1, PETSC_INT);
2947: PetscBinaryRead(fd, &m->numCorners, 1, PETSC_INT);
2948: PetscBinaryRead(fd, &m->numCellCorners, 1, PETSC_INT);
2949: PetscBinaryRead(fd, &m->numHoles, 1, PETSC_INT);
2951: PartitionSerialize(m, &m->part, viewer, store);
2952: PetscLogObjectParent(m, m->part);
2954: PartitionGetNumOverlapElements(m->part, &numFaces);
2955: PartitionGetNumOverlapNodes(m->part, &numNodes);
2956: numEdges = m->numEdges;
2957: numCorners = m->numCorners;
2958: numBd = m->numBd;
2959: numBdNodes = m->numBdNodes;
2960: numBdEdges = m->numBdEdges;
2961: PetscMalloc(numNodes*2 * sizeof(double), &tri->nodes);
2962: PetscMalloc(numNodes * sizeof(int), &tri->markers);
2963: PetscMalloc(numNodes * sizeof(int), &tri->degrees);
2964: PetscMalloc(numEdges*2 * sizeof(int), &tri->edges);
2965: PetscMalloc(numEdges * sizeof(int), &tri->edgemarkers);
2966: PetscMalloc(numFaces*numCorners * sizeof(int), &tri->faces);
2967: PetscMalloc(numFaces*3 * sizeof(int), &tri->neighbors);
2968: PetscMalloc(numBdNodes * sizeof(int), &tri->bdNodes);
2969: PetscMalloc(numBdEdges * sizeof(int), &tri->bdEdges);
2970: PetscMalloc(numBd * sizeof(int), &tri->bdMarkers);
2971: PetscMalloc((numBd+1) * sizeof(int), &tri->bdBegin);
2972: PetscMalloc((numBd+1) * sizeof(int), &tri->bdEdgeBegin);
2973: if (m->numHoles > 0) {
2974: PetscMalloc(m->numHoles*2 * sizeof(double), &m->holes);
2975: } else {
2976: m->holes = PETSC_NULL;
2977: }
2978: PetscLogObjectMemory(m, (numNodes*2 + m->numHoles*2) * sizeof(double) +
2979: (numNodes*2 + numEdges*3 + numFaces*(numCorners+1) + numBdNodes + numBdEdges + numBd*3+2)*sizeof(int));
2980: PetscBinaryRead(fd, tri->nodes, numNodes*2, PETSC_DOUBLE);
2981: PetscBinaryRead(fd, &old, 1, PETSC_INT);
2982: if (old) {
2983: PetscMalloc(numNodes*2 * sizeof(double), &tri->nodesOld);
2984: PetscLogObjectMemory(m, numNodes*2 * sizeof(double));
2985: PetscBinaryRead(fd, tri->nodesOld, numNodes*2, PETSC_DOUBLE);
2986: } else {
2987: tri->nodesOld = PETSC_NULL;
2988: }
2989: PetscBinaryRead(fd, tri->markers, numNodes, PETSC_INT);
2990: PetscBinaryRead(fd, tri->degrees, numNodes, PETSC_INT);
2991: PetscBinaryRead(fd, tri->edges, numEdges*2, PETSC_INT);
2992: PetscBinaryRead(fd, tri->edgemarkers, numEdges, PETSC_INT);
2993: PetscBinaryRead(fd, tri->faces, numFaces*numCorners, PETSC_INT);
2994: PetscBinaryRead(fd, tri->neighbors, numFaces*3, PETSC_INT);
2995: PetscBinaryRead(fd, tri->bdNodes, numBdNodes, PETSC_INT);
2996: PetscBinaryRead(fd, tri->bdEdges, numBdEdges, PETSC_INT);
2997: PetscBinaryRead(fd, tri->bdMarkers, numBd, PETSC_INT);
2998: PetscBinaryRead(fd, tri->bdBegin, numBd+1, PETSC_INT);
2999: PetscBinaryRead(fd, tri->bdEdgeBegin, numBd+1, PETSC_INT);
3000: PetscBinaryRead(fd, m->holes, m->numHoles*m->dim, PETSC_DOUBLE);
3002: PetscBinaryRead(fd, &m->maxAspectRatio, 1, PETSC_DOUBLE);
3004: #ifdef PETSC_USE_BOPT_g
3005: /* Check mesh integrity */
3006: MeshDebug_Triangular_2D(m, PETSC_TRUE);
3007: #endif
3008: *mesh = m;
3009: }
3011: return(0);
3012: }
3013: EXTERN_C_END
3015: EXTERN_C_BEGIN
3016: int MeshCreate_Triangular_2D(Mesh mesh) {
3017: Mesh_Triangular *tri;
3018: int (*f)(MeshBoundary2D *, int, Mesh);
3019: MPI_Comm comm;
3020: PetscTruth opt;
3021: int ierr;
3024: PetscNew(Mesh_Triangular, &tri);
3025: PetscMemcpy(mesh->ops, &MeshOps, sizeof(struct _MeshOps));
3026: mesh->data = (void *) tri;
3027: PetscStrallocpy(MESH_SER_TRIANGULAR_2D_BINARY, &mesh->serialize_name);
3028: mesh->dim = 2;
3030: /* Setup the periodic domain */
3031: if (mesh->isPeriodic == PETSC_TRUE) {
3032: MeshBoundary2DSetBoundingBox(mesh, mesh->bdCtx);
3033: }
3035: /* Setup default mechanisms */
3036: PetscObjectQueryFunction((PetscObject) mesh, "MeshTriangular2D_Create_C", (PetscVoidFunction) &f);
3037: if (f == PETSC_NULL) {
3038: #ifdef PETSC_HAVE_TRIANGLE
3039: if (mesh->isPeriodic == PETSC_TRUE) {
3040: PetscObjectComposeFunction((PetscObject) mesh, "MeshTriangular2D_Create_C", "MeshCreate_Periodic",
3041: (void (*)(void)) MeshCreate_Periodic);
3042:
3043: PetscObjectComposeFunction((PetscObject) mesh, "MeshTriangular2D_Refine_C", "MeshRefine_Periodic",
3044: (void (*)(void)) MeshRefine_Periodic);
3045:
3046: } else {
3047: PetscObjectComposeFunction((PetscObject) mesh, "MeshTriangular2D_Create_C", "MeshCreate_Triangle",
3048: (void (*)(void)) MeshCreate_Triangle);
3049:
3050: PetscObjectComposeFunction((PetscObject) mesh, "MeshTriangular2D_Refine_C", "MeshRefine_Triangle",
3051: (void (*)(void)) MeshRefine_Triangle);
3052: }
3053: #endif
3054: }
3056: MeshCheckBoundary_Triangular_2D(mesh);
3058: PetscObjectQueryFunction((PetscObject) mesh, "MeshTriangular2D_Create_C", (PetscVoidFunction) &f);
3059: (*f)(mesh->bdCtx, mesh->numCorners, mesh);
3061: /* Calculate maximum degree of vertices */
3062: MeshSetupSupport_Triangular_2D(mesh);
3064: /* Construct derived and boundary information */
3065: MeshSetupBoundary_Triangular_2D(mesh);
3067: /* Reorder nodes before distributing mesh */
3068: PetscOptionsHasName(mesh->prefix, "-mesh_reorder", &opt);
3069: if (opt == PETSC_TRUE) {
3070: /* MUST FIX: Since we have duplicated the whole mesh, we may impersonate a serial mesh */
3071: comm = mesh->comm;
3072: mesh->comm = PETSC_COMM_SELF;
3073: MeshGetOrdering(mesh, MESH_ORDER_TRIANGULAR_2D_RCM, &mesh->nodeOrdering);
3074: mesh->comm = comm;
3076: /* Permute arrays implicitly numbered by node numbers */
3077: AOApplicationToPetscPermuteReal(mesh->nodeOrdering, mesh->dim, tri->nodes);
3078: AOApplicationToPetscPermuteInt(mesh->nodeOrdering, 1, tri->markers);
3079: AOApplicationToPetscPermuteInt(mesh->nodeOrdering, 1, tri->degrees);
3080: /* Renumber arrays dependent on the canonical node numbering */
3081: AOApplicationToPetsc(mesh->nodeOrdering, mesh->numEdges*2, tri->edges);
3082: AOApplicationToPetsc(mesh->nodeOrdering, mesh->numFaces*mesh->numCorners, tri->faces);
3083: AOApplicationToPetsc(mesh->nodeOrdering, mesh->numBdNodes, tri->bdNodes);
3084: }
3086: /* Initialize partition */
3087: MeshPartition(mesh);
3089: /* Initialize save space */
3090: tri->nodesOld = PETSC_NULL;
3092: /* Reset the midnodes */
3093: if (mesh->isPeriodic == PETSC_TRUE) {
3094: MeshResetNodes(mesh, PETSC_TRUE);
3095: }
3097: return(0);
3098: }
3099: EXTERN_C_END