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