- Home
- About DUNE
- Download
- Documentation
- Community
- Development
00001 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*- 00002 // vi: set et ts=4 sw=4 sts=4: 00003 00004 #ifndef DUNE_GMSHREADER_HH 00005 #define DUNE_GMSHREADER_HH 00006 00007 #include <cstdarg> 00008 #include <cstdio> 00009 #include <cstring> 00010 #include <fstream> 00011 #include <iostream> 00012 #include <map> 00013 #include <string> 00014 #include <vector> 00015 00016 #include <stdio.h> 00017 00018 #include <dune/common/exceptions.hh> 00019 #include <dune/common/fvector.hh> 00020 #include <dune/common/geometrytype.hh> 00021 00022 #include <dune/grid/common/boundarysegment.hh> 00023 #include <dune/grid/common/gridfactory.hh> 00024 00025 namespace Dune 00026 { 00027 00033 00034 struct GmshReaderOptions 00035 { 00036 enum GeometryOrder { 00038 firstOrder, 00040 secondOrder 00041 }; 00042 }; 00043 00044 namespace { 00045 00046 // arbitrary dimension, implementation is in specialization 00047 template< int dimension, int dimWorld = dimension > 00048 class GmshReaderQuadraticBoundarySegment 00049 { 00050 }; 00051 00052 // quadratic boundary segments in 1d 00053 /* 00054 Note the points 00055 00056 (0) (alpha) (1) 00057 00058 are mapped to the points in global coordinates 00059 00060 p0 p2 p1 00061 00062 alpha is determined automatically from the given points. 00063 */ 00064 template< int dimWorld > 00065 struct GmshReaderQuadraticBoundarySegment< 2, dimWorld > 00066 : public Dune::BoundarySegment< 2, dimWorld > 00067 { 00068 typedef Dune::FieldVector< double, dimWorld > GlobalVector; 00069 00070 GmshReaderQuadraticBoundarySegment ( const GlobalVector &p0_, const GlobalVector &p1_, const GlobalVector &p2_) 00071 : p0(p0_), p1(p1_), p2(p2_) 00072 { 00073 GlobalVector d1 = p1; 00074 d1 -= p0; 00075 GlobalVector d2 = p2; 00076 d2 -= p1; 00077 00078 alpha=d1.two_norm()/(d1.two_norm()+d2.two_norm()); 00079 if (alpha<1E-6 || alpha>1-1E-6) 00080 DUNE_THROW(Dune::IOError, "ration in quadratic boundary segment bad"); 00081 } 00082 00083 virtual GlobalVector operator() ( const Dune::FieldVector<double,1> &local ) const 00084 { 00085 GlobalVector y; 00086 y = 0.0; 00087 y.axpy((local[0]-alpha)*(local[0]-1.0)/alpha,p0); 00088 y.axpy(local[0]*(local[0]-1.0)/(alpha*(alpha-1.0)),p1); 00089 y.axpy(local[0]*(local[0]-alpha)/(1.0-alpha),p2); 00090 return y; 00091 } 00092 00093 private: 00094 GlobalVector p0,p1,p2; 00095 double alpha; 00096 }; 00097 00098 00099 // quadratic boundary segments in 2d 00100 /* numbering of points corresponding to gmsh: 00101 00102 2 00103 00104 5 4 00105 00106 0 3 1 00107 00108 Note: The vertices 3, 4, 5 are not necessarily at the edge midpoints but can 00109 be placed with parameters alpha, beta , gamma at the following positions 00110 in local coordinates: 00111 00112 00113 2 = (0,1) 00114 00115 5 = (0,beta) 4 = (1-gamma/sqrt(2),gamma/sqrt(2)) 00116 00117 0 = (0,0) 3 = (alpha,0) 1 = (1,0) 00118 00119 The parameters alpha, beta, gamma are determined from the given vertices in 00120 global coordinates. 00121 */ 00122 template<> 00123 class GmshReaderQuadraticBoundarySegment< 3, 3 > 00124 : public Dune::BoundarySegment< 3 > 00125 { 00126 public: 00127 GmshReaderQuadraticBoundarySegment (Dune::FieldVector<double,3> p0_, Dune::FieldVector<double,3> p1_, 00128 Dune::FieldVector<double,3> p2_, Dune::FieldVector<double,3> p3_, 00129 Dune::FieldVector<double,3> p4_, Dune::FieldVector<double,3> p5_) 00130 : p0(p0_), p1(p1_), p2(p2_), p3(p3_), p4(p4_), p5(p5_) 00131 { 00132 sqrt2 = sqrt(2.0); 00133 Dune::FieldVector<double,3> d1,d2; 00134 00135 d1 = p3; d1 -= p0; 00136 d2 = p1; d2 -= p3; 00137 alpha=d1.two_norm()/(d1.two_norm()+d2.two_norm()); 00138 if (alpha<1E-6 || alpha>1-1E-6) 00139 DUNE_THROW(Dune::IOError, "alpha in quadratic boundary segment bad"); 00140 00141 d1 = p5; d1 -= p0; 00142 d2 = p2; d2 -= p5; 00143 beta=d1.two_norm()/(d1.two_norm()+d2.two_norm()); 00144 if (beta<1E-6 || beta>1-1E-6) 00145 DUNE_THROW(Dune::IOError, "beta in quadratic boundary segment bad"); 00146 00147 d1 = p4; d1 -= p1; 00148 d2 = p2; d2 -= p4; 00149 gamma=sqrt2*(d1.two_norm()/(d1.two_norm()+d2.two_norm())); 00150 if (gamma<1E-6 || gamma>1-1E-6) 00151 DUNE_THROW(Dune::IOError, "gamma in quadratic boundary segment bad"); 00152 } 00153 00154 virtual Dune::FieldVector<double,3> operator() (const Dune::FieldVector<double,2>& local) const 00155 { 00156 Dune::FieldVector<double,3> y; 00157 y = 0.0; 00158 y.axpy(phi0(local),p0); 00159 y.axpy(phi1(local),p1); 00160 y.axpy(phi2(local),p2); 00161 y.axpy(phi3(local),p3); 00162 y.axpy(phi4(local),p4); 00163 y.axpy(phi5(local),p5); 00164 return y; 00165 } 00166 00167 private: 00168 // The six Lagrange basis function on the reference element 00169 // for the points given above 00170 00171 double phi0 (const Dune::FieldVector<double,2>& local) const 00172 { 00173 return (alpha*beta-beta*local[0]-alpha*local[1])*(1-local[0]-local[1])/(alpha*beta); 00174 } 00175 double phi3 (const Dune::FieldVector<double,2>& local) const 00176 { 00177 return local[0]*(1-local[0]-local[1])/(alpha*(1-alpha)); 00178 } 00179 double phi1 (const Dune::FieldVector<double,2>& local) const 00180 { 00181 return local[0]*(gamma*local[0]-(sqrt2-gamma-sqrt2*alpha)*local[1]-alpha*gamma)/(gamma*(1-alpha)); 00182 } 00183 double phi5 (const Dune::FieldVector<double,2>& local) const 00184 { 00185 return local[1]*(1-local[0]-local[1])/(beta*(1-beta)); 00186 } 00187 double phi4 (const Dune::FieldVector<double,2>& local) const 00188 { 00189 return local[0]*local[1]/((1-gamma/sqrt2)*gamma/sqrt2); 00190 } 00191 double phi2 (const Dune::FieldVector<double,2>& local) const 00192 { 00193 return local[1]*(beta*(1-gamma/sqrt2)-local[0]*(beta-gamma/sqrt2)-local[1]*(1-gamma/sqrt2))/((1-gamma/sqrt2)*(beta-1)); 00194 } 00195 00196 Dune::FieldVector<double,3> p0,p1,p2,p3,p4,p5; 00197 double alpha,beta,gamma,sqrt2; 00198 }; 00199 00200 } // end empty namespace 00201 00203 template<typename GridType> 00204 class GmshReaderParser 00205 { 00206 protected: 00207 // private data 00208 Dune::GridFactory<GridType>& factory; 00209 bool verbose; 00210 bool insert_boundary_segments; 00211 unsigned int number_of_real_vertices; 00212 int boundary_element_count; 00213 int element_count; 00214 // read buffer 00215 char buf[512]; 00216 std::string fileName; 00217 // exported data 00218 std::vector<int> boundary_id_to_physical_entity; 00219 std::vector<int> element_index_to_physical_entity; 00220 00221 // static data 00222 static const int dim = GridType::dimension; 00223 static const int dimWorld = GridType::dimensionworld; 00224 dune_static_assert( (dimWorld <= 3), "GmshReader requires dimWorld <= 3." ); 00225 00226 // typedefs 00227 typedef FieldVector< double, dimWorld > GlobalVector; 00228 00229 // template<class T1, class T2, class T3, class T4, class T5, class T6, 00230 // class T7, class T8, class T9> 00231 // void readfile(FILE * file, int cnt, const char * format, 00232 // T1& t1, T2& t2, T3& t3, T4& t4, T5& t5, T6& t6, T7& t7, T8& t8, T9& t9) 00233 // { 00234 // readfile(file, cnt, &t1, &t2, &t3, &t4, &t5, &t6, &t7, &t8, &t9); 00235 // }; 00236 00237 // don't use something like 00238 // readfile(file, 1, "%s\n", buf); 00239 // to skip the rest of of the line -- that will only skip the next 00240 // whitespace-seperated word! Use skipline() instead. 00241 void readfile(FILE * file, int cnt, const char * format, 00242 void* t1, void* t2 = 0, void* t3 = 0, void* t4 = 0, 00243 void* t5 = 0, void* t6 = 0, void* t7 = 0, void* t8 = 0, 00244 void* t9 = 0, void* t10 = 0) 00245 { 00246 off_t pos = ftello(file); 00247 int c = fscanf(file, format, t1, t2, t3, t4, t5, t6, t7, t8, t9, t10); 00248 if (c != cnt) 00249 DUNE_THROW(Dune::IOError, "Error parsing " << fileName << " " 00250 "file pos " << pos 00251 << ": Expected '" << format << "', only read " << c << " entries instead of " << cnt << "."); 00252 } 00253 00254 // skip over the rest of the line, including the terminating newline 00255 void skipline(FILE * file) 00256 { 00257 int c; 00258 do { 00259 c = std::fgetc(file); 00260 } while(c != '\n' && c != EOF); 00261 } 00262 00263 public: 00264 00265 GmshReaderParser(Dune::GridFactory<GridType>& _factory, bool v, bool i) : 00266 factory(_factory), verbose(v), insert_boundary_segments(i) {} 00267 00268 std::vector<int> & boundaryIdMap() 00269 { 00270 return boundary_id_to_physical_entity; 00271 } 00272 00273 std::vector<int> & elementIndexMap() 00274 { 00275 return element_index_to_physical_entity; 00276 } 00277 00278 void read (const std::string& f) 00279 { 00280 if (verbose) std::cout << "Reading " << dim << "d Gmsh grid..." << std::endl; 00281 00282 // open file name, we use C I/O 00283 fileName = f; 00284 FILE* file = fopen(fileName.c_str(),"r"); 00285 if (file==0) 00286 DUNE_THROW(Dune::IOError, "Could not open " << fileName); 00287 00288 //========================================= 00289 // Header: Read vertices into vector 00290 // Check vertices that are needed 00291 //========================================= 00292 00293 number_of_real_vertices = 0; 00294 boundary_element_count = 0; 00295 element_count = 0; 00296 00297 // process header 00298 double version_number; 00299 int file_type, data_size; 00300 00301 readfile(file,1,"%s\n",buf); 00302 if (strcmp(buf,"$MeshFormat")!=0) 00303 DUNE_THROW(Dune::IOError, "expected $MeshFormat in first line"); 00304 readfile(file,3,"%lg %d %d\n",&version_number,&file_type,&data_size); 00305 if( (version_number < 2.0) || (version_number > 2.2) ) 00306 DUNE_THROW(Dune::IOError, "can only read Gmsh version 2 files"); 00307 if (verbose) std::cout << "version " << version_number << " Gmsh file detected" << std::endl; 00308 readfile(file,1,"%s\n",buf); 00309 if (strcmp(buf,"$EndMeshFormat")!=0) 00310 DUNE_THROW(Dune::IOError, "expected $EndMeshFormat"); 00311 00312 // node section 00313 int number_of_nodes; 00314 00315 readfile(file,1,"%s\n",buf); 00316 if (strcmp(buf,"$Nodes")!=0) 00317 DUNE_THROW(Dune::IOError, "expected $Nodes"); 00318 readfile(file,1,"%d\n",&number_of_nodes); 00319 if (verbose) std::cout << "file contains " << number_of_nodes << " nodes" << std::endl; 00320 00321 // read nodes 00322 std::vector< GlobalVector > nodes( number_of_nodes+1 ); // store positions 00323 { 00324 int id; 00325 double x[ 3 ]; 00326 for( int i = 1; i <= number_of_nodes; ++i ) 00327 { 00328 readfile(file,4, "%d %lg %lg %lg\n", &id, &x[ 0 ], &x[ 1 ], &x[ 2 ] ); 00329 if( id != i ) 00330 DUNE_THROW( Dune::IOError, "Expected id " << i << "(got id " << id << "." ); 00331 00332 // just store node position 00333 for( int j = 0; j < dimWorld; ++j ) 00334 nodes[ i ][ j ] = x[ j ]; 00335 } 00336 readfile(file,1,"%s\n",buf); 00337 if (strcmp(buf,"$EndNodes")!=0) 00338 DUNE_THROW(Dune::IOError, "expected $EndNodes"); 00339 } 00340 00341 // element section 00342 readfile(file,1,"%s\n",buf); 00343 if (strcmp(buf,"$Elements")!=0) 00344 DUNE_THROW(Dune::IOError, "expected $Elements"); 00345 int number_of_elements; 00346 readfile(file,1,"%d\n",&number_of_elements); 00347 if (verbose) std::cout << "file contains " << number_of_elements << " elements" << std::endl; 00348 00349 //========================================= 00350 // Pass 1: Renumber needed vertices 00351 //========================================= 00352 00353 long section_element_offset = ftell(file); 00354 std::map<int,unsigned int> renumber; 00355 for (int i=1; i<=number_of_elements; i++) 00356 { 00357 int id, elm_type, number_of_tags; 00358 readfile(file,3,"%d %d %d ",&id,&elm_type,&number_of_tags); 00359 int physical_entity, elementary_entity; 00360 std::vector<int> mesh_partitions; 00361 if ( version_number < 2.2 ) 00362 { 00363 mesh_partitions.resize(1); 00364 } 00365 for (int k=1; k<=number_of_tags; k++) 00366 { 00367 int blub; 00368 readfile(file,1,"%d ",&blub); 00369 if (k==1) physical_entity = blub; 00370 if (k==2) elementary_entity = blub; 00371 if ( version_number < 2.2 ) 00372 { 00373 if (k==3) mesh_partitions[0] = blub; 00374 } 00375 else 00376 { 00377 if (k > 3) 00378 mesh_partitions[k-4] = blub; 00379 else 00380 mesh_partitions.resize(blub); 00381 } 00382 } 00383 pass1HandleElement(file, elm_type, renumber, nodes); 00384 } 00385 if (verbose) std::cout << "number of real vertices = " << number_of_real_vertices << std::endl; 00386 if (verbose) std::cout << "number of boundary elements = " << boundary_element_count << std::endl; 00387 if (verbose) std::cout << "number of elements = " << element_count << std::endl; 00388 readfile(file,1,"%s\n",buf); 00389 if (strcmp(buf,"$EndElements")!=0) 00390 DUNE_THROW(Dune::IOError, "expected $EndElements"); 00391 boundary_id_to_physical_entity.resize(boundary_element_count); 00392 element_index_to_physical_entity.resize(element_count); 00393 00394 //============================================== 00395 // Pass 2: Insert boundary segments and elements 00396 //============================================== 00397 00398 fseek(file, section_element_offset, SEEK_SET); 00399 boundary_element_count = 0; 00400 element_count = 0; 00401 for (int i=1; i<=number_of_elements; i++) 00402 { 00403 int id, elm_type, number_of_tags; 00404 readfile(file,3,"%d %d %d ",&id,&elm_type,&number_of_tags); 00405 int physical_entity = -1, elementary_entity = -1; 00406 std::vector<int> mesh_partitions; 00407 if ( version_number < 2.2 ) 00408 { 00409 mesh_partitions.resize(1); 00410 } 00411 for (int k=1; k<=number_of_tags; k++) 00412 { 00413 int blub; 00414 readfile(file,1,"%d ",&blub); 00415 if (k==1) physical_entity = blub; 00416 if (k==2) elementary_entity = blub; 00417 if ( version_number < 2.2 ) 00418 { 00419 if (k==3) mesh_partitions[0] = blub; 00420 } 00421 else 00422 { 00423 if (k > 3) 00424 mesh_partitions[k-4] = blub; 00425 else 00426 mesh_partitions.resize(blub); 00427 } 00428 } 00429 pass2HandleElement(file, elm_type, renumber, nodes, physical_entity); 00430 } 00431 readfile(file,1,"%s\n",buf); 00432 if (strcmp(buf,"$EndElements")!=0) 00433 DUNE_THROW(Dune::IOError, "expected $EndElements"); 00434 00435 fclose(file); 00436 } 00437 00438 // dimension dependent routines 00439 void pass1HandleElement(FILE* file, const int elm_type, 00440 std::map<int,unsigned int> & renumber, 00441 const std::vector< GlobalVector > & nodes) 00442 { 00443 // some data about gmsh elements 00444 const int nDofs[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 3, 6, -1, 10}; 00445 const int nVertices[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 2, 3, -1, 4}; 00446 const int elementDim[12] = {-1, 1, 2, 2, 3, 3, 3, 3, 1, 2, -1, 3}; 00447 00448 // test whether we support the element type 00449 if ( not (elm_type >= 0 && elm_type < 12 // index in suitable range? 00450 && (elementDim[elm_type] == dim || elementDim[elm_type] == (dim-1) ) ) ) // real element or boundary element? 00451 { 00452 skipline(file); // skip rest of line if element is unknown 00453 return; 00454 } 00455 00456 // The format string for parsing is n times '%d' in a row 00457 std::string formatString = "%d"; 00458 for (int i=1; i<nDofs[elm_type]; i++) 00459 formatString += " %d"; 00460 formatString += "\n"; 00461 00462 // '10' is the largest number of dofs we may encounter in a .msh file 00463 std::vector<int> elementDofs(10); 00464 00465 readfile(file,nDofs[elm_type], formatString.c_str(), 00466 &(elementDofs[0]),&(elementDofs[1]),&(elementDofs[2]), 00467 &(elementDofs[3]),&(elementDofs[4]),&(elementDofs[5]), 00468 &(elementDofs[6]),&(elementDofs[7]),&(elementDofs[8]), 00469 &(elementDofs[9])); 00470 00471 // insert each vertex if it hasn't been inserted already 00472 for (int i=0; i<nVertices[elm_type]; i++) 00473 if (renumber.find(elementDofs[i])==renumber.end()) 00474 { 00475 renumber[elementDofs[i]] = number_of_real_vertices++; 00476 factory.insertVertex(nodes[elementDofs[i]]); 00477 } 00478 00479 // count elements and boundary elements 00480 if (elementDim[elm_type] == dim) 00481 element_count++; 00482 else 00483 boundary_element_count++; 00484 00485 } 00486 00487 virtual void pass2HandleElement(FILE* file, const int elm_type, 00488 std::map<int,unsigned int> & renumber, 00489 const std::vector< GlobalVector > & nodes, 00490 const int physical_entity) 00491 { 00492 // some data about gmsh elements 00493 const int nDofs[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 3, 6, -1, 10}; 00494 const int nVertices[12] = {-1, 2, 3, 4, 4, 8, 6, 5, 2, 3, -1, 4}; 00495 const int elementDim[12] = {-1, 1, 2, 2, 3, 3, 3, 3, 1, 2, -1, 3}; 00496 00497 // test whether we support the element type 00498 if ( not (elm_type >= 0 && elm_type < 12 // index in suitable range? 00499 && (elementDim[elm_type] == dim || elementDim[elm_type] == (dim-1) ) ) ) // real element or boundary element? 00500 { 00501 skipline(file); // skip rest of line if element is unknown 00502 return; 00503 } 00504 00505 // The format string for parsing is n times '%d' in a row 00506 std::string formatString = "%d"; 00507 for (int i=1; i<nDofs[elm_type]; i++) 00508 formatString += " %d"; 00509 formatString += "\n"; 00510 00511 // '10' is the largest number of dofs we may encounter in a .msh file 00512 std::vector<int> elementDofs(10); 00513 00514 readfile(file,nDofs[elm_type], formatString.c_str(), 00515 &(elementDofs[0]),&(elementDofs[1]),&(elementDofs[2]), 00516 &(elementDofs[3]),&(elementDofs[4]),&(elementDofs[5]), 00517 &(elementDofs[6]),&(elementDofs[7]),&(elementDofs[8]), 00518 &(elementDofs[9])); 00519 00520 // correct differences between gmsh and Dune in the local vertex numbering 00521 switch (elm_type) 00522 { 00523 case 3: // 4-node quadrilateral 00524 std::swap(elementDofs[2],elementDofs[3]); 00525 break; 00526 case 5: // 8-node hexahedron 00527 std::swap(elementDofs[2],elementDofs[3]); 00528 std::swap(elementDofs[6],elementDofs[7]); 00529 break; 00530 case 7: // 5-node pyramid 00531 std::swap(elementDofs[2],elementDofs[3]); 00532 break; 00533 } 00534 00535 // renumber corners to account for the explicitly given vertex 00536 // numbering in the file 00537 std::vector<unsigned int> vertices(nVertices[elm_type]); 00538 00539 for (int i=0; i<nVertices[elm_type]; i++) 00540 vertices[i] = renumber[elementDofs[i]]; 00541 00542 // If it is an element, insert it as such 00543 if (elementDim[elm_type] == dim) { 00544 00545 switch (elm_type) 00546 { 00547 case 1: // 2-node line 00548 factory.insertElement(Dune::GeometryType(Dune::GeometryType::simplex,dim),vertices); 00549 break; 00550 case 2: // 3-node triangle 00551 factory.insertElement(Dune::GeometryType(Dune::GeometryType::simplex,dim),vertices); 00552 break; 00553 case 3: // 4-node quadrilateral 00554 factory.insertElement(Dune::GeometryType(Dune::GeometryType::cube,dim),vertices); 00555 break; 00556 case 4: // 4-node tetrahedron 00557 factory.insertElement(Dune::GeometryType(Dune::GeometryType::simplex,dim),vertices); 00558 break; 00559 case 5: // 8-node hexahedron 00560 factory.insertElement(Dune::GeometryType(Dune::GeometryType::cube,dim),vertices); 00561 break; 00562 case 6: // 6-node prism 00563 factory.insertElement(Dune::GeometryType(Dune::GeometryType::prism,dim),vertices); 00564 break; 00565 case 7: // 5-node pyramid 00566 factory.insertElement(Dune::GeometryType(Dune::GeometryType::pyramid,dim),vertices); 00567 break; 00568 case 9: // 6-node triangle 00569 factory.insertElement(Dune::GeometryType(Dune::GeometryType::simplex,dim),vertices); 00570 break; 00571 case 11: // 10-node tetrahedron 00572 factory.insertElement(Dune::GeometryType(Dune::GeometryType::simplex,dim),vertices); 00573 break; 00574 } 00575 00576 } else { 00577 // it must be a boundary segment then 00578 if (insert_boundary_segments) { 00579 00580 switch (elm_type) 00581 { 00582 case 1: // 2-node line 00583 factory.insertBoundarySegment(vertices); 00584 break; 00585 00586 case 2: // 3-node triangle 00587 factory.insertBoundarySegment(vertices); 00588 break; 00589 00590 case 8: { // 3-node line 00591 array<FieldVector<double,dim>, 3> v; 00592 for (int i=0; i<dim; i++) { 00593 v[0][i] = nodes[elementDofs[0]][i]; 00594 v[1][i] = nodes[elementDofs[2]][i]; // yes, the renumbering is intended! 00595 v[2][i] = nodes[elementDofs[1]][i]; 00596 } 00597 BoundarySegment<dim,dimWorld>* newBoundarySegment 00598 = (BoundarySegment<dim,dimWorld>*)new GmshReaderQuadraticBoundarySegment< 2, dimWorld >(v[0], v[1], v[2]); 00599 factory.insertBoundarySegment(vertices, 00600 shared_ptr<BoundarySegment<dim,dimWorld> >(newBoundarySegment)); 00601 break; 00602 } 00603 case 9: { // 6-node triangle 00604 array<FieldVector<double,dim>, 6> v; 00605 for (int i=0; i<6; i++) 00606 for (int j=0; j<dim; j++) 00607 v[i][j] = nodes[elementDofs[i]][j]; 00608 00609 BoundarySegment<dim,dimWorld>* newBoundarySegment 00610 = (BoundarySegment<dim,dimWorld>*)new GmshReaderQuadraticBoundarySegment< 3, 3 >(v[0], v[1], v[2], 00611 v[3], v[4], v[5]); 00612 00613 factory.insertBoundarySegment(vertices, 00614 shared_ptr<BoundarySegment<dim,dimWorld> >(newBoundarySegment)); 00615 break; 00616 } 00617 00618 } 00619 00620 } 00621 } 00622 00623 // count elements and boundary elements 00624 if (elementDim[elm_type] == dim) { 00625 element_index_to_physical_entity[element_count] = physical_entity; 00626 element_count++; 00627 } else { 00628 boundary_id_to_physical_entity[boundary_element_count] = physical_entity; 00629 boundary_element_count++; 00630 } 00631 00632 } 00633 00634 }; 00635 00652 template<typename GridType> 00653 class GmshReader 00654 { 00655 public: 00657 static GridType* read (const std::string& fileName, bool verbose = true, bool insert_boundary_segments=true) 00658 { 00659 // make a grid factory 00660 Dune::GridFactory<GridType> factory; 00661 00662 // create parse object 00663 GmshReaderParser<GridType> parser(factory,verbose,insert_boundary_segments); 00664 parser.read(fileName); 00665 00666 return factory.createGrid(); 00667 } 00668 00670 static GridType* read (const std::string& fileName, 00671 std::vector<int>& boundary_id_to_physical_entity, 00672 std::vector<int>& element_index_to_physical_entity, 00673 bool verbose = true, bool insert_boundary_segments=true) 00674 { 00675 // make a grid factory 00676 Dune::GridFactory<GridType> factory; 00677 00678 // create parse object 00679 GmshReaderParser<GridType> parser(factory,verbose,insert_boundary_segments); 00680 parser.read(fileName); 00681 00682 boundary_id_to_physical_entity.swap(parser.boundaryIdMap()); 00683 element_index_to_physical_entity.swap(parser.elementIndexMap()); 00684 00685 return factory.createGrid(); 00686 } 00687 00689 static GridType* read (GridType& grid, const std::string& fileName, 00690 bool verbose = true, bool insert_boundary_segments=true) DUNE_DEPRECATED 00691 { 00692 // make a grid factory 00693 Dune::GridFactory<GridType> factory(&grid); 00694 00695 // create parse object 00696 GmshReaderParser<GridType> parser(factory,verbose,insert_boundary_segments); 00697 parser.read(fileName); 00698 00699 return factory.createGrid(); 00700 } 00701 00703 static GridType* read (GridType& grid, const std::string& fileName, 00704 std::vector<int>& boundary_id_to_physical_entity, 00705 std::vector<int>& element_index_to_physical_entity, 00706 bool verbose = true, bool insert_boundary_segments=true) DUNE_DEPRECATED 00707 { 00708 // make a grid factory 00709 Dune::GridFactory<GridType> factory(&grid); 00710 00711 // create parse object 00712 GmshReaderParser<GridType> parser(factory,verbose,insert_boundary_segments); 00713 parser.read(fileName); 00714 00715 boundary_id_to_physical_entity.swap(parser.boundaryIdMap()); 00716 element_index_to_physical_entity.swap(parser.elementIndexMap()); 00717 00718 return factory.createGrid(); 00719 } 00720 00721 00723 static void read (Dune::GridFactory<GridType>& factory, const std::string& fileName, 00724 bool verbose = true, bool insert_boundary_segments=true) 00725 { 00726 // create parse object 00727 GmshReaderParser<GridType> parser(factory,verbose,insert_boundary_segments); 00728 parser.read(fileName); 00729 } 00730 00732 static void read (Dune::GridFactory<GridType>& factory, 00733 const std::string& fileName, 00734 std::vector<int>& boundary_id_to_physical_entity, 00735 std::vector<int>& element_index_to_physical_entity, 00736 bool verbose = true, bool insert_boundary_segments=true) 00737 { 00738 // create parse object 00739 GmshReaderParser<GridType> parser(factory,verbose,insert_boundary_segments); 00740 parser.read(fileName); 00741 00742 boundary_id_to_physical_entity.swap(parser.boundaryIdMap()); 00743 element_index_to_physical_entity.swap(parser.elementIndexMap()); 00744 } 00745 }; 00746 00749 } // namespace Dune 00750 00751 #endif
Generated on Fri Apr 29 2011 with Doxygen (ver 1.7.1) [doxygen-log,error-log].