starcdreader.hh

00001 #ifndef DUNE_STARCD_READER_HH
00002 #define DUNE_STARCD_READER_HH
00003 
00004 #include <fstream>
00005 
00006 namespace Dune {
00007 
00041     template <class GridType>
00042     class StarCDReader {
00043 
00044     public:
00045 
00051         static void read(GridType& grid, const std::string& fileName, bool verbose = true)
00052         {
00053             // extract the grid dimension
00054             const int dim = GridType::dimension;
00055             
00056             // currently only dim = 3 is implemented
00057             if (dim != 3)
00058                 DUNE_THROW(Dune::NotImplemented, 
00059                            "Reading Star-CD format is not implemented for dimension " << dim);
00060         
00061             // set up the internal grid creation process
00062             grid.createBegin();
00063         
00064             // set the name of the vertex file 
00065             std::string vertexFileName = fileName + ".vrt";
00066         
00067             // set the vertex input stream 
00068             std::ifstream vertexFile(vertexFileName.c_str());
00069             if (!vertexFile)
00070                 DUNE_THROW(Dune::IOError, "Could not open " << vertexFileName);
00071             
00072             // read the vertices 
00073             int dummyIdx;
00074             int numberOfVertices = 0;
00075             while (vertexFile >> dummyIdx) {
00076                 numberOfVertices++;
00077                 
00078                 Dune::FieldVector<double,dim> position;
00079                 
00080                 for (int k = 0; k < dim; k++)
00081                     vertexFile >> position[k];
00082                 
00083                 grid.insertVertex(position);
00084             }
00085             if (verbose)
00086                 std::cout << numberOfVertices << " vertices read." << std::endl;
00087             
00088             // set the name of the element file 
00089             std::string elementFileName = fileName + ".cel";
00090             
00091             // set the element input stream 
00092             std::ifstream elementFile(elementFileName.c_str());
00093             if (!elementFile)
00094                 DUNE_THROW(Dune::IOError, "Could not open " << elementFileName);
00095             
00096             // read the elements 
00097             int numberOfElements = 0;
00098             int numberOfSimplices = 0;
00099             int numberOfPyramids = 0;
00100             int numberOfPrisms = 0;
00101             int numberOfCubes = 0;;
00102             int maxNumberOfVertices = (int)pow(2, dim);
00103             int isVolume = 1;
00104             while (elementFile >> dummyIdx) {
00105                 std::vector<unsigned int> vertices(maxNumberOfVertices); 
00106                 for (int k = 0; k < maxNumberOfVertices; k++)
00107                     elementFile >> vertices[k];
00108                 
00109                 int boundaryId;
00110                 elementFile >> boundaryId;
00111             
00112                 int volumeOrSurface[2];
00113                 elementFile >> volumeOrSurface[0] >> volumeOrSurface[1]; 
00114                 
00115                 if (volumeOrSurface[0] == isVolume) {
00116                     numberOfElements++;
00117                     
00118                     if (vertices[2] == vertices[3]) { // simplex or prism
00119                         if (vertices[4] == vertices[5]) { // simplex
00120                             numberOfSimplices++;
00121                             std::vector<unsigned int> simplexVertices(4);
00122                             for (int k = 0; k < 3; k++)
00123                                 simplexVertices[k] = vertices[k] - 1; 
00124                             simplexVertices[3] = vertices[4] - 1;
00125                             grid.insertElement(Dune::GeometryType(Dune::GeometryType::simplex,dim), simplexVertices);
00126                         }
00127                         else { // prism
00128                             numberOfPrisms++;
00129                             std::vector<unsigned int> prismVertices(6);
00130                             for (int k = 0; k < 3; k++)
00131                                 prismVertices[k] = vertices[k] - 1; 
00132                             for (int k = 3; k < 6; k++)
00133                                 prismVertices[k] = vertices[k+1] - 1; 
00134                             grid.insertElement(Dune::GeometryType(Dune::GeometryType::prism,dim), prismVertices);
00135                         }
00136                     }
00137                     else { // cube or pyramid 
00138                         if (vertices[4] == vertices[5]) { // pyramid
00139                             numberOfPyramids++;
00140                             std::vector<unsigned int> pyramidVertices(5);
00141                             for (int k = 0; k < 5; k++)
00142                                 pyramidVertices[k] = vertices[k] - 1; 
00143                             grid.insertElement(Dune::GeometryType(Dune::GeometryType::pyramid,dim), pyramidVertices);
00144                         }
00145                         else { // cube
00146                             numberOfCubes++;
00147                             std::vector<unsigned int> cubeVertices(8);
00148                             for (int k = 0; k < 8; k++)
00149                                 cubeVertices[k] = vertices[k] - 1; 
00150                             std::swap(cubeVertices[2], cubeVertices[3]);
00151                             std::swap(cubeVertices[6], cubeVertices[7]);
00152                             grid.insertElement(Dune::GeometryType(Dune::GeometryType::cube,dim), cubeVertices);
00153                         }
00154                     }
00155                 }
00156             }
00157             if (verbose) 
00158                 std::cout << numberOfElements << " elements read: " 
00159                           << numberOfSimplices << " simplices, " << numberOfPyramids << " pyramids, " 
00160                           << numberOfPrisms << " prisms, " << numberOfCubes << " cubes." << std::endl;
00161             
00162             // finish off the construction of the grid object
00163             if (verbose) 
00164                 std::cout << "Starting createEnd() ... " << std::flush;
00165             grid.createEnd();
00166             if (verbose) 
00167                 std::cout << "finished." << std::endl;
00168             
00169         }
00170         
00171     };
00172 
00173 }
00174 
00175 #endif

Generated on 6 Nov 2008 with Doxygen (ver 1.5.6) [logfile].