Actual source code: meshpflotran.c

  1: #include <petscmesh_formats.hh>   /*I      "petscmesh.h"   I*/

  3: #if defined(PETSC_HAVE_HDF5)
  4: #include <hdf5.h>
  5: #endif

  7: namespace ALE {
  8:   namespace PFLOTRAN {
  9:     //
 10:     // Builder methods
 11:     //
 12:     void Builder::readConnectivity(MPI_Comm comm, const std::string& filename,
 13:                                    int& corners, const bool useZeroBase,
 14:                                    int& numElements, int *vertices[]) {
 15:       PetscViewer    viewer;
 16:       PetscInt       numCells;
 17:       PetscInt      *verts;
 18:       PetscInt       commRank;
 20: #if defined(PETSC_HAVE_HDF5)
 21:       herr_t         status;
 22:       hid_t          file_id;
 23:       hid_t          group_id;
 24:       hid_t          dataset_id;
 25:       hid_t          prop_id;
 26:       hid_t          type_id;
 27:       hid_t          attribute_id;
 28:       hid_t          string_id;
 29:       H5T_class_t    class_type;
 30:       char           element_type[33];
 31: #endif

 33:       MPI_Comm_rank(comm, &commRank);

 35:       if (commRank != 0) return;

 37: #if defined(PETSC_HAVE_HDF5)
 38:       PetscViewerCreate(PETSC_COMM_SELF, &viewer);
 39:       PetscViewerSetType(viewer, PETSC_VIEWER_HDF5);
 40:       PetscViewerFileSetMode(viewer, FILE_MODE_READ);
 41:       PetscViewerFileSetName(viewer, filename.c_str());
 42:       if (ierr) {
 43:         ostringstream txt;
 44:         txt << "Could not open PFLOTRAN connectivity file: " << filename;
 45:         throw ALE::Exception(txt.str().c_str());
 46:       }
 47:       PetscViewerHDF5GetFileId(viewer, &file_id);

 49:       group_id = H5Gopen(file_id, "Connectivity");

 51: // read in number of cells
 52:       attribute_id = H5Aopen_name(group_id, "Number of Elements");
 53:       type_id = H5Aget_type(attribute_id);
 54:       class_type = H5Tget_class(type_id);
 55:       if (class_type != H5T_INTEGER) {
 56:         PetscPrintf(PETSC_COMM_WORLD,
 57:                     "ERROR: 'Number of Elements' attribute should be an int\n");
 58:       }
 59:       status = H5Tclose(type_id);
 60:       status = H5Aread(attribute_id, H5T_NATIVE_INT, &numCells);
 61:       status = H5Aclose(attribute_id);

 63: // read in cell type

 65:       attribute_id = H5Aopen_name(group_id, "Element Type");
 66:       type_id = H5Aget_type(attribute_id);
 67:       class_type = H5Tget_class(type_id);
 68:       if (class_type != H5T_STRING) {
 69:         // print error message
 70:       }
 71:       size_t strsize = H5Tget_size(type_id);
 72:       status = H5Tclose(type_id);
 73:       if (strsize != 32) {  // right now, 32 is arbitrary
 74:         PetscPrintf(PETSC_COMM_WORLD,
 75:                     "ERROR: Size of attribute string should be 32\n");
 76:       }
 77:       string_id = H5Tcopy(H5T_C_S1);
 78:       status = H5Tset_strpad(string_id, H5T_STR_NULLTERM);
 79:       status = H5Tset_size(string_id, strsize+1);
 80:       status = H5Aread(attribute_id, string_id, element_type);
 81:       status = H5Aclose(attribute_id);


 84:       if (!strncmp(element_type, "tri", 3)) {
 85:         corners = 3;
 86:       }
 87:       else if (!strncmp(element_type, "quad", 4) ||
 88:               !strncmp(element_type, "tet", 3)) {
 89:         corners = 4;
 90:       }
 91:       else if (!strncmp(element_type, "hex", 3)) {
 92:         corners = 8;
 93:       }
 94:       PetscMalloc(numCells*corners * sizeof(PetscInt), &verts);

 96:       dataset_id = H5Dopen(group_id, "Connectivity");
 97:       prop_id = H5Pcreate(H5P_DATASET_XFER);
 98: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
 99:       H5Pset_dxpl_mpio(prop_id, H5FD_MPIO_INDEPENDENT);
100: #endif
101:       H5Dread(dataset_id, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, prop_id, verts);
102:       H5Pclose(prop_id);
103:       H5Dclose(dataset_id);
104:       H5Gclose(group_id);

106:       if (!useZeroBase) {
107:         for (int i=0; i<numCells*corners; i++)
108:           verts[i] -= 1;
109:       }
110:       PetscViewerDestroy(viewer);
111:       numElements = numCells;
112:       *vertices = verts;
113:       PetscPrintf(PETSC_COMM_WORLD,"%d %s elements read.\n",numCells,
114:                   element_type);
115: #else
116:       SETERRABORT(comm,PETSC_ERR_SUP,"PETSc has not been compiled with hdf5 enabled.");
117: #endif
118:     };
119:     void Builder::readCoordinates(MPI_Comm comm, const std::string& filename,
120:                                   const int dim, int& numVertices,
121:                                   double *coordinates[]) {
122:       PetscViewer    viewer;
123:       PetscInt       numVerts;
124:       PetscScalar   *coords, *coord;
125:       PetscInt       c;
126:       PetscInt       commRank;
128: #if defined(PETSC_HAVE_HDF5)
129:       herr_t         status;
130:       hid_t          file_id;
131:       hid_t          group_id;
132:       hid_t          dataset_id;
133:       hid_t          prop_id;
134:       hid_t          type_id;
135:       hid_t          attribute_id;
136:       H5T_class_t    class_type;
137: #endif

139:       MPI_Comm_rank(comm, &commRank);

141:       if (commRank != 0) return;

143: #if defined(PETSC_HAVE_HDF5)
144:       PetscViewerCreate(PETSC_COMM_SELF, &viewer);
145:       PetscViewerSetType(viewer, PETSC_VIEWER_HDF5);
146:       PetscViewerFileSetMode(viewer, FILE_MODE_READ);
147:       PetscViewerFileSetName(viewer, filename.c_str());
148:       if (ierr) {
149:         ostringstream txt;
150:         txt << "Could not open PFLOTRAN connectivity file: " << filename;
151:         throw ALE::Exception(txt.str().c_str());
152:       }
153:       PetscViewerHDF5GetFileId(viewer, &file_id);

155:       group_id = H5Gopen(file_id, "Coordinates");

157: // read in number of vertices
158:       attribute_id = H5Aopen_name(group_id, "Number of Vertices");
159:       type_id = H5Aget_type(attribute_id);
160:       class_type = H5Tget_class(type_id);
161:       if (type_id != H5T_INTEGER) {
162:         // print error message
163:       }
164:       status = H5Tclose(type_id);
165:       status = H5Aread(attribute_id, H5T_NATIVE_INT, &numVerts);
166:       status = H5Aclose(attribute_id);

168:       PetscMalloc(numVerts*dim * sizeof(PetscScalar), &coords);
169:       PetscMalloc(numVerts * sizeof(PetscScalar), &coord);

171:       dataset_id = H5Dopen(group_id, "X-Coordinates");
172:       prop_id = H5Pcreate(H5P_DATASET_XFER);
173: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
174:       H5Pset_dxpl_mpio(prop_id, H5FD_MPIO_INDEPENDENT);
175: #endif
176:       H5Dread(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, prop_id, coord);
177:       H5Pclose(prop_id);
178:       H5Dclose(dataset_id);
179:       c = 0;
180:       for (int i=0; i<numVerts; i++) {
181:         coords[c] = coord[i];
182:         c += dim;
183:       }

185:       if (dim > 1) {
186:         dataset_id = H5Dopen(group_id, "Y-Coordinates");
187:         prop_id = H5Pcreate(H5P_DATASET_XFER);
188: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
189:         H5Pset_dxpl_mpio(prop_id, H5FD_MPIO_INDEPENDENT);
190: #endif
191:         H5Dread(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, prop_id,
192:                 coord);
193:         H5Pclose(prop_id);
194:         H5Dclose(dataset_id);
195:         c = 1;
196:         for (int i=0; i<numVerts; i++) {
197:           coords[c] = coord[i];
198:           c += dim;
199:         }
200:       }

202:       if (dim > 2) {
203:         dataset_id = H5Dopen(group_id, "Z-Coordinates");
204:         prop_id = H5Pcreate(H5P_DATASET_XFER);
205: #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
206:         H5Pset_dxpl_mpio(prop_id, H5FD_MPIO_INDEPENDENT);
207: #endif
208:         H5Dread(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, prop_id,
209:                 coord);
210:         H5Pclose(prop_id);
211:         H5Dclose(dataset_id);
212:         c = 2;
213:         for (int i=0; i<numVerts; i++) {
214:           coords[c] = coord[i];
215:           c += dim;
216:         }
217:       }
218:       H5Gclose(group_id);
219:       PetscFree(coord);

221:       PetscViewerDestroy(viewer);
222:       numVertices = numVerts;
223:       *coordinates = coords;
224:       PetscPrintf(PETSC_COMM_WORLD,"%d vertices read.\n",numVerts);
225: #else
226:       SETERRABORT(comm,PETSC_ERR_SUP,"PETSc has not been compiled with hdf5 enabled.");
227: #endif
228:     };
229:     Obj<ALE::Mesh> Builder::readMesh(MPI_Comm comm, const int dim,
230:                                      const std::string& basename,
231:                                      const bool useZeroBase = true,
232:                                      const bool interpolate = true,
233:                                      const int debug = 0) {
234:       return readMesh(comm, dim, basename+".h5", basename+".h5",
235:                       useZeroBase, interpolate, debug);
236:     };
237:     Obj<ALE::Mesh> Builder::readMesh(MPI_Comm comm, const int dim,
238:                                      const std::string& coordFilename,
239:                                      const std::string& adjFilename,
240:                                      const bool useZeroBase = true,
241:                                      const bool interpolate = true,
242:                                      const int debug = 0) {
243:       Obj<Mesh>          mesh     = new Mesh(comm, dim, debug);
244:       Obj<sieve_type>    sieve    = new sieve_type(comm, debug);
245:       int    *cells = NULL;
246:       double *coordinates = NULL;
247:       int     numCells = 0, numVertices = 0, numCorners = dim+1;

250:       ALE::PFLOTRAN::Builder::readConnectivity(comm, adjFilename, numCorners,
251:                                                useZeroBase, numCells, &cells);
252:       ALE::PFLOTRAN::Builder::readCoordinates(comm, coordFilename, dim,
253:                                               numVertices, &coordinates);
254:       ALE::SieveBuilder<ALE::Mesh>::buildTopology(sieve, dim, numCells, cells,
255:                                                   numVertices, interpolate,
256:                                                   numCorners, -1,
257:                                           mesh->getArrowSection("orientation"));
258:       mesh->setSieve(sieve);
259:       mesh->stratify();
260:       ALE::SieveBuilder<ALE::Mesh>::buildCoordinates(mesh, dim, coordinates);
261:       if (cells) {PetscFree(cells);}
262:       if (coordinates) {PetscFree(coordinates);}
263:       return mesh;
264:     };
265:     // Creates boundary sections:
266:     //   IBC[NBFS,2]:     ALL
267:     //     BL[NBFS,1]:
268:     //     BNVEC[NBFS,2]:
269:     //   BCFUNC[NBCF,NV]: ALL
270:     //   IBNDFS[NBN,2]:   STILL NEED 4-5
271:     //     BNNV[NBN,2]
272:     void Builder::readBoundary(const Obj<Mesh>& mesh, const std::string& bcFilename) {
273:       PetscViewer    viewer;
274:       FILE          *f;
275:       char           buf[2048];

278:       const Obj<Mesh::int_section_type>&  ibc    = mesh->getIntSection("IBC");
279:       const Obj<Mesh::int_section_type>&  ibndfs = mesh->getIntSection("IBNDFS");
280:       const Obj<Mesh::int_section_type>&  ibcnum = mesh->getIntSection("IBCNUM");
281:       const Obj<Mesh::int_section_type>&  ibfcon = mesh->getIntSection("IBFCON");
282:       const Obj<Mesh::real_section_type>& bl     = mesh->getRealSection("BL");
283:       const Obj<Mesh::real_section_type>& bnvec  = mesh->getRealSection("BNVEC");
284:       const Obj<Mesh::real_section_type>& bnnv   = mesh->getRealSection("BNNV");
285:       if (mesh->commRank() != 0) {
286: #if 0
287:         mesh->distributeBCValues();
288: #endif
289:         return;
290:       }
291:       PetscViewerCreate(PETSC_COMM_SELF, &viewer);
292:       PetscViewerSetType(viewer, PETSC_VIEWER_ASCII);
293:       PetscViewerFileSetMode(viewer, FILE_MODE_READ);
294:       PetscViewerFileSetName(viewer, bcFilename.c_str());
295:       PetscViewerASCIIGetPointer(viewer, &f);
296:       // Create IBC section
297:       int  numBdFaces = atoi(strtok(fgets(buf, 2048, f), " "));
298:       int *tmpIBC     = new int[numBdFaces*4];
299:       std::map<int,std::set<int> > elem2Idx;
300:       std::map<int,int> bfReorder;
301:       for(int bf = 0; bf < numBdFaces; bf++) {
302:         const char *x = strtok(fgets(buf, 2048, f), " ");

304:         // Ignore boundary face number
305:         x = strtok(NULL, " ");
306:         tmpIBC[bf*4+0] = atoi(x);
307:         x = strtok(NULL, " ");
308:         tmpIBC[bf*4+1] = atoi(x);
309:         x = strtok(NULL, " ");
310:         tmpIBC[bf*4+2] = atoi(x);
311:         x = strtok(NULL, " ");
312:         tmpIBC[bf*4+3] = atoi(x);
313:         const int elem = tmpIBC[bf*4+0]-1;

315:         ibc->addFiberDimension(elem, 4);
316:         ibcnum->addFiberDimension(elem, 1);
317:         ibfcon->addFiberDimension(elem, 2);
318:         bl->addFiberDimension(elem, 1);
319:         bnvec->addFiberDimension(elem, 2);
320:         elem2Idx[elem].insert(bf);
321:       }
322:       mesh->allocate(ibc);
323:       mesh->allocate(ibcnum);
324:       mesh->allocate(ibfcon);
325:       mesh->allocate(bl);
326:       mesh->allocate(bnvec);
327:       const Mesh::int_section_type::chart_type& chart = ibc->getChart();
328:       int num = 1;

330:       for(Mesh::int_section_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
331:         const int elem = *p_iter;
332:         int bfNum[2];
333:         int k = 0;

335:         for(std::set<int>::const_iterator i_iter = elem2Idx[elem].begin(); i_iter != elem2Idx[elem].end(); ++i_iter) {
336:           bfReorder[(*i_iter)+1] = num;
337:           bfNum[k++] = num;
338:           num++;
339:         }
340:         ibcnum->updatePoint(elem, bfNum);
341:       }
342:       for(int bf = 0; bf < numBdFaces; bf++) {
343:         const int elem = tmpIBC[bf*4]-1;

345:         if (elem2Idx[elem].size() > 1) {
346:           if (*elem2Idx[elem].begin() == bf) {
347:             int values[8];
348:             int k = 0;

350:             for(std::set<int>::const_iterator i_iter = elem2Idx[elem].begin(); i_iter != elem2Idx[elem].end(); ++i_iter) {
351:               for(int v = 0; v < 4; ++v) {
352:                 values[k*4+v] = tmpIBC[*i_iter*4+v];
353:               }
354:               k++;
355:             }
356:             ibc->updatePoint(elem, values);
357:           }
358:         } else {
359:           ibc->updatePoint(elem, &tmpIBC[bf*4]);
360:         }
361:       }
362:       delete [] tmpIBC;
363:       // Create BCFUNC section
364:       int numBcFunc = atoi(strtok(fgets(buf, 2048, f), " "));
365:       if (numBcFunc != 0) {throw ALE::Exception("Cannot handle BCFUNCS after rewrite");}
366:       for(int bc = 0; bc < numBcFunc; bc++) {
367: #if 0
368:         const char *x = strtok(fgets(buf, 2048, f), " ");
369:         Mesh::bc_value_type value;

371:         // Ignore function number
372:         x = strtok(NULL, " ");
373:         value.rho = atof(x);
374:         x = strtok(NULL, " ");
375:         value.u   = atof(x);
376:         x = strtok(NULL, " ");
377:         value.v   = atof(x);
378:         x = strtok(NULL, " ");
379:         value.p   = atof(x);
380:         mesh->setBCValue(bc+1, value);
381: #endif
382:       }
383: #if 0
384:       mesh->distributeBCValues();
385: #endif
386:       // Create IBNDFS section
387:       int       numBdVertices = atoi(strtok(fgets(buf, 2048, f), " "));
388:       const int numElements   = mesh->heightStratum(0)->size();
389:       int      *tmpIBNDFS     = new int[numBdVertices*3];

391:       for(int bv = 0; bv < numBdVertices; bv++) {
392:         const char *x = strtok(fgets(buf, 2048, f), " ");

394:         // Ignore boundary node number
395:         x = strtok(NULL, " ");
396:         tmpIBNDFS[bv*3+0] = atoi(x);
397:         x = strtok(NULL, " ");
398:         tmpIBNDFS[bv*3+1] = atoi(x);
399:         x = strtok(NULL, " ");
400:         tmpIBNDFS[bv*3+2] = atoi(x);
401:         ibndfs->setFiberDimension(tmpIBNDFS[bv*3+0]-1+numElements, 6);
402:       }
403:       mesh->allocate(ibndfs);
404:       for(int bv = 0; bv < numBdVertices; bv++) {
405:         int values[5];

407:         values[0] = tmpIBNDFS[bv*3+0];
408:         // Covert to new boundary face numbers
409:         values[1] = bfReorder[tmpIBNDFS[bv*3+1]];
410:         values[2] = bfReorder[tmpIBNDFS[bv*3+2]];
411:         values[3] = 0;
412:         values[4] = 0;
413:         ibndfs->updatePoint(values[0]-1+numElements, values);
414:       }
415:       PetscViewerDestroy(viewer);
416:       // Create BNNV[NBN,2]
417:       const int dim = mesh->getDimension();

419:       for(int bv = 0; bv < numBdVertices; bv++) {
420:         bnnv->setFiberDimension(tmpIBNDFS[bv*3+0]-1+numElements, dim);
421:       }
422:       mesh->allocate(bnnv);
423:       delete [] tmpIBNDFS;
424:     };
425:     void Builder::outputVerticesLocal(const Obj<Mesh>& mesh, int *numVertices, int *dim, double *coordinates[], const bool columnMajor) {
426:       const Obj<Mesh::real_section_type>& coordSec = mesh->getRealSection("coordinates");
427:       if (!coordSec->size()) {
428:         *numVertices = 0;
429:         *dim         = 0;
430:         *coordinates = NULL;
431:         return;
432:       }
433:       const Obj<Mesh::label_sequence>& vertices   = mesh->depthStratum(0);
434:       const Obj<Mesh::numbering_type>& vNumbering = mesh->getFactory()->getLocalNumbering(mesh, 0);
435:       int            size     = vertices->size();
436:       int            embedDim = coordSec->getFiberDimension(*vertices->begin());
437:       double        *coords;

440:       PetscMalloc(vertices->size()*embedDim * sizeof(double), &coords);
441:       for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
442:         const Mesh::real_section_type::value_type *array = coordSec->restrictPoint(*v_iter);
443:         const int                                  row   = vNumbering->getIndex(*v_iter);

445:         if (columnMajor) {
446:           for(int d = 0; d < embedDim; d++) {
447:             coords[d*size + row] = array[d];
448:           }
449:         } else {
450:           for(int d = 0; d < embedDim; d++) {
451:             coords[row*embedDim + d] = array[d];
452:           }
453:         }
454:       }
455:       *numVertices = size;
456:       *dim         = embedDim;
457:       *coordinates = coords;
458:     };
459:     void Builder::outputElementsLocal(const Obj<Mesh>& mesh, int *numElements, int *numCorners, int *vertices[], const bool columnMajor) {
460:       if (!mesh->heightStratum(0)->size()) {
461:         *numElements = 0;
462:         *numCorners  = 0;
463:         *vertices    = NULL;
464:         return;
465:       }
466:       const Obj<Mesh::sieve_type>&     sieve      = mesh->getSieve();
467:       const Obj<Mesh::label_sequence>& elements   = mesh->heightStratum(0);
468:       const Obj<Mesh::numbering_type>& eNumbering = mesh->getFactory()->getLocalNumbering(mesh, mesh->depth());
469:       const Obj<Mesh::numbering_type>& vNumbering = mesh->getFactory()->getLocalNumbering(mesh, 0);
470:       int            size         = elements->size();
471:       //int            corners      = sieve->nCone(*elements->begin(), topology->depth())->size();
472:       int            corners      = sieve->cone(*elements->begin())->size();
473:       int           *v;

476:       PetscMalloc(elements->size()*corners * sizeof(int), &v);
477:       for(Mesh::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
478:         const Obj<Mesh::sieve_type::traits::coneSequence> cone  = sieve->cone(*e_iter);
479:         Mesh::sieve_type::traits::coneSequence::iterator  begin = cone->begin();
480:         Mesh::sieve_type::traits::coneSequence::iterator  end   = cone->end();

482:         const int row = eNumbering->getIndex(*e_iter);
483:         int       c   = -1;
484:         if (columnMajor) {
485:           for(Mesh::sieve_type::traits::coneSequence::iterator c_iter = begin; c_iter != end; ++c_iter) {
486:             v[(++c)*size + row] = vNumbering->getIndex(*c_iter)+1;
487:           }
488:         } else {
489:           for(Mesh::sieve_type::traits::coneSequence::iterator c_iter = begin; c_iter != end; ++c_iter) {
490:             v[row*corners + ++c] = vNumbering->getIndex(*c_iter)+1;
491:           }
492:         }
493:       }
494:       *numElements = size;
495:       *numCorners  = corners;
496:       *vertices    = v;
497:     };
500:     PetscErrorCode Viewer::writeVertices(const ALE::Obj<Mesh>& mesh, PetscViewer viewer) {
501:       ALE::Obj<Mesh::real_section_type> coordinates = mesh->getRealSection("coordinates");
502: #if 0
503:       Mesh::field_type::patch_type patch;
504:       const double  *array = coordinates->restrict(patch);
505:       int            numVertices;

509:       //FIX:
510:       if (vertexBundle->getGlobalOffsets()) {
511:         numVertices = vertexBundle->getGlobalOffsets()[mesh->commSize()];
512:       } else {
513:         numVertices = mesh->getTopology()->depthStratum(0)->size();
514:       }
515:       PetscViewerASCIIPrintf(viewer, "%D\n", numVertices);
516:       if (mesh->commRank() == 0) {
517:         int numLocalVertices = mesh->getTopology()->depthStratum(0)->size();
518:         int embedDim = coordinates->getFiberDimension(patch, *mesh->getTopology()->depthStratum(0)->begin());
519:         int vertexCount = 1;

521:         for(int v = 0; v < numLocalVertices; v++) {
522:           PetscViewerASCIIPrintf(viewer, "%7D   ", vertexCount++);
523:           for(int d = 0; d < embedDim; d++) {
524:             if (d > 0) {
525:               PetscViewerASCIIPrintf(viewer, " ");
526:             }
527:             PetscViewerASCIIPrintf(viewer, "% 12.5E", array[v*embedDim+d]);
528:           }
529:           PetscViewerASCIIPrintf(viewer, "\n");
530:         }
531:         for(int p = 1; p < mesh->commSize(); p++) {
532:           double    *remoteCoords;
533:           MPI_Status status;

535:           MPI_Recv(&numLocalVertices, 1, MPI_INT, p, 1, mesh->comm(), &status);
536:           PetscMalloc(numLocalVertices*embedDim * sizeof(double), &remoteCoords);
537:           MPI_Recv(remoteCoords, numLocalVertices*embedDim, MPI_DOUBLE, p, 1, mesh->comm(), &status);
538:           for(int v = 0; v < numLocalVertices; v++) {
539:             PetscViewerASCIIPrintf(viewer,"%7D   ", vertexCount++);
540:             for(int d = 0; d < embedDim; d++) {
541:               if (d > 0) {
542:                 PetscViewerASCIIPrintf(viewer, " ");
543:               }
544:               PetscViewerASCIIPrintf(viewer, "% 12.5E", remoteCoords[v*embedDim+d]);
545:             }
546:             PetscViewerASCIIPrintf(viewer, "\n");
547:           }
548:         }
549:       } else {
550:         ALE::Obj<Mesh::bundle_type>                           globalOrder = coordinates->getGlobalOrder();
551:         ALE::Obj<Mesh::bundle_type::order_type::coneSequence> cone        = globalOrder->getPatch(patch);
552:         const int *offsets = coordinates->getGlobalOffsets();
553:         int        embedDim = coordinates->getFiberDimension(patch, *mesh->getTopology()->depthStratum(0)->begin());
554:         int        numLocalVertices = (offsets[mesh->commRank()+1] - offsets[mesh->commRank()])/embedDim;
555:         double    *localCoords;
556:         int        k = 0;

558:         PetscMalloc(numLocalVertices*embedDim * sizeof(double), &localCoords);
559:         for(Mesh::bundle_type::order_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
560:           int dim = globalOrder->getFiberDimension(patch, *p_iter);

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

565:             for(int i = offset; i < offset+dim; ++i) {
566:               localCoords[k++] = array[i];
567:             }
568:           }
569:         }
570:         if (k != numLocalVertices*embedDim) {
571:           SETERRQ2(PETSC_ERR_PLIB, "Invalid number of coordinates to send %d should be %d", k, numLocalVertices*embedDim);
572:         }
573:         MPI_Send(&numLocalVertices, 1, MPI_INT, 0, 1, mesh->comm());
574:         MPI_Send(localCoords, numLocalVertices*embedDim, MPI_DOUBLE, 0, 1, mesh->comm());
575:         PetscFree(localCoords);
576:       }
577: #endif
578:       return(0);
579:     };
582:     PetscErrorCode Viewer::writeElements(const ALE::Obj<Mesh>& mesh, PetscViewer viewer) {
583: #if 0
584:       ALE::Obj<Mesh::sieve_type::traits::heightSequence> elements = topology->heightStratum(0);
585:       ALE::Obj<Mesh::bundle_type> elementBundle = mesh->getBundle(topology->depth());
586:       ALE::Obj<Mesh::bundle_type> vertexBundle = mesh->getBundle(0);
587:       ALE::Obj<Mesh::bundle_type> globalVertex = vertexBundle->getGlobalOrder();
588:       ALE::Obj<Mesh::bundle_type> globalElement = elementBundle->getGlobalOrder();
589:       Mesh::bundle_type::patch_type patch;
590:       std::string    orderName("element");
591:       int            dim  = mesh->getDimension();
592:       int            corners = topology->nCone(*elements->begin(), topology->depth())->size();
593:       int            numElements;

597:       if (corners != dim+1) {
598:         SETERRQ(PETSC_ERR_SUP, "PFLOTRAN only supports simplicies");
599:       }
600:       if (!globalVertex) {
601:         globalVertex = vertexBundle;
602:       }
603:       if (elementBundle->getGlobalOffsets()) {
604:         numElements = elementBundle->getGlobalOffsets()[mesh->commSize()];
605:       } else {
606:         numElements = mesh->getTopology()->heightStratum(0)->size();
607:       }
608:       if (mesh->commRank() == 0) {
609:         int elementCount = 1;

611:         PetscViewerASCIIPrintf(viewer, "%d\n", numElements);
612:         for(Mesh::sieve_type::traits::heightSequence::iterator e_itor = elements->begin(); e_itor != elements->end(); ++e_itor) {
613:           ALE::Obj<Mesh::bundle_type::order_type::coneSequence> cone = vertexBundle->getPatch(orderName, *e_itor);

615:           PetscViewerASCIIPrintf(viewer, "%7d", elementCount++);
616:           for(Mesh::bundle_type::order_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); ++c_itor) {
617:             PetscViewerASCIIPrintf(viewer, " %7d", globalVertex->getIndex(patch, *c_itor).prefix);
618:           }
619:           PetscViewerASCIIPrintf(viewer, "\n");
620:         }
621:         for(int p = 1; p < mesh->commSize(); p++) {
622:           int        numLocalElements;
623:           int       *remoteVertices;
624:           MPI_Status status;

626:           MPI_Recv(&numLocalElements, 1, MPI_INT, p, 1, mesh->comm(), &status);
627:           PetscMalloc(numLocalElements*corners * sizeof(int), &remoteVertices);
628:           MPI_Recv(remoteVertices, numLocalElements*corners, MPI_INT, p, 1, mesh->comm(), &status);
629:           for(int e = 0; e < numLocalElements; e++) {
630:             PetscViewerASCIIPrintf(viewer, "%7d", elementCount++);
631:             for(int c = 0; c < corners; c++) {
632:               PetscViewerASCIIPrintf(viewer, " %7d", remoteVertices[e*corners+c]);
633:             }
634:             PetscViewerASCIIPrintf(viewer, "\n");
635:           }
636:           PetscFree(remoteVertices);
637:         }
638:       } else {
639:         const int *offsets = elementBundle->getGlobalOffsets();
640:         int        numLocalElements = offsets[mesh->commRank()+1] - offsets[mesh->commRank()];
641:         int       *localVertices;
642:         int        k = 0;

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

648:           if (globalElement->getFiberDimension(patch, *e_itor) > 0) {
649:             for(Mesh::bundle_type::order_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); ++c_itor) {
650:               localVertices[k++] = globalVertex->getIndex(patch, *c_itor).prefix;
651:             }
652:           }
653:         }
654:         if (k != numLocalElements*corners) {
655:           SETERRQ2(PETSC_ERR_PLIB, "Invalid number of vertices to send %d should be %d", k, numLocalElements*corners);
656:         }
657:         MPI_Send(&numLocalElements, 1, MPI_INT, 0, 1, mesh->comm());
658:         MPI_Send(localVertices, numLocalElements*corners, MPI_INT, 0, 1, mesh->comm());
659:         PetscFree(localVertices);
660:       }
661: #endif
662:       return(0);
663:     };
666:     PetscErrorCode Viewer::writeVerticesLocal(const Obj<Mesh>& mesh, PetscViewer viewer) {
667:       Obj<Mesh::real_section_type>     coordinates = mesh->getRealSection("coordinates");
668:       const Obj<Mesh::label_sequence>& vertices    = mesh->depthStratum(0);
669:       const Obj<Mesh::numbering_type>& vNumbering  = mesh->getFactory()->getLocalNumbering(mesh, 0);
670:       int            embedDim = coordinates->getFiberDimension(*vertices->begin());

674:       PetscViewerASCIIPrintf(viewer, "%D\n", vertices->size());
675:       for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
676:         const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);

