00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020 #include <GmshFormatReader.hpp>
00021
00022 #include <SurfaceMeshOfQuadrangles.hpp>
00023 #include <MeshOfHexahedra.hpp>
00024
00025 #include <SurfaceMeshOfTriangles.hpp>
00026 #include <MeshOfTetrahedra.hpp>
00027
00028 #include <EndianConverter.hpp>
00029
00030 void
00031 GmshFormatReader::__readVertices()
00032 {
00033 const int numberOfVerices = this->__getInteger();
00034 ffout(3) << "- Number Of Vertices: " << numberOfVerices << '\n';
00035 if (numberOfVerices<0) {
00036 throw ErrorHandler(__FILE__,__LINE__,
00037 "reading file '"+this->__fileName
00038 +"': number of vertices is negative",
00039 ErrorHandler::normal);
00040 }
00041
00042 __verticesNumbers.resize(numberOfVerices);
00043 __vertices = new VerticesSet(numberOfVerices);
00044 VerticesSet& vertices = *__vertices;
00045
00046 if (not __binary) {
00047 for (int i=0; i<numberOfVerices; ++i) {
00048 __verticesNumbers[i] = this->__getInteger();
00049 const real_t x = this->__getReal();
00050 const real_t y = this->__getReal();
00051 const real_t z = this->__getReal();
00052 vertices[i] = TinyVector<3, real_t>(x,y,z);
00053 }
00054 } else {
00055 fseek(__ifh,1L,SEEK_CUR);
00056 for (int i=0; i<numberOfVerices; ++i) {
00057 int number = 0;
00058 fread(reinterpret_cast<char*>(&number),sizeof(int),1,__ifh);
00059 __verticesNumbers[i] = number;
00060 TinyVector<3,double> x;
00061 fread(reinterpret_cast<char*>(&(x[0])),sizeof(double),3,__ifh);
00062 for (size_t j=0; j<3; ++j) {
00063 vertices[i][j] = x[j];
00064 }
00065 }
00066
00067 if (__convertEndian) {
00068 for (int i=0; i<numberOfVerices; ++i) {
00069 invertEndianess(__verticesNumbers[i]);
00070 for (size_t j=0; j<3; ++j) {
00071 invertEndianess(vertices[i][j]);
00072 }
00073 }
00074 }
00075 }
00076 }
00077
00078 void
00079 GmshFormatReader::__readElements1()
00080 {
00081 const int numberOfElements = this->__getInteger();
00082 ffout(3) << "- Number Of Elements: " << numberOfElements << '\n';
00083 if (numberOfElements<0) {
00084 throw ErrorHandler(__FILE__,__LINE__,
00085 "reading file '"+__fileName
00086 +"': number of elements is negative",
00087 ErrorHandler::normal);
00088 }
00089
00090 __elementType.resize(numberOfElements);
00091 __references.resize(numberOfElements);
00092 __elementVertices.resize(numberOfElements);
00093
00094 for (int i=0; i<numberOfElements; ++i) {
00095 this->__getInteger();
00096 const int elementType = this->__getInteger()-1;
00097
00098 if ((elementType < 0) or (elementType > 14)) {
00099 throw ErrorHandler(__FILE__,__LINE__,
00100 "reading file '"+__fileName
00101 +"': unknown element type '"+stringify(elementType)+"'",
00102 ErrorHandler::normal);
00103 }
00104 __elementType[i] = elementType;
00105 __references[i] = this->__getInteger();
00106 this->__getInteger();
00107
00108 const int numberOfVertices = this->__getInteger();
00109 if (numberOfVertices != __numberOfPrimitiveNodes[elementType]) {
00110 std::stringstream errorMsg;
00111 errorMsg << "reading file '" <<__fileName
00112 << "':number of vertices (" << numberOfVertices
00113 << ") does not match expected ("
00114 << __numberOfPrimitiveNodes[elementType] << ")" << std::ends;
00115
00116 throw ErrorHandler(__FILE__,__LINE__,
00117 errorMsg.str(),
00118 ErrorHandler::normal);
00119 }
00120 __elementVertices[i].resize(numberOfVertices);
00121 for (int j=0; j<numberOfVertices; ++j) {
00122 __elementVertices[i][j] = this->__getInteger();
00123 }
00124 }
00125 }
00126
00127 void
00128 GmshFormatReader::__readElements2()
00129 {
00130 const int numberOfElements = this->__getInteger();
00131 ffout(3) << "- Number Of Elements: " << numberOfElements << '\n';
00132 if (numberOfElements<0) {
00133 throw ErrorHandler(__FILE__,__LINE__,
00134 "reading file '"+__fileName
00135 +"': number of elements is negative",
00136 ErrorHandler::normal);
00137 }
00138
00139 __elementType.resize(numberOfElements);
00140 __references.resize(numberOfElements);
00141 __elementVertices.resize(numberOfElements);
00142
00143 if (not __binary) {
00144 for (int i=0; i<numberOfElements; ++i) {
00145 this->__getInteger();
00146 const int elementType = this->__getInteger()-1;
00147
00148 if ((elementType < 0) or (elementType > 14)) {
00149 throw ErrorHandler(__FILE__,__LINE__,
00150 "reading file '"+__fileName
00151 +"': unknown element type '"+stringify(elementType)+"'",
00152 ErrorHandler::normal);
00153 }
00154 __elementType[i] = elementType;
00155 const int numberOfTags = this->__getInteger();
00156 __references[i] = this->__getInteger();
00157 for (int tag=1; tag<numberOfTags; ++tag) {
00158 this->__getInteger();
00159 }
00160
00161 const int numberOfVertices = __numberOfPrimitiveNodes[elementType];
00162 __elementVertices[i].resize(numberOfVertices);
00163 for (int j=0; j<numberOfVertices; ++j) {
00164 __elementVertices[i][j] = this->__getInteger();
00165 }
00166 }
00167 } else {
00168 fseek(__ifh,1L,SEEK_CUR);
00169 int i=0;
00170 for (;i<numberOfElements;) {
00171 int elementType = 0;
00172 fread(reinterpret_cast<char*>(&elementType),sizeof(int),1,__ifh);
00173 if (__convertEndian) {
00174 invertEndianess(elementType);
00175 }
00176 --elementType;
00177
00178 if ((elementType < 0) or (elementType > 14)) {
00179 throw ErrorHandler(__FILE__,__LINE__,
00180 "reading file '"+__fileName
00181 +"': unknown element type '"+stringify(elementType)+"'",
00182 ErrorHandler::normal);
00183 }
00184
00185 int elementTypeNumber = 0;
00186 fread(reinterpret_cast<char*>(&elementTypeNumber),sizeof(int),1,__ifh);
00187 if (__convertEndian) {
00188 invertEndianess(elementTypeNumber);
00189 }
00190
00191 int numberOfTags = 0;
00192 fread(reinterpret_cast<char*>(&numberOfTags),sizeof(int),1,__ifh);
00193 if (__convertEndian) {
00194 invertEndianess(numberOfTags);
00195 }
00196
00197 for (int k=0; k<elementTypeNumber; ++k) {
00198 __elementType[i] = elementType;
00199
00200 int elementNumber = 0;
00201 fread(reinterpret_cast<char*>(&elementNumber),sizeof(int),1,__ifh);
00202
00203 int reference = 0;
00204 fread(reinterpret_cast<char*>(&reference),sizeof(int),1,__ifh);
00205
00206 __references[i] = reference;
00207 for (int tag=1; tag<numberOfTags; ++tag) {
00208 int t;
00209 fread(reinterpret_cast<char*>(&t),sizeof(int),1,__ifh);
00210 }
00211
00212 const int numberOfVertices = __numberOfPrimitiveNodes[elementType];
00213 __elementVertices[i].resize(numberOfVertices);
00214 fread(reinterpret_cast<char*>(&(__elementVertices[i][0])),sizeof(int),numberOfVertices,__ifh);
00215 ++i;
00216 }
00217 }
00218
00219 if (__convertEndian) {
00220 for (size_t i=0; i<__references.size(); ++i) {
00221 invertEndianess(__references[i]);
00222 }
00223 for (size_t i=0; i<__elementVertices.size(); ++i) {
00224 for (size_t j=0; j<__elementVertices[i].size(); ++i) {
00225 invertEndianess(__elementVertices[i][j]);
00226 }
00227 }
00228 }
00229 }
00230 }
00231
00232 void
00233 GmshFormatReader::__proceedData()
00234 {
00235 TinyVector<15,int> elementNumber = 0;
00236 std::vector<std::set<size_t> > elementReferences(15);
00237
00238 for (size_t i=0; i<__elementType.size(); ++i) {
00239 const int elementType = __elementType[i];
00240 elementNumber[elementType]++;
00241 elementReferences[elementType].insert(__references[i]);
00242 }
00243
00244 for (size_t i=0; i<elementNumber.size(); ++i) {
00245 if (elementNumber[i]>0) {
00246 ffout(3) << " - Number of "
00247 << __primitivesNames[i]
00248 << ": " << elementNumber[i];
00249 if (not(__supportedPrimitives[i])) {
00250 ffout(3) << " [IGNORED]";
00251 } else {
00252 ffout(3) << " references: ";
00253 for(std::set<size_t>::const_iterator iref
00254 = elementReferences[i].begin();
00255 iref != elementReferences[i].end(); ++iref) {
00256 if (iref != elementReferences[i].begin()) ffout(3) << ',';
00257 ffout(3) << *iref;
00258 }
00259
00260 switch (i) {
00261
00262 case 1: {
00263 __triangles = new Vector<Triangle>(elementNumber[i]);
00264 break;
00265 }
00266 case 2: {
00267 __quadrilaterals = new Vector<Quadrangle>(elementNumber[i]);
00268 break;
00269 }
00270 case 3: {
00271 __tetrahedra = new Vector<Tetrahedron>(elementNumber[i]);
00272 break;
00273 }
00274 case 4: {
00275 __hexahedra = new Vector<Hexahedron>(elementNumber[i]);
00276 }
00277
00278 case 0:
00279 case 14: {
00280 break;
00281 }
00282
00283 case 5:
00284 case 6:
00285 case 7:
00286 case 8:
00287 case 9:
00288 case 10:
00289 case 11:
00290 case 12:
00291 case 13:
00292 default: {
00293 throw ErrorHandler(__FILE__,__LINE__,
00294 "reading file '"+__fileName
00295 +"': ff3d cannot read this kind of element",
00296 ErrorHandler::normal);
00297 }
00298 }
00299 }
00300 ffout(3) << '\n';
00301 }
00302 }
00303
00304 ffout(3) << "- Building correspondance table\n";
00305 int minNumber = std::numeric_limits<int>::max();
00306 int maxNumber = std::numeric_limits<int>::min();
00307 for (size_t i=0; i<__verticesNumbers.size(); ++i) {
00308 const int& vertexNumber = __verticesNumbers[i];
00309 minNumber = std::min(minNumber,vertexNumber);
00310 maxNumber = std::max(maxNumber,vertexNumber);
00311 }
00312
00313 if (minNumber<0) {
00314 throw ErrorHandler(__FILE__,__LINE__,
00315 "reading file '"+__fileName
00316 +"': ff3d does not implement negative vertices number",
00317 ErrorHandler::normal);
00318 }
00319
00320
00321 __verticesCorrepondance.resize(maxNumber+1,-1);
00322
00323 for (size_t i=0; i<__verticesNumbers.size(); ++i) {
00324
00325 __verticesCorrepondance[__verticesNumbers[i]] = i;
00326 }
00327
00328
00329 elementNumber = 0;
00330
00331 VerticesSet& vertices = *__vertices;
00332
00333
00334 for (size_t i=0; i<__elementType.size(); ++i) {
00335 std::vector<int>& elementVertices = __elementVertices[i];
00336 switch (__elementType[i]) {
00337
00338 case 1: {
00339 int& triangleNumber = elementNumber[1];
00340
00341 const size_t a = __verticesCorrepondance[elementVertices[0]];
00342 const size_t b = __verticesCorrepondance[elementVertices[1]];
00343 const size_t c = __verticesCorrepondance[elementVertices[2]];
00344 if ((a<0) or (b<0) or (c<0)) {
00345 throw ErrorHandler(__FILE__,__LINE__,
00346 "reading file '"+__fileName
00347 +"': error reading element "+stringify(i)
00348 +" [bad vertices definition]",
00349 ErrorHandler::normal);
00350 }
00351 (*__triangles)[triangleNumber]
00352 = Triangle(vertices[a], vertices[b], vertices[c],
00353 __references[i]);
00354 triangleNumber++;
00355 break;
00356 }
00357 case 2: {
00358 int& quadrilateralNumber = elementNumber[2];
00359
00360 const size_t a = __verticesCorrepondance[elementVertices[0]];
00361 const size_t b = __verticesCorrepondance[elementVertices[1]];
00362 const size_t c = __verticesCorrepondance[elementVertices[2]];
00363 const size_t d = __verticesCorrepondance[elementVertices[3]];
00364 if ((a<0) or (b<0) or (c<0) or (d<0)) {
00365 throw ErrorHandler(__FILE__,__LINE__,
00366 "reading file '"+__fileName
00367 +"': error reading element "+stringify(i)
00368 +" [bad vertices definition]",
00369 ErrorHandler::normal);
00370 }
00371 (*__quadrilaterals)[quadrilateralNumber]
00372 = Quadrangle(vertices[a], vertices[b], vertices[c], vertices[d],
00373 __references[i]);
00374 quadrilateralNumber++;
00375 break;
00376 }
00377 case 3: {
00378 int& tetrahedronNumber = elementNumber[3];
00379
00380 const int a = __verticesCorrepondance[elementVertices[0]];
00381 const int b = __verticesCorrepondance[elementVertices[1]];
00382 const int c = __verticesCorrepondance[elementVertices[2]];
00383 const int d = __verticesCorrepondance[elementVertices[3]];
00384 if ((a<0) or (b<0) or (c<0) or (d<0)) {
00385 throw ErrorHandler(__FILE__,__LINE__,
00386 "reading file '"+__fileName
00387 +"': error reading element "+stringify(i)
00388 +" [bad vertices definition]",
00389 ErrorHandler::normal);
00390 }
00391 (*__tetrahedra)[tetrahedronNumber]
00392 = Tetrahedron(vertices[a], vertices[b], vertices[c], vertices[d],
00393 __references[i]);
00394 tetrahedronNumber++;
00395 break;
00396 }
00397 case 4: {
00398 int& hexahedronNumber = elementNumber[4];
00399
00400 const int a = __verticesCorrepondance[elementVertices[0]];
00401 const int b = __verticesCorrepondance[elementVertices[1]];
00402 const int c = __verticesCorrepondance[elementVertices[2]];
00403 const int d = __verticesCorrepondance[elementVertices[3]];
00404 const int e = __verticesCorrepondance[elementVertices[4]];
00405 const int f = __verticesCorrepondance[elementVertices[5]];
00406 const int g = __verticesCorrepondance[elementVertices[6]];
00407 const int h = __verticesCorrepondance[elementVertices[7]];
00408 if ((a<0) or (b<0) or (c<0) or (d<0) or (e<0) or (f<0) or (g<0) or (h<0)) {
00409 throw ErrorHandler(__FILE__,__LINE__,
00410 "reading file '"+__fileName
00411 +"': error reading element "+stringify(i)
00412 +" [bad vertices definition]",
00413 ErrorHandler::normal);
00414 }
00415 (*__hexahedra)[hexahedronNumber]
00416 = Hexahedron(vertices[a], vertices[b], vertices[c], vertices[d],
00417 vertices[e], vertices[f], vertices[g], vertices[h],
00418 __references[i]);
00419 hexahedronNumber++;
00420 break; }
00421
00422 case 0:
00423 case 14:{
00424 break;
00425 }
00426 case 5:
00427 case 6:
00428 case 7:
00429 case 8:
00430 case 9:
00431 case 10:
00432 case 11:
00433 case 12:
00434 case 13:{
00435 }
00436 default: {
00437 throw ErrorHandler(__FILE__,__LINE__,
00438 "reading file '"+__fileName
00439 +"': ff3d cannot read this kind of element",
00440 ErrorHandler::normal);
00441 }
00442 }
00443 }
00444 }
00445
00446 GmshFormatReader::Keyword
00447 GmshFormatReader::__nextKeyword()
00448 {
00449 GmshFormatReader::Keyword kw;
00450
00451 char pKeyword[256];
00452 int retval = fscanf(__ifh,"%255s",pKeyword);
00453 const std::string aKeyword = pKeyword;
00454
00455 if (retval < 0) {
00456 kw.second = EndOfFile;
00457 return kw;
00458 } else if (retval == 0) {
00459 kw.second = Unknown;
00460 return kw;
00461 }
00462
00463 KeywordList::iterator i = __keywordList.find(aKeyword.c_str());
00464
00465 if (i != __keywordList.end()) {
00466 kw.first = (*i).first;
00467 kw.second = (*i).second;
00468 return kw;
00469 }
00470
00471 throw ErrorHandler(__FILE__,__LINE__,
00472 "reading file '"+__fileName
00473 +"': unknown keyword '"+aKeyword+"'",
00474 ErrorHandler::normal);
00475
00476 kw.first = aKeyword;
00477 kw.second = Unknown;
00478 return kw;
00479 }
00480
00481 GmshFormatReader::GmshFormatReader(const std::string & s)
00482 : MeshReader(s)
00483 {
00484
00485 __keywordList["$NOD"] = NOD;
00486 __keywordList["$ENDNOD"] = ENDNOD;
00487 __keywordList["$ELM"] = ELM;
00488 __keywordList["$ENDELM"] = ENDELM;
00489
00490
00491 __keywordList["$MeshFormat"] = MESHFORMAT;
00492 __keywordList["$EndMeshFormat"] = ENDMESHFORMAT;
00493 __keywordList["$Nodes"] = NODES;
00494 __keywordList["$EndNodes"] = ENDNODES;
00495 __keywordList["$Elements"] = ELEMENTS;
00496 __keywordList["$EndElements"] = ENDELEMENTS;
00497
00498 __numberOfPrimitiveNodes.resize(16);
00499 __numberOfPrimitiveNodes[ 0] = 2;
00500 __numberOfPrimitiveNodes[ 1] = 3;
00501 __numberOfPrimitiveNodes[ 2] = 4;
00502 __numberOfPrimitiveNodes[ 3] = 4;
00503 __numberOfPrimitiveNodes[ 4] = 8;
00504 __numberOfPrimitiveNodes[ 5] = 6;
00505 __numberOfPrimitiveNodes[ 6] = 5;
00506 __numberOfPrimitiveNodes[ 7] = 3;
00507 __numberOfPrimitiveNodes[ 8] = 6;
00508 __numberOfPrimitiveNodes[ 9] = 9;
00509 __numberOfPrimitiveNodes[10] = 10;
00510 __numberOfPrimitiveNodes[11] = 27;
00511 __numberOfPrimitiveNodes[12] = 18;
00512 __numberOfPrimitiveNodes[13] = 14;
00513 __numberOfPrimitiveNodes[14] = 1;
00514
00515 __primitivesNames[0] = "edges";
00516 __supportedPrimitives[0] = false;
00517 __primitivesNames[1] = "triangles";
00518 __supportedPrimitives[1] = true;
00519 __primitivesNames[2] = "quadrangles";
00520 __supportedPrimitives[2] = true;
00521 __primitivesNames[3] = "tetrahedra";
00522 __supportedPrimitives[3] = true;
00523 __primitivesNames[4] = "hexahedra";
00524 __supportedPrimitives[4] = true;
00525 __primitivesNames[5] = "prisms";
00526 __supportedPrimitives[5] = false;
00527 __primitivesNames[6] = "pyramids";
00528 __supportedPrimitives[6] = false;
00529 __primitivesNames[7] = "second order edges";
00530 __supportedPrimitives[7] = false;
00531 __primitivesNames[8] = "second order triangles";
00532 __supportedPrimitives[8] = false;
00533 __primitivesNames[9] = "second order quadrangles";
00534 __supportedPrimitives[9] = false;
00535 __primitivesNames[10] = "second order tetrahedra";
00536 __supportedPrimitives[10]= false;
00537 __primitivesNames[11] = "second order hexahedra";
00538 __supportedPrimitives[11]= false;
00539 __primitivesNames[12] = "second order prisms";
00540 __supportedPrimitives[12]= false;
00541 __primitivesNames[13] = "second order pyramids";
00542 __supportedPrimitives[13]= false;
00543 __primitivesNames[14] = "point";
00544 __supportedPrimitives[14]= false;
00545
00546 ffout(3) << "\n";
00547 ffout(3) << "Reading file '" << s << "'\n";
00548
00549
00550 GmshFormatReader::Keyword kw = this->__nextKeyword();
00551 switch(kw.second) {
00552 case NOD: {
00553 this->__readGmshFormat1();
00554 break;
00555 }
00556 case MESHFORMAT: {
00557 real_t fileVersion = this->__getReal();
00558 if (fileVersion != 2) {
00559 throw ErrorHandler(__FILE__,__LINE__,
00560 "Cannot read Gmsh format '"+stringify(fileVersion)+"'",
00561 ErrorHandler::normal);
00562 }
00563 int fileType = this->__getInteger();
00564 __binary = (fileType == 1);
00565 if ((fileType < 0) or (fileType > 1)) {
00566 throw ErrorHandler(__FILE__,__LINE__,
00567 "Cannot read Gmsh file type '"+stringify(fileType)+"'",
00568 ErrorHandler::normal);
00569 }
00570
00571 int dataSize = this->__getInteger();
00572 if (dataSize != sizeof(double)) {
00573 throw ErrorHandler(__FILE__,__LINE__,
00574 "Data size not supported '"+stringify(dataSize)+"'",
00575 ErrorHandler::normal);
00576 }
00577
00578 if (__binary) {
00579 int one=0;
00580
00581 fseek(__ifh,1L,SEEK_CUR);
00582 fread(reinterpret_cast<char*>(&one),sizeof(int),1,__ifh);
00583
00584 if (one == 1) {
00585 __convertEndian = false;
00586 } else {
00587 invertEndianess(one);
00588
00589 if (one == 1) {
00590 __convertEndian = true;
00591 } else {
00592 throw ErrorHandler(__FILE__,__LINE__,
00593 "Cannot determine endianess",
00594 ErrorHandler::normal);
00595 }
00596 }
00597
00598 ffout(3) << "- Binary format: ";
00599 #ifdef WORDS_BIGENDIAN
00600 if (not __convertEndian) {
00601 ffout(3) << "big endian\n";
00602 } else {
00603 ffout(3) << "little endian\n";
00604 }
00605 #else // WORDS_BIGENDIAN
00606 if (not __convertEndian) {
00607 ffout(3) << "little endian\n";
00608 } else {
00609 ffout(3) << "big endian\n";
00610 }
00611 #endif // WORDS_BIGENDIAN
00612 }
00613
00614 kw = this->__nextKeyword();
00615 if (kw.second != ENDMESHFORMAT) {
00616 throw ErrorHandler(__FILE__,__LINE__,
00617 "reading file '"+__fileName
00618 +"': expecting $EndMeshFormat, '"+kw.first+"' was found",
00619 ErrorHandler::normal);
00620 }
00621
00622 this->__readGmshFormat2();
00623
00624 break;
00625 }
00626 default: {
00627 throw ErrorHandler(__FILE__,__LINE__,
00628 "cannot determine format version of '"+__fileName+"'",
00629 ErrorHandler::normal);
00630 }
00631 }
00632
00633 this->__proceedData();
00634 this->__createMesh();
00635 }
00636
00637 void GmshFormatReader::
00638 __readGmshFormat1()
00639 {
00640 ffout(3) << "- Reading Gmsh format 1.0\n";
00641 this->__readVertices();
00642
00643 GmshFormatReader::Keyword kw;
00644 kw = this->__nextKeyword();
00645 if (kw.second != ENDNOD) {
00646 throw ErrorHandler(__FILE__,__LINE__,
00647 "reading file '"+__fileName
00648 +"': expecting $ENDNOD, '"+kw.first+"' was found",
00649 ErrorHandler::normal);
00650 }
00651
00652
00653 kw = this->__nextKeyword();
00654 if (kw.second != ELM) {
00655 throw ErrorHandler(__FILE__,__LINE__,
00656 "reading file '"+__fileName
00657 +"': expecting $ELM, '"+kw.first+"' was found",
00658 ErrorHandler::normal);
00659 }
00660
00661 this->__readElements1();
00662
00663 kw = this->__nextKeyword();
00664 if (kw.second != ENDELM) {
00665 throw ErrorHandler(__FILE__,__LINE__,
00666 "reading file '"+__fileName
00667 +"': expecting $ENDELM, '"+kw.first+"' was found",
00668 ErrorHandler::normal);
00669 }
00670 }
00671
00672 void GmshFormatReader::
00673 __readGmshFormat2()
00674 {
00675 ffout(3) << "- Reading Gmsh format 2.0\n";
00676 GmshFormatReader::Keyword kw;
00677 kw = this->__nextKeyword();
00678 if (kw.second != NODES) {
00679 throw ErrorHandler(__FILE__,__LINE__,
00680 "reading file '"+__fileName
00681 +"': expecting $Nodes, '"+kw.first+"' was found",
00682 ErrorHandler::normal);
00683 }
00684
00685 this->__readVertices();
00686
00687 kw = this->__nextKeyword();
00688 if (kw.second != ENDNODES) {
00689 throw ErrorHandler(__FILE__,__LINE__,
00690 "reading file '"+__fileName
00691 +"': expecting $EndNodes, '"+kw.first+"' was found",
00692 ErrorHandler::normal);
00693 }
00694
00695
00696 kw = this->__nextKeyword();
00697 if (kw.second != ELEMENTS) {
00698 throw ErrorHandler(__FILE__,__LINE__,
00699 "reading file '"+__fileName
00700 +"': expecting $Elements, '"+kw.first+"' was found",
00701 ErrorHandler::normal);
00702 }
00703
00704 this->__readElements2();
00705
00706 kw = this->__nextKeyword();
00707 if (kw.second != ENDELEMENTS) {
00708 throw ErrorHandler(__FILE__,__LINE__,
00709 "reading file '"+__fileName
00710 +"': expecting $EndElements, '"+kw.first+"' was found",
00711 ErrorHandler::normal);
00712 }
00713 }