Actual source code: tri2dView.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: tri2dView.c,v 1.4 2000/10/17 13:48:56 knepley Exp $";
3: #endif
5: /* Viewers for 2d triangular grids */
6: #include "src/mesh/impls/triangular/2d/2dimpl.h" /*I "mesh.h" I*/
7: #include "src/sys/src/viewer/impls/silo/vsilo.h" /* For Silo viewer */
8: #include "tri2dView.h"
10: #undef __FUNCT__
12: static int MeshView_Triangular_2D_File(Mesh mesh, PetscViewer viewer)
13: {
14: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
15: Partition p = mesh->part;
16: FILE *fd;
17: int i, j;
18: int ierr;
21: PetscViewerASCIIPrintf(viewer, "Mesh Object:n");
22: if (mesh->isPeriodic == PETSC_FALSE) {
23: PetscViewerASCIIPrintf(viewer, " Triangular 2D grid withn");
24: } else {
25: PetscViewerASCIIPrintf(viewer, " Periodic Triangular 2D grid withn");
26: }
27: if (mesh->numBd == 1) {
28: PetscViewerASCIIPrintf(viewer, " %d closed boundaryn", mesh->numBd);
29: } else {
30: PetscViewerASCIIPrintf(viewer, " %d closed boundariesn", mesh->numBd);
31: }
32: for(i = 0; i < mesh->numBd; i++) {
33: PetscViewerASCIIPrintf(viewer, " Boundary %d: %d nodesn", tri->bdMarkers[i], tri->bdBegin[i+1] - tri->bdBegin[i]);
34: }
35: PetscViewerASCIIPrintf(viewer, " Total boundary nodes: %d edges: %dn", mesh->numBdNodes, mesh->numBdEdges);
36: PetscViewerASCIIGetPointer(viewer, &fd);
37: PetscSynchronizedFPrintf(mesh->comm, fd, " Local graph %d: %d nodes %d bdEdgesn", p->rank, mesh->numNodes, mesh->numBdEdges);
38: for(i = 0; i < mesh->numNodes; i++) {
39: PetscSynchronizedFPrintf(mesh->comm, fd, " %d %g %g %dn", i, tri->nodes[i*2], tri->nodes[i*2+1], tri->markers[i]);
40: }
41: PetscSynchronizedFlush(mesh->comm);
42: PetscSynchronizedFPrintf(mesh->comm, fd, " Local graph %d: %d edgesn", p->rank, mesh->numEdges);
43: for(i = 0; i < mesh->numEdges; i++) {
44: PetscSynchronizedFPrintf(mesh->comm, fd, " %d %d %d %dn", i, tri->edges[i*2], tri->edges[i*2+1], tri->edgemarkers[i]);
45: }
46: PetscSynchronizedFlush(mesh->comm);
47: PetscSynchronizedFPrintf(mesh->comm, fd, " Local graph %d: %d faces with %d nodes per facen",
48: p->rank, mesh->numFaces, mesh->numCorners);
49: for(i = 0; i < mesh->numFaces; i++) {
50: PetscSynchronizedFPrintf(mesh->comm, fd, " %d", i);
51: for(j = 0; j < mesh->numCorners; j++) PetscSynchronizedFPrintf(mesh->comm, fd, " %d", tri->faces[i*mesh->numCorners+j]);
52: PetscSynchronizedFPrintf(mesh->comm, fd, "n");
53: }
54: for(i = 0; i < mesh->numFaces; i++) {
55: PetscSynchronizedFPrintf(mesh->comm, fd, " %d %d %d %dn", i, tri->neighbors[i*3], tri->neighbors[i*3+1],
56: tri->neighbors[i*3+2]);
57: }
58: PetscSynchronizedFlush(mesh->comm);
59:
60: if (tri->areas != PETSC_NULL) {
61: PetscSynchronizedFPrintf(mesh->comm, fd, " Local graph %d: element areasn", p->rank);
62: for(i = 0; i < mesh->numFaces; i++)
63: PetscSynchronizedFPrintf(mesh->comm, fd, " %d %gn", i, tri->areas[i]);
64: PetscSynchronizedFlush(mesh->comm);
65: }
66: if (tri->aspectRatios != PETSC_NULL) {
67: PetscSynchronizedFPrintf(mesh->comm, fd, " Local graph %d: element aspect ratiosn", p->rank);
68: for(i = 0; i < mesh->numFaces; i++)
69: PetscSynchronizedFPrintf(mesh->comm, fd, " %d %gn", i, tri->aspectRatios[i]);
70: PetscSynchronizedFlush(mesh->comm);
71: }
73: if (mesh->partitioned) {
74: PartitionView(mesh->part, viewer);
75: }
76: return(0);
77: }
79: #undef __FUNCT__
81: static int MeshViewLocalizeMesh_Private(Mesh mesh, double **nodes, int **faces) {
82: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
83: Partition_Triangular_2D *q;
84: Partition part;
85: MPI_Comm comm;
86: int *numRecvNodes;
87: int *numRecvFaces;
88: int *nodeOffsets;
89: int *faceOffsets;
90: int *locFaces;
91: double *locNodes;
92: int dim, numLocNodes, numNodes, numLocFaces, numFaces, numCorners;
93: int numProcs, rank;
94: int proc, node;
95: int ierr;
98: PetscObjectGetComm((PetscObject) mesh, &comm);
99: MPI_Comm_size(comm, &numProcs);
100: MPI_Comm_rank(comm, &rank);
101: MeshGetPartition(mesh, &part);
102: MeshGetDimension(mesh, &dim);
103: MeshGetNumCorners(mesh, &numCorners);
104: PartitionGetNumNodes(part, &numLocNodes);
105: PartitionGetTotalNodes(part, &numNodes);
106: PartitionGetNumElements(part, &numLocFaces);
107: PartitionGetTotalElements(part, &numFaces);
108: if (numProcs > 1) {
109: q = (Partition_Triangular_2D *) part->data;
110: /* Allocate global arrays */
111: PetscMalloc(numNodes*dim * sizeof(double), nodes);
112: PetscMalloc(numFaces*numCorners * sizeof(int), faces);
113: PetscMalloc(numLocFaces*numCorners * sizeof(int), &locFaces);
114: locNodes = tri->nodes;
116: /* Calculate offsets */
117: PetscMalloc(numProcs * sizeof(int), &numRecvNodes);
118: PetscMalloc(numProcs * sizeof(int), &numRecvFaces);
119: PetscMalloc(numProcs * sizeof(int), &nodeOffsets);
120: PetscMalloc(numProcs * sizeof(int), &faceOffsets);
121: for(proc = 0; proc < numProcs; proc++) {
122: numRecvNodes[proc] = (q->firstNode[proc+1] - q->firstNode[proc])*dim;
123: numRecvFaces[proc] = (part->firstElement[proc+1] - part->firstElement[proc])*numCorners;
124: }
125: nodeOffsets[0] = 0;
126: faceOffsets[0] = 0;
127: for(proc = 1; proc < numProcs; proc++) {
128: nodeOffsets[proc] = numRecvNodes[proc-1] + nodeOffsets[proc-1];
129: faceOffsets[proc] = numRecvFaces[proc-1] + faceOffsets[proc-1];
130: }
132: /* Local to global node number conversion */
133: for(node = 0; node < numLocFaces*numCorners; node++) {
134: PartitionLocalToGlobalNodeIndex(part, tri->faces[node], &locFaces[node]);
135: }
137: /* Collect global arrays */
138: MPI_Gatherv(locNodes, numLocNodes*dim, MPI_DOUBLE, *nodes, numRecvNodes, nodeOffsets, MPI_DOUBLE, 0, comm);
139: MPI_Gatherv(locFaces, numLocFaces*numCorners, MPI_INT, *faces, numRecvFaces, faceOffsets, MPI_INT, 0, comm);
141: /* Cleanup */
142: PetscFree(locFaces);
143: PetscFree(numRecvNodes);
144: PetscFree(numRecvFaces);
145: PetscFree(nodeOffsets);
146: PetscFree(faceOffsets);
148: if (rank == 0) {
149: /* We must globally renumber and permute so that midnodes come after vertices */
150: AOPetscToApplication(q->nodeOrdering, numFaces*numCorners, *faces);
151: AOPetscToApplicationPermuteReal(q->nodeOrdering, dim, *nodes);
152: if (mesh->nodeOrdering) {
153: AOPetscToApplication(mesh->nodeOrdering, numFaces*numCorners, *faces);
154: AOPetscToApplicationPermuteReal(mesh->nodeOrdering, dim, *nodes);
155: }
156: }
157: } else {
158: *nodes = tri->nodes;
159: *faces = tri->faces;
160: }
161: return(0);
162: }
164: #undef __FUNCT__
166: static int MeshViewLocalizeDestroy_Private(Mesh mesh, double *nodes, int *faces) {
167: MPI_Comm comm;
168: int numProcs;
169: int ierr;
172: PetscObjectGetComm((PetscObject) mesh, &comm);
173: MPI_Comm_size(comm, &numProcs);
174: if (numProcs > 1) {
175: PetscFree(nodes);
176: PetscFree(faces);
177: }
178: return(0);
179: }
181: #undef __FUNCT__
183: static int MeshView_Triangular_2D_VU(Mesh mesh, PetscViewer viewer) {
184: MPI_Comm comm;
185: Partition part;
186: FILE *fp;
187: double *nodes;
188: int *faces;
189: int dim, numNodes, numFaces, numCorners;
190: int d, node, elem;
191: int ierr;
194: PetscObjectGetComm((PetscObject) mesh, &comm);
195: MeshGetPartition(mesh, &part);
196: MeshGetDimension(mesh, &dim);
197: MeshGetNumCorners(mesh, &numCorners);
198: PartitionGetTotalNodes(part, &numNodes);
199: PartitionGetTotalElements(part, &numFaces);
200: PetscViewerVUGetPointer(viewer, &fp);
201: MeshViewLocalizeMesh_Private(mesh, &nodes, &faces);
202: /* Write the node coordinates */
203: PetscFPrintf(comm, fp, "// Node coordinatesnFIELD Coordinates() = n{n");
204: for(node = 0; node < numNodes; node++) {
205: PetscFPrintf(comm, fp, " ");
206: for(d = 0; d < dim-1; d++) {
207: PetscFPrintf(comm, fp, "%g ", nodes[node*dim+d]);
208: }
209: PetscFPrintf(comm, fp, "%gn", nodes[node*dim+(dim-1)]);
210: }
211: PetscFPrintf(comm, fp, "};nn");
212: /* Write the node connectivity */
213: if (numCorners == 6) {
214: PetscFPrintf(comm, fp, "// Full mesh connectivitynFIELD<int> FullConnectivity() =n{n");
215: for (elem = 0; elem < numFaces; elem++) {
216: PetscFPrintf(comm, fp, " ");
217: for (node = 0; node < numCorners-1; node++) {
218: PetscFPrintf(comm, fp, " %d ", faces[elem*numCorners+node]+1);
219: }
220: PetscFPrintf(comm, fp, "%dn", faces[elem*numCorners+numCorners-1]+1);
221: }
222: PetscFPrintf(comm, fp, "};nn");
223: }
224: PetscFPrintf(comm, fp, "// Vertex mesh connectivitynFIELD<int> Connectivity() =n{n");
225: for (elem = 0; elem < numFaces; elem++) {
226: PetscFPrintf(comm, fp, " %d %d %dn", faces[elem*numCorners]+1, faces[elem*numCorners+1]+1,
227: faces[elem*numCorners+2]+1);
228: }
229: PetscFPrintf(comm, fp, "};nn");
230: /* Define mesh itself */
231: if (numCorners == 6) {
232: PetscFPrintf(comm, fp, "MESH PetscMesh() =n{n ZONE Zone1(LagrTrian06, Coordinates, FullConnectivity);n};nn");
233: } else {
234: PetscFPrintf(comm, fp, "MESH PetscMesh() =n{n ZONE Zone1(LagrTrian03, Coordinates, Connectivity);n};nn");
235: }
237: MeshViewLocalizeDestroy_Private(mesh, nodes, faces);
238: return(0);
239: }
241: #if 0
242: #undef __FUNCT__
244: static int MeshView_Triangular_2D_Poum(Mesh mesh, PetscViewer viewer)
245: {
246: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
247: Partition p = mesh->part;
248: Partition_Triangular_2D *q = (Partition_Triangular_2D *) p->data;
249: int numCorners = tri->numCorners;
250: int *faces;
251: double *nodes;
252: FILE *fp;
253: int proc, elem, node;
254: int ierr;
257: MeshViewLocalizeMesh_Private(mesh, &nodes, &faces);
258: ViewerPoumGetMeshPointer(viewer, &fp);
259: /* Writes the nodes description */
260: if (p->rank == 0) {
261: PetscFPrintf(PETSC_COMM_SELF, fp, "Nodes Nmattn");
262: for (node = 0; node < tri->numVertices; node++)
263: PetscFPrintf(PETSC_COMM_SELF, fp, " %dt %ft %f 0.0n", node + 1, nodes[node*2], nodes[node*2+1]);
264: }
266: /* Writes the element description */
267: if (p->rank == 0) {
268: PetscFPrintf(PETSC_COMM_SELF, fp, "Elements Ematt using Nmattn");
269: for (elem = 0; elem < p->numElements; elem++)
270: PetscFPrintf(PETSC_COMM_SELF, fp, " %dt 4t %dt %dt %dn", elem + 1, faces[elem*numCorners]+1,
271: faces[elem*numCorners+1]+1, faces[elem*numCorners+2]+1);
272: }
274: ViewerPoumGetPartitionPointer(viewer, &fp);
275: if (p->rank == 0) {
276: /* Print the partition */
277: PetscFPrintf(PETSC_COMM_SELF, fp, "Decomposition Dmatt using Emattn %dn", p->numProcs);
278: for(proc = 0; proc< p->numProcs; proc++) {
279: PetscFPrintf(PETSC_COMM_SELF, fp, "%dn", p->firstElement[proc+1] - p->firstElement[proc]);
280: for (elem = 0; elem < p->firstElement[proc+1] - p->firstElement[proc]; elem++)
281: PetscFPrintf(PETSC_COMM_SELF, fp, " %dn", elem + p->firstElement[proc] + 1);
282: }
283: }
284: MeshViewLocalizeDestroy_Private(mesh, nodes, faces);
286: return(0);
287: }
288: #endif
290: #ifdef PETSC_HAVE_SILO
291: #undef __FUNCT__
293: static int ViewerSiloNodeValues_Private(PetscViewer viewer, Mesh mesh, Vec vec, char *fieldName)
294: {
295: Viewer_Silo *vsilo = (Viewer_Silo *) viewer->data;
296: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
297: int dim = mesh->dim;
298: int numNodes = tri->numNodes;
299: DBfile *fp;
300: DBoptlist *optlist;
301: Scalar *array;
302: float **newArray;
303: char **subNames;
304: char name[256];
305: char buf[1024];
306: int size, nameLen, len;
307: int node, c;
308: int ierr;
311: ViewerSiloCheckMesh(viewer, mesh);
312: VecGetSize(vec, &size);
313: ViewerSiloGetFilePointer(viewer, &fp);
314: VecGetArray(vec, &array);
315: /* Name vector */
316: if (vsilo->objName != PETSC_NULL) {
317: PetscStrncpy(name, vsilo->objName, 240);
318: } else {
319: PetscStrcpy(name, "");
320: }
321: PetscStrcat(name, fieldName);
322: /* Allocate storage compatible with SILO calls */
323: PetscMalloc(dim * sizeof(char *), &subNames);
324: PetscMalloc(dim * sizeof(float *), &newArray);
325: PetscStrlen(name, &nameLen);
326: nameLen += (int) log10((double) dim) + 6;
327: for(c = 0; c < dim; c++) {
328: PetscMalloc(nameLen * sizeof(char), &subNames[c]);
329: PetscMalloc(numNodes * sizeof(float), &newArray[c]);
330: sprintf(subNames[c], "%scomp%d", name, c);
331: }
332: /* Convert vector to SILO format */
333: for(node = 0; node < numNodes; node++) {
334: for(c = 0; c < dim; c++) {
335: newArray[c][node] = array[node*dim+c];
336: }
337: }
338: /* Add each component */
339: for(c = 0; c < dim; c++) {
340: optlist = DBMakeOptlist(3);
341: DBAddOption(optlist, DBOPT_LABEL, name);
342: if (vsilo->meshName == PETSC_NULL) {
343: DBPutUcdvar1(fp, subNames[c], "PetscMesh", newArray[c], numNodes, PETSC_NULL, 0, DB_FLOAT, DB_NODECENT, optlist);
344:
345: } else {
346: DBPutUcdvar1(fp, subNames[c], vsilo->meshName, newArray[c], numNodes, PETSC_NULL, 0, DB_FLOAT, DB_NODECENT, optlist);
347:
348: }
349: DBFreeOptlist(optlist);
350: }
352: if (dim > 1) {
353: PetscMemzero(buf, 1024 * sizeof(char));
354: len = DBGetVarLength(fp, "_meshtv_defvars");
355: if (len > 0) {
356: if (DBGetVarType(fp, "_meshtv_defvars") != DB_CHAR) {
357: SETERRQ(PETSC_ERR_FILE_READ, "Invalid type for variable _meshtv_defvars");
358: }
359: if (len > 1024) SETERRQ(PETSC_ERR_SUP, "Need to do dyanmic allocation here");
360: DBReadVar(fp, "_meshtv_defvars", buf);
361: PetscStrcat(buf, ";vec");
362: } else {
363: PetscStrcpy(buf, "vec");
364: }
365: PetscStrcat(buf, name);
366: PetscStrcat(buf, " vector {");
367: for(c = 0; c < dim-1; c++) {
368: PetscStrcat(buf, subNames[c]);
369: PetscStrcat(buf, ",");
370: }
371: PetscStrcat(buf, subNames[c]);
372: PetscStrcat(buf, "}");
373: PetscStrlen(buf, &len);
374: if (len > 1024) SETERRQ(PETSC_ERR_SUP, "Need to do dyanmic allocation here");
375: len = 1024;
376: DBWrite(fp, "_meshtv_defvars", buf, &len, 1, DB_CHAR);
377: }
379: for(c = 0; c < dim; c++) {
380: PetscFree(subNames[c]);
381: PetscFree(newArray[c]);
382: }
383: PetscFree(subNames);
384: PetscFree(newArray);
385: VecRestoreArray(vec, &array);
386: return(0);
387: }
388: #endif
390: #undef __FUNCT__
392: static int MeshView_Triangular_2D_Silo(Mesh mesh, PetscViewer viewer)
393: {
394: #ifdef HAVE_SILO
395: Viewer_Silo *vsilo = (Viewer_Silo *) viewer->data;
396: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
397: int numCorners = tri->numCorners;
398: int numElements = tri->numFaces;
399: int *faces = tri->faces;
400: int *edges = tri->edges;
401: int numNodes = tri->numNodes;
402: double *nodes = tri->nodes;
403: int numBdEdges = tri->numBdEdges;
404: int *bdEdges = tri->bdEdges;
405: double sizeX = mesh->sizeX;
406: double sizeY = mesh->sizeY;
407: PetscTruth isPeriodic = mesh->isPeriodic;
408: PetscTruth xPer = mesh->isPeriodicDim[0];
409: PetscTruth yPer = mesh->isPeriodicDim[1];
410: int pointsPerEdge = 2; /* The number of vertices on an edge - should always be 2 for 2D mesh */
411: int pointsPerFace = 3; /* The number of vertices on a face - should always be 3 for triangular mesh */
412: int numShapes = 1; /* The number of different shapes - we only have triangles */
413: char edgeName[256]; /* The name of the edge list */
414: char zoneName[256]; /* The name of the zone list */
415: char meshName[256]; /* The name of the mesh */
416: DBfile *fp;
417: DBoptlist *optlist;
418: int *facelist, *zonelist;
419: float *xcoords, *ycoords;
420: float *coords[2];
421: int numFaces, numZones;
422: int node, face, zone;
423: int ierr;
424:
425:
427: /* Setup base names */
428: if (vsilo->objName != PETSC_NULL) {
429: PetscStrncpy(edgeName, vsilo->objName, 251);
430: PetscStrncpy(zoneName, vsilo->objName, 251);
431: PetscStrncpy(meshName, vsilo->objName, 251);
432: } else {
433: PetscStrcpy(edgeName, "Petsc");
434: PetscStrcpy(zoneName, "Petsc");
435: PetscStrcpy(meshName, "Petsc");
436: }
437: /* Add suffixes */
438: PetscStrcat(edgeName, "Face");
439: PetscStrcat(zoneName, "Zone");
440: PetscStrcat(meshName, "Mesh");
442: /* Allocate space for the new arrays that have to be created */
443: PetscMalloc(numNodes * sizeof(float), &xcoords);
444: PetscMalloc(numNodes * sizeof(float), &ycoords);
445: PetscMalloc(numBdEdges*pointsPerEdge * sizeof(int), &facelist);
446: PetscMalloc(numElements*pointsPerFace * sizeof(int), &zonelist);
448: if (isPeriodic == PETSC_TRUE) {
449: for(face = 0, numFaces = 0; face < numBdEdges; face++) {
450: if (((xPer == PETSC_TRUE) && (PetscAbsScalar(nodes[edges[bdEdges[face]*2]*2] - nodes[edges[bdEdges[face]*2+1]*2]) > 0.5*sizeX)) ||
451: ((yPer == PETSC_TRUE) && (PetscAbsScalar(nodes[edges[bdEdges[face]*2]*2+1] - nodes[edges[bdEdges[face]*2+1]*2+1]) > 0.5*sizeY)))
452: continue;
453: facelist[numFaces*2] = edges[bdEdges[face]*2];
454: facelist[numFaces*2+1] = edges[bdEdges[face]*2+1];
455: numFaces++;
456: }
457: /* DBPutZoneList only wants the points on the corners of the triangle, so we remove the other nodes */
458: for(zone = 0, numZones = 0; zone < numElements; zone++) {
459: if (((xPer == PETSC_TRUE) && (PetscAbsScalar(nodes[faces[zone*numCorners]*2] - nodes[faces[zone*numCorners+1]*2]) > 0.5*sizeX)) ||
460: ((yPer == PETSC_TRUE) && (PetscAbsScalar(nodes[faces[zone*numCorners]*2+1] - nodes[faces[zone*numCorners+1]*2+1]) > 0.5*sizeY)))
461: continue;
462: if (((xPer == PETSC_TRUE) && (PetscAbsScalar(nodes[faces[zone*numCorners]*2] - nodes[faces[zone*numCorners+2]*2]) > 0.5*sizeX)) ||
463: ((yPer == PETSC_TRUE) && (PetscAbsScalar(nodes[faces[zone*numCorners]*2+1] - nodes[faces[zone*numCorners+2]*2+1]) > 0.5*sizeY)))
464: continue;
465: zonelist[numZones*3] = faces[zone*numCorners];
466: zonelist[numZones*3+1] = faces[zone*numCorners+1];
467: zonelist[numZones*3+2] = faces[zone*numCorners+2];
468: numZones++;
469: }
470: } else {
471: numFaces = numBdEdges;
472: for(face = 0; face < numBdEdges; face++) {
473: facelist[face*2] = edges[bdEdges[face]*2];
474: facelist[face*2+1] = edges[bdEdges[face]*2+1];
475: }
476: /* DBPutZoneList only wants the points on the corners of the triangle, so we remove the other nodes */
477: numZones = numElements;
478: for(zone = 0; zone < numElements; zone++) {
479: zonelist[zone*3] = faces[zone*numCorners];
480: zonelist[zone*3+1] = faces[zone*numCorners+1];
481: zonelist[zone*3+2] = faces[zone*numCorners+2];
482: }
483: }
484:
485: /* DBPutUcdMesh expects x- and y- coordinates to be in separate arrays, so we split them up */
486: for(node = 0; node < numNodes; node++) {
487: xcoords[node] = nodes[node*2];
488: ycoords[node] = nodes[node*2+1];
489: }
490: coords[0] = xcoords;
491: coords[1] = ycoords;
492:
493: ViewerSiloGetFilePointer(viewer, &fp);
494: DBPutFacelist(fp, edgeName, numFaces, 2, facelist, numFaces*pointsPerEdge, 0, PETSC_NULL,
495: &pointsPerEdge, &numFaces, numShapes, PETSC_NULL, PETSC_NULL, 0);
496:
497: DBPutZonelist(fp, zoneName, numZones, 2, zonelist, numZones*pointsPerFace, 0, &pointsPerFace, &numZones, 1);
498:
499: PetscFree(facelist);
500: PetscFree(zonelist);
501:
502: /* Assorted options to make the graph prettier */
503: optlist = DBMakeOptlist(5);
504: /* DBAddOption(optlist, DBOPT_COORDSYS, DB_CARTESIAN); */
505: DBAddOption(optlist, DBOPT_XLABEL, "X") ;
506: DBAddOption(optlist, DBOPT_YLABEL, "Y");
507: DBAddOption(optlist, DBOPT_XUNITS, "cm");
508: DBAddOption(optlist, DBOPT_YUNITS, "cm");
509: DBPutUcdmesh(fp, meshName, 2, PETSC_NULL, coords, numNodes, numZones, zoneName, edgeName, DB_FLOAT, optlist);
510:
511: DBFreeOptlist(optlist);
513: /* Check for moving mesh and visualize velocity/acceleration */
514: if (mesh->movingMesh == PETSC_TRUE) {
515: if (mesh->nodeVelocities != PETSC_NULL) {
516: ViewerSiloSetName(viewer, "Mesh");
517: ViewerSiloNodeValues_Private(viewer, mesh, mesh->nodeVelocities, "Velocity");
518: }
519: if (mesh->nodeAccelerations != PETSC_NULL) {
520: ViewerSiloSetName(viewer, "Mesh");
521: ViewerSiloNodeValues_Private(viewer, mesh, mesh->nodeAccelerations, "Acceleration");
522: }
523: ViewerSiloClearName(viewer);
524: }
526: #ifdef NOT_WORKING
527: /* Set mesh name */
528: if (vsilo->objName != PETSC_NULL) {
529: PetscMalloc((PetscStrlen(meshName)+1) * sizeof(char), &newMeshName);
530: PetscStrcpy(newMeshName, meshName);
531: ViewerSiloSetMeshName(viewer, newMeshName);
532: }
533: #endif
534: return(0);
535: #else
536: SETERRQ(PETSC_ERR_SUP, "Need LLNL's Silo package");
537: #endif
538: }
540: #undef __FUNCT__
542: static int MeshView_Triangular_2D_Draw_Edge(Mesh mesh, PetscDraw draw, int edge, int color) {
543: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
544: int numCorners = mesh->numCorners;
545: int *edges = tri->edges;
546: int *elements = tri->faces;
547: double *nodes = tri->nodes;
548: int numOverlapElements;
549: int elem, corner, corner2, midCorner, midNode;
550: int ierr;
553: if (numCorners == 3) {
554: MeshDrawLine(mesh, draw, nodes[edges[edge*2]*2], nodes[edges[edge*2]*2+1], nodes[edges[edge*2+1]*2],
555: nodes[edges[edge*2+1]*2+1], color);
556:
557: } else if (numCorners == 6) {
558: /* Find midnode */
559: PartitionGetNumOverlapElements(mesh->part, &numOverlapElements);
560: for(elem = 0, midNode = -1; elem < numOverlapElements; elem++) {
561: for(corner = 0; corner < 3; corner++) {
562: if (elements[elem*numCorners+corner] == edges[edge*2]) break;
563: }
564: if (corner < 3) {
565: for(corner2 = 0; corner2 < 3; corner2++) {
566: if (elements[elem*numCorners+corner2] == edges[edge*2+1]) break;
567: }
568: if (corner2 < 3) {
569: midCorner = ((corner+corner2)*2)%3 + 3;
570: midNode = elements[elem*numCorners+midCorner];
571: break;
572: }
573: }
574: }
575: if (midNode < 0) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid mesh connectivity");
576: MeshDrawLine(mesh, draw, nodes[edges[edge*2]*2], nodes[edges[edge*2]*2+1], nodes[midNode*2],
577: nodes[midNode*2+1], color);
578:
579: MeshDrawLine(mesh, draw, nodes[midNode*2], nodes[midNode*2+1], nodes[edges[edge*2+1]*2],
580: nodes[edges[edge*2+1]*2+1], color);
581:
582: } else {
583: SETERRQ1(PETSC_ERR_SUP, "Invalid number of corners %d", numCorners);
584: }
585: return(0);
586: }
588: #undef __FUNCT__
590: static int MeshView_Triangular_2D_Draw_Element(Mesh mesh, PetscDraw draw, int elem, int color, PetscTruth drawEdges) {
591: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
592: int numCorners = mesh->numCorners;
593: int *elements = tri->faces;
594: int *neighbors = tri->faces;
595: double *nodes = tri->nodes;
596: int corner;
597: int ierr;
600: if (numCorners == 3) {
601: MeshDrawTriangle(mesh, draw, nodes[elements[elem*numCorners]*2], nodes[elements[elem*numCorners]*2+1],
602: nodes[elements[elem*numCorners+1]*2], nodes[elements[elem*numCorners+1]*2+1],
603: nodes[elements[elem*numCorners+2]*2], nodes[elements[elem*numCorners+2]*2+1],
604: color, color, color);
605:
606: } else if (numCorners == 6) {
607: /* Draw 4 interior triangles */
608: MeshDrawTriangle(mesh, draw, nodes[elements[elem*numCorners]*2], nodes[elements[elem*numCorners]*2+1],
609: nodes[elements[elem*numCorners+5]*2], nodes[elements[elem*numCorners+5]*2+1],
610: nodes[elements[elem*numCorners+4]*2], nodes[elements[elem*numCorners+4]*2+1],
611: color, color, color);
612:
613: MeshDrawTriangle(mesh, draw, nodes[elements[elem*numCorners+1]*2], nodes[elements[elem*numCorners+1]*2+1],
614: nodes[elements[elem*numCorners+3]*2], nodes[elements[elem*numCorners+3]*2+1],
615: nodes[elements[elem*numCorners+5]*2], nodes[elements[elem*numCorners+5]*2+1],
616: color, color, color);
617:
618: MeshDrawTriangle(mesh, draw, nodes[elements[elem*numCorners+2]*2], nodes[elements[elem*numCorners+2]*2+1],
619: nodes[elements[elem*numCorners+4]*2], nodes[elements[elem*numCorners+4]*2+1],
620: nodes[elements[elem*numCorners+3]*2], nodes[elements[elem*numCorners+3]*2+1],
621: color, color, color);
622:
623: MeshDrawTriangle(mesh, draw, nodes[elements[elem*numCorners+3]*2], nodes[elements[elem*numCorners+3]*2+1],
624: nodes[elements[elem*numCorners+4]*2], nodes[elements[elem*numCorners+4]*2+1],
625: nodes[elements[elem*numCorners+5]*2], nodes[elements[elem*numCorners+5]*2+1],
626: color, color, color);
627:
628: } else {
629: SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid number of corners %d", numCorners);
630: }
632: if (drawEdges == PETSC_TRUE) {
633: for(corner = 0; corner < 3; corner++) {
634: if (neighbors[elem*3+corner] < 0) {
635: color = PETSC_DRAW_RED;
636: } else {
637: color = PETSC_DRAW_BLACK;
638: }
639: if (numCorners == 3) {
640: MeshDrawLine(mesh, draw, nodes[elements[elem*numCorners+corner]*2], nodes[elements[elem*numCorners+corner]*2+1],
641: nodes[elements[elem*numCorners+((corner+1)%3)]*2], nodes[elements[elem*numCorners+((corner+1)%3)]*2+1], color);
642:
643: } else if (numCorners == 6) {
644: MeshDrawLine(mesh, draw, nodes[elements[elem*numCorners+corner]*2], nodes[elements[elem*numCorners+corner]*2+1],
645: nodes[elements[elem*numCorners+corner+3]*2], nodes[elements[elem*numCorners+corner+3]*2+1], color);
646:
647: MeshDrawLine(mesh, draw, nodes[elements[elem*numCorners+corner+3]*2], nodes[elements[elem*numCorners+corner+3]*2+1],
648: nodes[elements[elem*numCorners+((corner+1)%3)]*2], nodes[elements[elem*numCorners+((corner+1)%3)]*2+1], color);
649:
650: }
651: }
652: }
653: return(0);
654: }
656: #undef __FUNCT__
658: static int MeshView_Triangular_2D_Draw_Mesh(Mesh mesh, PetscDraw draw) {
659: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
660: int *edgemarkers = tri->edgemarkers;
661: int numElements = mesh->numFaces;
662: int numEdges = mesh->numEdges;
663: int hElement = mesh->highlightElement;
664: int color, rank;
665: int elem, edge;
666: int ierr;
669: MPI_Comm_rank(mesh->comm, &rank);
670: for(elem = 0; elem < numElements; elem++) {
671: if (hElement == elem) {
672: color = PETSC_DRAW_WHITE;
673: } else {
674: color = PETSC_DRAW_RED + rank + 1;
675: }
676: MeshView_Triangular_2D_Draw_Element(mesh, draw, elem, color, PETSC_FALSE);
677: }
678: /* Needed so that triangles do not overwrite edges */
679: PetscDrawSynchronizedFlush(draw);
681: for(edge = 0; edge < numEdges; edge++) {
682: if (edgemarkers[edge]) {
683: color = PETSC_DRAW_RED;
684: } else {
685: color = PETSC_DRAW_BLACK;
686: }
687: MeshView_Triangular_2D_Draw_Edge(mesh, draw, edge, color);
688: }
689: PetscDrawSynchronizedFlush(draw);
690: return(0);
691: }
693: #undef __FUNCT__
695: /*
696: Remember that this function sets the coordinates for the window.
697: */
698: static int MeshView_DrawStatusLine_Private(Mesh mesh, PetscDraw draw, double startX, double startY, double endX,
699: double endY, double offX, double offY) {
700: Partition part = mesh->part;
701: MPI_Comm comm;
702: char statusLine[1024];
703: PetscReal textWidth, textHeight;
704: int hElement, rank;
705: int ierr;
708: PetscObjectGetComm((PetscObject) mesh, &comm);
709: MPI_Comm_rank(comm, &rank);
710: MeshGetHighlightElement(mesh, &hElement);
711: PartitionLocalToGlobalElementIndex(part, hElement, &hElement);
712: sprintf(statusLine, "Element %d Proc: %d", hElement, rank);
713: PetscDrawStringGetSize(draw, &textWidth, &textHeight);
714: PetscDrawSetCoordinates(draw, startX-offX, startY-offY-textHeight, endX+offX, endY+offY);
715: if (hElement >= 0) {
716: PetscDrawString(draw, startX-offX, startY-offY-textHeight, PETSC_DRAW_BLACK, statusLine);
717: }
718: return(0);
719: }
721: #undef __FUNCT__
723: static int MeshViewElement_Triangular_2D_Draw(Mesh mesh, int elem, PetscViewer v) {
724: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
725: int dim = mesh->dim;
726: int numCorners = mesh->numCorners;
727: int *elements = tri->faces;
728: double *nodes = tri->nodes;
729: double startX, endX;
730: double startY, endY;
731: double sizeX, sizeY;
732: double offX, offY;
733: PetscDraw draw;
734: PetscTruth isnull;
735: PetscDrawButton button;
736: int pause, color;
737: PetscReal xc, yc, scale;
738: int xsize, ysize;
739: int corner;
740: int ierr;
743: startX = endX = nodes[elements[elem*numCorners]*dim+0];
744: startY = endY = nodes[elements[elem*numCorners]*dim+1];
745: for(corner = 1; corner < numCorners; corner++) {
746: if (nodes[elements[elem*numCorners+corner]*dim+0] < startX) startX = nodes[elements[elem*numCorners+corner]*dim+0];
747: if (nodes[elements[elem*numCorners+corner]*dim+0] > endX) endX = nodes[elements[elem*numCorners+corner]*dim+0];
748: if (nodes[elements[elem*numCorners+corner]*dim+1] < startY) startY = nodes[elements[elem*numCorners+corner]*dim+1];
749: if (nodes[elements[elem*numCorners+corner]*dim+1] > endY) endY = nodes[elements[elem*numCorners+corner]*dim+1];
750: }
751: sizeX = endX - startX;
752: sizeY = endY - startY;
753: color = PETSC_DRAW_WHITE;
754: PetscViewerDrawGetDraw(v, 0, &draw);
755: PetscDrawIsNull(draw, &isnull);
756: if (isnull) return(0);
757: if (sizeX > sizeY) {
758: ysize = 300;
759: if (sizeY > 0.0) {
760: xsize = ysize * (int) (sizeX/sizeY);
761: } else {
762: return(0);
763: }
764: } else {
765: xsize = 300;
766: if (sizeX > 0.0) {
767: ysize = xsize * (int) (sizeY/sizeX);
768: } else {
769: return(0);
770: }
771: }
772: if ((xsize == 0) || (ysize == 0)) return(0);
773: offX = 0.02*sizeX;
774: offY = 0.02*sizeY;
775: PetscDrawResizeWindow(draw, xsize, ysize);
776: MeshView_DrawStatusLine_Private(mesh, draw, startX, startY, endX, endY, offX, offY);
777: MeshView_Triangular_2D_Draw_Element(mesh, draw, elem, color, PETSC_TRUE);
778: PetscDrawGetPause(draw, &pause);
779: if (pause >= 0) {
780: PetscSleep(pause);
781: return(0);
782: }
784: /* Allow the mesh to zoom or shrink */
785: PetscDrawCheckResizedWindow(draw);
786: PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
787: while (button != BUTTON_RIGHT) {
788: PetscDrawSynchronizedClear(draw);
789: switch (button)
790: {
791: case BUTTON_LEFT:
792: scale = 0.5;
793: break;
794: case BUTTON_CENTER:
795: scale = 2.0;
796: break;
797: default:
798: scale = 1.0;
799: break;
800: }
801: if (scale != 1.0) {
802: startX = scale*(startX + sizeX - xc) + xc - sizeX*scale;
803: endX = scale*(endX - sizeX - xc) + xc + sizeX*scale;
804: startY = scale*(startY + sizeY - yc) + yc - sizeY*scale;
805: endY = scale*(endY - sizeY - yc) + yc + sizeY*scale;
806: sizeX *= scale;
807: sizeY *= scale;
808: offX *= scale;
809: offY *= scale;
810: }
812: MeshView_DrawStatusLine_Private(mesh, draw, startX, startY, endX, endY, offX, offY);
813: MeshView_Triangular_2D_Draw_Element(mesh, draw, elem, color, PETSC_TRUE);
814: PetscDrawCheckResizedWindow(draw);
815: PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
816: }
818: return(0);
819: }
821: #undef __FUNCT__
823: static int MeshView_Triangular_2D_Draw(Mesh mesh, PetscViewer v) {
824: double startX = mesh->startX;
825: double endX = mesh->endX;
826: double startY = mesh->startY;
827: double endY = mesh->endY;
828: double sizeX = mesh->sizeX;
829: double sizeY = mesh->sizeY;
830: double offX, offY;
831: PetscDraw draw;
832: PetscTruth isnull;
833: PetscDrawButton button;
834: int pause, hElement;
835: PetscReal xc, yc, scale;
836: int xsize, ysize;
837: int ierr;
840: PetscViewerDrawGetDraw(v, 0, &draw);
841: PetscDrawIsNull(draw, &isnull);
842: if (isnull) return(0);
843: if (sizeX > sizeY) {
844: ysize = 300;
845: if (sizeY > 0.0) {
846: xsize = ysize * (int) (sizeX/sizeY);
847: } else {
848: return(0);
849: }
850: } else {
851: xsize = 300;
852: if (sizeX > 0.0) {
853: ysize = xsize * (int) (sizeY/sizeX);
854: } else {
855: return(0);
856: }
857: }
858: if ((xsize == 0) || (ysize == 0)) return(0);
859: offX = 0.02*sizeX;
860: offY = 0.02*sizeY;
861: PetscDrawResizeWindow(draw, xsize, ysize);
862: MeshView_DrawStatusLine_Private(mesh, draw, startX, startY, endX, endY, offX, offY);
863: MeshView_Triangular_2D_Draw_Mesh(mesh, draw);
864: PetscDrawGetPause(draw, &pause);
865: if (pause >= 0) {
866: PetscSleep(pause);
867: return(0);
868: }
870: /* Allow the mesh to zoom or shrink */
871: PetscDrawCheckResizedWindow(draw);
872: PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
873: while (button != BUTTON_RIGHT) {
874: PetscDrawSynchronizedClear(draw);
875: switch (button)
876: {
877: case BUTTON_LEFT:
878: scale = 0.5;
879: break;
880: case BUTTON_CENTER:
881: scale = 2.0;
882: break;
883: case BUTTON_LEFT_SHIFT:
884: scale = 1.0;
885: MeshLocatePoint(mesh, xc, yc, 0.0, &hElement);
886: MeshSetHighlightElement(mesh, hElement);
887: break;
888: case BUTTON_CENTER_SHIFT:
889: scale = 1.0;
890: MeshLocatePoint(mesh, xc, yc, 0.0, &hElement);
891: MeshViewElement_Triangular_2D_Draw(mesh, hElement, v);
892: break;
893: default:
894: scale = 1.0;
895: break;
896: }
897: if (scale != 1.0) {
898: startX = scale*(startX + sizeX - xc) + xc - sizeX*scale;
899: endX = scale*(endX - sizeX - xc) + xc + sizeX*scale;
900: startY = scale*(startY + sizeY - yc) + yc - sizeY*scale;
901: endY = scale*(endY - sizeY - yc) + yc + sizeY*scale;
902: sizeX *= scale;
903: sizeY *= scale;
904: offX *= scale;
905: offY *= scale;
906: }
908: MeshView_DrawStatusLine_Private(mesh, draw, startX, startY, endX, endY, offX, offY);
909: MeshView_Triangular_2D_Draw_Mesh(mesh, draw);
910: PetscDrawCheckResizedWindow(draw);
911: PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
912: }
914: return(0);
915: }
917: #undef __FUNCT__
919: int MeshView_Triangular_2D(Mesh mesh, PetscViewer viewer)
920: {
921: PetscTruth isascii, issilo, isdraw, isvu, ismathematica;
922: int ierr;
925: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &isascii);
926: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_SILO, &issilo);
927: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_DRAW, &isdraw);
928: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_VU, &isvu);
929: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_MATHEMATICA, &ismathematica);
930: if (isascii == PETSC_TRUE) {
931: MeshView_Triangular_2D_File(mesh, viewer);
932: } else if (issilo == PETSC_TRUE) {
933: MeshView_Triangular_2D_Silo(mesh, viewer);
934: } else if (isdraw == PETSC_TRUE) {
935: MeshView_Triangular_2D_Draw(mesh, viewer);
936: } else if (isvu == PETSC_TRUE) {
937: MeshView_Triangular_2D_VU(mesh, viewer);
938: } else if (ismathematica == PETSC_TRUE) {
939: #ifdef PETSC_HAVE_MATHEMATICA
940: ViewerMathematica_Mesh_Triangular_2D(viewer, mesh);
941: #endif
942: }
944: return(0);
945: }