1 // Gmsh - Copyright (C) 1997-2021 C. Geuzaine, J.-F. Remacle 2 // 3 // See the LICENSE.txt file in the Gmsh root directory for license information. 4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues. 5 6 #ifndef GMODEL_H 7 #define GMODEL_H 8 9 #include <algorithm> 10 #include <vector> 11 #include <set> 12 #include <map> 13 #include <string> 14 #include <unordered_map> 15 #include <functional> 16 #include "GVertex.h" 17 #include "GEdge.h" 18 #include "GFace.h" 19 #include "GRegion.h" 20 #include "SPoint3.h" 21 #include "SBoundingBox3d.h" 22 #include "MFaceHash.h" 23 #include "MEdgeHash.h" 24 25 #define hashmapMFace \ 26 std::unordered_map<MFace, std::size_t, MFaceHash, MFaceEqual> 27 #define hashmapMEdge \ 28 std::unordered_map<MEdge, std::size_t, MEdgeHash, MEdgeEqual> 29 30 template <class scalar> class simpleFunction; 31 32 class GEO_Internals; 33 class OCC_Internals; 34 class ACIS_Internals; 35 class Parasolid_Internals; 36 class smooth_normals; 37 class FieldManager; 38 class gLevelset; 39 class discreteFace; 40 class discreteRegion; 41 class MElementOctree; 42 43 // A geometric model. The model is a "not yet" non-manifold B-Rep. 44 class GModel { 45 private: 46 std::multimap<std::pair<const std::vector<int>, const std::vector<int> >, 47 std::pair<const std::string, const std::vector<int> > > 48 _homologyRequests; 49 std::set<GRegion *, GEntityPtrLessThan> _chainRegions; 50 std::set<GFace *, GEntityPtrLessThan> _chainFaces; 51 std::set<GEdge *, GEntityPtrLessThan> _chainEdges; 52 std::set<GVertex *, GEntityPtrLessThan> _chainVertices; 53 hashmapMEdge _mapEdgeNum; 54 hashmapMFace _mapFaceNum; 55 // the maximum vertex and element id number in the mesh 56 std::size_t _maxVertexNum, _maxElementNum; 57 std::size_t _checkPointedMaxVertexNum, _checkPointedMaxElementNum; 58 59 private: 60 int _readMSH2(const std::string &name); 61 int _writeMSH2(const std::string &name, double version, bool binary, 62 bool saveAll, bool saveParametric, double scalingFactor, 63 int elementStartNum, int saveSinglePartition, bool append, 64 bool renumberVertices); 65 int _writePartitionedMSH2(const std::string &baseName, bool binary, 66 bool saveAll, bool saveParametric, 67 double scalingFactor); 68 int _readMSH3(const std::string &name); 69 int _writeMSH3(const std::string &name, double version, bool binary, 70 bool saveAll, bool saveParametric, double scalingFactor, 71 int elementStartNum, int saveSinglePartition, bool append); 72 int _writePartitionedMSH3(const std::string &baseName, double version, 73 bool binary, bool saveAll, bool saveParametric, 74 double scalingFactor); 75 int _readMSH4(const std::string &name); 76 int _writeMSH4(const std::string &name, double version, bool binary, 77 bool saveAll, bool saveParametric, double scalingFactor, 78 bool append, int partitionToSave = 0); 79 int _writePartitionedMSH4(const std::string &baseName, double version, 80 bool binary, bool saveAll, bool saveParametric, 81 double scalingFactor); 82 83 protected: 84 // the name of the model 85 std::string _name; 86 87 // the name of the file the model was read from 88 std::string _fileName; 89 std::set<std::string> _fileNames; 90 91 // the visibility flag 92 char _visible; 93 94 // vertex and element caches to speed-up direct access by tag (mostly 95 // used for post-processing I/O) 96 std::vector<MVertex *> _vertexVectorCache; 97 std::map<int, MVertex *> _vertexMapCache; 98 std::vector<std::pair<MElement *, int> > _elementVectorCache; 99 std::map<int, std::pair<MElement *, int> > _elementMapCache; 100 std::map<int, int> _elementIndexCache; 101 102 // ghost cell information (stores partitions for each element acting 103 // as a ghost cell) 104 // /!\ Use only for compatibility with mesh format msh2 and msh3 105 std::multimap<MElement *, short> _ghostCells; 106 107 // an octree for fast mesh element lookup 108 MElementOctree *_elementOctree; 109 110 // global cache storage of discrete curvatures 111 std::map<MVertex *, std::pair<SVector3, SVector3> > _curvatures; 112 113 // GEO (Gmsh native) model internal data 114 GEO_Internals *_geo_internals; 115 // OpenCASCADE model internal data 116 OCC_Internals *_occ_internals; 117 // ACIS model internal data 118 ACIS_Internals *_acis_internals; 119 // Parasolid model internal data 120 Parasolid_Internals *_parasolid_internals; 121 122 // characteristic length (mesh size) fields 123 FieldManager *_fields; 124 125 // entity that is currently being meshed (used for error reporting) 126 GEntity *_currentMeshEntity; 127 128 // last entities/vertices where a meshing error has been reported 129 std::vector<GEntity *> _lastMeshEntityError; 130 std::vector<MVertex *> _lastMeshVertexError; 131 132 // index of the current model (in the static list of all loaded 133 // models) 134 static int _current; 135 136 // the sets of geometrical regions, faces, edges and vertices in the 137 // model 138 std::set<GRegion *, GEntityPtrLessThan> regions; 139 std::set<GFace *, GEntityPtrLessThan> faces; 140 std::set<GEdge *, GEntityPtrLessThan> edges; 141 std::set<GVertex *, GEntityPtrLessThan> vertices; 142 143 // map between the pair <dimension, elementary or physical number> 144 // and an optional associated name 145 std::map<std::pair<int, int>, std::string> _physicalNames, _elementaryNames; 146 147 // the set of all used mesh partition numbers 148 std::size_t _numPartitions; 149 150 protected: 151 // store the elements given in the map (indexed by elementary region 152 // number) into the model, creating discrete geometrical entities on 153 // the fly if needed 154 void _storeElementsInEntities(std::map<int, std::vector<MElement *> > &map); 155 156 // store the parent's pointer back into MSubElements (replacing numeric id) 157 void _storeParentsInSubElements(std::map<int, std::vector<MElement *> > &map); 158 159 // loop over all vertices connected to elements and associate 160 // geometrical entity 161 void _associateEntityWithMeshVertices(); 162 163 // store the vertices in the geometrical entity they are associated 164 // with, and delete those that are not associated with any entity 165 void _storeVerticesInEntities(std::map<int, MVertex *> &vertices); 166 void _storeVerticesInEntities(std::vector<MVertex *> &vertices); 167 168 // store the physical tags in the geometrical entities 169 void 170 _storePhysicalTagsInEntities(int dim, 171 std::map<int, std::map<int, std::string> > &map); 172 173 public: 174 // region, face, edge and vertex iterators 175 typedef std::set<GRegion *, GEntityPtrLessThan>::iterator riter; 176 typedef std::set<GFace *, GEntityPtrLessThan>::iterator fiter; 177 typedef std::set<GEdge *, GEntityPtrLessThan>::iterator eiter; 178 typedef std::set<GVertex *, GEntityPtrLessThan>::iterator viter; 179 180 typedef std::set<GRegion *, GEntityPtrLessThan>::const_iterator const_riter; 181 typedef std::set<GFace *, GEntityPtrLessThan>::const_iterator const_fiter; 182 typedef std::set<GEdge *, GEntityPtrLessThan>::const_iterator const_eiter; 183 typedef std::set<GVertex *, GEntityPtrLessThan>::const_iterator const_viter; 184 185 // elementary/physical name iterator 186 typedef std::map<std::pair<int, int>, std::string>::iterator piter; 187 188 public: 189 GModel(const std::string &name = ""); 190 virtual ~GModel(); 191 // Required for python bindings on Windows 192 #ifndef SWIG 193 // the static list of all loaded models 194 static std::vector<GModel *> list; 195 #endif 196 197 // return the current model, and sets the current model index if 198 // index >= 0 199 static GModel *current(int index = -1); 200 201 // sets a model to current 202 static int setCurrent(GModel *m); setAsCurrent()203 int setAsCurrent() { return setCurrent(this); } 204 205 // find a model by name; if fileName is given, return model only if it does 206 // *not* have a link to the fileName 207 static GModel *findByName(const std::string &name, 208 const std::string &fileName = ""); 209 210 // delete everything in a GModel (optionally keep name and fileName) 211 void destroy(bool keepName = false); 212 213 // get/set global vertex/element num getMaxVertexNumber()214 std::size_t getMaxVertexNumber() const { return _maxVertexNum; } getMaxElementNumber()215 std::size_t getMaxElementNumber() const { return _maxElementNum; } setMaxVertexNumber(std::size_t num)216 void setMaxVertexNumber(std::size_t num) 217 { 218 #pragma omp atomic write 219 _maxVertexNum = _maxVertexNum > num ? _maxVertexNum : num; 220 } setMaxElementNumber(std::size_t num)221 void setMaxElementNumber(std::size_t num) 222 { 223 #pragma omp atomic write 224 _maxElementNum = _maxElementNum > num ? _maxElementNum : num; 225 } 226 227 // increment and get global vertex/element num incrementAndGetMaxVertexNumber()228 std::size_t incrementAndGetMaxVertexNumber() 229 { 230 std::size_t _myVertexNum; 231 #pragma omp atomic capture 232 { 233 ++_maxVertexNum; 234 _myVertexNum = _maxVertexNum; 235 } 236 return _myVertexNum; 237 } incrementAndGetMaxElementNumber()238 std::size_t incrementAndGetMaxElementNumber() 239 { 240 std::size_t _myElementNum; 241 #pragma omp atomic capture 242 { 243 ++_maxElementNum; 244 _myElementNum = _maxElementNum; 245 } 246 return _myElementNum; 247 } 248 249 // decrement global vertex num decrementMaxVertexNumber()250 void decrementMaxVertexNumber() 251 { 252 #pragma omp atomic update 253 --_maxVertexNum; 254 } 255 checkPointMaxNumbers()256 void checkPointMaxNumbers() 257 { 258 _checkPointedMaxVertexNum = _maxVertexNum; 259 _checkPointedMaxVertexNum = _maxVertexNum; 260 } getCheckPointedMaxNumbers(std::size_t & maxv,std::size_t & maxe)261 void getCheckPointedMaxNumbers(std::size_t &maxv, std::size_t &maxe) const 262 { 263 maxv = _checkPointedMaxVertexNum; 264 maxe = _checkPointedMaxElementNum; 265 } 266 267 // add a mesh edge or face in the global edge or face map with number "num", 268 // or number it (starting at 1) if num == 0 269 std::size_t addMEdge(MEdge &edge, std::size_t num = 0); 270 std::size_t addMFace(MFace &face, std::size_t num = 0); 271 // get the edge of face and its global number given mesh nodes (return 0 if 272 // the edge or face does not exist in the edge or face map) 273 std::size_t getMEdge(MVertex *v0, MVertex *v1, MEdge &edge); 274 std::size_t getMFace(MVertex *v0, MVertex *v1, MVertex *v2, MVertex *v3, 275 MFace &face); 276 // iterate on edges and faces firstMEdge()277 hashmapMEdge::const_iterator firstMEdge() { return _mapEdgeNum.begin(); } lastMEdge()278 hashmapMEdge::const_iterator lastMEdge() { return _mapEdgeNum.end(); } firstMFace()279 hashmapMFace::const_iterator firstMFace() { return _mapFaceNum.begin(); } lastMFace()280 hashmapMFace::const_iterator lastMFace() { return _mapFaceNum.end(); } 281 282 // renumber mesh vertices and elements in a continuous sequence (this 283 // invalidates the mesh caches) 284 void renumberMeshVertices(); 285 void renumberMeshElements(); 286 287 // delete all the mesh-related caches (this must be called when the 288 // mesh is changed) 289 void destroyMeshCaches(); 290 // delete the mesh stored in entities and call destroMeshCaches 291 void deleteMesh(); 292 void deleteMesh(const std::vector<GEntity *> &entities); 293 // delete the vertex arrays used for efficient mesh drawing 294 void deleteVertexArrays(); 295 296 // remove all mesh vertex associations to geometrical entities and remove 297 // vertices from geometrical entities, then _associateEntityWithMeshVertices 298 // and _storeVerticesInEntities are called to rebuild the associations 299 void pruneMeshVertexAssociations(); 300 301 // access internal CAD representations 302 void createGEOInternals(); 303 void createOCCInternals(); 304 void createACISInternals(); 305 void createParasolidInternals(); getOCCInternals()306 OCC_Internals *getOCCInternals() { return _occ_internals; } getGEOInternals()307 GEO_Internals *getGEOInternals() { return _geo_internals; } getACISInternals()308 ACIS_Internals *getACISInternals() { return _acis_internals; } getParasolidInternals()309 Parasolid_Internals *getParasolidInternals() { return _parasolid_internals; } 310 void deleteGEOInternals(); 311 void deleteOCCInternals(); 312 void resetOCCInternals(); 313 void deleteACISInternals(); 314 void deleteParasolidInternals(); 315 316 // access characteristic length (mesh size) fields getFields()317 FieldManager *getFields() { return _fields; } 318 319 // get/set the model name setName(const std::string & name)320 void setName(const std::string &name) { _name = name; } getName()321 std::string getName() const { return _name; } 322 323 // get/set the model file name 324 void setFileName(const std::string &fileName); getFileName()325 std::string getFileName() const { return _fileName; } hasFileName(const std::string & name)326 bool hasFileName(const std::string &name) const 327 { 328 return _fileNames.find(name) != _fileNames.end(); 329 } 330 331 // get/set the visibility flag getVisibility()332 char getVisibility() const { return _visible; } setVisibility(char val)333 void setVisibility(char val) { _visible = val; } 334 335 // get the number of entities in this model getNumRegions()336 std::size_t getNumRegions() const { return regions.size(); } getNumFaces()337 std::size_t getNumFaces() const { return faces.size(); } getNumEdges()338 std::size_t getNumEdges() const { return edges.size(); } getNumVertices()339 std::size_t getNumVertices() const { return vertices.size(); } 340 341 // quickly check if the model is empty (i.e., if it contains no 342 // entities) 343 bool empty() const; 344 345 // get an iterator initialized to the first/last entity in this model firstRegion()346 riter firstRegion() { return regions.begin(); } firstFace()347 fiter firstFace() { return faces.begin(); } firstEdge()348 eiter firstEdge() { return edges.begin(); } firstVertex()349 viter firstVertex() { return vertices.begin(); } lastRegion()350 riter lastRegion() { return regions.end(); } lastFace()351 fiter lastFace() { return faces.end(); } lastEdge()352 eiter lastEdge() { return edges.end(); } lastVertex()353 viter lastVertex() { return vertices.end(); } firstRegion()354 const_riter firstRegion() const { return regions.begin(); } firstFace()355 const_fiter firstFace() const { return faces.begin(); } firstEdge()356 const_eiter firstEdge() const { return edges.begin(); } firstVertex()357 const_viter firstVertex() const { return vertices.begin(); } lastRegion()358 const_riter lastRegion() const { return regions.end(); } lastFace()359 const_fiter lastFace() const { return faces.end(); } lastEdge()360 const_eiter lastEdge() const { return edges.end(); } lastVertex()361 const_viter lastVertex() const { return vertices.end(); } 362 363 // get the set of entities getRegions()364 std::set<GRegion *, GEntityPtrLessThan> getRegions() const 365 { 366 return regions; 367 }; getFaces()368 std::set<GFace *, GEntityPtrLessThan> getFaces() const { return faces; }; getEdges()369 std::set<GEdge *, GEntityPtrLessThan> getEdges() const { return edges; }; getVertices()370 std::set<GVertex *, GEntityPtrLessThan> getVertices() const 371 { 372 return vertices; 373 }; 374 375 // find the entity with the given tag 376 GRegion *getRegionByTag(int n) const; 377 GFace *getFaceByTag(int n) const; 378 GEdge *getEdgeByTag(int n) const; 379 GVertex *getVertexByTag(int n) const; 380 GEntity *getEntityByTag(int dim, int n) const; 381 382 // change entity tag (modifies the model entity sets) 383 bool changeEntityTag(int dim, int tag, int newTag); 384 385 // add/remove an entity in the model add(GRegion * r)386 bool add(GRegion *r) { return regions.insert(r).second; } add(GFace * f)387 bool add(GFace *f) { return faces.insert(f).second; } add(GEdge * e)388 bool add(GEdge *e) { return edges.insert(e).second; } add(GVertex * v)389 bool add(GVertex *v) { return vertices.insert(v).second; } 390 bool remove(GRegion *r); 391 bool remove(GFace *f); 392 bool remove(GEdge *e); 393 bool remove(GVertex *v); 394 void remove(int dim, int tag, std::vector<GEntity*> &removed, 395 bool recursive = false); 396 void remove(const std::vector<std::pair<int, int> > &dimTags, 397 std::vector<GEntity*> &removed, 398 bool recursive = false); 399 void remove(); 400 401 // snap vertices on model edges by using geometry tolerance 402 void snapVertices(); 403 404 // fill a vector containing all the entities in the model 405 void getEntities(std::vector<GEntity *> &entities, int dim = -1) const; 406 407 // fill a vector containing all the entities in a given bounding box 408 void getEntitiesInBox(std::vector<GEntity *> &entities, 409 const SBoundingBox3d &box, int dim = -1) const; 410 411 // get tags of entities of the boundary of the given input entities 412 bool getBoundaryTags(const std::vector<std::pair<int, int> > &inDimTags, 413 std::vector<std::pair<int, int> > &outDimTags, 414 bool combined, bool oriented = true, 415 bool recursive = false); 416 417 // return the highest number associated with an elementary entity of 418 // a given dimension (or the highest overall if dim < 0) 419 int getMaxElementaryNumber(int dim); 420 421 // check if there are no physical entities in the model 422 bool noPhysicalGroups(); 423 424 // return all physical groups (one map per dimension: 0-D to 3-D) 425 void 426 getPhysicalGroups(std::map<int, std::vector<GEntity *> > groups[4]) const; 427 void getPhysicalGroups(int dim, 428 std::map<int, std::vector<GEntity *> > &groups) const; getPhysicalNames()429 const std::map<std::pair<int, int>, std::string> &getPhysicalNames() const 430 { 431 return _physicalNames; 432 } setPhysicalNames(const std::map<std::pair<int,int>,std::string> & names)433 void setPhysicalNames(const std::map<std::pair<int, int>, std::string> &names) 434 { 435 _physicalNames = names; 436 } 437 438 // add a physical group (made of elementary entities "tags") 439 void addPhysicalGroup(int dim, int tag, const std::vector<int> &tags); 440 441 // remove physical groups 442 void removePhysicalGroups(); 443 void removePhysicalGroup(int dim, int num); 444 445 // return the highest number associated with a physical entity of a 446 // given dimension (or highest for all dimenions if dim < 0) 447 int getMaxPhysicalNumber(int dim); 448 449 // get an iterator on the elementary/physical names firstPhysicalName()450 piter firstPhysicalName() { return _physicalNames.begin(); } lastPhysicalName()451 piter lastPhysicalName() { return _physicalNames.end(); } firstElementaryName()452 piter firstElementaryName() { return _elementaryNames.begin(); } lastElementaryName()453 piter lastElementaryName() { return _elementaryNames.end(); } 454 455 // get the number of physical names numPhysicalNames()456 int numPhysicalNames() const { return (int)_physicalNames.size(); } 457 458 // get iterators to the last physical name of each dimension 459 void getInnerPhysicalNamesIterators(std::vector<piter> &iterators); 460 461 // associate a name with a physical entity of dimension "dim" and 462 // number "num" (returns a new number id if "num"==0) 463 int setPhysicalName(const std::string &name, int dim, int num = 0); 464 piter setPhysicalName(piter pos, const std::string &name, int dim, 465 int num = 0); 466 467 // get the name (if any) of a given physical group of dimension 468 // "dim" and id number "num" 469 std::string getPhysicalName(int dim, int num) const; 470 471 // remove physical name(s) 472 void removePhysicalName(const std::string &name); 473 474 // get the number of a given physical group of dimension 475 // "dim" and name "name". return -1 if not found 476 int getPhysicalNumber(const int &dim, const std::string &name); 477 478 // get all tags of elementary entities associated with a given physical group 479 // name 480 std::vector<int> getTagsForPhysicalName(int dim, const std::string &name); 481 482 // set physical tags to entities in a given bounding box 483 void setPhysicalNumToEntitiesInBox(int EntityDimension, int PhysicalNumber, 484 std::vector<double> p1, 485 std::vector<double> p2); 486 void setPhysicalNumToEntitiesInBox(int EntityDimension, int PhysicalNumber, 487 const SBoundingBox3d &box); 488 489 // get the name (if any) of a given elementary entity of dimension 490 // "dim" and id number "num" 491 std::string getElementaryName(int dim, int tag); setElementaryName(int dim,int tag,const std::string & name)492 void setElementaryName(int dim, int tag, const std::string &name) 493 { 494 _elementaryNames[std::make_pair(dim, tag)] = name; 495 } 496 497 // remove elememtary name(s) 498 void removeElementaryName(const std::string &name); 499 500 // get the highest dimension of the GModel 501 int getDim() const; 502 503 // get the highest dimension of the mesh in the GModel 504 int getMeshDim() const; 505 506 // set the selection flag on all entities 507 void setSelection(int val); 508 509 // the bounding box 510 SBoundingBox3d bounds(bool aroundVisible = false); 511 512 // return the mesh status for the entire model 513 int getMeshStatus(bool countDiscrete = true); 514 515 // return the total number of elements in the mesh 516 std::size_t getNumMeshElements(int dim = -1) const; 517 std::size_t getNumMeshParentElements() const; 518 519 // get the number of each type of element in the mesh at the largest 520 // dimension and return the dimension 521 std::size_t getNumMeshElements(unsigned c[6]); 522 523 // access a mesh element by coordinates (using an octree search) 524 MElement *getMeshElementByCoord(SPoint3 &p, SPoint3 ¶m, int dim = -1, 525 bool strict = true); 526 std::vector<MElement *> getMeshElementsByCoord(SPoint3 &p, int dim = -1, 527 bool strict = true); 528 529 // access a mesh element by tag, using the element cache getMeshElementByTag(int n)530 MElement *getMeshElementByTag(int n) 531 { 532 int tag; 533 return getMeshElementByTag(n, tag); 534 } 535 MElement *getMeshElementByTag(int n, int &entityTag); 536 537 // access temporary mesh element index 538 int getMeshElementIndex(MElement *e); 539 void setMeshElementIndex(MElement *e, int index); 540 541 // return the total number of vertices in the mesh 542 std::size_t getNumMeshVertices(int dim = -1) const; 543 544 // recompute _vertexVectorCache if there is a dense vertex numbering or 545 // _vertexMapCache if not. 546 void rebuildMeshVertexCache(bool onlyIfNecessary = false); 547 548 // recompute _elementVectorCache if there is a dense element numbering or 549 // _elementMapCache if not. 550 void rebuildMeshElementCache(bool onlyIfNecessary = false); 551 552 // access a mesh vertex by tag, using the vertex cache 553 MVertex *getMeshVertexByTag(int n); 554 555 // add a mesh vertex to the global mesh vertex cache 556 void addMVertexToVertexCache(MVertex* v); 557 558 // get all the mesh vertices associated with the physical group 559 // of dimension "dim" and id number "num" 560 void getMeshVerticesForPhysicalGroup(int dim, int num, 561 std::vector<MVertex *> &); 562 563 // index all the (used) mesh vertices in a continuous sequence, 564 // starting at 1 565 std::size_t indexMeshVertices(bool all, int singlePartition = 0, 566 bool renumber = true); 567 568 // scale the mesh by the given factor 569 void scaleMesh(double factor); 570 571 // set/get entity that is currently being meshed (for error reporting) 572 void setCurrentMeshEntity(GEntity *e); getCurrentMeshEntity()573 GEntity *getCurrentMeshEntity() { return _currentMeshEntity; } 574 575 // set/get entities/vertices linked meshing errors clearLastMeshEntityError()576 void clearLastMeshEntityError() { _lastMeshEntityError.clear(); } addLastMeshEntityError(GEntity * e)577 void addLastMeshEntityError(GEntity *e) { _lastMeshEntityError.push_back(e); } getLastMeshEntityError()578 std::vector<GEntity *> getLastMeshEntityError() 579 { 580 return _lastMeshEntityError; 581 } clearLastMeshVertexError()582 void clearLastMeshVertexError() { _lastMeshVertexError.clear(); } addLastMeshVertexError(MVertex * v)583 void addLastMeshVertexError(MVertex *v) { _lastMeshVertexError.push_back(v); } getLastMeshVertexError()584 std::vector<MVertex *> getLastMeshVertexError() 585 { 586 return _lastMeshVertexError; 587 } 588 589 // delete or reverse all invisble mesh elements 590 std::size_t removeInvisibleElements(); 591 std::size_t reverseInvisibleElements(); 592 593 // the list of partitions getNumPartitions()594 std::size_t getNumPartitions() const { return _numPartitions; } setNumPartitions(std::size_t npart)595 void setNumPartitions(std::size_t npart) { _numPartitions = npart; } 596 597 // partition the mesh 598 int partitionMesh(int num, 599 std::vector<std::pair<MElement *, int> > elementPartition = 600 std::vector<std::pair<MElement *, int> >()); 601 // unpartition the mesh 602 int unpartitionMesh(); 603 // import a mesh partitionned by a tag given by element (i.e. the old way we 604 // stored partitions) and create the new topology-based partition entitiesx 605 int convertOldPartitioningToNewOne(); 606 // write the partitioned topology file 607 int writePartitionedTopology(std::string &name); 608 609 // /!\ Use only for compatibility with mesh format msh2 and msh3 getGhostCells()610 std::multimap<MElement *, short> &getGhostCells() { return _ghostCells; } addGhostCells(MElement * elm,short partition)611 void addGhostCells(MElement *elm, short partition) 612 { 613 _ghostCells.insert(std::make_pair(elm, partition)); 614 } 615 616 // perform various coherence tests on the mesh 617 void checkMeshCoherence(double tolerance); 618 619 // remove duplicate mesh vertices 620 int removeDuplicateMeshVertices(double tolerance); 621 622 // create a geometry (i.e. a parametrization for curves and surfaces) for the 623 // given discrete entities (or all of them if dimTags is empty) 624 void createGeometryOfDiscreteEntities( 625 const std::vector<std::pair<int, int> > &dimTags = 626 std::vector<std::pair<int, int> >()); 627 628 // make discrete entities simply connected 629 void makeDiscreteRegionsSimplyConnected(); 630 void makeDiscreteFacesSimplyConnected(); 631 632 // create topology from mesh 633 void createTopologyFromMesh(); 634 635 // align periodic boundaries 636 void alignPeriodicBoundaries(); 637 638 // a container for smooth normals 639 smooth_normals *normals; 640 641 // mesh the model 642 int mesh(int dimension); 643 644 // adapt 3d mesh 645 int adaptMesh(); 646 647 // adapt the mesh anisotropically using metrics that are computed from a set 648 // of functions f(x,y,z). The algorithm first generate a mesh if no one is 649 // available; see the cpp for parameter documentation 650 int adaptMesh(std::vector<int> technique, 651 std::vector<simpleFunction<double> *> f, 652 std::vector<std::vector<double> > parameters, int niter, 653 bool meshAll = false); 654 655 // ensure that the Jacobian of all volume elements is positive 656 bool setAllVolumesPositive(); 657 void setAllVolumesPositiveTopology(); 658 659 // make the mesh a high order mesh at order N (linear is 1 if the high order 660 // points are not placed on the geometry of the model; incomplete is 1 if 661 // incomplete basis are used) 662 int setOrderN(int order, int linear, int incomplete); 663 664 // refine the mesh by splitting all elements 665 int refineMesh(int linear, bool splitIntoQuads = false, 666 bool splitIntoHexas = false, bool barycentric = false); 667 668 // optimize the mesh 669 int optimizeMesh(const std::string &how, bool force = false, int niter = 1); 670 671 // recombine the mesh 672 int recombineMesh(); 673 674 // fill the vertex arrays, given the current option and data 675 bool fillVertexArrays(); 676 677 // reclassify a surface mesh, using an angle threshold to tag edges and faces 678 void classifySurfaces(double angleThreshold, bool includeBoundary, 679 bool forReparametrization, double curveAngleThreshold); 680 681 // build a new GModel by cutting the elements crossed by the levelset ls 682 // if cutElem is set to false, split the model without cutting the elements 683 GModel *buildCutGModel(gLevelset *ls, bool cutElem = true, 684 bool saveTri = false); 685 686 // store mesh elements of a chain in a new elementary and physical entity 687 void storeChain(int dim, std::map<int, std::vector<MElement *> > &entityMap, 688 std::map<int, std::map<int, std::string> > &physicalMap); 689 690 // request homology computation 691 void addHomologyRequest(const std::string &type, 692 const std::vector<int> &domain, 693 const std::vector<int> &subdomain, 694 const std::vector<int> &dim); 695 void computeHomology(); 696 697 // mesh size callback 698 std::function<double(int, int, double, double, double, double)> lcCallback; 699 700 // compute automatic sizing field from curvature 701 void computeSizeField(); 702 703 // access global cache of discrete curvatures getCurvatures()704 std::map<MVertex *, std::pair<SVector3, SVector3> > &getCurvatures() 705 { 706 return _curvatures; 707 } 708 709 // "automatic" IO based on Gmsh global functions 710 void load(const std::string &fileName); 711 void save(const std::string &fileName); 712 713 // Gmsh native CAD format (readGEO is static, since it can create 714 // multiple models) 715 static int readGEO(const std::string &name); 716 int writeGEO(const std::string &name, bool printLabels = true, 717 bool onlyPhysicals = false); 718 int exportDiscreteGEOInternals(); 719 720 // OCC model 721 int readOCCBREP(const std::string &name); 722 int readOCCSTEP(const std::string &name); 723 int readOCCIGES(const std::string &name); 724 int writeOCCSTEP(const std::string &name); 725 int writeOCCBREP(const std::string &name); 726 int importOCCShape(const void *shape); 727 GVertex *getVertexForOCCShape(const void *shape); 728 GEdge *getEdgeForOCCShape(const void *shape); 729 GFace *getFaceForOCCShape(const void *shape); 730 GRegion *getRegionForOCCShape(const void *shape); 731 732 // ACIS Model 733 int readACISSAT(const std::string &name); 734 735 // Parasolid Model 736 int readParasolidXMT(const std::string &name); 737 int writeParasolidXMT(const std::string &name); 738 int readParasolidSTEP(const std::string &name); 739 int writeParasolidSTEP(const std::string &name); 740 741 // Gmsh mesh file format 742 int readMSH(const std::string &name); 743 int writeMSH(const std::string &name, double version = 2.2, 744 bool binary = false, bool saveAll = false, 745 bool saveParametric = false, double scalingFactor = 1.0, 746 int elementStartNum = 0, int saveSinglePartition = 0, 747 bool append = false); 748 int writePartitionedMSH(const std::string &baseName, double version = 2.2, 749 bool binary = false, bool saveAll = false, 750 bool saveParametric = false, 751 double scalingFactor = 1.0); 752 753 // Iridium file format 754 int writeIR3(const std::string &name, int elementTagType, bool saveAll, 755 double scalingFactor); 756 757 // mesh statistics (saved as a Gmsh post-processing view) 758 int writePOS(const std::string &name, bool printElementary, 759 bool printElementNumber, bool printSICN, bool printSIGE, 760 bool printGamma, bool printDisto, bool saveAll = false, 761 double scalingFactor = 1.0); 762 763 // Stereo lithography format 764 int readSTL(const std::string &name, double tolerance = 1.e-3); 765 int writeSTL(const std::string &name, bool binary = false, 766 bool saveAll = false, double scalingFactor = 1.0, 767 int oneSolidPerSurface = 0); 768 769 // X3D (only output from OCCT's triangulation) 770 int writeX3D(const std::string &name, bool saveAll = false, 771 double scalingFactor = 1.0, int x3dsurfaces = 1, 772 int x3dedges = 0, int x3dvertices = 0); 773 774 // PLY(2) format (ascii text format) 775 int readPLY(const std::string &name); 776 int readPLY2(const std::string &name); 777 int writePLY2(const std::string &name); 778 779 // Inventor/VRML format 780 int readVRML(const std::string &name); 781 int writeVRML(const std::string &name, bool saveAll = false, 782 double scalingFactor = 1.0); 783 784 // I-deas universal mesh format 785 int readUNV(const std::string &name, bool readGroupsOfElements = true); 786 int writeUNV(const std::string &name, bool saveAll = false, 787 int saveGroupsOfElements = 0, int saveGroupsOfNodes = 0, 788 double scalingFactor = 1.0); 789 790 // Medit (INRIA) mesh format 791 int readMESH(const std::string &name); 792 int writeMESH(const std::string &name, int elementTagType = 1, 793 bool saveAll = false, double scalingFactor = 1.0); 794 795 // Object file format (OFF) 796 int readOFF(const std::string &name); 797 int writeOFF(const std::string &name, bool saveAll = false, 798 double scalingFactor = 1.0); 799 800 // Nastran Bulk Data File format 801 int readBDF(const std::string &name); 802 int writeBDF(const std::string &name, int format = 0, int elementTagType = 1, 803 bool saveAll = false, double scalingFactor = 1.0); 804 805 // Actran mesh 806 int readACTRAN(const std::string &name); 807 808 // Sameced mesh 809 int readSAMCEF(const std::string &name); 810 811 // Plot3D structured mesh format 812 int readP3D(const std::string &name); 813 int writeP3D(const std::string &name, bool saveAll = false, 814 double scalingFactor = 1.0); 815 816 // CFD General Notation System files 817 int readCGNS(const std::string &name, 818 std::vector<std::vector<MVertex *> > &vertPerZone, 819 std::vector<std::vector<MElement *> > &eltPerZone); 820 int writeCGNS(const std::string &name, bool saveAll = false, 821 double scalingFactor = 1.0, bool structured = false); 822 823 // Med "Modele d'Echange de Donnees" file format (the static routine 824 // is allowed to load multiple models/meshes) 825 static int readMED(const std::string &name); 826 int readMED(const std::string &name, int meshIndex); 827 int writeMED(const std::string &name, bool saveAll = false, 828 double scalingFactor = 1.0); 829 830 // VTK format 831 int readVTK(const std::string &name, bool bigEndian = false); 832 int writeVTK(const std::string &name, bool binary = false, 833 bool saveAll = false, double scalingFactor = 1.0, 834 bool bigEndian = false); 835 836 // Matlab format 837 int writeMATLAB(const std::string &name, bool binary = false, 838 bool saveAll = false, double scalingFactor = 1.0, 839 int filetype = 1); 840 841 // Tochnog format 842 int writeTOCHNOG(const std::string &name, bool saveAll = false, 843 int saveGroupsOfNodes = 0, double scalingFactor = 1.0); 844 845 // DIFFPACK format 846 int readDIFF(const std::string &name); 847 int writeDIFF(const std::string &name, bool binary = false, 848 bool saveAll = false, double scalingFactor = 1.0); 849 850 // Abaqus 851 int writeINP(const std::string &name, bool saveAll = false, 852 int saveGroupsOfElements = 0, int saveGroupsOfNodes = 0, 853 double scalingFactor = 1.0); 854 855 // LSDYNA 856 int writeKEY(const std::string &name, int saveAll = 0, 857 int saveGroupsOfNodes = 0, double scalingFactor = 1.0); 858 859 // CELUM 860 int writeCELUM(const std::string &name, bool saveAll = false, 861 double scalingFactor = 1.0); 862 863 // Geomview mesh 864 int readGEOM(const std::string &name); 865 866 // CEA triangulation 867 int writeMAIL(const std::string &name, bool saveAll, double scalingFactor); 868 869 // SU2 mesh file 870 int writeSU2(const std::string &name, bool saveAll, double scalingFactor); 871 872 // GAMBIT neutral mesh file (.neu) 873 int writeNEU(const std::string &name, bool saveAll, double scalingFactor); 874 }; 875 876 #endif 877