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_TRUE)
 21:     PetscViewerASCIIPrintf(viewer, "  Triangular 2D grid withn");
 22:   else
 23:     PetscViewerASCIIPrintf(viewer, "  Periodic Triangular 2D grid withn");
 24:   if (mesh->numBd == 1)
 25:     PetscViewerASCIIPrintf(viewer, "    %d closed boundaryn", mesh->numBd);
 26:   else
 27:     PetscViewerASCIIPrintf(viewer, "    %d closed boundariesn", mesh->numBd);
 28:   for(i = 0; i < mesh->numBd; i++)
 29:     PetscViewerASCIIPrintf(viewer, "      Boundary %d: %d nodesn", tri->bdMarkers[i], tri->bdBegin[i+1] - tri->bdBegin[i]);
 30:   PetscViewerASCIIPrintf(viewer, "      Total boundary nodes: %d edges: %dn", mesh->numBdNodes, mesh->numBdEdges);
 31:   PetscViewerASCIIGetPointer(viewer, &fd);
 32:   PetscSynchronizedFPrintf(mesh->comm, fd, "  Local graph %d: %d nodes %d bdEdgesn", p->rank, mesh->numNodes, mesh->numBdEdges);
 33:   for(i = 0; i < mesh->numNodes; i++)
 34:     PetscSynchronizedFPrintf(mesh->comm, fd, "    %d %g %g %dn", i, tri->nodes[i*2], tri->nodes[i*2+1], tri->markers[i]);
 35:   PetscSynchronizedFlush(mesh->comm);
 36:   PetscSynchronizedFPrintf(mesh->comm, fd, "  Local graph %d: %d edgesn", p->rank, mesh->numEdges);
 37:   for(i = 0; i < mesh->numEdges; i++)
 38:     PetscSynchronizedFPrintf(mesh->comm, fd, "    %d %d %d %dn", i, tri->edges[i*2], tri->edges[i*2+1], tri->edgemarkers[i]);
 39:   PetscSynchronizedFlush(mesh->comm);
 40:   PetscSynchronizedFPrintf(mesh->comm, fd, "  Local graph %d: %d faces with %d nodes per facen",
 41:                            p->rank, mesh->numFaces, mesh->numCorners);
 42:   for(i = 0; i < mesh->numFaces; i++) {
 43:     PetscSynchronizedFPrintf(mesh->comm, fd, "    %d", i);
 44:     for(j = 0; j < mesh->numCorners; j++)
 45:       PetscSynchronizedFPrintf(mesh->comm, fd, " %d", tri->faces[i*mesh->numCorners+j]);
 46:     PetscSynchronizedFPrintf(mesh->comm, fd, "n");
 47:   }
 48:   for(i = 0; i < mesh->numFaces; i++)
 49:     PetscSynchronizedFPrintf(mesh->comm, fd, "    %d %d %d %dn", i, tri->neighbors[i*3], tri->neighbors[i*3+1],
 50:                              tri->neighbors[i*3+2]);
 51:   PetscSynchronizedFlush(mesh->comm);
 52: 
 53:   if (tri->areas != PETSC_NULL) {
 54:     PetscSynchronizedFPrintf(mesh->comm, fd, "  Local graph %d: element areasn", p->rank);
 55:     for(i = 0; i < mesh->numFaces; i++)
 56:       PetscSynchronizedFPrintf(mesh->comm, fd, "    %d %gn", i, tri->areas[i]);
 57:     PetscSynchronizedFlush(mesh->comm);
 58:   }
 59:   if (tri->aspectRatios != PETSC_NULL) {
 60:     PetscSynchronizedFPrintf(mesh->comm, fd, "  Local graph %d: element aspect ratiosn", p->rank);
 61:     for(i = 0; i < mesh->numFaces; i++)
 62:       PetscSynchronizedFPrintf(mesh->comm, fd, "    %d %gn", i, tri->aspectRatios[i]);
 63:     PetscSynchronizedFlush(mesh->comm);
 64:   }

 66:   if (mesh->partitioned) {
 67:     PartitionView(mesh->part, viewer);
 68:   }
 69:   return(0);
 70: }

 72: static int MeshViewLocalizeMesh_Private(Mesh mesh, double **nodes, int **faces) {
 73:   Mesh_Triangular         *tri = (Mesh_Triangular *) mesh->data;
 74:   Partition_Triangular_2D *q;
 75:   Partition                part;
 76:   MPI_Comm                 comm;
 77:   int                     *numRecvNodes;
 78:   int                     *numRecvFaces;
 79:   int                     *nodeOffsets;
 80:   int                     *faceOffsets;
 81:   int                     *locFaces;
 82:   double                  *locNodes;
 83:   int                      dim, numLocNodes, numNodes, numLocFaces, numFaces, numCorners;
 84:   int                      numProcs, rank;
 85:   int                      proc, node;
 86:   int                      ierr;

 89:   PetscObjectGetComm((PetscObject) mesh, &comm);
 90:   MPI_Comm_size(comm, &numProcs);
 91:   MPI_Comm_rank(comm, &rank);
 92:   MeshGetPartition(mesh, &part);
 93:   MeshGetDimension(mesh, &dim);
 94:   MeshGetNumCorners(mesh, &numCorners);
 95:   PartitionGetNumNodes(part, &numLocNodes);
 96:   PartitionGetTotalNodes(part, &numNodes);
 97:   PartitionGetNumElements(part, &numLocFaces);
 98:   PartitionGetTotalElements(part, &numFaces);
 99:   if (numProcs > 1) {
100:     q = (Partition_Triangular_2D *) part->data;
101:     /* Allocate global arrays */
102:     PetscMalloc(numNodes*dim           * sizeof(double), nodes);
103:     PetscMalloc(numFaces*numCorners    * sizeof(int),    faces);
104:     PetscMalloc(numLocFaces*numCorners * sizeof(int),    &locFaces);
105:     locNodes = tri->nodes;

107:     /* Calculate offsets */
108:     PetscMalloc(numProcs * sizeof(int), &numRecvNodes);
109:     PetscMalloc(numProcs * sizeof(int), &numRecvFaces);
110:     PetscMalloc(numProcs * sizeof(int), &nodeOffsets);
111:     PetscMalloc(numProcs * sizeof(int), &faceOffsets);
112:     for(proc = 0; proc < numProcs; proc++) {
113:       numRecvNodes[proc] = (q->firstNode[proc+1]       - q->firstNode[proc])*dim;
114:       numRecvFaces[proc] = (part->firstElement[proc+1] - part->firstElement[proc])*numCorners;
115:     }
116:     nodeOffsets[0] = 0;
117:     faceOffsets[0] = 0;
118:     for(proc = 1; proc < numProcs; proc++) {
119:       nodeOffsets[proc] = numRecvNodes[proc-1] + nodeOffsets[proc-1];
120:       faceOffsets[proc] = numRecvFaces[proc-1] + faceOffsets[proc-1];
121:     }

123:     /* Local to global node number conversion */
124:     for(node = 0; node < numLocFaces*numCorners; node++) {
125:       PartitionLocalToGlobalNodeIndex(part, tri->faces[node], &locFaces[node]);
126:     }

128:     /* Collect global arrays */
129:     MPI_Gatherv(locNodes, numLocNodes*dim,        MPI_DOUBLE, *nodes, numRecvNodes, nodeOffsets, MPI_DOUBLE, 0, comm);
130:     MPI_Gatherv(locFaces, numLocFaces*numCorners, MPI_INT,    *faces, numRecvFaces, faceOffsets, MPI_INT,    0, comm);

132:     /* Cleanup */
133:     PetscFree(locFaces);
134:     PetscFree(numRecvNodes);
135:     PetscFree(numRecvFaces);
136:     PetscFree(nodeOffsets);
137:     PetscFree(faceOffsets);

139:     if (rank == 0) {
140:       /* We must globally renumber and permute so that midnodes come after vertices */
141:       AOPetscToApplication(q->nodeOrdering, numFaces*numCorners, *faces);
142:       AOPetscToApplicationPermuteReal(q->nodeOrdering, dim, *nodes);
143:       if (mesh->nodeOrdering) {
144:         AOPetscToApplication(mesh->nodeOrdering, numFaces*numCorners, *faces);
145:         AOPetscToApplicationPermuteReal(mesh->nodeOrdering, dim, *nodes);
146:       }
147:     }
148:   } else {
149:     *nodes = tri->nodes;
150:     *faces = tri->faces;
151:   }
152:   return(0);
153: }