678:         PetscViewerASCIIPrintf(viewer, "%7D   ", vNumbering->getIndex(*v_iter)+1);
679:         for(int d = 0; d < embedDim; d++) {
680:           if (d > 0) {
681:             PetscViewerASCIIPrintf(viewer, " ");
682:           }
683:           PetscViewerASCIIPrintf(viewer, "% 12.5E", array[d]);
684:         }
685:         PetscViewerASCIIPrintf(viewer, "\n");
686:       }
687:       return(0);
688:     };
691:     PetscErrorCode Viewer::writeRestart(const Obj<Mesh>& mesh, PetscViewer viewer) {
692:       const Obj<Mesh::real_section_type>&   velocity    = mesh->getRealSection("VELN");
693:       const Obj<Mesh::real_section_type>&   pressure    = mesh->getRealSection("PN");
694:       const Obj<Mesh::real_section_type>&   temperature = mesh->getRealSection("TN");
695:       const Obj<Mesh::numbering_type>& cNumbering  = mesh->getFactory()->getNumbering(mesh, mesh->depth());
696:       const Obj<Mesh::numbering_type>& vNumbering  = mesh->getFactory()->getNumbering(mesh, 0);
697:       const int                        numCells    = cNumbering->getGlobalSize();

701:       int          blen[2];
702:       MPI_Aint     indices[2];
703:       MPI_Datatype oldtypes[2], newtype;
704:       blen[0] = 1; indices[0] = 0;           oldtypes[0] = MPI_INT;
705:       blen[1] = 4; indices[1] = sizeof(int); oldtypes[1] = MPI_DOUBLE;
706:       MPI_Type_struct(2, blen, indices, oldtypes, &newtype);
707:       MPI_Type_commit(&newtype);

709:       if (mesh->commRank() == 0) {
710:         const Obj<Mesh::label_sequence>& vertices = mesh->depthStratum(0);

712:         for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
713:           if (vNumbering->isLocal(*v_iter)) {
714:             const Mesh::real_section_type::value_type *veln = velocity->restrictPoint(*v_iter);
715:             const Mesh::real_section_type::value_type *pn   = pressure->restrictPoint(*v_iter);
716:             const Mesh::real_section_type::value_type *tn   = temperature->restrictPoint(*v_iter);

718:             PetscViewerASCIIPrintf(viewer, "%6d% 16.8E% 16.8E% 16.8E% 16.8E\n", *v_iter-numCells+1, veln[0], veln[1], pn[0], tn[0]);
719:           }
720:         }
721:         for(int p = 1; p < mesh->commSize(); p++) {
722:           RestartType *remoteValues;
723:           int          numLocalElements;
724:           MPI_Status   status;

726:           MPI_Recv(&numLocalElements, 1, MPI_INT, p, 1, mesh->comm(), &status);
727:           PetscMalloc(numLocalElements * sizeof(RestartType), &remoteValues);
728:           MPI_Recv(remoteValues, numLocalElements, newtype, p, 1, mesh->comm(), &status);
729:           for(int e = 0; e < numLocalElements; e++) {
730:             PetscViewerASCIIPrintf(viewer, "%6d% 16.8E% 16.8E% 16.8E% 16.8E\n", remoteValues[e].vertex-numCells+1, remoteValues[e].veln_x, remoteValues[e].veln_y, remoteValues[e].pn, remoteValues[e].tn);
731:           }
732:         }
733:       } else {
734:         const Obj<Mesh::label_sequence>& vertices = mesh->depthStratum(0);
735:         RestartType *localValues;
736:         int numLocalElements = vNumbering->getLocalSize();
737:         int k = 0;

739:         PetscMalloc(numLocalElements * sizeof(RestartType), &localValues);
740:         for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
741:           if (vNumbering->isLocal(*v_iter)) {
742:             const Mesh::real_section_type::value_type *veln = velocity->restrictPoint(*v_iter);
743:             const Mesh::real_section_type::value_type *pn   = pressure->restrictPoint(*v_iter);
744:             const Mesh::real_section_type::value_type *tn   = temperature->restrictPoint(*v_iter);

746:             localValues[k].vertex = *v_iter;
747:             localValues[k].veln_x = veln[0];
748:             localValues[k].veln_y = veln[1];
749:             localValues[k].pn     = pn[0];
750:             localValues[k].tn     = tn[0];
751:             k++;
752:           }
753:         }
754:         if (k != numLocalElements) {
755:           SETERRQ2(PETSC_ERR_PLIB, "Invalid number of values to send for field, %d should be %d", k, numLocalElements);
756:         }
757:         MPI_Send(&numLocalElements, 1, MPI_INT, 0, 1, mesh->comm());
758:         MPI_Send(localValues, numLocalElements, newtype, 0, 1, mesh->comm());
759:         PetscFree(localValues);
760:       }
761:       MPI_Type_free(&newtype);
762:       return(0);
763:     };

