Actual source code: meshpylith.c

  1: #include "src/dm/mesh/meshpylith.h"   /*I      "petscmesh.h"   I*/

  3: #include<list>

  5: namespace ALE {
  6:   namespace PyLith {
  7:     using ALE::Mesh;
  8:     //
  9:     // Builder methods
 10:     //
 11:     inline void Builder::ignoreComments(char *buf, PetscInt bufSize, FILE *f) {
 12:       while((fgets(buf, bufSize, f) != NULL) && ((buf[0] == '#') || (buf[0] == '\0'))) {}
 13:     };
 14:     void Builder::readConnectivity(MPI_Comm comm, const std::string& filename, int& corners, const bool useZeroBase, int& numElements, int *vertices[], int *materials[]) {
 15:       PetscViewer    viewer;
 16:       FILE          *f;
 17:       PetscInt       maxCells = 1024, cellCount = 0;
 18:       PetscInt      *verts;
 19:       PetscInt      *mats;
 20:       char           buf[2048];
 21:       PetscInt       c;
 22:       PetscInt       commRank;

 25:       MPI_Comm_rank(comm, &commRank);
 26:       if (commRank != 0) return;
 27:       PetscViewerCreate(PETSC_COMM_SELF, &viewer);
 28:       PetscViewerSetType(viewer, PETSC_VIEWER_ASCII);
 29:       PetscViewerFileSetMode(viewer, FILE_MODE_READ);
 30:       PetscViewerFileSetName(viewer, filename.c_str());
 31:       PetscViewerASCIIGetPointer(viewer, &f);
 32:       /* Ignore comments */
 33:       ignoreComments(buf, 2048, f);
 34:       do {
 35:         const char *v = strtok(buf, " ");
 36:         int         elementType;

 38:         if (cellCount == maxCells) {
 39:           PetscInt *vtmp, *mtmp;

 41:           vtmp = verts;
 42:           mtmp = mats;
 43:           PetscMalloc2(maxCells*2*corners,PetscInt,&verts,maxCells*2,PetscInt,&mats);
 44:           PetscMemcpy(verts, vtmp, maxCells*corners * sizeof(PetscInt));
 45:           PetscMemcpy(mats,  mtmp, maxCells         * sizeof(PetscInt));
 46:           PetscFree2(vtmp,mtmp);
 47:           maxCells *= 2;
 48:         }
 49:         /* Ignore cell number */
 50:         v = strtok(NULL, " ");
 51:         /* Get element type */
 52:         elementType = atoi(v);
 53:         if (elementType == 1) {
 54:           corners = 8;
 55:         } else if (elementType == 5) {
 56:           corners = 4;
 57:         } else {
 58:           ostringstream msg;

 60:           msg << "We do not accept element type " << elementType << " right now";
 61:           throw ALE::Exception(msg.str().c_str());
 62:         }
 63:         if (cellCount == 0) {
 64:           PetscMalloc2(maxCells*corners,PetscInt,&verts,maxCells,PetscInt,&mats);
 65:         }
 66:         v = strtok(NULL, " ");
 67:         /* Store material type */
 68:         mats[cellCount] = atoi(v);
 69:         v = strtok(NULL, " ");
 70:         /* Ignore infinite domain element code */
 71:         v = strtok(NULL, " ");
 72:         for(c = 0; c < corners; c++) {
 73:           int vertex = atoi(v);
 74: 
 75:           if (!useZeroBase) vertex -= 1;
 76:           verts[cellCount*corners+c] = vertex;
 77:           v = strtok(NULL, " ");
 78:         }
 79:         cellCount++;
 80:       } while(fgets(buf, 2048, f) != NULL);
 81:       PetscViewerDestroy(viewer);
 82:       numElements = cellCount;
 83:       *vertices   = verts;
 84:       *materials  = mats;
 85:     };
 86:     void Builder::readCoordinates(MPI_Comm comm, const std::string& filename, const int dim, int& numVertices, double *coordinates[]) {
 87:       PetscViewer    viewer;
 88:       FILE          *f;
 89:       PetscInt       maxVerts = 1024, vertexCount = 0;
 90:       PetscScalar   *coords;
 91:       double         scaleFactor = 1.0;
 92:       char           buf[2048];
 93:       PetscInt       c;
 94:       PetscInt       commRank;

 97:       MPI_Comm_rank(comm, &commRank);
 98:       if (commRank == 0) {
 99:         PetscViewerCreate(PETSC_COMM_SELF, &viewer);
100:         PetscViewerSetType(viewer, PETSC_VIEWER_ASCII);
101:         PetscViewerFileSetMode(viewer, FILE_MODE_READ);
102:         PetscViewerFileSetName(viewer, filename.c_str());
103:         PetscViewerASCIIGetPointer(viewer, &f);
104:         /* Ignore comments */
105:         ignoreComments(buf, 2048, f);
106:         PetscMalloc(maxVerts*dim * sizeof(PetscScalar), &coords);
107:         /* Read units */
108:         const char *units = strtok(buf, " ");
109:         if (strcmp(units, "coord_units")) {
110:           throw ALE::Exception("Invalid coordinate units line");
111:         }
112:         units = strtok(NULL, " ");
113:         if (strcmp(units, "=")) {
114:           throw ALE::Exception("Invalid coordinate units line");
115:         }
116:         units = strtok(NULL, " ");
117:         if (!strcmp(units, "km")) {
118:           /* Should use Pythia to do units conversion */
119:           scaleFactor = 1000.0;
120:         }
121:         /* Ignore comments */
122:         ignoreComments(buf, 2048, f);
123:         do {
124:           const char *x = strtok(buf, " ");

126:           if (vertexCount == maxVerts) {
127:             PetscScalar *ctmp;

129:             ctmp = coords;
130:             PetscMalloc(maxVerts*2*dim * sizeof(PetscScalar), &coords);
131:             PetscMemcpy(coords, ctmp, maxVerts*dim * sizeof(PetscScalar));
132:             PetscFree(ctmp);
133:             maxVerts *= 2;
134:           }
135:           /* Ignore vertex number */
136:           x = strtok(NULL, " ");
137:           for(c = 0; c < dim; c++) {
138:             coords[vertexCount*dim+c] = atof(x)*scaleFactor;
139:             x = strtok(NULL, " ");
140:           }
141:           vertexCount++;
142:         } while(fgets(buf, 2048, f) != NULL);
143:         PetscViewerDestroy(viewer);
144:         numVertices = vertexCount;
145:         *coordinates = coords;
146:       }
147:     };
148:     // numSplit is the number of split nodes
149:     // splitInd[] is an array of numSplit pairs, <element, vertex>
150:     // splitValues[] is an array of numSplit*dim displacements
151:     void Builder::readSplit(MPI_Comm comm, const std::string& filename, const int dim, const bool useZeroBase, int& numSplit, int *splitInd[], double *splitValues[]) {
152:       PetscViewer    viewer;
153:       FILE          *f;
154:       PetscInt       maxSplit = 1024, splitCount = 0;
155:       PetscInt      *splitId;
156:       PetscScalar   *splitVal;
157:       char           buf[2048];
158:       PetscInt       c;
159:       PetscInt       commRank;

162:       MPI_Comm_rank(comm, &commRank);
163:       if (dim != 3) {
164:         throw ALE::Exception("PyLith only works in 3D");
165:       }
166:       if (commRank != 0) return;
167:       PetscViewerCreate(PETSC_COMM_SELF, &viewer);
168:       PetscViewerSetType(viewer, PETSC_VIEWER_ASCII);
169:       PetscViewerFileSetMode(viewer, FILE_MODE_READ);
170:       PetscExceptionTry1(PetscViewerFileSetName(viewer, filename.c_str()), PETSC_ERR_FILE_OPEN);
171:       if (PetscExceptionValue(ierr)) {
172:         // this means that a caller above me has also tryed this exception so I don't handle it here, pass it up
173:       } else if (PetscExceptionCaught(ierr,PETSC_ERR_FILE_OPEN)) {
174:         // File does not exist
175:         return;
176:       }
177:       PetscViewerASCIIGetPointer(viewer, &f);
178:       /* Ignore comments */
179:       ignoreComments(buf, 2048, f);
180:       PetscMalloc2(maxSplit*2,PetscInt,&splitId,maxSplit*dim,PetscScalar,&splitVal);
181:       do {
182:         const char *s = strtok(buf, " ");

184:         if (splitCount == maxSplit) {
185:           PetscInt    *sitmp;
186:           PetscScalar *svtmp;

188:           sitmp = splitId;
189:           svtmp = splitVal;
190:           PetscMalloc2(maxSplit*2*2,PetscInt,&splitId,maxSplit*dim*2,PetscScalar,&splitVal);
191:           PetscMemcpy(splitId,  sitmp, maxSplit*2   * sizeof(PetscInt));
192:           PetscMemcpy(splitVal, svtmp, maxSplit*dim * sizeof(PetscScalar));
193:           PetscFree2(sitmp,svtmp);
194:           maxSplit *= 2;
195:         }
196:         /* Get element number */
197:         int elem = atoi(s);
198:         if (!useZeroBase) elem -= 1;
199:         splitId[splitCount*2+0] = elem;
200:         s = strtok(NULL, " ");
201:         /* Get node number */
202:         int node = atoi(s);
203:         if (!useZeroBase) node -= 1;
204:         splitId[splitCount*2+1] = node;
205:         s = strtok(NULL, " ");
206:         /* Ignore load history number */
207:         s = strtok(NULL, " ");
208:         /* Get split values */
209:         for(c = 0; c < dim; c++) {
210:           splitVal[splitCount*dim+c] = atof(s);
211:           s = strtok(NULL, " ");
212:         }
213:         splitCount++;
214:       } while(fgets(buf, 2048, f) != NULL);
215:       PetscViewerDestroy(viewer);
216:       numSplit     = splitCount;
217:       *splitInd    = splitId;
218:       *splitValues = splitVal;
219:     };
220:     void Builder::buildSplit(const Obj<split_section_type>& splitField, int numCells, int numSplit, int splitInd[], double splitVals[]) {
221:       const Obj<split_section_type::atlas_type>& atlas = splitField->getAtlas();
222:       const split_section_type::patch_type       patch = 0;
223:       split_section_type::value_type             values[3];
224:       std::map<split_section_type::point_type, std::set<int> > elem2index;

226:       for(int e = 0; e < numSplit; e++) {
227:         atlas->addFiberDimension(patch, splitInd[e*2+0], 1);
228:         elem2index[splitInd[e*2+0]].insert(e);
229:       }
230:       atlas->orderPatches();
231:       splitField->allocate();
232:       for(std::map<split_section_type::point_type, std::set<int> >::const_iterator e_iter = elem2index.begin(); e_iter != elem2index.end(); ++e_iter) {
233:         const split_section_type::point_type& e = e_iter->first;
234:         int                                   k = 0;

236:         for(std::set<int>::const_iterator i_iter = e_iter->second.begin(); i_iter != e_iter->second.end(); ++i_iter, ++k) {
237:           const int& i = *i_iter;

239:           values[k].first    = splitInd[i*2+1] + numCells;
240:           values[k].second.x = splitVals[i*3+0];
241:           values[k].second.y = splitVals[i*3+1];
242:           values[k].second.z = splitVals[i*3+2];
243:         }
244:         splitField->update(patch, e, values);
245:       }
246:     };
247:     void Builder::buildCoordinates(const Obj<section_type>& coords, const int embedDim, const double coordinates[]) {
248:       const section_type::patch_type patch = 0;
249:       const Obj<topology_type::label_sequence>& vertices = coords->getAtlas()->getTopology()->depthStratum(patch, 0);
250:       const int numCells = coords->getAtlas()->getTopology()->heightStratum(patch, 0)->size();

252:       coords->getAtlas()->setFiberDimensionByDepth(patch, 0, embedDim);
253:       coords->getAtlas()->orderPatches();
254:       coords->allocate();
255:       for(topology_type::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
256:         coords->update(patch, *v_iter, &(coordinates[(*v_iter - numCells)*embedDim]));
257:       }
258:     };
259:     void Builder::buildMaterials(const Obj<Mesh::section_type>& matField, const int materials[]) {
260:       const Mesh::section_type::patch_type patch = 0;
261:       const Obj<Mesh::section_type::topology_type::label_sequence>& elements = matField->getAtlas()->getTopology()->heightStratum(patch, 0);

263:       matField->getAtlas()->setFiberDimensionByHeight(patch, 0, 1);
264:       matField->getAtlas()->orderPatches();
265:       matField->allocate();
266:       for(Mesh::section_type::topology_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
267:         double mat = (double) materials[*e_iter];
268:         matField->update(patch, *e_iter, &mat);
269:       }
270:     };
271:     Obj<Mesh> Builder::readMesh(MPI_Comm comm, const int dim, const std::string& basename, const bool useZeroBase = false, const bool interpolate = false, const int debug = 0) {
272:       Obj<Mesh>          mesh     = Mesh(comm, dim, debug);
273:       Obj<sieve_type>    sieve    = new sieve_type(comm, debug);
274:       Obj<topology_type> topology = new topology_type(comm, debug);
275:       int    *cells, *materials;
276:       double *coordinates;
277:       int    *splitInd;
278:       double *splitValues;
279:       int     numCells = 0, numVertices = 0, numCorners = dim+1, numSplit = 0, hasSplit;

281:       ALE::PyLith::Builder::readConnectivity(comm, basename+".connect", numCorners, useZeroBase, numCells, &cells, &materials);
282:       ALE::PyLith::Builder::readCoordinates(comm, basename+".coord", dim, numVertices, &coordinates);
283:       ALE::PyLith::Builder::readSplit(comm, basename+".split", dim, useZeroBase, numSplit, &splitInd, &splitValues);
284:       ALE::New::SieveBuilder<sieve_type>::buildTopology(sieve, dim, numCells, cells, numVertices, interpolate, numCorners);
285:       sieve->stratify();
286:       topology->setPatch(0, sieve);
287:       topology->stratify();
288:       mesh->setTopologyNew(topology);
289:       buildCoordinates(mesh->getSection("coordinates"), dim, coordinates);
290:       buildMaterials(mesh->getSection("material"), materials);
291:       MPI_Allreduce(&numSplit, &hasSplit, 1, MPI_INT, MPI_MAX, comm);
292:       if (hasSplit) {
293:         Obj<split_section_type> splitField = new split_section_type(comm, debug);

295:         splitField->getAtlas()->setTopology(topology);
296:         buildSplit(splitField, numCells, numSplit, splitInd, splitValues);
297:         mesh->setSplitSection(splitField);
298:       }
299:       return mesh;
300:     };
301:     //
302:     // Viewer methods
303:     //
306:     PetscErrorCode Viewer::writeVertices(const Obj<Mesh>& mesh, PetscViewer viewer) {
307:       Obj<Mesh::section_type> coordinates  = mesh->getSection("coordinates");
308:       //Mesh::section_type::patch_type patch;
309:       //const double  *array = coordinates->restrict(Mesh::section_type::patch_type());
310:       //int            dim = mesh->getDimension();
311:       //int            numVertices;

315: #if 0
316:       //FIX:
317:       if (vertexBundle->getGlobalOffsets()) {
318:         numVertices = vertexBundle->getGlobalOffsets()[mesh->commSize()];
319:       } else {
320:         numVertices = mesh->getTopology()->depthStratum(0)->size();
321:       }
322:       PetscViewerASCIIPrintf(viewer,"#\n");
323:       PetscViewerASCIIPrintf(viewer,"coord_units = m\n");
324:       PetscViewerASCIIPrintf(viewer,"#\n");
325:       PetscViewerASCIIPrintf(viewer,"#\n");
326:       PetscViewerASCIIPrintf(viewer,"#  Node      X-coord           Y-coord           Z-coord\n");
327:       PetscViewerASCIIPrintf(viewer,"#\n");
328:       if (mesh->commRank() == 0) {
329:         int numLocalVertices = mesh->getTopology()->depthStratum(0)->size();
330:         int vertexCount = 1;

332:         for(int v = 0; v < numLocalVertices; v++) {
333:           PetscViewerASCIIPrintf(viewer,"%7D ", vertexCount++);
334:           for(int d = 0; d < dim; d++) {
335:             if (d > 0) {
336:               PetscViewerASCIIPrintf(viewer," ");
337:             }
338:             PetscViewerASCIIPrintf(viewer,"% 16.8E", array[v*dim+d]);
339:           }
340:           PetscViewerASCIIPrintf(viewer,"\n");
341:         }
342:         for(int p = 1; p < mesh->commSize(); p++) {
343:           double    *remoteCoords;
344:           MPI_Status status;

346:           MPI_Recv(&numLocalVertices, 1, MPI_INT, p, 1, mesh->comm(), &status);
347:           PetscMalloc(numLocalVertices*dim * sizeof(double), &remoteCoords);
348:           MPI_Recv(remoteCoords, numLocalVertices*dim, MPI_DOUBLE, p, 1, mesh->comm(), &status);
349:           for(int v = 0; v < numLocalVertices; v++) {
350:             PetscViewerASCIIPrintf(viewer,"%7D   ", vertexCount++);
351:             for(int d = 0; d < dim; d++) {
352:               if (d > 0) {
353:                 PetscViewerASCIIPrintf(viewer, " ");
354:               }
355:               PetscViewerASCIIPrintf(viewer, "% 16.8E", remoteCoords[v*dim+d]);
356:             }
357:             PetscViewerASCIIPrintf(viewer, "\n");
358:           }
359:         }
360:       } else {
361:         Obj<Mesh::bundle_type> globalOrder = coordinates->getGlobalOrder();
362:         Obj<Mesh::field_type::order_type::coneSequence> cone = globalOrder->getPatch(patch);
363:         const int *offsets = coordinates->getGlobalOffsets();
364:         int        numLocalVertices = (offsets[mesh->commRank()+1] - offsets[mesh->commRank()])/dim;
365:         double    *localCoords;
366:         int        k = 0;

368:         PetscMalloc(numLocalVertices*dim * sizeof(double), &localCoords);
369:         for(Mesh::field_type::order_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
370:           int dim = globalOrder->getFiberDimension(patch, *p_iter);

372:           if (dim > 0) {
373:             int offset = coordinates->getFiberOffset(patch, *p_iter);

375:             for(int i = offset; i < offset+dim; ++i) {
376:               localCoords[k++] = array[i];
377:             }
378:           }
379:         }
380:         if (k != numLocalVertices*dim) {
381:           SETERRQ2(PETSC_ERR_PLIB, "Invalid number of coordinates to send %d should be %d", k, numLocalVertices*dim);
382:         }
383:         MPI_Send(&numLocalVertices, 1, MPI_INT, 0, 1, mesh->comm());
384:         MPI_Send(localCoords, numLocalVertices*dim, MPI_DOUBLE, 0, 1, mesh->comm());
385:         PetscFree(localCoords);
386:       }
387: #endif
388:       return(0);
389:     };
392:     PetscErrorCode Viewer::writeElements(const Obj<Mesh>& mesh, PetscViewer viewer) {
393:       Obj<Mesh::topology_type> topology = mesh->getTopologyNew();
394: #if 0
395:       Obj<Mesh::sieve_type::traits::heightSequence> elements = topology->heightStratum(0);
396:       Obj<Mesh::bundle_type> elementBundle = mesh->getBundle(topology->depth());
397:       Obj<Mesh::bundle_type> vertexBundle = mesh->getBundle(0);
398:       Obj<Mesh::bundle_type> globalVertex = vertexBundle->getGlobalOrder();
399:       Obj<Mesh::bundle_type> globalElement = elementBundle->getGlobalOrder();
400:       Obj<Mesh::field_type> material;
401:       Mesh::bundle_type::patch_type patch;
402:       std::string    orderName("element");
403:       bool           hasMaterial = mesh->hasField("material");
404:       int            dim  = mesh->getDimension();
405:       int            corners = topology->nCone(*elements->begin(), topology->depth())->size();
406:       int            elementType = -1;

410:       if (dim != 3) {
411:         SETERRQ(PETSC_ERR_SUP, "PyLith only supports 3D meshes.");
412:       }
413:       if (corners == 4) {
414:         // Linear tetrahedron
415:         elementType = 5;
416:       } else if (corners == 8) {
417:         // Linear hexahedron
418:         elementType = 1;
419:       } else {
420:         SETERRQ1(PETSC_ERR_SUP, "PyLith Error: Unsupported number of elements vertices: %d", corners);
421:       }
422:       if (hasMaterial) {
423:         material = mesh->getField("material");
424:       }
425:       PetscViewerASCIIPrintf(viewer,"#\n");
426:       PetscViewerASCIIPrintf(viewer,"#     N ETP MAT INF     N1     N2     N3     N4     N5     N6     N7     N8\n");
427:       PetscViewerASCIIPrintf(viewer,"#\n");
428:       if (mesh->commRank() == 0) {
429:         int elementCount = 1;

431:         for(Mesh::sieve_type::traits::heightSequence::iterator e_itor = elements->begin(); e_itor != elements->end(); ++e_itor) {
432:           Obj<Mesh::bundle_type::order_type::coneSequence> cone = vertexBundle->getPatch(orderName, *e_itor);

434:           PetscViewerASCIIPrintf(viewer, "%7d %3d", elementCount++, elementType);
435:           if (hasMaterial) {
436:             // No infinite elements
437:             PetscViewerASCIIPrintf(viewer, " %3d %3d", (int) material->restrict(patch, *e_itor)[0], 0);
438:           } else {
439:             // No infinite elements
440:             PetscViewerASCIIPrintf(viewer, " %3d %3d", 1, 0);
441:           }
442:           for(Mesh::bundle_type::order_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); ++c_itor) {
443:             PetscViewerASCIIPrintf(viewer, " %6d", globalVertex->getIndex(patch, *c_itor).prefix+1);
444:           }
445:           PetscViewerASCIIPrintf(viewer, "\n");
446:         }
447:         for(int p = 1; p < mesh->commSize(); p++) {
448:           int         numLocalElements;
449:           int        *remoteVertices;
450:           MPI_Status  status;

452:           MPI_Recv(&numLocalElements, 1, MPI_INT, p, 1, mesh->comm(), &status);
453:           PetscMalloc(numLocalElements*(corners+1) * sizeof(int), &remoteVertices);
454:           MPI_Recv(remoteVertices, numLocalElements*(corners+1), MPI_INT, p, 1, mesh->comm(), &status);
455:           for(int e = 0; e < numLocalElements; e++) {
456:             // Only linear tetrahedra, material, no infinite elements
457:             int mat = remoteVertices[e*(corners+1)+corners];

459:             PetscViewerASCIIPrintf(viewer, "%7d %3d %3d %3d", elementCount++, elementType, mat, 0);
460:             for(int c = 0; c < corners; c++) {
461:               PetscViewerASCIIPrintf(viewer, " %6d", remoteVertices[e*(corners+1)+c]);
462:             }
463:             PetscViewerASCIIPrintf(viewer, "\n");
464:           }
465:           PetscFree(remoteVertices);
466:         }
467:       } else {
468:         const int *offsets = elementBundle->getGlobalOffsets();
469:         int        numLocalElements = offsets[mesh->commRank()+1] - offsets[mesh->commRank()];
470:         int       *localVertices;
471:         int        k = 0;

473:         PetscMalloc(numLocalElements*(corners+1) * sizeof(int), &localVertices);
474:         for(Mesh::sieve_type::traits::heightSequence::iterator e_itor = elements->begin(); e_itor != elements->end(); ++e_itor) {
475:           Obj<Mesh::bundle_type::order_type::coneSequence> cone = vertexBundle->getPatch(orderName, *e_itor);

477:           if (globalElement->getFiberDimension(patch, *e_itor) > 0) {
478:             for(Mesh::bundle_type::order_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); ++c_itor) {
479:               localVertices[k++] = globalVertex->getIndex(patch, *c_itor).prefix;
480:             }
481:             if (hasMaterial) {
482:               localVertices[k++] = (int) material->restrict(patch, *e_itor)[0];
483:             } else {
484:               localVertices[k++] = 1;
485:             }
486:           }
487:         }
488:         if (k != numLocalElements*corners) {
489:           SETERRQ2(PETSC_ERR_PLIB, "Invalid number of vertices to send %d should be %d", k, numLocalElements*corners);
490:         }
491:         MPI_Send(&numLocalElements, 1, MPI_INT, 0, 1, mesh->comm());
492:         MPI_Send(localVertices, numLocalElements*(corners+1), MPI_INT, 0, 1, mesh->comm());
493:         PetscFree(localVertices);
494:       }
495: #endif
496:       return(0);
497:     };
500:     PetscErrorCode Viewer::writeVerticesLocal(const Obj<Mesh>& mesh, PetscViewer viewer) {
501:       const Mesh::section_type::patch_type            patch       = 0;
502:       Obj<Mesh::section_type>                         coordinates = mesh->getSection("coordinates");
503:       const Obj<Mesh::topology_type>&                 topology    = mesh->getTopologyNew();
504:       const Obj<Mesh::topology_type::label_sequence>& vertices    = topology->depthStratum(patch, 0);
505:       const Obj<Mesh::numbering_type>&                vNumbering  = mesh->getLocalNumbering(0);
506:       int            embedDim = coordinates->getAtlas()->getFiberDimension(patch, *vertices->begin());

510:       PetscViewerASCIIPrintf(viewer,"#\n");
511:       PetscViewerASCIIPrintf(viewer,"coord_units = m\n");
512:       PetscViewerASCIIPrintf(viewer,"#\n");
513:       PetscViewerASCIIPrintf(viewer,"#\n");
514:       PetscViewerASCIIPrintf(viewer,"#  Node      X-coord           Y-coord           Z-coord\n");
515:       PetscViewerASCIIPrintf(viewer,"#\n");

517:       for(Mesh::topology_type::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
518:         const Mesh::section_type::value_type *array = coordinates->restrict(patch, *v_iter);

520:         PetscViewerASCIIPrintf(viewer, "%7D ", vNumbering->getIndex(*v_iter)+1);
521:         for(int d = 0; d < embedDim; d++) {
522:           if (d > 0) {
523:             PetscViewerASCIIPrintf(viewer, " ");
524:           }
525:           PetscViewerASCIIPrintf(viewer, "% 16.8E", array[d]);
526:         }
527:         PetscViewerASCIIPrintf(viewer, "\n");
528:       }
529:       return(0);
530:     };
533:     PetscErrorCode Viewer::writeElementsLocal(const Obj<Mesh>& mesh, PetscViewer viewer) {
534:       const Mesh::topology_type::patch_type           patch      = 0;
535:       const Obj<Mesh::topology_type>&                 topology   = mesh->getTopologyNew();
536:       const Obj<Mesh::sieve_type>&                    sieve      = topology->getPatch(patch);
537:       const Obj<Mesh::topology_type::label_sequence>& elements   = topology->heightStratum(patch, 0);
538:       const Obj<Mesh::numbering_type>&                eNumbering = mesh->getLocalNumbering(topology->depth());
539:       const Obj<Mesh::numbering_type>&                vNumbering = mesh->getLocalNumbering(0);
540:       Obj<Mesh::section_type>                         material;
541:       int            dim          = mesh->getDimension();
542:       int            corners      = sieve->nCone(*elements->begin(), topology->depth())->size();
543:       bool           hasMaterial  = mesh->hasSection("material");
544:       int            elementType  = -1;

548:       if (dim != 3) {
549:         SETERRQ(PETSC_ERR_SUP, "PyLith only supports 3D meshes.");
550:       }
551:       if (corners == 4) {
552:         // Linear tetrahedron
553:         elementType = 5;
554:       } else if (corners == 8) {
555:         // Linear hexahedron
556:         elementType = 1;
557:       } else {
558:         SETERRQ1(PETSC_ERR_SUP, "PyLith Error: Unsupported number of elements vertices: %d", corners);
559:       }
560:       if (hasMaterial) {
561:         material = mesh->getSection("material");
562:       }
563:       PetscViewerASCIIPrintf(viewer,"#\n");
564:       PetscViewerASCIIPrintf(viewer,"#     N ETP MAT INF     N1     N2     N3     N4     N5     N6     N7     N8\n");
565:       PetscViewerASCIIPrintf(viewer,"#\n");
566:       for(Mesh::topology_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
567:         Obj<Mesh::sieve_type::traits::coneSequence> cone = sieve->cone(*e_iter);

569:         PetscViewerASCIIPrintf(viewer, "%7d %3d", eNumbering->getIndex(*e_iter)+1, elementType);
570:         if (hasMaterial) {
571:           // No infinite elements
572:           PetscViewerASCIIPrintf(viewer, " %3d %3d", (int) material->restrict(patch, *e_iter)[0], 0);
573:         } else {
574:           // No infinite elements
575:           PetscViewerASCIIPrintf(viewer, " %3d %3d", 1, 0);
576:         }
577:         for(Mesh::sieve_type::traits::coneSequence::iterator c_iter = cone->begin(); c_iter != cone->end(); ++c_iter) {
578:           //FIX: Need a global ordering here
579:           PetscViewerASCIIPrintf(viewer, " %6d", vNumbering->getIndex(*c_iter)+1);
580:         }
581:         PetscViewerASCIIPrintf(viewer, "\n");
582:       }
583:       return(0);
584:     };
587:     // The elements seem to be implicitly numbered by appearance, which makes it impossible to
588:     //   number here by bundle, but we can fix it by traversing the elements like the vertices
589:     PetscErrorCode Viewer::writeSplitLocal(const Obj<Mesh>& mesh, const Obj<Builder::split_section_type>& splitField, PetscViewer viewer) {
590:       const Obj<Mesh::numbering_type>&        eNumbering = mesh->getLocalNumbering(mesh->getTopologyNew()->depth());
591:       const Obj<Mesh::numbering_type>&        vNumbering = mesh->getLocalNumbering(0);
592:       Builder::split_section_type::patch_type patch      = 0;

596:       const Builder::split_section_type::atlas_type::chart_type& chart = splitField->getAtlas()->getChart(patch);

598:       for(Mesh::atlas_type::chart_type::const_iterator c_iter = chart.begin(); c_iter != chart.end(); ++c_iter) {
599:         const Builder::split_section_type::point_type& e      = c_iter->first;
600:         const Builder::split_section_type::value_type *values = splitField->restrict(patch, e);
601:         const int                                      size   = splitField->getAtlas()->getFiberDimension(patch, e);

603:         for(int i = 0; i < size; i++) {
604:           const Builder::split_section_type::point_type& v      = values[i].first;
605:           const Builder::split_value&                    split  = values[i].second;

607:           // No time history
608:           PetscViewerASCIIPrintf(viewer, "%6d %6d 0 %15.9g %15.9g %15.9g\n", eNumbering->getIndex(e)+1, vNumbering->getIndex(v)+1, split.x, split.y, split.z);
609:         }
610:       }
611:       return(0);
612:     };
613:   };
614: };