155: static int MeshViewLocalizeDestroy_Private(Mesh mesh, double *nodes, int *faces) {
156:   MPI_Comm comm;
157:   int      numProcs;
158:   int      ierr;

161:   PetscObjectGetComm((PetscObject) mesh, &comm);
162:   MPI_Comm_size(comm, &numProcs);
163:   if (numProcs > 1) {
164:     PetscFree(nodes);
165:     PetscFree(faces);
166:   }
167:   return(0);
168: }

170: static int MeshView_Triangular_2D_VU(Mesh mesh, PetscViewer viewer) {
171:   MPI_Comm  comm;
172:   Partition part;
173:   FILE     *fp;
174:   double   *nodes;
175:   int      *faces;
176:   int       dim, numNodes, numFaces, numCorners;
177:   int       d, node, elem;
178:   int       ierr;

181:   PetscObjectGetComm((PetscObject) mesh, &comm);
182:   MeshGetPartition(mesh, &part);
183:   MeshGetDimension(mesh, &dim);
184:   MeshGetNumCorners(mesh, &numCorners);
185:   PartitionGetTotalNodes(part, &numNodes);
186:   PartitionGetTotalElements(part, &numFaces);
187:   PetscViewerVUGetPointer(viewer, &fp);
188:   MeshViewLocalizeMesh_Private(mesh, &nodes, &faces);
189:   /* Write the node coordinates */
190:   PetscFPrintf(comm, fp, "// Node coordinatesnFIELD Coordinates() = n{n");
191:   for(node = 0; node < numNodes; node++) {
192:     PetscFPrintf(comm, fp, "  ");
193:     for(d = 0; d < dim-1; d++) {
194:       PetscFPrintf(comm, fp, "%g ", nodes[node*dim+d]);
195:     }
196:     PetscFPrintf(comm, fp, "%gn", nodes[node*dim+(dim-1)]);
197:   }
198:   PetscFPrintf(comm, fp, "};nn");
199:   /* Write the node connectivity */
200:   if (numCorners == 6) {
201:     PetscFPrintf(comm, fp, "// Full mesh connectivitynFIELD<int> FullConnectivity() =n{n");
202:     for (elem = 0; elem < numFaces; elem++) {
203:       PetscFPrintf(comm, fp, "  ");
204:       for (node = 0; node < numCorners-1; node++) {
205:         PetscFPrintf(comm, fp, " %d ", faces[elem*numCorners+node]+1);
206:       }
207:       PetscFPrintf(comm, fp, "%dn", faces[elem*numCorners+numCorners-1]+1);
208:     }
209:     PetscFPrintf(comm, fp, "};nn");
210:   }
211:   PetscFPrintf(comm, fp, "// Vertex mesh connectivitynFIELD<int> Connectivity() =n{n");
212:   for (elem = 0; elem < numFaces; elem++) {
213:     PetscFPrintf(comm, fp, "  %d %d %dn", faces[elem*numCorners]+1, faces[elem*numCorners+1]+1,
214:                  faces[elem*numCorners+2]+1);
215:   }
216:   PetscFPrintf(comm, fp, "};nn");
217:   /* Define mesh itself */
218:   if (numCorners == 6) {
219:     PetscFPrintf(comm, fp, "MESH PetscMesh() =n{n  ZONE Zone1(LagrTrian06, Coordinates, FullConnectivity);n};nn");
220:   } else {
221:     PetscFPrintf(comm, fp, "MESH PetscMesh() =n{n  ZONE Zone1(LagrTrian03, Coordinates, Connectivity);n};nn");
222:   }

224:   MeshViewLocalizeDestroy_Private(mesh, nodes, faces);
225:   return(0);
226: }