765:     //   This class reconstructs the local pieces of the boundary that distributed PFLOTRAN needs.
766:     // The boundary along with the boundary conditions is encoded in a collection of sections
767:     // over the PFLOTRAN mesh.  These sections contain a description of the boundary topology
768:     // using elements' global names.  This is unacceptable for PFLOTRAN, since it interprets
769:     // elements of the connectivity data arrays as local offsets into (some other of) these arrays.
770:     //   This subroutine performs the renumbering based on the local numbering on the distributed
771:     // mesh.  Essentially, we replace each global element id by its local number.
772:     //
773:     // Note: Any vertex or element number from PFLOTRAN is 1-based, but in Sieve we are 0-based. Thus
774:     //       we add and subtract 1 during conversion. Also, Sieve vertices are numbered after cells.
775:     void fuseBoundary(const ALE::Obj<ALE::Mesh>& mesh) {
776:       // Extract PFLOTRAN boundary sections
777:       ALE::Obj<ALE::Mesh::int_section_type> IBCsec    = mesh->getIntSection("IBC");
778:       ALE::Obj<ALE::Mesh::int_section_type> IBNDFSsec = mesh->getIntSection("IBNDFS");
779:       ALE::Obj<ALE::Mesh::int_section_type> IBCNUMsec = mesh->getIntSection("IBCNUM");

781:       // Look at the sections, if debugging
782:       if (mesh->debug()) {
783:         IBCsec->view("IBCsec", mesh->comm());
784:         IBNDFSsec->view("IBNDFSsec", mesh->comm());
785:       }

787:       // Extract the local numberings
788:       ALE::Obj<ALE::Mesh::numbering_type> vertexNum  = mesh->getFactory()->getLocalNumbering(mesh, 0);
789:       ALE::Obj<ALE::Mesh::numbering_type> elementNum = mesh->getFactory()->getLocalNumbering(mesh, mesh->depth());
790:       const int numElements = mesh->getFactory()->getNumbering(mesh, mesh->depth())->getGlobalSize();
791:       std::map<int,int> bfMap;
792:       // Declare points to the extracted numbering data
793:       const ALE::Mesh::numbering_type::value_type *vNum, *eNum;

795:       // Create map from serial bdFace numbers to local bdFace numbers
796:       {
797:         const ALE::Mesh::int_section_type::chart_type     chart = IBCNUMsec->getChart();
798:         ALE::Mesh::int_section_type::chart_type::iterator begin = chart.begin();
799:         ALE::Mesh::int_section_type::chart_type::iterator end   = chart.end();
800:         int num = 0;

802:         for(ALE::Mesh::int_section_type::chart_type::iterator p_iter = begin; p_iter != end; ++p_iter) {
803:           const int  fiberDim  = IBCNUMsec->getFiberDimension(*p_iter);
804:           const int *globalNum = IBCNUMsec->restrictPoint(*p_iter);

806:           for(int n = 0; n < fiberDim; ++n) {
807:             bfMap[globalNum[n]] = ++num;
808:           }
809:         }
810:       }
811:       // Renumber vertex section IBC
812:       {
813:         const ALE::Mesh::int_section_type::chart_type     IBCchart = IBCsec->getChart();
814:         ALE::Mesh::int_section_type::chart_type::iterator begin    = IBCchart.begin();
815:         ALE::Mesh::int_section_type::chart_type::iterator end      = IBCchart.end();
816:         for(ALE::Mesh::int_section_type::chart_type::iterator p_iter = begin; p_iter != end; ++p_iter) {
817:           ALE::Mesh::point_type p = *p_iter;
818:           const ALE::Mesh::int_section_type::value_type *ibc_in = IBCsec->restrictPoint(p);
819:           int fiberDimension = IBCsec->getFiberDimension(p);
820:           ALE::Mesh::int_section_type::value_type ibc_out[8];
821:           // k controls the update of edge connectivity for each containing element;
822:           // if fiberDimension is 4, only one boundary face is connected to the element, and that edge's data
823:           // are contained in entries 0 - 3 of the section over the element p;
824:           // if fiberDimension is 8, two boundary faces are connected to the element, so the second edge's data
825:           // are contained in entries 4 - 7
826:           for(int k = 0; k < fiberDimension/4; k++) {
827:             // Extract IBC entry 1 (entry kk*4) for edge kk connected to element p.
828:             // This is the entry that needs renumbering for renumbering (2,3 & 4 are invariant under distribution),
829:             // see IBC's description.
830:             // Here we assume that elementNum's domain contains all boundary elements found in IBC,
831:             // so we can restrict to the extracted entry.
832:             eNum = elementNum->restrictPoint((ALE::Mesh::numbering_type::point_type) ibc_in[k*4]-1);
833:             ibc_out[k*4+0] = eNum[0]+1;
834:             // Copy the other entries right over
835:             ibc_out[k*4+1] = ibc_in[k*4+1];
836:             ibc_out[k*4+2] = ibc_in[k*4+2];
837:             ibc_out[k*4+3] = ibc_in[k*4+3];
838:           }
839:           // Update IBC
840:           IBCsec->updatePoint(p, ibc_out);
841:         }
842:       }
843:       {
844:         // Renumber vertex section IBNDFS
845:         const ALE::Mesh::int_section_type::chart_type     IBNDFSchart = IBNDFSsec->getChart();
846:         ALE::Mesh::int_section_type::chart_type::iterator begin       = IBNDFSchart.begin();
847:         ALE::Mesh::int_section_type::chart_type::iterator end         = IBNDFSchart.end();
848:         for(ALE::Mesh::int_section_type::chart_type::iterator p_iter = begin; p_iter != end; ++p_iter) {
849:           ALE::Mesh::point_type p = *p_iter;
850:           const ALE::Mesh::int_section_type::value_type *ibndfs_in = IBNDFSsec->restrictPoint(p);
851:           // Here we assume the fiber dimension is 5
852:           ALE::Mesh::int_section_type::value_type ibndfs_out[5];
853:           // Renumber entries 1,2 & 3 (4 & 5 are invariant under distribution), see IBNDFS's description
854:           // Here we assume that vertexNum's domain contains all boundary verticies found in IBNDFS, so we can restrict to the first extracted entry
855:           vNum= vertexNum->restrictPoint((ALE::Mesh::numbering_type::point_type)ibndfs_in[0]-1+numElements);
856:           ibndfs_out[0] = vNum[0]+1;
857:           // Map serial bdFace numbers to local bdFace numbers
858:           ibndfs_out[1] = bfMap[ibndfs_in[1]];
859:           ibndfs_out[2] = bfMap[ibndfs_in[2]];
860:           // Copy the other entries right over
861:           ibndfs_out[3] = ibndfs_in[3];
862:           ibndfs_out[4] = ibndfs_in[4];
863:           // Update IBNDFS
864:           IBNDFSsec->updatePoint(p,ibndfs_out);
865:         }
866:       }
867:       if (mesh->debug()) {
868:         IBCsec->view("Renumbered IBCsec", mesh->comm());
869:         IBNDFSsec->view("Renumbered IBNDFSsec", mesh->comm());
870:       }
871:       // It's not clear whether IBFCON needs to be renumbered (what does it mean that its entries are not "GLOBAL NODE NUMBER(s)" -- see IBFCON's descriptions
872:     };
873:   };
874: };