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: PetscInt commRank;
17: #if defined(PETSC_HAVE_HDF5)
18: PetscViewer viewer;
19: PetscInt numCells;
20: PetscInt *verts;
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: PetscInt commRank;
124: #if defined(PETSC_HAVE_HDF5)
125: PetscViewer viewer;
126: PetscInt numVerts;
127: PetscScalar *coords, *coord;
128: PetscInt c;
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<PETSC_MESH_TYPE> 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: #ifdef PETSC_OPT_SIEVE
238: Obj<PETSC_MESH_TYPE> Builder::readMesh(MPI_Comm comm, const int dim, const std::string& coordFilename, const std::string& adjFilename, const bool useZeroBase = true, const bool interpolate = true, const int debug = 0) {
239: throw ALE::Exception("Not implemented for optimized sieves");
240: };
241: void Builder::readBoundary(const Obj<Mesh>& mesh, const std::string& bcFilename) {
242: throw ALE::Exception("Not implemented for optimized sieves");
243: };
244: void Builder::outputVerticesLocal(const Obj<Mesh>& mesh, int *numVertices, int *dim, double *coordinates[], const bool columnMajor) {
245: throw ALE::Exception("Not implemented for optimized sieves");
246: };
247: void Builder::outputElementsLocal(const Obj<Mesh>& mesh, int *numElements, int *numCorners, int *vertices[], const bool columnMajor) {
248: throw ALE::Exception("Not implemented for optimized sieves");
249: };
250: PetscErrorCode Viewer::writeVertices(const ALE::Obj<Mesh>& mesh, PetscViewer viewer) {
251: throw ALE::Exception("Not implemented for optimized sieves");
252: };
253: PetscErrorCode Viewer::writeElements(const ALE::Obj<Mesh>& mesh, PetscViewer viewer) {
254: throw ALE::Exception("Not implemented for optimized sieves");
255: };
256: PetscErrorCode Viewer::writeVerticesLocal(const Obj<Mesh>& mesh, PetscViewer viewer) {
257: throw ALE::Exception("Not implemented for optimized sieves");
258: };
259: PetscErrorCode Viewer::writeRestart(const Obj<Mesh>& mesh, PetscViewer viewer) {
260: throw ALE::Exception("Not implemented for optimized sieves");
261: };
262: void fuseBoundary(const ALE::Obj<PETSC_MESH_TYPE>& mesh) {
263: throw ALE::Exception("Not implemented for optimized sieves");
264: };
265: #else
266: Obj<PETSC_MESH_TYPE> Builder::readMesh(MPI_Comm comm, const int dim,
267: const std::string& coordFilename,
268: const std::string& adjFilename,
269: const bool useZeroBase = true,
270: const bool interpolate = true,
271: const int debug = 0) {
272: Obj<Mesh> mesh = new Mesh(comm, dim, debug);
273: Obj<sieve_type> sieve = new sieve_type(comm, debug);
274: int *cells = NULL;
275: double *coordinates = NULL;
276: int numCells = 0, numVertices = 0, numCorners = dim+1;
279: ALE::PFLOTRAN::Builder::readConnectivity(comm, adjFilename, numCorners,
280: useZeroBase, numCells, &cells);
281: ALE::PFLOTRAN::Builder::readCoordinates(comm, coordFilename, dim,
282: numVertices, &coordinates);
283: ALE::SieveBuilder<PETSC_MESH_TYPE>::buildTopology(sieve, dim, numCells, cells,
284: numVertices, interpolate,
285: numCorners, -1,
286: mesh->getArrowSection("orientation"));
287: mesh->setSieve(sieve);
288: mesh->stratify();
289: ALE::SieveBuilder<PETSC_MESH_TYPE>::buildCoordinates(mesh, dim, coordinates);
290: if (cells) {PetscFree(cells);}
291: if (coordinates) {PetscFree(coordinates);}
292: return mesh;
293: };
294: // Creates boundary sections:
295: // IBC[NBFS,2]: ALL
296: // BL[NBFS,1]:
297: // BNVEC[NBFS,2]:
298: // BCFUNC[NBCF,NV]: ALL
299: // IBNDFS[NBN,2]: STILL NEED 4-5
300: // BNNV[NBN,2]
301: void Builder::readBoundary(const Obj<Mesh>& mesh, const std::string& bcFilename) {
302: PetscViewer viewer;
303: FILE *f;
304: char buf[2048];
307: const Obj<Mesh::int_section_type>& ibc = mesh->getIntSection("IBC");
308: const Obj<Mesh::int_section_type>& ibndfs = mesh->getIntSection("IBNDFS");
309: const Obj<Mesh::int_section_type>& ibcnum = mesh->getIntSection("IBCNUM");
310: const Obj<Mesh::int_section_type>& ibfcon = mesh->getIntSection("IBFCON");
311: const Obj<Mesh::real_section_type>& bl = mesh->getRealSection("BL");
312: const Obj<Mesh::real_section_type>& bnvec = mesh->getRealSection("BNVEC");
313: const Obj<Mesh::real_section_type>& bnnv = mesh->getRealSection("BNNV");
314: if (mesh->commRank() != 0) {
315: #if 0
316: mesh->distributeBCValues();
317: #endif
318: return;
319: }
320: PetscViewerCreate(PETSC_COMM_SELF, &viewer);
321: PetscViewerSetType(viewer, PETSC_VIEWER_ASCII);
322: PetscViewerFileSetMode(viewer, FILE_MODE_READ);
323: PetscViewerFileSetName(viewer, bcFilename.c_str());
324: PetscViewerASCIIGetPointer(viewer, &f);
325: // Create IBC section
326: int numBdFaces = atoi(strtok(fgets(buf, 2048, f), " "));
327: int *tmpIBC = new int[numBdFaces*4];
328: std::map<int,std::set<int> > elem2Idx;
329: std::map<int,int> bfReorder;
330: for(int bf = 0; bf < numBdFaces; bf++) {
331: const char *x = strtok(fgets(buf, 2048, f), " ");
333: // Ignore boundary face number
334: x = strtok(NULL, " ");
335: tmpIBC[bf*4+0] = atoi(x);
336: x = strtok(NULL, " ");
337: tmpIBC[bf*4+1] = atoi(x);
338: x = strtok(NULL, " ");
339: tmpIBC[bf*4+2] = atoi(x);
340: x = strtok(NULL, " ");
341: tmpIBC[bf*4+3] = atoi(x);
342: const int elem = tmpIBC[bf*4+0]-1;
344: ibc->addFiberDimension(elem, 4);
345: ibcnum->addFiberDimension(elem, 1);
346: ibfcon->addFiberDimension(elem, 2);
347: bl->addFiberDimension(elem, 1);
348: bnvec->addFiberDimension(elem, 2);
349: elem2Idx[elem].insert(bf);
350: }
351: mesh->allocate(ibc);
352: mesh->allocate(ibcnum);
353: mesh->allocate(ibfcon);
354: mesh->allocate(bl);
355: mesh->allocate(bnvec);
356: const Mesh::int_section_type::chart_type& chart = ibc->getChart();
357: int num = 1;
359: for(Mesh::int_section_type::chart_type::const_iterator p_iter = chart.begin(); p_iter != chart.end(); ++p_iter) {
360: const int elem = *p_iter;
361: int bfNum[2];
362: int k = 0;
364: for(std::set<int>::const_iterator i_iter = elem2Idx[elem].begin(); i_iter != elem2Idx[elem].end(); ++i_iter) {
365: bfReorder[(*i_iter)+1] = num;
366: bfNum[k++] = num;
367: num++;
368: }
369: ibcnum->updatePoint(elem, bfNum);
370: }
371: for(int bf = 0; bf < numBdFaces; bf++) {
372: const int elem = tmpIBC[bf*4]-1;
374: if (elem2Idx[elem].size() > 1) {
375: if (*elem2Idx[elem].begin() == bf) {
376: int values[8];
377: int k = 0;
379: for(std::set<int>::const_iterator i_iter = elem2Idx[elem].begin(); i_iter != elem2Idx[elem].end(); ++i_iter) {
380: for(int v = 0; v < 4; ++v) {
381: values[k*4+v] = tmpIBC[*i_iter*4+v];
382: }
383: k++;
384: }
385: ibc->updatePoint(elem, values);
386: }
387: } else {
388: ibc->updatePoint(elem, &tmpIBC[bf*4]);
389: }
390: }
391: delete [] tmpIBC;
392: // Create BCFUNC section
393: int numBcFunc = atoi(strtok(fgets(buf, 2048, f), " "));
394: if (numBcFunc != 0) {throw ALE::Exception("Cannot handle BCFUNCS after rewrite");}
395: for(int bc = 0; bc < numBcFunc; bc++) {
396: #if 0
397: const char *x = strtok(fgets(buf, 2048, f), " ");
398: Mesh::bc_value_type value;
400: // Ignore function number
401: x = strtok(NULL, " ");
402: value.rho = atof(x);
403: x = strtok(NULL, " ");
404: value.u = atof(x);
405: x = strtok(NULL, " ");
406: value.v = atof(x);
407: x = strtok(NULL, " ");
408: value.p = atof(x);
409: mesh->setBCValue(bc+1, value);
410: #endif
411: }
412: #if 0
413: mesh->distributeBCValues();
414: #endif
415: // Create IBNDFS section
416: int numBdVertices = atoi(strtok(fgets(buf, 2048, f), " "));
417: const int numElements = mesh->heightStratum(0)->size();
418: int *tmpIBNDFS = new int[numBdVertices*3];
420: for(int bv = 0; bv < numBdVertices; bv++) {
421: const char *x = strtok(fgets(buf, 2048, f), " ");
423: // Ignore boundary node number
424: x = strtok(NULL, " ");
425: tmpIBNDFS[bv*3+0] = atoi(x);
426: x = strtok(NULL, " ");
427: tmpIBNDFS[bv*3+1] = atoi(x);
428: x = strtok(NULL, " ");
429: tmpIBNDFS[bv*3+2] = atoi(x);
430: ibndfs->setFiberDimension(tmpIBNDFS[bv*3+0]-1+numElements, 6);
431: }
432: mesh->allocate(ibndfs);
433: for(int bv = 0; bv < numBdVertices; bv++) {
434: int values[5];
436: values[0] = tmpIBNDFS[bv*3+0];
437: // Covert to new boundary face numbers
438: values[1] = bfReorder[tmpIBNDFS[bv*3+1]];
439: values[2] = bfReorder[tmpIBNDFS[bv*3+2]];
440: values[3] = 0;
441: values[4] = 0;
442: ibndfs->updatePoint(values[0]-1+numElements, values);
443: }
444: PetscViewerDestroy(viewer);
445: // Create BNNV[NBN,2]
446: const int dim = mesh->getDimension();
448: for(int bv = 0; bv < numBdVertices; bv++) {
449: bnnv->setFiberDimension(tmpIBNDFS[bv*3+0]-1+numElements, dim);
450: }
451: mesh->allocate(bnnv);
452: delete [] tmpIBNDFS;
453: };
454: void Builder::outputVerticesLocal(const Obj<Mesh>& mesh, int *numVertices, int *dim, double *coordinates[], const bool columnMajor) {
455: const Obj<Mesh::real_section_type>& coordSec = mesh->getRealSection("coordinates");
456: if (!coordSec->size()) {
457: *numVertices = 0;
458: *dim = 0;
459: *coordinates = NULL;
460: return;
461: }
462: const Obj<Mesh::label_sequence>& vertices = mesh->depthStratum(0);
463: const Obj<Mesh::numbering_type>& vNumbering = mesh->getFactory()->getLocalNumbering(mesh, 0);
464: int size = vertices->size();
465: int embedDim = coordSec->getFiberDimension(*vertices->begin());
466: double *coords;
469: PetscMalloc(vertices->size()*embedDim * sizeof(double), &coords);
470: for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
471: const Mesh::real_section_type::value_type *array = coordSec->restrictPoint(*v_iter);
472: const int row = vNumbering->getIndex(*v_iter);
474: if (columnMajor) {
475: for(int d = 0; d < embedDim; d++) {
476: coords[d*size + row] = array[d];
477: }
478: } else {
479: for(int d = 0; d < embedDim; d++) {
480: coords[row*embedDim + d] = array[d];
481: }
482: }
483: }
484: *numVertices = size;
485: *dim = embedDim;
486: *coordinates = coords;
487: };
488: void Builder::outputElementsLocal(const Obj<Mesh>& mesh, int *numElements, int *numCorners, int *vertices[], const bool columnMajor) {
489: if (!mesh->heightStratum(0)->size()) {
490: *numElements = 0;
491: *numCorners = 0;
492: *vertices = NULL;
493: return;
494: }
495: const Obj<Mesh::sieve_type>& sieve = mesh->getSieve();
496: const Obj<Mesh::label_sequence>& elements = mesh->heightStratum(0);
497: const Obj<Mesh::numbering_type>& eNumbering = mesh->getFactory()->getLocalNumbering(mesh, mesh->depth());
498: const Obj<Mesh::numbering_type>& vNumbering = mesh->getFactory()->getLocalNumbering(mesh, 0);
499: int size = elements->size();
500: //int corners = sieve->nCone(*elements->begin(), topology->depth())->size();
501: int corners = sieve->cone(*elements->begin())->size();
502: int *v;
505: PetscMalloc(elements->size()*corners * sizeof(int), &v);
506: for(Mesh::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
507: const Obj<Mesh::sieve_type::traits::coneSequence> cone = sieve->cone(*e_iter);
508: Mesh::sieve_type::traits::coneSequence::iterator begin = cone->begin();
509: Mesh::sieve_type::traits::coneSequence::iterator end = cone->end();
511: const int row = eNumbering->getIndex(*e_iter);
512: int c = -1;
513: if (columnMajor) {
514: for(Mesh::sieve_type::traits::coneSequence::iterator c_iter = begin; c_iter != end; ++c_iter) {
515: v[(++c)*size + row] = vNumbering->getIndex(*c_iter)+1;
516: }
517: } else {
518: for(Mesh::sieve_type::traits::coneSequence::iterator c_iter = begin; c_iter != end; ++c_iter) {
519: v[row*corners + ++c] = vNumbering->getIndex(*c_iter)+1;
520: }
521: }
522: }
523: *numElements = size;
524: *numCorners = corners;
525: *vertices = v;
526: };
529: PetscErrorCode Viewer::writeVertices(const ALE::Obj<Mesh>& mesh, PetscViewer viewer) {
530: ALE::Obj<Mesh::real_section_type> coordinates = mesh->getRealSection("coordinates");
531: #if 0
532: Mesh::field_type::patch_type patch;
533: const double *array = coordinates->restrict(patch);
534: int numVertices;
538: //FIX:
539: if (vertexBundle->getGlobalOffsets()) {
540: numVertices = vertexBundle->getGlobalOffsets()[mesh->commSize()];
541: } else {
542: numVertices = mesh->getTopology()->depthStratum(0)->size();
543: }
544: PetscViewerASCIIPrintf(viewer, "%D\n", numVertices);
545: if (mesh->commRank() == 0) {
546: int numLocalVertices = mesh->getTopology()->depthStratum(0)->size();
547: int embedDim = coordinates->getFiberDimension(patch, *mesh->getTopology()->depthStratum(0)->begin());
548: int vertexCount = 1;
550: for(int v = 0; v < numLocalVertices; v++) {
551: PetscViewerASCIIPrintf(viewer, "%7D ", vertexCount++);
552: for(int d = 0; d < embedDim; d++) {
553: if (d > 0) {
554: PetscViewerASCIIPrintf(viewer, " ");
555: }
556: PetscViewerASCIIPrintf(viewer, "% 12.5E", array[v*embedDim+d]);
557: }
558: PetscViewerASCIIPrintf(viewer, "\n");
559: }
560: for(int p = 1; p < mesh->commSize(); p++) {
561: double *remoteCoords;
562: MPI_Status status;
564: MPI_Recv(&numLocalVertices, 1, MPI_INT, p, 1, mesh->comm(), &status);
565: PetscMalloc(numLocalVertices*embedDim * sizeof(double), &remoteCoords);
566: MPI_Recv(remoteCoords, numLocalVertices*embedDim, MPI_DOUBLE, p, 1, mesh->comm(), &status);
567: for(int v = 0; v < numLocalVertices; v++) {
568: PetscViewerASCIIPrintf(viewer,"%7D ", vertexCount++);
569: for(int d = 0; d < embedDim; d++) {
570: if (d > 0) {
571: PetscViewerASCIIPrintf(viewer, " ");
572: }
573: PetscViewerASCIIPrintf(viewer, "% 12.5E", remoteCoords[v*embedDim+d]);
574: }
575: PetscViewerASCIIPrintf(viewer, "\n");
576: }
577: }
578: } else {
579: ALE::Obj<Mesh::bundle_type> globalOrder = coordinates->getGlobalOrder();
580: ALE::Obj<Mesh::bundle_type::order_type::coneSequence> cone = globalOrder->getPatch(patch);
581: const int *offsets = coordinates->getGlobalOffsets();
582: int embedDim = coordinates->getFiberDimension(patch, *mesh->getTopology()->depthStratum(0)->begin());
583: int numLocalVertices = (offsets[mesh->commRank()+1] - offsets[mesh->commRank()])/embedDim;
584: double *localCoords;
585: int k = 0;
587: PetscMalloc(numLocalVertices*embedDim * sizeof(double), &localCoords);
588: for(Mesh::bundle_type::order_type::coneSequence::iterator p_iter = cone->begin(); p_iter != cone->end(); ++p_iter) {
589: int dim = globalOrder->getFiberDimension(patch, *p_iter);
591: if (dim > 0) {
592: int offset = coordinates->getFiberOffset(patch, *p_iter);
594: for(int i = offset; i < offset+dim; ++i) {
595: localCoords[k++] = array[i];
596: }
597: }
598: }
599: if (k != numLocalVertices*embedDim) {
600: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of coordinates to send %d should be %d", k, numLocalVertices*embedDim);
601: }
602: MPI_Send(&numLocalVertices, 1, MPI_INT, 0, 1, mesh->comm());
603: MPI_Send(localCoords, numLocalVertices*embedDim, MPI_DOUBLE, 0, 1, mesh->comm());
604: PetscFree(localCoords);
605: }
606: #endif
607: return(0);
608: };
611: PetscErrorCode Viewer::writeElements(const ALE::Obj<Mesh>& mesh, PetscViewer viewer) {
612: #if 0
613: ALE::Obj<Mesh::sieve_type::traits::heightSequence> elements = topology->heightStratum(0);
614: ALE::Obj<Mesh::bundle_type> elementBundle = mesh->getBundle(topology->depth());
615: ALE::Obj<Mesh::bundle_type> vertexBundle = mesh->getBundle(0);
616: ALE::Obj<Mesh::bundle_type> globalVertex = vertexBundle->getGlobalOrder();
617: ALE::Obj<Mesh::bundle_type> globalElement = elementBundle->getGlobalOrder();
618: Mesh::bundle_type::patch_type patch;
619: std::string orderName("element");
620: int dim = mesh->getDimension();
621: int corners = topology->nCone(*elements->begin(), topology->depth())->size();
622: int numElements;
626: if (corners != dim+1) {
627: SETERRQ(PETSC_ERR_SUP, "PFLOTRAN only supports simplicies");
628: }
629: if (!globalVertex) {
630: globalVertex = vertexBundle;
631: }
632: if (elementBundle->getGlobalOffsets()) {
633: numElements = elementBundle->getGlobalOffsets()[mesh->commSize()];
634: } else {
635: numElements = mesh->getTopology()->heightStratum(0)->size();
636: }
637: if (mesh->commRank() == 0) {
638: int elementCount = 1;
640: PetscViewerASCIIPrintf(viewer, "%d\n", numElements);
641: for(Mesh::sieve_type::traits::heightSequence::iterator e_itor = elements->begin(); e_itor != elements->end(); ++e_itor) {
642: ALE::Obj<Mesh::bundle_type::order_type::coneSequence> cone = vertexBundle->getPatch(orderName, *e_itor);
644: PetscViewerASCIIPrintf(viewer, "%7d", elementCount++);
645: for(Mesh::bundle_type::order_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); ++c_itor) {
646: PetscViewerASCIIPrintf(viewer, " %7d", globalVertex->getIndex(patch, *c_itor).prefix);
647: }
648: PetscViewerASCIIPrintf(viewer, "\n");
649: }
650: for(int p = 1; p < mesh->commSize(); p++) {
651: int numLocalElements;
652: int *remoteVertices;
653: MPI_Status status;
655: MPI_Recv(&numLocalElements, 1, MPI_INT, p, 1, mesh->comm(), &status);
656: PetscMalloc(numLocalElements*corners * sizeof(int), &remoteVertices);
657: MPI_Recv(remoteVertices, numLocalElements*corners, MPI_INT, p, 1, mesh->comm(), &status);
658: for(int e = 0; e < numLocalElements; e++) {
659: PetscViewerASCIIPrintf(viewer, "%7d", elementCount++);
660: for(int c = 0; c < corners; c++) {
661: PetscViewerASCIIPrintf(viewer, " %7d", remoteVertices[e*corners+c]);
662: }
663: PetscViewerASCIIPrintf(viewer, "\n");
664: }
665: PetscFree(remoteVertices);
666: }
667: } else {
668: const int *offsets = elementBundle->getGlobalOffsets();
669: int numLocalElements = offsets[mesh->commRank()+1] - offsets[mesh->commRank()];
670: int *localVertices;
671: int k = 0;
673: PetscMalloc(numLocalElements*corners * sizeof(int), &localVertices);
674: for(Mesh::sieve_type::traits::heightSequence::iterator e_itor = elements->begin(); e_itor != elements->end(); ++e_itor) {
675: ALE::Obj<Mesh::bundle_type::order_type::coneSequence> cone = vertexBundle->getPatch(orderName, *e_itor);
677: if (globalElement->getFiberDimension(patch, *e_itor) > 0) {
678: for(Mesh::bundle_type::order_type::coneSequence::iterator c_itor = cone->begin(); c_itor != cone->end(); ++c_itor) {
679: localVertices[k++] = globalVertex->getIndex(patch, *c_itor).prefix;
680: }
681: }
682: }
683: if (k != numLocalElements*corners) {
684: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of vertices to send %d should be %d", k, numLocalElements*corners);
685: }
686: MPI_Send(&numLocalElements, 1, MPI_INT, 0, 1, mesh->comm());
687: MPI_Send(localVertices, numLocalElements*corners, MPI_INT, 0, 1, mesh->comm());
688: PetscFree(localVertices);
689: }
690: #endif
691: return(0);
692: };
695: PetscErrorCode Viewer::writeVerticesLocal(const Obj<Mesh>& mesh, PetscViewer viewer) {
696: Obj<Mesh::real_section_type> coordinates = mesh->getRealSection("coordinates");
697: const Obj<Mesh::label_sequence>& vertices = mesh->depthStratum(0);
698: const Obj<Mesh::numbering_type>& vNumbering = mesh->getFactory()->getLocalNumbering(mesh, 0);
699: int embedDim = coordinates->getFiberDimension(*vertices->begin());
703: PetscViewerASCIIPrintf(viewer, "%D\n", vertices->size());
704: for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
705: const Mesh::real_section_type::value_type *array = coordinates->restrictPoint(*v_iter);
707: PetscViewerASCIIPrintf(viewer, "%7D ", vNumbering->getIndex(*v_iter)+1);
708: for(int d = 0; d < embedDim; d++) {
709: if (d > 0) {
710: PetscViewerASCIIPrintf(viewer, " ");
711: }
712: PetscViewerASCIIPrintf(viewer, "% 12.5E", array[d]);
713: }
714: PetscViewerASCIIPrintf(viewer, "\n");
715: }
716: return(0);
717: };
720: PetscErrorCode Viewer::writeRestart(const Obj<Mesh>& mesh, PetscViewer viewer) {
721: const Obj<Mesh::real_section_type>& velocity = mesh->getRealSection("VELN");
722: const Obj<Mesh::real_section_type>& pressure = mesh->getRealSection("PN");
723: const Obj<Mesh::real_section_type>& temperature = mesh->getRealSection("TN");
724: const Obj<Mesh::numbering_type>& cNumbering = mesh->getFactory()->getNumbering(mesh, mesh->depth());
725: const Obj<Mesh::numbering_type>& vNumbering = mesh->getFactory()->getNumbering(mesh, 0);
726: const int numCells = cNumbering->getGlobalSize();
730: int blen[2];
731: MPI_Aint indices[2];
732: MPI_Datatype oldtypes[2], newtype;
733: blen[0] = 1; indices[0] = 0; oldtypes[0] = MPI_INT;
734: blen[1] = 4; indices[1] = sizeof(int); oldtypes[1] = MPI_DOUBLE;
735: MPI_Type_struct(2, blen, indices, oldtypes, &newtype);
736: MPI_Type_commit(&newtype);
738: if (mesh->commRank() == 0) {
739: const Obj<Mesh::label_sequence>& vertices = mesh->depthStratum(0);
741: for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
742: if (vNumbering->isLocal(*v_iter)) {
743: const Mesh::real_section_type::value_type *veln = velocity->restrictPoint(*v_iter);
744: const Mesh::real_section_type::value_type *pn = pressure->restrictPoint(*v_iter);
745: const Mesh::real_section_type::value_type *tn = temperature->restrictPoint(*v_iter);
747: PetscViewerASCIIPrintf(viewer, "%6d% 16.8E% 16.8E% 16.8E% 16.8E\n", *v_iter-numCells+1, veln[0], veln[1], pn[0], tn[0]);
748: }
749: }
750: for(int p = 1; p < mesh->commSize(); p++) {
751: RestartType *remoteValues;
752: int numLocalElements;
753: MPI_Status status;
755: MPI_Recv(&numLocalElements, 1, MPI_INT, p, 1, mesh->comm(), &status);
756: PetscMalloc(numLocalElements * sizeof(RestartType), &remoteValues);
757: MPI_Recv(remoteValues, numLocalElements, newtype, p, 1, mesh->comm(), &status);
758: for(int e = 0; e < numLocalElements; e++) {
759: 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);
760: }
761: }
762: } else {
763: const Obj<Mesh::label_sequence>& vertices = mesh->depthStratum(0);
764: RestartType *localValues;
765: int numLocalElements = vNumbering->getLocalSize();
766: int k = 0;
768: PetscMalloc(numLocalElements * sizeof(RestartType), &localValues);
769: for(Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
770: if (vNumbering->isLocal(*v_iter)) {
771: const Mesh::real_section_type::value_type *veln = velocity->restrictPoint(*v_iter);
772: const Mesh::real_section_type::value_type *pn = pressure->restrictPoint(*v_iter);
773: const Mesh::real_section_type::value_type *tn = temperature->restrictPoint(*v_iter);
775: localValues[k].vertex = *v_iter;
776: localValues[k].veln_x = veln[0];
777: localValues[k].veln_y = veln[1];
778: localValues[k].pn = pn[0];
779: localValues[k].tn = tn[0];
780: k++;
781: }
782: }
783: if (k != numLocalElements) {
784: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of values to send for field, %d should be %d", k, numLocalElements);
785: }
786: MPI_Send(&numLocalElements, 1, MPI_INT, 0, 1, mesh->comm());
787: MPI_Send(localValues, numLocalElements, newtype, 0, 1, mesh->comm());
788: PetscFree(localValues);
789: }
790: MPI_Type_free(&newtype);
791: return(0);
792: };
794: // This class reconstructs the local pieces of the boundary that distributed PFLOTRAN needs.
795: // The boundary along with the boundary conditions is encoded in a collection of sections
796: // over the PFLOTRAN mesh. These sections contain a description of the boundary topology
797: // using elements' global names. This is unacceptable for PFLOTRAN, since it interprets
798: // elements of the connectivity data arrays as local offsets into (some other of) these arrays.
799: // This subroutine performs the renumbering based on the local numbering on the distributed
800: // mesh. Essentially, we replace each global element id by its local number.
801: //
802: // Note: Any vertex or element number from PFLOTRAN is 1-based, but in Sieve we are 0-based. Thus
803: // we add and subtract 1 during conversion. Also, Sieve vertices are numbered after cells.
804: void fuseBoundary(const ALE::Obj<PETSC_MESH_TYPE>& mesh) {
805: // Extract PFLOTRAN boundary sections
806: ALE::Obj<PETSC_MESH_TYPE::int_section_type> IBCsec = mesh->getIntSection("IBC");
807: ALE::Obj<PETSC_MESH_TYPE::int_section_type> IBNDFSsec = mesh->getIntSection("IBNDFS");
808: ALE::Obj<PETSC_MESH_TYPE::int_section_type> IBCNUMsec = mesh->getIntSection("IBCNUM");
810: // Look at the sections, if debugging
811: if (mesh->debug()) {
812: IBCsec->view("IBCsec", mesh->comm());
813: IBNDFSsec->view("IBNDFSsec", mesh->comm());
814: }
816: // Extract the local numberings
817: ALE::Obj<PETSC_MESH_TYPE::numbering_type> vertexNum = mesh->getFactory()->getLocalNumbering(mesh, 0);
818: ALE::Obj<PETSC_MESH_TYPE::numbering_type> elementNum = mesh->getFactory()->getLocalNumbering(mesh, mesh->depth());
819: const int numElements = mesh->getFactory()->getNumbering(mesh, mesh->depth())->getGlobalSize();
820: std::map<int,int> bfMap;
821: // Declare points to the extracted numbering data
822: const PETSC_MESH_TYPE::numbering_type::value_type *vNum, *eNum;
824: // Create map from serial bdFace numbers to local bdFace numbers
825: {
826: const PETSC_MESH_TYPE::int_section_type::chart_type chart = IBCNUMsec->getChart();
827: PETSC_MESH_TYPE::int_section_type::chart_type::const_iterator begin = chart.begin();
828: PETSC_MESH_TYPE::int_section_type::chart_type::const_iterator end = chart.end();
829: int num = 0;
831: for(PETSC_MESH_TYPE::int_section_type::chart_type::const_iterator p_iter = begin; p_iter != end; ++p_iter) {
832: const int fiberDim = IBCNUMsec->getFiberDimension(*p_iter);
833: const int *globalNum = IBCNUMsec->restrictPoint(*p_iter);
835: for(int n = 0; n < fiberDim; ++n) {
836: bfMap[globalNum[n]] = ++num;
837: }
838: }
839: }
840: // Renumber vertex section IBC
841: {
842: const PETSC_MESH_TYPE::int_section_type::chart_type IBCchart = IBCsec->getChart();
843: PETSC_MESH_TYPE::int_section_type::chart_type::const_iterator begin = IBCchart.begin();
844: PETSC_MESH_TYPE::int_section_type::chart_type::const_iterator end = IBCchart.end();
845: for(PETSC_MESH_TYPE::int_section_type::chart_type::const_iterator p_iter = begin; p_iter != end; ++p_iter) {
846: PETSC_MESH_TYPE::point_type p = *p_iter;
847: const PETSC_MESH_TYPE::int_section_type::value_type *ibc_in = IBCsec->restrictPoint(p);
848: int fiberDimension = IBCsec->getFiberDimension(p);
849: PETSC_MESH_TYPE::int_section_type::value_type ibc_out[8];
850: // k controls the update of edge connectivity for each containing element;
851: // if fiberDimension is 4, only one boundary face is connected to the element, and that edge's data
852: // are contained in entries 0 - 3 of the section over the element p;
853: // if fiberDimension is 8, two boundary faces are connected to the element, so the second edge's data
854: // are contained in entries 4 - 7
855: for(int k = 0; k < fiberDimension/4; k++) {
856: // Extract IBC entry 1 (entry kk*4) for edge kk connected to element p.
857: // This is the entry that needs renumbering for renumbering (2,3 & 4 are invariant under distribution),
858: // see IBC's description.
859: // Here we assume that elementNum's domain contains all boundary elements found in IBC,
860: // so we can restrict to the extracted entry.
861: eNum = elementNum->restrictPoint((PETSC_MESH_TYPE::numbering_type::point_type) ibc_in[k*4]-1);
862: ibc_out[k*4+0] = eNum[0]+1;
863: // Copy the other entries right over
864: ibc_out[k*4+1] = ibc_in[k*4+1];
865: ibc_out[k*4+2] = ibc_in[k*4+2];
866: ibc_out[k*4+3] = ibc_in[k*4+3];
867: }
868: // Update IBC
869: IBCsec->updatePoint(p, ibc_out);
870: }
871: }
872: {
873: // Renumber vertex section IBNDFS
874: const PETSC_MESH_TYPE::int_section_type::chart_type IBNDFSchart = IBNDFSsec->getChart();
875: PETSC_MESH_TYPE::int_section_type::chart_type::const_iterator begin = IBNDFSchart.begin();
876: PETSC_MESH_TYPE::int_section_type::chart_type::const_iterator end = IBNDFSchart.end();
877: for(PETSC_MESH_TYPE::int_section_type::chart_type::const_iterator p_iter = begin; p_iter != end; ++p_iter) {
878: PETSC_MESH_TYPE::point_type p = *p_iter;
879: const PETSC_MESH_TYPE::int_section_type::value_type *ibndfs_in = IBNDFSsec->restrictPoint(p);
880: // Here we assume the fiber dimension is 5
881: PETSC_MESH_TYPE::int_section_type::value_type ibndfs_out[5];
882: // Renumber entries 1,2 & 3 (4 & 5 are invariant under distribution), see IBNDFS's description
883: // Here we assume that vertexNum's domain contains all boundary verticies found in IBNDFS, so we can restrict to the first extracted entry
884: vNum= vertexNum->restrictPoint((PETSC_MESH_TYPE::numbering_type::point_type)ibndfs_in[0]-1+numElements);
885: ibndfs_out[0] = vNum[0]+1;
886: // Map serial bdFace numbers to local bdFace numbers
887: ibndfs_out[1] = bfMap[ibndfs_in[1]];
888: ibndfs_out[2] = bfMap[ibndfs_in[2]];
889: // Copy the other entries right over
890: ibndfs_out[3] = ibndfs_in[3];
891: ibndfs_out[4] = ibndfs_in[4];
892: // Update IBNDFS
893: IBNDFSsec->updatePoint(p,ibndfs_out);
894: }
895: }
896: if (mesh->debug()) {
897: IBCsec->view("Renumbered IBCsec", mesh->comm());
898: IBNDFSsec->view("Renumbered IBNDFSsec", mesh->comm());
899: }
900: // 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
901: };
902: #endif
903: };
904: };