228: #if 0
229: static int MeshView_Triangular_2D_Poum(Mesh mesh, PetscViewer viewer)
230: {
231:   Mesh_Triangular         *tri = (Mesh_Triangular *) mesh->data;
232:   Partition                p   = mesh->part;
233:   Partition_Triangular_2D *q   = (Partition_Triangular_2D *) p->data;
234:   int                      numCorners = tri->numCorners;
235:   int                     *faces;
236:   double                  *nodes;
237:   FILE                    *fp;
238:   int                      proc, elem, node;
239:   int                      ierr;

242:   MeshViewLocalizeMesh_Private(mesh, &nodes, &faces);
243:   ViewerPoumGetMeshPointer(viewer, &fp);
244:   /* Writes the nodes description */
245:   if (p->rank == 0) {
246:     PetscFPrintf(PETSC_COMM_SELF, fp, "Nodes Nmattn");
247:     for (node = 0; node < tri->numVertices; node++)
248:       PetscFPrintf(PETSC_COMM_SELF, fp, " %dt %ft %f 0.0n", node + 1, nodes[node*2], nodes[node*2+1]);
249:   }

251:   /* Writes the element description */
252:   if (p->rank == 0) {
253:     PetscFPrintf(PETSC_COMM_SELF, fp, "Elements Ematt using Nmattn");
254:     for (elem = 0; elem < p->numElements; elem++)
255:       PetscFPrintf(PETSC_COMM_SELF, fp, " %dt 4t %dt %dt %dn", elem + 1, faces[elem*numCorners]+1,
256:                    faces[elem*numCorners+1]+1, faces[elem*numCorners+2]+1);
257:   }

259:   ViewerPoumGetPartitionPointer(viewer, &fp);
260:   if (p->rank == 0) {
261:     /* Print the partition */
262:     PetscFPrintf(PETSC_COMM_SELF, fp, "Decomposition Dmatt using Emattn %dn", p->numProcs);
263:     for(proc = 0; proc< p->numProcs; proc++) {
264:       PetscFPrintf(PETSC_COMM_SELF, fp, "%dn", p->firstElement[proc+1] - p->firstElement[proc]);
265:       for (elem = 0; elem < p->firstElement[proc+1] - p->firstElement[proc]; elem++)
266:         PetscFPrintf(PETSC_COMM_SELF, fp, " %dn", elem + p->firstElement[proc] + 1);
267:     }
268:   }
269:   MeshViewLocalizeDestroy_Private(mesh, nodes, faces);

271:   return(0);
272: }
273: #endif

275: #ifdef PETSC_HAVE_SILO
276: static int ViewerSiloNodeValues_Private(PetscViewer viewer, Mesh mesh, Vec vec, char *fieldName)
277: {
278:   Viewer_Silo     *vsilo    = (Viewer_Silo *) viewer->data;
279:   Mesh_Triangular *tri      = (Mesh_Triangular *) mesh->data;
280:   int              dim      = mesh->dim;
281:   int              numNodes = tri->numNodes;
282:   DBfile          *fp;
283:   DBoptlist       *optlist;
284:   Scalar          *array;
285:   float          **newArray;
286:   char           **subNames;
287:   char             name[256];
288:   char             buf[1024];
289:   int              size, nameLen, len;
290:   int              node, c;
291:   int              ierr;

294:   ViewerSiloCheckMesh(viewer, mesh);
295:   VecGetSize(vec, &size);
296:   ViewerSiloGetFilePointer(viewer, &fp);
297:   VecGetArray(vec, &array);
298:   /* Name vector */
299:   if (vsilo->objName != PETSC_NULL) {
300:     PetscStrncpy(name, vsilo->objName, 240);
301:   } else {
302:     PetscStrcpy(name, "");
303:   }
304:   PetscStrcat(name, fieldName);
305:   /* Allocate storage compatible with SILO calls */
306:   PetscMalloc(dim * sizeof(char *),  &subNames);
307:   PetscMalloc(dim * sizeof(float *), &newArray);
308:   PetscStrlen(name, &nameLen);
309:   nameLen += (int) log10((double) dim) + 6;
310:   for(c = 0; c < dim; c++) {
311:     PetscMalloc(nameLen  * sizeof(char),  &subNames[c]);
312:     PetscMalloc(numNodes * sizeof(float), &newArray[c]);
313:     sprintf(subNames[c], "%scomp%d", name, c);
314:   }
315:   /* Convert vector to SILO format */
316:   for(node = 0; node < numNodes; node++) {
317:     for(c = 0; c < dim; c++) {
318:       newArray[c][node] = array[node*dim+c];
319:     }
320:   }
321:   /* Add each component */
322:   for(c = 0; c < dim; c++) {
323:     optlist = DBMakeOptlist(3);
324:     DBAddOption(optlist, DBOPT_LABEL, name);
325:     if (vsilo->meshName == PETSC_NULL) {
326:       DBPutUcdvar1(fp, subNames[c], "PetscMesh",     newArray[c], numNodes, PETSC_NULL, 0, DB_FLOAT, DB_NODECENT, optlist);
327: 
328:     } else {
329:       DBPutUcdvar1(fp, subNames[c], vsilo->meshName, newArray[c], numNodes, PETSC_NULL, 0, DB_FLOAT, DB_NODECENT, optlist);
330: 
331:     }
332:     DBFreeOptlist(optlist);
333:   }

335:   if (dim > 1) {
336:     PetscMemzero(buf, 1024 * sizeof(char));
337:     len  = DBGetVarLength(fp, "_meshtv_defvars");
338:     if (len > 0) {
339:       if (DBGetVarType(fp, "_meshtv_defvars") != DB_CHAR) {
340:         SETERRQ(PETSC_ERR_FILE_READ, "Invalid type for variable _meshtv_defvars");
341:       }
342:       if (len > 1024) SETERRQ(PETSC_ERR_SUP, "Need to do dyanmic allocation here");
343:       DBReadVar(fp, "_meshtv_defvars", buf);
344:       PetscStrcat(buf, ";vec");
345:     } else {
346:       PetscStrcpy(buf, "vec");
347:     }
348:     PetscStrcat(buf, name);
349:     PetscStrcat(buf, " vector {");
350:     for(c = 0; c < dim-1; c++) {
351:       PetscStrcat(buf, subNames[c]);
352:       PetscStrcat(buf, ",");
353:     }
354:     PetscStrcat(buf, subNames[c]);
355:     PetscStrcat(buf, "}");
356:     PetscStrlen(buf, &len);
357:     if (len > 1024) SETERRQ(PETSC_ERR_SUP, "Need to do dyanmic allocation here");
358:     len  = 1024;
359:     DBWrite(fp, "_meshtv_defvars", buf, &len, 1, DB_CHAR);
360:   }

362:   for(c = 0; c < dim; c++) {
363:     PetscFree(subNames[c]);
364:     PetscFree(newArray[c]);
365:   }
366:   PetscFree(subNames);
367:   PetscFree(newArray);
368:   VecRestoreArray(vec, &array);
369:   return(0);
370: }
371: #endif

