1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585
|
// This file is part of ff3d - http://www.freefem.org/ff3d
// Copyright (C) 2001, 2002, 2003 Stphane Del Pino
// This program is free software; you can redistribute it and/or modify
// it under the terms of the GNU General Public License as published by
// the Free Software Foundation; either version 2, or (at your option)
// any later version.
// This program is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
// GNU General Public License for more details.
// You should have received a copy of the GNU General Public License
// along with this program; if not, write to the Free Software Foundation,
// Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
// $Id: GmshFormatReader.cpp,v 1.10 2007/03/18 19:29:05 delpinux Exp $
#include <GmshFormatReader.hpp>
#include <SurfaceMeshOfQuadrangles.hpp>
#include <MeshOfHexahedra.hpp>
#include <SurfaceMeshOfTriangles.hpp>
#include <MeshOfTetrahedra.hpp>
void
GmshFormatReader::__readVertices()
{
const int numberOfVerices = this->__getInteger();
ffout(3) << "- Number Of Vertices: " << numberOfVerices << '\n';
if (numberOfVerices<0) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+this->__fileName
+"': number of vertices is negative",
ErrorHandler::normal);
}
__verticesNumbers.resize(numberOfVerices);
__vertices = new VerticesSet(numberOfVerices);
VerticesSet& vertices = *__vertices;
for (int i=0; i<numberOfVerices; ++i) {
__verticesNumbers[i] = this->__getInteger();
const real_t x = this->__getReal();
const real_t y = this->__getReal();
const real_t z = this->__getReal();
vertices[i] = TinyVector<3, real_t>(x,y,z);
}
}
void
GmshFormatReader::__readElements1()
{
const int numberOfElements = this->__getInteger();
ffout(3) << "- Number Of Elements: " << numberOfElements << '\n';
if (numberOfElements<0) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': number of elements is negative",
ErrorHandler::normal);
}
__elementType.resize(numberOfElements);
__references.resize(numberOfElements);
__elementVertices.resize(numberOfElements);
for (int i=0; i<numberOfElements; ++i) {
this->__getInteger(); // drops element number
const int elementType = this->__getInteger()-1;
if ((elementType < 0) or (elementType > 14)) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': unknown element type '"+stringify(elementType)+"'",
ErrorHandler::normal);
}
__elementType[i] = elementType;
__references[i] = this->__getInteger(); // physical reference
this->__getInteger(); // drops entity reference
const int numberOfVertices = this->__getInteger();
if (numberOfVertices != __numberOfPrimitiveNodes[elementType]) {
std::stringstream errorMsg;
errorMsg << "reading file '" <<__fileName
<< "':number of vertices (" << numberOfVertices
<< ") does not match expected ("
<< __numberOfPrimitiveNodes[elementType] << ")" << std::ends;
throw ErrorHandler(__FILE__,__LINE__,
errorMsg.str(),
ErrorHandler::normal);
}
__elementVertices[i].resize(numberOfVertices);
for (int j=0; j<numberOfVertices; ++j) {
__elementVertices[i][j] = this->__getInteger();
}
}
}
void
GmshFormatReader::__readElements2()
{
const int numberOfElements = this->__getInteger();
ffout(3) << "- Number Of Elements: " << numberOfElements << '\n';
if (numberOfElements<0) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': number of elements is negative",
ErrorHandler::normal);
}
__elementType.resize(numberOfElements);
__references.resize(numberOfElements);
__elementVertices.resize(numberOfElements);
for (int i=0; i<numberOfElements; ++i) {
this->__getInteger(); // drops element number
const int elementType = this->__getInteger()-1;
if ((elementType < 0) or (elementType > 14)) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': unknown element type '"+stringify(elementType)+"'",
ErrorHandler::normal);
}
__elementType[i] = elementType;
const int numberOfTags = this->__getInteger();
__references[i] = this->__getInteger(); // physical reference
for (int tag=1; tag<numberOfTags; ++tag) {
this->__getInteger(); // drops remaining tags
}
const int numberOfVertices = __numberOfPrimitiveNodes[elementType];
__elementVertices[i].resize(numberOfVertices);
for (int j=0; j<numberOfVertices; ++j) {
__elementVertices[i][j] = this->__getInteger();
}
}
}
void
GmshFormatReader::__proceedData()
{
TinyVector<15,int> elementNumber = 0;
std::vector<std::set<size_t> > elementReferences(15);
for (size_t i=0; i<__elementType.size(); ++i) {
const int elementType = __elementType[i];
elementNumber[elementType]++;
elementReferences[elementType].insert(__references[i]);
}
for (size_t i=0; i<elementNumber.size(); ++i) {
if (elementNumber[i]>0) {
ffout(3) << " - Number of "
<< __primitivesNames[i]
<< ": " << elementNumber[i];
if (not(__supportedPrimitives[i])) {
ffout(3) << " [IGNORED]";
} else {
ffout(3) << " references: ";
for(std::set<size_t>::const_iterator iref
= elementReferences[i].begin();
iref != elementReferences[i].end(); ++iref) {
if (iref != elementReferences[i].begin()) ffout(3) << ',';
ffout(3) << *iref;
}
switch (i) {
// Supported entities
case 1: { // triangles
__triangles = new Vector<Triangle>(elementNumber[i]);
break;
}
case 2: { // quadrangles
__quadrilaterals = new Vector<Quadrangle>(elementNumber[i]);
break;
}
case 3: { // tetrahedra
__tetrahedra = new Vector<Tetrahedron>(elementNumber[i]);
break;
}
case 4: {// hexahedra
__hexahedra = new Vector<Hexahedron>(elementNumber[i]);
}
// Ignored entities
case 0: // edge
case 14: {// point
break;
}
// Unsupported entities
case 5: // prism
case 6: // pyramid
case 7: // second order edge
case 8: // second order triangle
case 9: // second order quadrangle
case 10: // second order tetrahedron
case 11: // second order hexahedron
case 12: // second order prism
case 13: // second order pyramid
default: {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': ff3d cannot read this kind of element",
ErrorHandler::normal);
}
}
}
ffout(3) << '\n';
}
}
ffout(3) << "- Building correspondance table\n";
int minNumber = std::numeric_limits<int>::max();
int maxNumber = std::numeric_limits<int>::min();
for (size_t i=0; i<__verticesNumbers.size(); ++i) {
const int& vertexNumber = __verticesNumbers[i];
minNumber = std::min(minNumber,vertexNumber);
maxNumber = std::max(maxNumber,vertexNumber);
}
if (minNumber<0) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': ff3d does not implement negative vertices number",
ErrorHandler::normal);
}
// A value of -1 means that the vertex is unknown
__verticesCorrepondance.resize(maxNumber+1,-1);
for (size_t i=0; i<__verticesNumbers.size(); ++i) {
__verticesCorrepondance[__verticesNumbers[i]] = i;
}
// reset element number to count them while filling structures
elementNumber = 0;
VerticesSet& vertices = *__vertices;
// Element structures filling
for (size_t i=0; i<__elementType.size(); ++i) {
std::vector<int>& elementVertices = __elementVertices[i];
switch (__elementType[i]) {
// Supported entities
case 1: { // triangles
int& triangleNumber = elementNumber[1];
const size_t a = __verticesCorrepondance[elementVertices[0]];
const size_t b = __verticesCorrepondance[elementVertices[1]];
const size_t c = __verticesCorrepondance[elementVertices[2]];
if ((a<0) or (b<0) or (c<0)) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': error reading element "+stringify(i)
+" [bad vertices definition]",
ErrorHandler::normal);
}
(*__triangles)[triangleNumber]
= Triangle(vertices[a], vertices[b], vertices[c],
__references[i]);
triangleNumber++; // one more triangle
break;
}
case 2: { // quadrangle
int& quadrilateralNumber = elementNumber[2];
const size_t a = __verticesCorrepondance[elementVertices[0]];
const size_t b = __verticesCorrepondance[elementVertices[1]];
const size_t c = __verticesCorrepondance[elementVertices[2]];
const size_t d = __verticesCorrepondance[elementVertices[3]];
if ((a<0) or (b<0) or (c<0) or (d<0)) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': error reading element "+stringify(i)
+" [bad vertices definition]",
ErrorHandler::normal);
}
(*__quadrilaterals)[quadrilateralNumber]
= Quadrangle(vertices[a], vertices[b], vertices[c], vertices[d],
__references[i]);
quadrilateralNumber++; // one more quadrangle
break;
}
case 3: { // tetrahedra
int& tetrahedronNumber = elementNumber[3];
const int a = __verticesCorrepondance[elementVertices[0]];
const int b = __verticesCorrepondance[elementVertices[1]];
const int c = __verticesCorrepondance[elementVertices[2]];
const int d = __verticesCorrepondance[elementVertices[3]];
if ((a<0) or (b<0) or (c<0) or (d<0)) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': error reading element "+stringify(i)
+" [bad vertices definition]",
ErrorHandler::normal);
}
(*__tetrahedra)[tetrahedronNumber]
= Tetrahedron(vertices[a], vertices[b], vertices[c], vertices[d],
__references[i]);
tetrahedronNumber++; // one more tetrahedron
break;
}
case 4: { // hexaredron
int& hexahedronNumber = elementNumber[4];
const int a = __verticesCorrepondance[elementVertices[0]];
const int b = __verticesCorrepondance[elementVertices[1]];
const int c = __verticesCorrepondance[elementVertices[2]];
const int d = __verticesCorrepondance[elementVertices[3]];
const int e = __verticesCorrepondance[elementVertices[4]];
const int f = __verticesCorrepondance[elementVertices[5]];
const int g = __verticesCorrepondance[elementVertices[6]];
const int h = __verticesCorrepondance[elementVertices[7]];
if ((a<0) or (b<0) or (c<0) or (d<0) or (e<0) or (f<0) or (g<0) or (h<0)) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': error reading element "+stringify(i)
+" [bad vertices definition]",
ErrorHandler::normal);
}
(*__hexahedra)[hexahedronNumber]
= Hexahedron(vertices[a], vertices[b], vertices[c], vertices[d],
vertices[e], vertices[f], vertices[g], vertices[h],
__references[i]);
hexahedronNumber++; // one more hexahedron
break; }
// Unsupported entities
case 0: // edge
case 14:{// point
break;
}
case 5: // prism
case 6: // pyramid
case 7: // second order edge
case 8: // second order triangle
case 9: // second order quadrangle
case 10: // second order tetrahedron
case 11: // second order hexahedron
case 12: // second order prism
case 13:{// second order pyramid
}
default: {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': ff3d cannot read this kind of element",
ErrorHandler::normal);
}
}
}
}
GmshFormatReader::Keyword
GmshFormatReader::__nextKeyword()
{
GmshFormatReader::Keyword kw;
char pKeyword[256];
int retval = fscanf(__ifh,"%255s",pKeyword);
const std::string aKeyword = pKeyword;
if (retval < 0) {
kw.second = EndOfFile;
return kw;
} else if (retval == 0) {
kw.second = Unknown;
return kw;
}
KeywordList::iterator i = __keywordList.find(aKeyword.c_str());
if (i != __keywordList.end()) {
kw.first = (*i).first;
kw.second = (*i).second;
return kw;
}
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': unknown keyword "+aKeyword,
ErrorHandler::normal);
kw.first = aKeyword;
kw.second = Unknown;
return kw;
}
GmshFormatReader::GmshFormatReader(const std::string & s)
: MeshReader(s)
{
// gmsh 1.0 format keywords
__keywordList["$NOD"] = NOD;
__keywordList["$ENDNOD"] = ENDNOD;
__keywordList["$ELM"] = ELM;
__keywordList["$ENDELM"] = ENDELM;
// gmsh 2.0 format keywords
__keywordList["$MeshFormat"] = MESHFORMAT;
__keywordList["$EndMeshFormat"] = ENDMESHFORMAT;
__keywordList["$Nodes"] = NODES;
__keywordList["$EndNodes"] = ENDNODES;
__keywordList["$Elements"] = ELEMENTS;
__keywordList["$EndElements"] = ENDELEMENTS;
__numberOfPrimitiveNodes.resize(16);
__numberOfPrimitiveNodes[ 0] = 2; // edge
__numberOfPrimitiveNodes[ 1] = 3; // triangle
__numberOfPrimitiveNodes[ 2] = 4; // quadrangle
__numberOfPrimitiveNodes[ 3] = 4; // Tetrahedron
__numberOfPrimitiveNodes[ 4] = 8; // Hexaredron
__numberOfPrimitiveNodes[ 5] = 6; // Prism
__numberOfPrimitiveNodes[ 6] = 5; // Pyramid
__numberOfPrimitiveNodes[ 7] = 3; // second order edge
__numberOfPrimitiveNodes[ 8] = 6; // second order triangle
__numberOfPrimitiveNodes[ 9] = 9; // second order quadrangle
__numberOfPrimitiveNodes[10] = 10; // second order tetrahedron
__numberOfPrimitiveNodes[11] = 27; // second order hexahedron
__numberOfPrimitiveNodes[12] = 18; // second order prism
__numberOfPrimitiveNodes[13] = 14; // second order pyramid
__numberOfPrimitiveNodes[14] = 1; // point
__primitivesNames[0] = "edges";
__supportedPrimitives[0] = false;
__primitivesNames[1] = "triangles";
__supportedPrimitives[1] = true;
__primitivesNames[2] = "quadrangles";
__supportedPrimitives[2] = true;
__primitivesNames[3] = "tetrahedra";
__supportedPrimitives[3] = true;
__primitivesNames[4] = "hexahedra";
__supportedPrimitives[4] = true;
__primitivesNames[5] = "prisms";
__supportedPrimitives[5] = false;
__primitivesNames[6] = "pyramids";
__supportedPrimitives[6] = false;
__primitivesNames[7] = "second order edges";
__supportedPrimitives[7] = false;
__primitivesNames[8] = "second order triangles";
__supportedPrimitives[8] = false;
__primitivesNames[9] = "second order quadrangles";
__supportedPrimitives[9] = false;
__primitivesNames[10] = "second order tetrahedra";
__supportedPrimitives[10]= false;
__primitivesNames[11] = "second order hexahedra";
__supportedPrimitives[11]= false;
__primitivesNames[12] = "second order prisms";
__supportedPrimitives[12]= false;
__primitivesNames[13] = "second order pyramids";
__supportedPrimitives[13]= false;
__primitivesNames[14] = "point";
__supportedPrimitives[14]= false;
ffout(3) << "\n";
ffout(3) << "Reading file '" << s << "'\n";
// Getting vertices list
GmshFormatReader::Keyword kw = this->__nextKeyword();
switch(kw.second) {
case NOD: {
this->__readGmshFormat1();
break;
}
case MESHFORMAT: {
real_t fileVersion = this->__getReal();
if (fileVersion != 2) {
throw ErrorHandler(__FILE__,__LINE__,
"Cannot read Gmsh format '"+stringify(fileVersion)+"'",
ErrorHandler::normal);
}
int fileType = this->__getInteger();
if (fileType != 0) {
throw ErrorHandler(__FILE__,__LINE__,
"Cannot read Gmsh file type '"+stringify(fileType)+"'",
ErrorHandler::normal);
}
int dataSize = this->__getInteger();
if (dataSize != sizeof(double)) {
throw ErrorHandler(__FILE__,__LINE__,
"Data size not supported '"+stringify(dataSize)+"'",
ErrorHandler::normal);
}
kw = this->__nextKeyword();
if (kw.second != ENDMESHFORMAT) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': expecting $EndMeshFormat, '"+kw.first+"' was found",
ErrorHandler::normal);
}
this->__readGmshFormat2();
break;
}
default: {
throw ErrorHandler(__FILE__,__LINE__,
"cannot determine format version of '"+__fileName+"'",
ErrorHandler::normal);
}
}
this->__proceedData();
this->__createMesh();
}
void GmshFormatReader::
__readGmshFormat1()
{
ffout(3) << "- Reading Gmsh format 1.0\n";
this->__readVertices();
GmshFormatReader::Keyword kw;
kw = this->__nextKeyword();
if (kw.second != ENDNOD) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': expecting $ENDNOD, '"+kw.first+"' was found",
ErrorHandler::normal);
}
// Getting elements list
kw = this->__nextKeyword();
if (kw.second != ELM) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': expecting $ELM, '"+kw.first+"' was found",
ErrorHandler::normal);
}
this->__readElements1();
kw = this->__nextKeyword();
if (kw.second != ENDELM) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': expecting $ENDELM, '"+kw.first+"' was found",
ErrorHandler::normal);
}
}
void GmshFormatReader::
__readGmshFormat2()
{
ffout(3) << "- Reading Gmsh format 2.0\n";
GmshFormatReader::Keyword kw;
kw = this->__nextKeyword();
if (kw.second != NODES) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': expecting $Nodes, '"+kw.first+"' was found",
ErrorHandler::normal);
}
this->__readVertices();
kw = this->__nextKeyword();
if (kw.second != ENDNODES) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': expecting $EndNodes, '"+kw.first+"' was found",
ErrorHandler::normal);
}
// Getting elements list
kw = this->__nextKeyword();
if (kw.second != ELEMENTS) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': expecting $Elements, '"+kw.first+"' was found",
ErrorHandler::normal);
}
this->__readElements2();
kw = this->__nextKeyword();
if (kw.second != ENDELEMENTS) {
throw ErrorHandler(__FILE__,__LINE__,
"reading file '"+__fileName
+"': expecting $EndElements, '"+kw.first+"' was found",
ErrorHandler::normal);
}
}
|