Actual source code: meshlagrit.c

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

  3: namespace ALE {
  4:   namespace LaGriT {
  5:     void Builder::readInpFile(MPI_Comm comm, const std::string& filename, const int dim, const int numCorners, int& numElements, int *vertices[], int& numVertices, double *coordinates[]) {
  6:       PetscViewer    viewer;
  7:       FILE          *f;
  8:       PetscScalar   *coords;
  9:       PetscInt      *verts;
 10:       PetscInt       commRank;
 11:       char           buf[2048];

 14:       MPI_Comm_rank(comm, &commRank);
 15:       if (commRank != 0) return;
 16:       PetscViewerCreate(PETSC_COMM_SELF, &viewer);
 17:       PetscViewerSetType(viewer, PETSC_VIEWER_ASCII);
 18:       PetscViewerFileSetMode(viewer, FILE_MODE_READ);
 19:       PetscViewerFileSetName(viewer, filename.c_str());
 20:       PetscViewerASCIIGetPointer(viewer, &f);
 21:       // Read header
 22:       if (fgets(buf, 2048, f) == NULL) throw ALE::Exception("File ended prematurely");
 23:       // Number of vertices
 24:       const char *x = strtok(buf, " ");
 25:       numVertices = atoi(x);
 26:       // Number of elements
 27:       x = strtok(NULL, " ");
 28:       numElements = atoi(x);
 29:       // Element type
 30:       x = strtok(NULL, " ");
 31:       // ???
 32:       x = strtok(NULL, " ");
 33:       // ???
 34:       x = strtok(NULL, " ");
 35:       PetscMalloc(numVertices*dim * sizeof(PetscScalar), &coords);
 36:       for(PetscInt i = 0; i < numVertices; ++i) {
 37:         if (fgets(buf, 2048, f) == NULL) throw ALE::Exception("File ended prematurely");
 38:         x = strtok(buf, " ");
 39:         // Ignore vertex number
 40:         x = strtok(NULL, " ");
 41:         for(int c = 0; c < dim; c++) {
 42:           coords[i*dim+c] = atof(x);
 43:           x = strtok(NULL, " ");
 44:         }
 45:       }
 46:       *coordinates = coords;
 47:       PetscMalloc(numElements*numCorners * sizeof(PetscInt), &verts);
 48:       for(PetscInt i = 0; i < numElements; ++i) {
 49:         if (fgets(buf, 2048, f) == NULL) throw ALE::Exception("File ended prematurely");
 50:         x = strtok(buf, " ");
 51:         // Ignore element number
 52:         x = strtok(NULL, " ");
 53:         // Ignore ??? (material)
 54:         x = strtok(NULL, " ");
 55:         // Ignore element type
 56:         x = strtok(NULL, " ");
 57:         for(int c = 0; c < numCorners; c++) {
 58:           verts[i*numCorners+c] = atoi(x) - 1;
 59:           x = strtok(NULL, " ");
 60:         }
 61:       }
 62:       *vertices = verts;
 63:       PetscViewerDestroy(viewer);
 64:     };
 65:     Obj<Builder::Mesh> Builder::readMesh(MPI_Comm comm, const int dim, const std::string& filename, const bool interpolate = false, const int debug = 0) {
 66:       Obj<Mesh>       mesh  = new Mesh(comm, dim, debug);
 67:       Obj<sieve_type> sieve = new sieve_type(comm, debug);
 68:       int    *cells;
 69:       double *coordinates;
 70:       int     numCells = 0, numVertices = 0, numCorners = dim+1;

 72:       Builder::readInpFile(comm, filename, dim, numCorners, numCells, &cells, numVertices, &coordinates);
 73:       ALE::SieveBuilder<Mesh>::buildTopology(sieve, dim, numCells, cells, numVertices, interpolate, numCorners);
 74:       mesh->setSieve(sieve);
 75:       mesh->stratify();
 76:       ALE::SieveBuilder<Mesh>::buildCoordinates(mesh, dim, coordinates);
 77:       return mesh;
 78:     };
 79:     void Builder::readFault(Obj<Builder::Mesh> mesh, const std::string& filename) {
 80:       const int      numCells = mesh->heightStratum(0)->size();
 81:       PetscViewer    viewer;
 82:       FILE          *f;
 83:       char           buf[2048];
 84:       PetscInt       numPsets;

 87:       if (mesh->commRank() != 0) return;
 88:       PetscViewerCreate(PETSC_COMM_SELF, &viewer);
 89:       PetscViewerSetType(viewer, PETSC_VIEWER_ASCII);
 90:       PetscViewerFileSetMode(viewer, FILE_MODE_READ);
 91:       PetscViewerFileSetName(viewer, filename.c_str());
 92:       PetscViewerASCIIGetPointer(viewer, &f);
 93:       // Read header
 94:       if (fgets(buf, 2048, f) == NULL) throw ALE::Exception("File ended prematurely");
 95:       // Check file type
 96:       const char *x = strtok(buf, " ");
 97:       std::string fileType("pset");
 98:       if (fileType != x) throw ALE::Exception("Invalid file type");
 99:       // Ignore set type
100:       x = strtok(NULL, " ");
101:       // Number of psets
102:       x = strtok(NULL, " ");
103:       numPsets = atoi(x);
104:       for(PetscInt p = 0; p < numPsets; ++p) {
105:         if (fgets(buf, 2048, f) == NULL) throw ALE::Exception("File ended prematurely");
106:         // Read name
107:         x = strtok(buf, " ");
108:         const Obj<Mesh::int_section_type>& fault = mesh->getIntSection(x);
109:         // Vertices per line
110:         x = strtok(NULL, " ");
111:         const PetscInt vertsPerLine = atoi(x);
112:         // Total vertices
113:         x = strtok(NULL, " ");
114:         const PetscInt totalVerts = atoi(x);

116:         for(PetscInt v = 0; v < totalVerts; ++v) {
117:           if (v%vertsPerLine == 0) {
118:             if (fgets(buf, 2048, f) == NULL) throw ALE::Exception("File ended prematurely");
119:             x = strtok(buf, " ");
120:           } else {
121:             x = strtok(NULL, " ");
122:           }
123:           const int v = atoi(x) + numCells - 1;

125:           fault->setFiberDimension(v, 1);
126:         }
127:         fault->allocatePoint();
128:       }
129:       PetscViewerDestroy(viewer);
130:     };
131:   }
132: }