373: static int MeshView_Triangular_2D_Silo(Mesh mesh, PetscViewer viewer)
374: {
375: #ifdef HAVE_SILO
376:   Viewer_Silo     *vsilo         = (Viewer_Silo *) viewer->data;
377:   Mesh_Triangular *tri           = (Mesh_Triangular *) mesh->data;
378:   int              numCorners    = tri->numCorners;
379:   int              numElements   = tri->numFaces;
380:   int             *faces         = tri->faces;
381:   int             *edges         = tri->edges;
382:   int              numNodes      = tri->numNodes;
383:   double          *nodes         = tri->nodes;
384:   int              numBdEdges    = tri->numBdEdges;
385:   int             *bdEdges       = tri->bdEdges;
386:   double           sizeX         = mesh->sizeX;
387:   double           sizeY         = mesh->sizeY;
388:   PetscTruth       isPeriodic    = mesh->isPeriodic;
389:   PetscTruth       xPer          = mesh->isPeriodicDim[0];
390:   PetscTruth       yPer          = mesh->isPeriodicDim[1];
391:   int              pointsPerEdge = 2;   /* The number of vertices on an edge - should always be 2 for 2D mesh */
392:   int              pointsPerFace = 3;   /* The number of vertices on a face - should always be 3 for triangular mesh */
393:   int              numShapes     = 1;   /* The number of different shapes - we only have triangles */
394:   char             edgeName[256];       /* The name of the edge list */
395:   char             zoneName[256];       /* The name of the zone list */
396:   char             meshName[256];       /* The name of the mesh */
397:   DBfile          *fp;
398:   DBoptlist       *optlist;
399:   int             *facelist, *zonelist;
400:   float           *xcoords, *ycoords;
401:   float           *coords[2];
402:   int              numFaces, numZones;
403:   int              node, face, zone;
404:   int              ierr;
405: 
406: 
408:   /* Setup base names */
409:   if (vsilo->objName != PETSC_NULL) {
410:     PetscStrncpy(edgeName, vsilo->objName, 251);
411:     PetscStrncpy(zoneName, vsilo->objName, 251);
412:     PetscStrncpy(meshName, vsilo->objName, 251);
413:   } else {
414:     PetscStrcpy(edgeName, "Petsc");
415:     PetscStrcpy(zoneName, "Petsc");
416:     PetscStrcpy(meshName, "Petsc");
417:   }
418:   /* Add suffixes */
419:   PetscStrcat(edgeName, "Face");
420:   PetscStrcat(zoneName, "Zone");
421:   PetscStrcat(meshName, "Mesh");

423:   /* Allocate space for the new arrays that have to be created */
424:   PetscMalloc(numNodes                  * sizeof(float), &xcoords);
425:   PetscMalloc(numNodes                  * sizeof(float), &ycoords);
426:   PetscMalloc(numBdEdges*pointsPerEdge  * sizeof(int),   &facelist);
427:   PetscMalloc(numElements*pointsPerFace * sizeof(int),   &zonelist);

429:   if (isPeriodic == PETSC_TRUE) {
430:     for(face = 0, numFaces = 0; face < numBdEdges; face++) {
431:       if (((xPer == PETSC_TRUE) && (PetscAbsScalar(nodes[edges[bdEdges[face]*2]*2]   - nodes[edges[bdEdges[face]*2+1]*2])   > 0.5*sizeX)) ||
432:           ((yPer == PETSC_TRUE) && (PetscAbsScalar(nodes[edges[bdEdges[face]*2]*2+1] - nodes[edges[bdEdges[face]*2+1]*2+1]) > 0.5*sizeY)))
433:         continue;
434:       facelist[numFaces*2]   = edges[bdEdges[face]*2];
435:       facelist[numFaces*2+1] = edges[bdEdges[face]*2+1];
436:       numFaces++;
437:     }
438:     /* DBPutZoneList only wants the points on the corners of the triangle, so we remove the other nodes */
439:     for(zone = 0, numZones = 0; zone < numElements; zone++) {
440:       if (((xPer == PETSC_TRUE) && (PetscAbsScalar(nodes[faces[zone*numCorners]*2]   - nodes[faces[zone*numCorners+1]*2])   > 0.5*sizeX)) ||
441:           ((yPer == PETSC_TRUE) && (PetscAbsScalar(nodes[faces[zone*numCorners]*2+1] - nodes[faces[zone*numCorners+1]*2+1]) > 0.5*sizeY)))
442:         continue;
443:       if (((xPer == PETSC_TRUE) && (PetscAbsScalar(nodes[faces[zone*numCorners]*2]   - nodes[faces[zone*numCorners+2]*2])   > 0.5*sizeX)) ||
444:           ((yPer == PETSC_TRUE) && (PetscAbsScalar(nodes[faces[zone*numCorners]*2+1] - nodes[faces[zone*numCorners+2]*2+1]) > 0.5*sizeY)))
445:         continue;
446:       zonelist[numZones*3]   = faces[zone*numCorners];
447:       zonelist[numZones*3+1] = faces[zone*numCorners+1];
448:       zonelist[numZones*3+2] = faces[zone*numCorners+2];
449:       numZones++;
450:     }
451:   } else {
452:     numFaces = numBdEdges;
453:     for(face = 0; face < numBdEdges; face++) {
454:       facelist[face*2]   = edges[bdEdges[face]*2];
455:       facelist[face*2+1] = edges[bdEdges[face]*2+1];
456:     }
457:     /* DBPutZoneList only wants the points on the corners of the triangle, so we remove the other nodes */
458:     numZones = numElements;
459:     for(zone = 0; zone < numElements; zone++) {
460:       zonelist[zone*3]   = faces[zone*numCorners];
461:       zonelist[zone*3+1] = faces[zone*numCorners+1];
462:       zonelist[zone*3+2] = faces[zone*numCorners+2];
463:     }
464:   }
465: 
466:   /* DBPutUcdMesh expects x- and y- coordinates to be in separate arrays, so we split them up */
467:   for(node = 0; node < numNodes; node++) {
468:     xcoords[node] = nodes[node*2];
469:     ycoords[node] = nodes[node*2+1];
470:   }
471:   coords[0] = xcoords;
472:   coords[1] = ycoords;
473: 
474:   ViewerSiloGetFilePointer(viewer, &fp);
475:   DBPutFacelist(fp, edgeName, numFaces, 2, facelist, numFaces*pointsPerEdge, 0, PETSC_NULL,
476:                        &pointsPerEdge, &numFaces, numShapes, PETSC_NULL, PETSC_NULL, 0);
477: 
478:   DBPutZonelist(fp, zoneName, numZones, 2, zonelist, numZones*pointsPerFace, 0, &pointsPerFace, &numZones, 1);
479: 
480:   PetscFree(facelist);
481:   PetscFree(zonelist);
482: 
483:   /* Assorted options to make the graph prettier */
484:   optlist = DBMakeOptlist(5);
485:   /* DBAddOption(optlist, DBOPT_COORDSYS, DB_CARTESIAN);                                            */
486:   DBAddOption(optlist, DBOPT_XLABEL, "X")    ;
487:   DBAddOption(optlist, DBOPT_YLABEL, "Y");
488:   DBAddOption(optlist, DBOPT_XUNITS, "cm");
489:   DBAddOption(optlist, DBOPT_YUNITS, "cm");
490:   DBPutUcdmesh(fp, meshName, 2, PETSC_NULL, coords, numNodes, numZones, zoneName, edgeName, DB_FLOAT, optlist);
491: 
492:   DBFreeOptlist(optlist);

494:   /* Check for moving mesh and visualize velocity/acceleration */
495:   if (mesh->movingMesh == PETSC_TRUE) {
496:     if (mesh->nodeVelocities    != PETSC_NULL) {
497:       ViewerSiloSetName(viewer, "Mesh");
498:       ViewerSiloNodeValues_Private(viewer, mesh, mesh->nodeVelocities, "Velocity");
499:     }
500:     if (mesh->nodeAccelerations != PETSC_NULL) {
501:       ViewerSiloSetName(viewer, "Mesh");
502:       ViewerSiloNodeValues_Private(viewer, mesh, mesh->nodeAccelerations, "Acceleration");
503:     }
504:     ViewerSiloClearName(viewer);
505:   }

507: #ifdef NOT_WORKING
508:   /* Set mesh name */
509:   if (vsilo->objName != PETSC_NULL) {
510:     PetscMalloc((PetscStrlen(meshName)+1) * sizeof(char), &newMeshName);
511:     PetscStrcpy(newMeshName, meshName);
512:     ViewerSiloSetMeshName(viewer, newMeshName);
513:   }
514: #endif
515:   return(0);
516: #else
517:   SETERRQ(PETSC_ERR_SUP, "Need LLNL's Silo package");
518: #endif
519: }

521: static int MeshView_Triangular_2D_Draw_Edge(Mesh mesh, PetscDraw draw, int edge, int color) {
522:   Mesh_Triangular *tri        = (Mesh_Triangular *) mesh->data;
523:   int              numCorners = mesh->numCorners;
524:   int             *edges      = tri->edges;
525:   int             *elements   = tri->faces;
526:   double          *nodes      = tri->nodes;
527:   int              numOverlapElements;
528:   int              elem, corner, corner2, midCorner, midNode;
529:   int              ierr;

532:   if (numCorners == 3) {
533:     MeshDrawLine(mesh, draw, nodes[edges[edge*2]*2], nodes[edges[edge*2]*2+1], nodes[edges[edge*2+1]*2],
534:                         nodes[edges[edge*2+1]*2+1], color);
535: 
536:   } else if (numCorners == 6) {
537:     /* Find midnode */
538:     PartitionGetNumOverlapElements(mesh->part, &numOverlapElements);
539:     for(elem = 0, midNode = -1; elem < numOverlapElements; elem++) {
540:       for(corner = 0; corner < 3; corner++) {
541:         if (elements[elem*numCorners+corner] == edges[edge*2]) break;
542:       }
543:       if (corner < 3) {
544:         for(corner2 = 0; corner2 < 3; corner2++) {
545:           if (elements[elem*numCorners+corner2] == edges[edge*2+1]) break;
546:         }
547:         if (corner2 < 3) {
548:           midCorner = ((corner+corner2)*2)%3 + 3;
549:           midNode   = elements[elem*numCorners+midCorner];
550:           break;
551:         }
552:       }
553:     }
554:     if (midNode < 0) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid mesh connectivity");
555:     MeshDrawLine(mesh, draw, nodes[edges[edge*2]*2], nodes[edges[edge*2]*2+1], nodes[midNode*2],
556:                         nodes[midNode*2+1],      color);
557: 
558:     MeshDrawLine(mesh, draw, nodes[midNode*2],    nodes[midNode*2+1],    nodes[edges[edge*2+1]*2],
559:                         nodes[edges[edge*2+1]*2+1], color);
560: 
561:   } else {
562:     SETERRQ1(PETSC_ERR_SUP, "Invalid number of corners %d", numCorners);
563:   }
564:   return(0);
565: }

567: static int MeshView_Triangular_2D_Draw_Element(Mesh mesh, PetscDraw draw, int elem, int color, PetscTruth drawEdges) {
568:   Mesh_Triangular *tri        = (Mesh_Triangular *) mesh->data;
569:   int              numCorners = mesh->numCorners;
570:   int             *elements   = tri->faces;
571:   int             *neighbors  = tri->faces;
572:   double          *nodes      = tri->nodes;
573:   int              corner;
574:   int              ierr;

577:   if (numCorners == 3) {
578:     MeshDrawTriangle(mesh, draw, nodes[elements[elem*numCorners]*2],   nodes[elements[elem*numCorners]*2+1],
579:                             nodes[elements[elem*numCorners+1]*2], nodes[elements[elem*numCorners+1]*2+1],
580:                             nodes[elements[elem*numCorners+2]*2], nodes[elements[elem*numCorners+2]*2+1],
581:                             color, color, color);
582: 
583:   } else if (numCorners == 6) {
584:     /* Draw 4 interior triangles */
585:     MeshDrawTriangle(mesh, draw, nodes[elements[elem*numCorners]*2],   nodes[elements[elem*numCorners]*2+1],
586:                             nodes[elements[elem*numCorners+5]*2], nodes[elements[elem*numCorners+5]*2+1],
587:                             nodes[elements[elem*numCorners+4]*2], nodes[elements[elem*numCorners+4]*2+1],
588:                             color, color, color);
589: 
590:     MeshDrawTriangle(mesh, draw, nodes[elements[elem*numCorners+1]*2], nodes[elements[elem*numCorners+1]*2+1],
591:                             nodes[elements[elem*numCorners+3]*2], nodes[elements[elem*numCorners+3]*2+1],
592:                             nodes[elements[elem*numCorners+5]*2], nodes[elements[elem*numCorners+5]*2+1],
593:                             color, color, color);
594: 
595:     MeshDrawTriangle(mesh, draw, nodes[elements[elem*numCorners+2]*2], nodes[elements[elem*numCorners+2]*2+1],
596:                             nodes[elements[elem*numCorners+4]*2], nodes[elements[elem*numCorners+4]*2+1],
597:                             nodes[elements[elem*numCorners+3]*2], nodes[elements[elem*numCorners+3]*2+1],
598:                             color, color, color);
599: 
600:     MeshDrawTriangle(mesh, draw, nodes[elements[elem*numCorners+3]*2], nodes[elements[elem*numCorners+3]*2+1],
601:                             nodes[elements[elem*numCorners+4]*2], nodes[elements[elem*numCorners+4]*2+1],
602:                             nodes[elements[elem*numCorners+5]*2], nodes[elements[elem*numCorners+5]*2+1],
603:                             color, color, color);
604: 
605:   } else {
606:     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid number of corners %d", numCorners);
607:   }

609:   if (drawEdges == PETSC_TRUE) {
610:     for(corner = 0; corner < 3; corner++) {
611:       if (neighbors[elem*3+corner] < 0) {
612:         color = PETSC_DRAW_RED;
613:       } else {
614:         color = PETSC_DRAW_BLACK;
615:       }
616:       if (numCorners == 3) {
617:         MeshDrawLine(mesh, draw, nodes[elements[elem*numCorners+corner]*2], nodes[elements[elem*numCorners+corner]*2+1],
618:                             nodes[elements[elem*numCorners+((corner+1)%3)]*2], nodes[elements[elem*numCorners+((corner+1)%3)]*2+1], color);
619: 
620:       } else if (numCorners == 6) {
621:         MeshDrawLine(mesh, draw, nodes[elements[elem*numCorners+corner]*2], nodes[elements[elem*numCorners+corner]*2+1],
622:                             nodes[elements[elem*numCorners+corner+3]*2], nodes[elements[elem*numCorners+corner+3]*2+1], color);
623: 
624:         MeshDrawLine(mesh, draw, nodes[elements[elem*numCorners+corner+3]*2], nodes[elements[elem*numCorners+corner+3]*2+1],
625:                             nodes[elements[elem*numCorners+((corner+1)%3)]*2], nodes[elements[elem*numCorners+((corner+1)%3)]*2+1], color);
626: 
627:       }
628:     }
629:   }
630:   return(0);
631: }

633: static int MeshView_Triangular_2D_Draw_Mesh(Mesh mesh, PetscDraw draw) {
634:   Mesh_Triangular *tri         = (Mesh_Triangular *) mesh->data;
635:   int             *edgemarkers = tri->edgemarkers;
636:   int              numElements = mesh->numFaces;
637:   int              numEdges    = mesh->numEdges;
638:   int              hElement    = mesh->highlightElement;
639:   int              color, rank;
640:   int              elem, edge;
641:   int              ierr;

644:   MPI_Comm_rank(mesh->comm, &rank);
645:   for(elem = 0; elem < numElements; elem++) {
646:     if (hElement == elem) {
647:       color = PETSC_DRAW_WHITE;
648:     } else {
649:       color = PETSC_DRAW_RED + rank + 1;
650:     }
651:     MeshView_Triangular_2D_Draw_Element(mesh, draw, elem, color, PETSC_FALSE);
652:   }
653:   /* Needed so that triangles do not overwrite edges */
654:   PetscDrawSynchronizedFlush(draw);

656:   for(edge = 0; edge < numEdges; edge++) {
657:     if (edgemarkers[edge]) {
658:       color = PETSC_DRAW_RED;
659:     } else {
660:       color = PETSC_DRAW_BLACK;
661:     }
662:     MeshView_Triangular_2D_Draw_Edge(mesh, draw, edge, color);
663:   }
664:   PetscDrawSynchronizedFlush(draw);
665:   return(0);
666: }

668: /*
669:   Remember that this function sets the coordinates for the window.
670: */
671: static int MeshView_DrawStatusLine_Private(Mesh mesh, PetscDraw draw, double startX, double startY, double endX,
672:                                            double endY, double offX, double offY) {
673:   Partition part = mesh->part;
674:   MPI_Comm  comm;
675:   char      statusLine[1024];
676:   PetscReal textWidth, textHeight;
677:   int       hElement, rank;
678:   int       ierr;

681:   PetscObjectGetComm((PetscObject) mesh, &comm);
682:   MPI_Comm_rank(comm, &rank);
683:   MeshGetHighlightElement(mesh, &hElement);
684:   PartitionLocalToGlobalElementIndex(part, hElement, &hElement);
685:   sprintf(statusLine, "Element %d Proc: %d", hElement, rank);
686:   PetscDrawStringGetSize(draw, &textWidth, &textHeight);
687:   PetscDrawSetCoordinates(draw, startX-offX, startY-offY-textHeight, endX+offX, endY+offY);
688:   if (hElement >= 0) {
689:     PetscDrawString(draw, startX-offX, startY-offY-textHeight, PETSC_DRAW_BLACK, statusLine);
690:   }
691:   return(0);
692: }

694: static int MeshViewElement_Triangular_2D_Draw(Mesh mesh, int elem, PetscViewer v) {
695:   Mesh_Triangular *tri       = (Mesh_Triangular *) mesh->data;
696:   int             dim        = mesh->dim;
697:   int             numCorners = mesh->numCorners;
698:   int            *elements   = tri->faces;
699:   double         *nodes      = tri->nodes;
700:   double          startX, endX;
701:   double          startY, endY;
702:   double          sizeX,  sizeY;
703:   double          offX,   offY;
704:   PetscDraw       draw;
705:   PetscTruth      isnull;
706:   PetscDrawButton button;
707:   int             pause, color;
708:   PetscReal       xc, yc, scale;
709:   int             xsize, ysize;
710:   int             corner;
711:   int             ierr;

714:   startX = endX = nodes[elements[elem*numCorners]*dim+0];
715:   startY = endY = nodes[elements[elem*numCorners]*dim+1];
716:   for(corner = 1; corner < numCorners; corner++) {
717:     if (nodes[elements[elem*numCorners+corner]*dim+0] < startX) startX = nodes[elements[elem*numCorners+corner]*dim+0];
718:     if (nodes[elements[elem*numCorners+corner]*dim+0] > endX)   endX   = nodes[elements[elem*numCorners+corner]*dim+0];
719:     if (nodes[elements[elem*numCorners+corner]*dim+1] < startY) startY = nodes[elements[elem*numCorners+corner]*dim+1];
720:     if (nodes[elements[elem*numCorners+corner]*dim+1] > endY)   endY   = nodes[elements[elem*numCorners+corner]*dim+1];
721:   }
722:   sizeX = endX - startX;
723:   sizeY = endY - startY;
724:   color = PETSC_DRAW_WHITE;
725:   PetscViewerDrawGetDraw(v, 0, &draw);
726:   PetscDrawIsNull(draw, &isnull);
727:   if (isnull) return(0);
728:   if (sizeX > sizeY) {
729:     ysize = 300;
730:     if (sizeY > 0.0) {
731:       xsize = ysize * (int) (sizeX/sizeY);
732:     } else {
733:       return(0);
734:     }
735:   } else {
736:     xsize = 300;
737:     if (sizeX > 0.0) {
738:       ysize = xsize * (int) (sizeY/sizeX);
739:     } else {
740:       return(0);
741:     }
742:   }
743:   if ((xsize == 0) || (ysize == 0)) return(0);
744:   offX = 0.02*sizeX;
745:   offY = 0.02*sizeY;
746:   PetscDrawResizeWindow(draw, xsize, ysize);
747:   MeshView_DrawStatusLine_Private(mesh, draw, startX, startY, endX, endY, offX, offY);
748:   MeshView_Triangular_2D_Draw_Element(mesh, draw, elem, color, PETSC_TRUE);
749:   PetscDrawGetPause(draw, &pause);
750:   if (pause >= 0) {
751:     PetscSleep(pause);
752:     return(0);
753:   }

755:   /* Allow the mesh to zoom or shrink */
756:   PetscDrawCheckResizedWindow(draw);
757:   PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
758:   while (button != BUTTON_RIGHT) {
759:     PetscDrawSynchronizedClear(draw);
760:     switch (button)
761:     {
762:     case BUTTON_LEFT:
763:       scale = 0.5;
764:       break;
765:     case BUTTON_CENTER:
766:       scale = 2.0;
767:       break;
768:     default:
769:       scale = 1.0;
770:       break;
771:     }
772:     if (scale != 1.0) {
773:       startX = scale*(startX + sizeX - xc) + xc - sizeX*scale;
774:       endX   = scale*(endX   - sizeX - xc) + xc + sizeX*scale;
775:       startY = scale*(startY + sizeY - yc) + yc - sizeY*scale;
776:       endY   = scale*(endY   - sizeY - yc) + yc + sizeY*scale;
777:       sizeX *= scale;
778:       sizeY *= scale;
779:       offX  *= scale;
780:       offY  *= scale;
781:     }

783:     MeshView_DrawStatusLine_Private(mesh, draw, startX, startY, endX, endY, offX, offY);
784:     MeshView_Triangular_2D_Draw_Element(mesh, draw, elem, color, PETSC_TRUE);
785:     PetscDrawCheckResizedWindow(draw);
786:     PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
787:   }

789:   return(0);
790: }

792: static int MeshView_Triangular_2D_Draw(Mesh mesh, PetscViewer v) {
793:   double          startX = mesh->startX;
794:   double          endX   = mesh->endX;
795:   double          startY = mesh->startY;
796:   double          endY   = mesh->endY;
797:   double          sizeX  = mesh->sizeX;
798:   double          sizeY  = mesh->sizeY;
799:   double          offX, offY;
800:   PetscDraw       draw;
801:   PetscTruth      isnull;
802:   PetscDrawButton button;
803:   int             pause, hElement;
804:   PetscReal       xc, yc, scale;
805:   int             xsize, ysize;
806:   int             ierr;

809:   PetscViewerDrawGetDraw(v, 0, &draw);
810:   PetscDrawIsNull(draw, &isnull);
811:   if (isnull) return(0);
812:   if (sizeX > sizeY) {
813:     ysize = 300;
814:     if (sizeY > 0.0) {
815:       xsize = ysize * (int) (sizeX/sizeY);
816:     } else {
817:       return(0);
818:     }
819:   } else {
820:     xsize = 300;
821:     if (sizeX > 0.0) {
822:       ysize = xsize * (int) (sizeY/sizeX);
823:     } else {
824:       return(0);
825:     }
826:   }
827:   if ((xsize == 0) || (ysize == 0)) return(0);
828:   offX = 0.02*sizeX;
829:   offY = 0.02*sizeY;
830:   PetscDrawResizeWindow(draw, xsize, ysize);
831:   MeshView_DrawStatusLine_Private(mesh, draw, startX, startY, endX, endY, offX, offY);
832:   MeshView_Triangular_2D_Draw_Mesh(mesh, draw);
833:   PetscDrawGetPause(draw, &pause);
834:   if (pause >= 0) {
835:     PetscSleep(pause);
836:     return(0);
837:   }

839:   /* Allow the mesh to zoom or shrink */
840:   PetscDrawCheckResizedWindow(draw);
841:   PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
842:   while (button != BUTTON_RIGHT) {
843:     PetscDrawSynchronizedClear(draw);
844:     switch (button)
845:     {
846:     case BUTTON_LEFT:
847:       scale = 0.5;
848:       break;
849:     case BUTTON_CENTER:
850:       scale = 2.0;
851:       break;
852:     case BUTTON_LEFT_SHIFT:
853:       scale = 1.0;
854:       MeshLocatePoint(mesh, xc, yc, 0.0, &hElement);
855:       MeshSetHighlightElement(mesh, hElement);
856:       break;
857:     case BUTTON_CENTER_SHIFT:
858:       scale = 1.0;
859:       MeshLocatePoint(mesh, xc, yc, 0.0, &hElement);
860:       MeshViewElement_Triangular_2D_Draw(mesh, hElement, v);
861:       break;
862:     default:
863:       scale = 1.0;
864:       break;
865:     }
866:     if (scale != 1.0) {
867:       startX = scale*(startX + sizeX - xc) + xc - sizeX*scale;
868:       endX   = scale*(endX   - sizeX - xc) + xc + sizeX*scale;
869:       startY = scale*(startY + sizeY - yc) + yc - sizeY*scale;
870:       endY   = scale*(endY   - sizeY - yc) + yc + sizeY*scale;
871:       sizeX *= scale;
872:       sizeY *= scale;
873:       offX  *= scale;
874:       offY  *= scale;
875:     }

877:     MeshView_DrawStatusLine_Private(mesh, draw, startX, startY, endX, endY, offX, offY);
878:     MeshView_Triangular_2D_Draw_Mesh(mesh, draw);
879:     PetscDrawCheckResizedWindow(draw);
880:     PetscDrawSynchronizedGetMouseButton(draw, &button, &xc, &yc, 0, 0);
881:   }

883:   return(0);
884: }

886: int MeshView_Triangular_2D(Mesh mesh, PetscViewer viewer)
887: {
888:   PetscTruth isascii, issilo, isdraw, isvu, ismathematica;
889:   int        ierr;

892:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII,       &isascii);
893:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_SILO,        &issilo);
894:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_DRAW,        &isdraw);
895:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_VU,          &isvu);
896:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_MATHEMATICA, &ismathematica);
897:   if (isascii == PETSC_TRUE) {
898:     MeshView_Triangular_2D_File(mesh, viewer);
899:   } else if (issilo == PETSC_TRUE) {
900:     MeshView_Triangular_2D_Silo(mesh, viewer);
901:   } else if (isdraw == PETSC_TRUE) {
902:     MeshView_Triangular_2D_Draw(mesh, viewer);
903:   } else if (isvu == PETSC_TRUE) {
904:     MeshView_Triangular_2D_VU(mesh, viewer);
905:   } else if (ismathematica == PETSC_TRUE) {
906: #ifdef PETSC_HAVE_MATHEMATICA
907:     ViewerMathematica_Mesh_Triangular_2D(viewer, mesh);
908: #endif
909:   }

911:   return(0);
912: }