1 // Copyright (C) 2012-2015 ALNEOS
2 // Copyright (C) 2016-2022 EDF R&D
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
18 // See http://www.alneos.com/ or email : contact@alneos.fr
19 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 #include "GMSHPlugin_Mesher.hxx"
22 #include "GMSHPlugin_Hypothesis_2D.hxx"
24 #include <SMDS_MeshElement.hxx>
25 #include <SMDS_MeshNode.hxx>
26 #include <SMESHDS_Mesh.hxx>
27 #include <SMESH_Comment.hxx>
28 #include <SMESH_ComputeError.hxx>
29 #include <SMESH_Gen_i.hxx>
30 #include <SMESH_Mesh.hxx>
31 #include <SMESH_MesherHelper.hxx>
32 #include <SMESH_subMesh.hxx>
33 #include <StdMeshers_FaceSide.hxx>
34 #include <utilities.h>
39 #include <TopExp_Explorer.hxx>
43 #include <MTriangle.h>
44 #include <MQuadrangle.h>
45 #if GMSH_MAJOR_VERSION >=4
46 #include <GmshGlobal.h>
47 #include <gmsh/Context.h>
50 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
60 SBoundingBox3d _bounds;
64 //================================================================================
66 * \brief Retrieve ShapeBounds from a compound GEdge
68 //================================================================================
70 bool getBoundsOfShapes( GEdge* gEdge,
71 std::vector< ShapeBounds > & topoEdges )
74 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
75 for ( size_t i = 0; i < gEdge->compound.size(); ++i )
77 GEdge* gE = static_cast< GEdge* >( gEdge->compound[ i ]);
78 topoEdges.push_back( ShapeBounds{ gE->bounds(), *((TopoDS_Edge*)gE->getNativePtr()) });
81 for ( size_t i = 0; i < gEdge->_compound.size(); ++i )
83 GEdge* gE = static_cast< GEdge* >( gEdge->_compound[ i ]);
84 topoEdges.push_back( ShapeBounds{ gE->bounds(), *((TopoDS_Edge*)gE->getNativePtr()) });
87 return topoEdges.size();
90 //================================================================================
92 * \brief Retrieve ShapeBounds from a compound GFace
94 //================================================================================
96 bool getBoundsOfShapes( GFace* gFace,
97 std::vector< ShapeBounds > & topoFaces )
100 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
101 for ( size_t i = 0; i < gFace->compound.size(); ++i )
103 GFace* gF = static_cast< GFace* >( gFace->compound[ i ]);
104 topoFaces.push_back( ShapeBounds{ gF->bounds(), *((TopoDS_Face*)gF->getNativePtr()) });
107 for ( size_t i = 0; i < gFace->_compound.size(); ++i )
109 GFace* gF = static_cast< GFace* >( gFace->_compound[ i ]);
110 topoFaces.push_back( ShapeBounds{ gF->bounds(), *((TopoDS_Face*)gF->getNativePtr()) });
113 return topoFaces.size();
115 //================================================================================
117 * \brief Find a shape whose bounding box includes a given point
119 //================================================================================
121 TopoDS_Shape getShapeAtPoint( const SPoint3& point, const std::vector< ShapeBounds > & shapes )
124 float distmin = std::numeric_limits<float>::max();
125 for ( size_t i = 0; i < shapes.size(); ++i )
127 float dist = GMSHPlugin_Mesher::DistBoundingBox( shapes[i]._bounds, point );
130 shape = shapes[i]._shape;
139 double segmentSize( const UVPtStructVec& nodeParam, size_t i )
141 double l1 = SMESH_NodeXYZ( nodeParam[i].node ).Distance( nodeParam[i-1].node );
142 double l2 = SMESH_NodeXYZ( nodeParam[i].node ).Distance( nodeParam[i+1].node );
143 return 0.5 * ( l1 + l2 );
147 //=============================================================================
151 //=============================================================================
153 GMSHPlugin_Mesher::GMSHPlugin_Mesher (SMESH_Mesh* mesh,
154 const TopoDS_Shape& aShape,
160 // il faudra peut être mettre un truc par defaut si l'utilisateur ne rentre rien en para
161 //defaultParameters();
164 //void GMSHPlugin_Mesher::defaultParameters(){}
166 void GMSHPlugin_Mesher::SetParameters(const GMSHPlugin_Hypothesis* hyp)
170 _algo2d = hyp->Get2DAlgo();
171 _algo3d = hyp->Get3DAlgo();
172 _recomb2DAlgo = hyp->GetRecomb2DAlgo();
173 _recombineAll = hyp->GetRecombineAll();
174 _subdivAlgo = hyp->GetSubdivAlgo();
175 _remeshAlgo = hyp->GetRemeshAlgo();
176 _remeshPara = hyp->GetRemeshPara();
177 _smouthSteps = hyp->GetSmouthSteps();
178 _sizeFactor = hyp->GetSizeFactor();
179 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=10
180 _meshCurvatureSize = hyp->GetMeshCurvatureSize();
182 _minSize = hyp->GetMinSize();
183 _maxSize = hyp->GetMaxSize();
184 _secondOrder = hyp->GetSecondOrder();
185 _useIncomplElem = hyp->GetUseIncomplElem();
186 _compounds = hyp->GetCompoundOnEntries();
193 _recombineAll = false;
199 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=10
200 _meshCurvatureSize = 0;
204 _secondOrder = false;
205 _useIncomplElem = true;
211 //================================================================================
213 * \brief Set maximum number of threads to be used by Gmsh
215 //================================================================================
217 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
218 void GMSHPlugin_Mesher::SetMaxThreadsGmsh()
220 MESSAGE("GMSHPlugin_Mesher::SetMaxThreadsGmsh");
221 // compound meshing (_compounds.size() > 0) and quad meshing (_algo2d >= 5) will
222 // not be multi-threaded
223 if (_compounds.size() > 0 || _algo2d >= 5){
227 _maxThreads = omp_get_max_threads();
231 //================================================================================
233 * \brief Set Gmsh Options
235 //================================================================================
237 void GMSHPlugin_Mesher::SetGmshOptions()
239 MESSAGE("GMSHPlugin_Mesher::SetGmshOptions");
241 printf("We chose _algo2d %d \n", _algo2d );
242 printf("We chose _algo3d %d \n", _algo3d );
243 printf("We chose _recomb2DAlgo %d \n", _recomb2DAlgo );
244 printf("We chose _recombineAll %d \n", (_recombineAll)?1:0);
245 printf("We chose _subdivAlgo %d \n", _subdivAlgo );
246 printf("We chose _remeshAlgo %d \n", _remeshAlgo );
247 printf("We chose _remeshPara %d \n", _remeshPara );
248 printf("We chose _smouthSteps %e \n", _smouthSteps );
249 printf("We chose _sizeFactor %e \n", _sizeFactor );
250 printf("We chose _minSize %e \n", _minSize );
251 printf("We chose _maxSize %e \n", _maxSize );
252 printf("We chose _secondOrder %d \n", (_secondOrder)?1:0);
253 printf("We chose _useIncomplElem %d \n", (_useIncomplElem)?1:0);
254 printf("We are in dimension %d \n", (_is2d)?2:3);
257 std::map <int,double> mapAlgo2d;
258 mapAlgo2d[0]=2; // Automatic
259 mapAlgo2d[1]=1; // MeshAdapt
260 mapAlgo2d[2]=5; // Delaunay
261 mapAlgo2d[3]=6; // Frontal-Delaunay
262 mapAlgo2d[4]=8; // DelQuad (Frontal-Delaunay for Quads)
263 mapAlgo2d[5]=9; // Packing of parallelograms
264 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=10
265 mapAlgo2d[6]=11;// Quasistructured quads with cross-fields
268 std::map <int,double> mapAlgo3d;
269 mapAlgo3d[0]=1; // Delaunay
270 mapAlgo3d[1]=4; // Frontal
271 mapAlgo3d[2]=7; // MMG3D
272 mapAlgo3d[3]=9; // R-tree
273 mapAlgo3d[4]=10;// HXT
276 ok = GmshSetOption("Mesh", "Algorithm" , mapAlgo2d[_algo2d]) ;
280 ok = GmshSetOption("Mesh", "Algorithm3D" , mapAlgo3d[_algo3d]) ;
283 ok = GmshSetOption("Mesh", "RecombinationAlgorithm" , (double)_recomb2DAlgo) ;
285 ok = GmshSetOption("Mesh", "RecombineAll" , (_recombineAll)?1.:0.) ;
287 ok = GmshSetOption("Mesh", "SubdivisionAlgorithm" , (double)_subdivAlgo) ;
289 ok = GmshSetOption("Mesh", "RemeshAlgorithm" , (double)_remeshAlgo) ;
291 ok = GmshSetOption("Mesh", "RemeshParametrization" , (double)_remeshPara) ;
293 ok = GmshSetOption("Mesh", "Smoothing" , (double)_smouthSteps) ;
295 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=10
296 ok = GmshSetOption("Mesh", "MeshSizeFromCurvature" , _meshCurvatureSize) ;
299 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
300 ok = GmshSetOption("Mesh", "MeshSizeFactor" , _sizeFactor) ;
302 ok = GmshSetOption("Mesh", "MeshSizeMin" , _minSize) ;
304 ok = GmshSetOption("Mesh", "MeshSizeMax" , _maxSize) ;
307 ok = GmshSetOption("Mesh", "CharacteristicLengthFactor" , _sizeFactor) ;
309 ok = GmshSetOption("Mesh", "CharacteristicLengthMin" , _minSize) ;
311 ok = GmshSetOption("Mesh", "CharacteristicLengthMax" , _maxSize) ;
314 ok = GmshSetOption("Mesh", "ElementOrder" , (_secondOrder)?2.:1.) ;
318 ok = GmshSetOption("Mesh", "SecondOrderIncomplete" ,(_useIncomplElem)?1.:0.);
322 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
323 /*ok = GmshSetOption("Mesh", "MaxNumThreads1D" , 0. ) ; // Coarse-grain algo threads
325 ok = GmshSetOption("Mesh", "MaxNumThreads2D" , 0. ) ; // Coarse-grain algo threads
327 ok = GmshSetOption("Mesh", "MaxNumThreads3D" , 0. ) ; // Fine-grain algo threads HXT
330 ok = GmshSetOption("General", "NumThreads" , _maxThreads ) ; // system default i.e. OMP_NUM_THREADS
333 if ( GMSHPlugin_Hypothesis::Algo3D::hxt == _algo3d ){
334 MESSAGE("GMSHPlugin_Mesher::SetGmshOptions: HXT algorithm is being used. Setting number of threads to 1.");
335 ok = GmshSetOption("Mesh", "MaxNumThreads3D" , 1. );
342 //================================================================================
344 * \brief Create and add Compounds into GModel _gModel.
346 //================================================================================
348 void GMSHPlugin_Mesher::CreateGmshCompounds()
350 MESSAGE("GMSHPlugin_Mesher::CreateGmshCompounds");
352 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
354 OCC_Internals* occgeo = _gModel->getOCCInternals();
355 bool toSynchronize = false;
357 for(std::set<std::string>::const_iterator its = _compounds.begin();its != _compounds.end(); ++its )
359 GEOM::GEOM_Object_var aGeomObj;
360 TopoDS_Shape geomShape = TopoDS_Shape();
361 SALOMEDS::SObject_var aSObj = SMESH_Gen_i::GetSMESHGen()->getStudyServant()->FindObjectID( (*its).c_str() );
362 SALOMEDS::GenericAttribute_var anAttr;
363 if (!aSObj->_is_nil() && aSObj->FindAttribute(anAttr, "AttributeIOR"))
365 SALOMEDS::AttributeIOR_var anIOR = SALOMEDS::AttributeIOR::_narrow(anAttr);
366 CORBA::String_var aVal = anIOR->Value();
367 CORBA::Object_var obj = SMESH_Gen_i::GetSMESHGen()->getStudyServant()->ConvertIORToObject(aVal);
368 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
370 geomShape = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
371 if ( geomShape.IsNull() )
374 TopAbs_ShapeEnum geomType = geomShape.ShapeType();
375 if ( geomType == TopAbs_COMPOUND)// voir s'il ne faut pas mettre une erreur dans le cas contraire
377 MESSAGE("shapeType == TopAbs_COMPOUND");
378 TopoDS_Iterator it(geomShape);
381 TopAbs_ShapeEnum shapeType = it.Value().ShapeType();
382 std::vector< std::pair< int, int > > dimTags;
383 for ( ; it.More(); it.Next())
385 const TopoDS_Shape& topoShape = it.Value();
386 ASSERT(topoShape.ShapeType() == shapeType);
387 if ( _mesh->GetMeshDS()->ShapeToIndex( topoShape ) > 0 )
388 occgeo->importShapes( &topoShape, false, dimTags );
391 TopoDS_Shape face = TopExp_Explorer( _shape, shapeType ).Current();
392 SMESH_subMesh* sm = _mesh->GetSubMesh( face );
393 sm->GetComputeError() =
394 SMESH_ComputeError::New
395 ( COMPERR_WARNING, "Compound shape does not belong to the main geometry. Ignored");
398 std::vector<int> tags;
399 int dim = ( shapeType == TopAbs_EDGE ) ? 1 : 2;
400 for ( size_t i = 0; i < dimTags.size(); ++i )
402 if ( dimTags[i].first == dim )
403 tags.push_back( dimTags[i].second );
407 _gModel->getGEOInternals()->setCompoundMesh( dim, tags );
408 toSynchronize = true;
411 _gModel->getGEOInternals()->synchronize(_gModel);
416 //================================================================================
418 * \brief For a compound mesh set the mesh components to be transmitted to SMESH
420 //================================================================================
422 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
423 void GMSHPlugin_Mesher::SetCompoundMeshVisibility()
426 // Loop over all faces, if the face belongs to a compound entry then
427 // for all (boundary) edges whithin the face visibility is set to 0,
428 // if the face doesn't belong to a compound entry then visibility is
429 // set to 1 for all its (boundary) edges. Later, in FillSMesh() func
430 // getVisibility() (returns either 1 or 0) is used to decide weather
431 // the mesh of edge should be transmitted to SMESH or not.
433 for ( GModel::fiter itF = _gModel->firstFace(); itF != _gModel->lastFace(); ++itF )
435 std::vector< GEdge *> faceEdges = (*itF)->edges();
437 for ( auto itE = faceEdges.begin(); itE != faceEdges.end(); ++itE )
439 if ( ((*itF)->compound.size()) )
440 (*itE)->setVisibility(0);
442 (*itE)->setVisibility(1);
447 // Loop over all edges, if the edge belongs to a compound entry then
448 // for all (boundary) vertices whithin the edge visibility is set to
449 // 0, if the edge doesn't belong to a compound entry then visibility
450 // is set to 1 for all its (boundary) vertices. Later, in FillSMesh()
451 // func getVisibility() (returns either 1 or 0) is used to decide we-
452 // ather the mesh of vertices should be transmitted to SMESH or not.
454 for ( GModel::eiter itE = _gModel->firstEdge(); itE != _gModel->lastEdge(); ++itE )
456 std::vector<GVertex *> bndVerticies = (*itE)->vertices();
458 for( auto itV = bndVerticies.begin(); itV != bndVerticies.end(); ++itV )
460 if(((*itE)->compound.size()))
461 (*itV)->setVisibility(0);
463 (*itV)->setVisibility(1);
470 //================================================================================
472 * \brief Get a node by a GMSH mesh vertex
474 //================================================================================
476 const SMDS_MeshNode* GMSHPlugin_Mesher::Node( const MVertex* v )
478 std::map< const MVertex *, const SMDS_MeshNode* >::iterator v2n = _nodeMap.find( v );
479 if ( v2n != _nodeMap.end() )
485 //================================================================================
487 * \brief Return a corresponding sub-mesh if a shape is meshed
489 //================================================================================
491 SMESHDS_SubMesh* GMSHPlugin_Mesher::HasSubMesh( const TopoDS_Shape& s )
493 if ( SMESHDS_SubMesh* sm = _mesh->GetMeshDS()->MeshElements( s ))
495 if ( s.ShapeType() == TopAbs_VERTEX )
496 return ( sm->NbNodes() > 0 ) ? sm : nullptr;
498 return ( sm->NbElements() > 0 ) ? sm : nullptr;
503 //================================================================================
505 * \brief Write mesh from GModel instance to SMESH instance
507 //================================================================================
509 void GMSHPlugin_Mesher::FillSMesh()
511 SMESHDS_Mesh* meshDS = _mesh->GetMeshDS();
514 for ( GModel::viter it = _gModel->firstVertex(); it != _gModel->lastVertex(); ++it)
516 GVertex *gVertex = *it;
518 // GET topoVertex CORRESPONDING TO gVertex
519 TopoDS_Vertex topoVertex = *((TopoDS_Vertex*)gVertex->getNativePtr());
521 if (gVertex->getVisibility() == 0) // belongs to a compound
523 SMESH_subMesh* sm = _mesh->GetSubMesh(topoVertex);
524 sm->SetIsAlwaysComputed(true); // prevent from displaying errors
528 // FILL SMESH FOR topoVertex
530 for( size_t i = 0; i < gVertex->mesh_vertices.size(); i++)
532 MVertex *v = gVertex->mesh_vertices[i];
533 if(v->getIndex() >= 0)
535 if ( SMESHDS_SubMesh* sm = HasSubMesh( topoVertex ))
537 const SMDS_MeshNode *node = sm->GetNodes()->next();
538 _nodeMap.insert({ v, node });
542 SMDS_MeshNode *node = meshDS->AddNode( v->x(),v->y(),v->z() );
543 meshDS->SetNodeOnVertex( node, topoVertex );
544 _nodeMap.insert({ v, node });
548 // WE DON'T ADD 0D ELEMENTS because it does not follow the salome meshers philosophy
550 // for(unsigned int i = 0; i < gVertex->getNumMeshElements(); i++)
552 // MElement *e = gVertex->getMeshElement(i);
553 // std::vector<MVertex*> verts;
554 // e->getVertices(verts);
555 // ASSERT(verts.size()==1);
556 // SMDS_Mesh0DElement* zeroDElement = 0;
557 // zeroDElement = meshDS->Add0DElementWithID(verts[0]->getNum(),e->getNum());
558 // meshDS->SetMeshElementOnShape(zeroDElement, topoVertex);
563 for(GModel::eiter it = _gModel->firstEdge(); it != _gModel->lastEdge(); ++it)
567 // GET topoEdge CORRESPONDING TO gEdge
568 TopoDS_Edge topoEdge;
569 std::vector< ShapeBounds > topoEdges;
570 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
571 if(gEdge->haveParametrization())
573 if ( gEdge->geomType() != GEntity::CompoundCurve )
576 topoEdge = *((TopoDS_Edge*)gEdge->getNativePtr());
577 if (gEdge->getVisibility() == 0) // belongs to a compound
579 SMESH_subMesh* sm = _mesh->GetSubMesh(topoEdge);
580 sm->SetIsAlwaysComputed(true); // prevent from displaying errors
583 if ( HasSubMesh( topoEdge ))
584 continue; // a meshed sub-mesh
586 bool isCompound = getBoundsOfShapes( gEdge, topoEdges );
588 // FILL SMESH FOR topoEdge
590 for ( size_t i = 0; i < gEdge->mesh_vertices.size(); i++ )
592 MVertex *v = gEdge->mesh_vertices[i];
593 if ( v->getIndex() >= 0 )
595 SMDS_MeshNode *node = meshDS->AddNode( v->x(),v->y(),v->z() );
598 topoEdge = TopoDS::Edge( getShapeAtPoint( v->point(), topoEdges ));
600 meshDS->SetNodeOnEdge( node, topoEdge );
601 _nodeMap.insert({ v, node });
606 for ( GModel::eiter it = _gModel->firstEdge(); it != _gModel->lastEdge(); ++it )
609 if ( gEdge->getVisibility() == 0) // belongs to a compound
612 TopoDS_Edge topoEdge;
613 std::vector< ShapeBounds > topoEdges;
614 bool isCompound = getBoundsOfShapes( gEdge, topoEdges );
616 topoEdge = *((TopoDS_Edge*)gEdge->getNativePtr());
618 if ( HasSubMesh( topoEdge ))
619 continue; // a meshed sub-mesh
622 std::vector<MVertex*> verts(3);
623 for ( size_t i = 0; i < gEdge->getNumMeshElements(); i++ )
625 MElement *e = gEdge->getMeshElement(i);
627 e->getVertices(verts);
629 // if a node wasn't set, it is assigned here
630 for ( size_t j = 0; j < verts.size(); j++ )
632 if ( verts[j]->onWhat()->getVisibility() == 0 )
634 SMDS_MeshNode *node = meshDS->AddNode(verts[j]->x(),verts[j]->y(),verts[j]->z() );
635 meshDS->SetNodeOnEdge( node, topoEdge );
636 verts[j]->setEntity(gEdge);
637 _nodeMap.insert({ verts[j], node });
641 SMDS_MeshEdge* edge = 0;
642 switch (verts.size())
645 edge = meshDS->AddEdge(Node( verts[0]),
649 edge = meshDS->AddEdge(Node( verts[0]),
658 topoEdge = TopoDS::Edge( getShapeAtPoint( e->barycenter(), topoEdges ));
660 meshDS->SetMeshElementOnShape( edge, topoEdge );
665 for ( GModel::fiter it = _gModel->firstFace(); it != _gModel->lastFace(); ++it)
669 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
670 // Gmsh since version 4.3 is now producing extra surface and mesh when
671 // compounds are involved. Since in Gmsh meshing procedure needs acess
672 // to each of the original topology and the meshed topology. Hence we
673 // bypass the additional mesh in case of compounds. Note, similar cri-
674 // teria also occurs in the following 'for' loop.
675 if ( _compounds.size() && gFace->geomType() == GEntity::DiscreteSurface )
679 // GET topoFace CORRESPONDING TO gFace
680 TopoDS_Face topoFace;
681 std::vector< ShapeBounds > topoFaces;
683 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
684 if(gFace->haveParametrization())
686 if ( gFace->geomType() != GEntity::CompoundSurface )
689 topoFace = *((TopoDS_Face*)gFace->getNativePtr());
690 if (gFace->getVisibility() == 0) // belongs to a compound
692 SMESH_subMesh* sm = _mesh->GetSubMesh(topoFace);
693 sm->SetIsAlwaysComputed(true); // prevent from displaying errors
696 if ( HasSubMesh( topoFace ))
697 continue; // a meshed sub-mesh
699 bool isCompound = getBoundsOfShapes( gFace, topoFaces );
701 // FILL SMESH FOR topoFace
703 for ( size_t i = 0; i < gFace->mesh_vertices.size(); i++ )
705 MVertex *v = gFace->mesh_vertices[i];
706 if ( v->getIndex() >= 0 )
708 SMDS_MeshNode *node = meshDS->AddNode( v->x(),v->y(),v->z() );
711 topoFace = TopoDS::Face( getShapeAtPoint( v->point(), topoFaces ));
713 meshDS->SetNodeOnFace( node, topoFace );
714 _nodeMap.insert({ v, node });
719 for ( GModel::fiter it = _gModel->firstFace(); it != _gModel->lastFace(); ++it)
723 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
724 if ( _compounds.size() && gFace->geomType() == GEntity::DiscreteSurface )
727 bool isCompound = (!gFace->haveParametrization());
729 bool isCompound = ( gFace->geomType() == GEntity::CompoundSurface );
732 if ( !isCompound && gFace->getVisibility() == 0 )
733 continue; // belongs to a compound
735 TopoDS_Face topoFace;
736 std::vector< ShapeBounds > topoFaces;
738 getBoundsOfShapes( gFace, topoFaces );
740 topoFace = *((TopoDS_Face*)gFace->getNativePtr());
742 if ( HasSubMesh( topoFace ))
743 continue; // a meshed sub-mesh
746 std::vector<MVertex*> verts;
747 for ( size_t i = 0; i < gFace->getNumMeshElements(); i++ )
749 MElement *e = gFace->getMeshElement(i);
751 e->getVertices(verts);
752 SMDS_MeshFace* face = 0;
754 // if a node wasn't set, it is assigned here
755 for ( size_t j = 0; j < verts.size(); j++)
757 if(verts[j]->onWhat()->getVisibility() == 0)
759 SMDS_MeshNode *node = meshDS->AddNode(verts[j]->x(),verts[j]->y(),verts[j]->z());
760 meshDS->SetNodeOnFace( node, topoFace );
761 _nodeMap.insert({ verts[j], node });
762 verts[j]->setEntity(gFace);
765 switch (verts.size())
768 face = meshDS->AddFace(Node( verts[0]),
773 face = meshDS->AddFace(Node( verts[0]),
779 face = meshDS->AddFace(Node( verts[0]),
787 face = meshDS->AddFace(Node( verts[0]),
797 face = meshDS->AddFace(Node( verts[0]),
813 topoFace = TopoDS::Face( getShapeAtPoint( e->barycenter(), topoFaces ));
815 meshDS->SetMeshElementOnShape(face, topoFace);
820 for ( GModel::riter it = _gModel->firstRegion(); it != _gModel->lastRegion(); ++it)
822 GRegion *gRegion = *it;
823 if (gRegion->getVisibility() == 0)
826 // GET topoSolid CORRESPONDING TO gRegion
827 TopoDS_Solid topoSolid = *((TopoDS_Solid*)gRegion->getNativePtr());
829 // FILL SMESH FOR topoSolid
832 for( size_t i = 0; i < gRegion->mesh_vertices.size(); i++)
834 MVertex *v = gRegion->mesh_vertices[i];
835 if(v->getIndex() >= 0)
837 SMDS_MeshNode *node = meshDS->AddNode( v->x(),v->y(),v->z() );
838 meshDS->SetNodeInVolume( node, topoSolid );
839 _nodeMap.insert({ v, node });
844 std::vector<MVertex*> verts;
845 for( size_t i = 0; i < gRegion->getNumMeshElements(); i++)
847 MElement *e = gRegion->getMeshElement(i);
849 e->getVertices(verts);
850 SMDS_MeshVolume* volume = 0;
851 switch (verts.size()){
853 volume = meshDS->AddVolume(Node( verts[0]),
859 volume = meshDS->AddVolume(Node( verts[0]),
866 volume = meshDS->AddVolume(Node( verts[0]),
874 volume = meshDS->AddVolume(Node( verts[0]),
884 volume = meshDS->AddVolume(Node( verts[0]),
896 volume = meshDS->AddVolume(Node( verts[0] ),
910 case 14: // same as case 13, because no pyra14 in smesh
911 volume = meshDS->AddVolume(Node( verts[0] ),
926 volume = meshDS->AddVolume(Node( verts[0] ),
942 case 18: // same as case 15, because no penta18 in smesh
943 volume = meshDS->AddVolume(Node( verts[0] ),
960 volume = meshDS->AddVolume(Node( verts[0] ),
982 volume = meshDS->AddVolume(Node( verts[0] ),
1014 meshDS->SetMeshElementOnShape(volume, topoSolid);
1021 //================================================================================
1023 * \brief Find if SPoint point is in SBoundingBox3d bounds
1025 //================================================================================
1027 float GMSHPlugin_Mesher::DistBoundingBox(const SBoundingBox3d& bounds, const SPoint3& point)
1029 SPoint3 min = bounds.min();
1030 SPoint3 max = bounds.max();
1034 if (point.x() < min.x())
1035 x = min.x()-point.x();
1036 else if (point.x() > max.x())
1037 x = point.x()-max.x();
1041 if (point.y() < min.y())
1042 y = min.y()-point.y();
1043 else if (point.y() > max.y())
1044 y = point.y()-max.y();
1048 if (point.z() < min.z())
1049 z = min.z()-point.z();
1050 else if (point.z() > max.z())
1051 z = point.z()-max.z();
1057 //================================================================================
1059 * \brief Reimplemented GmshMessage call. Actions done if errors occurs
1060 * during gmsh meshing. We define here what to display and what to do.
1062 //================================================================================
1063 void GMSHPlugin_Mesher::mymsg::operator()(std::string level, std::string msg)
1065 //MESSAGE("level="<< level.c_str() << ", msg=" << msg.c_str()<< "\n");
1066 printf("level=%s msg=%s\n", level.c_str(), msg.c_str());
1068 if(level == "Fatal" || level == "Error")
1070 std::ostringstream oss;
1071 if (level == "Fatal")
1072 oss << "Fatal error during Generation of Gmsh Mesh\n";
1074 oss << "Error during Generation of Gmsh Mesh\n";
1075 oss << " " << msg.c_str() << "\n";
1076 GEntity *e = _gModel->getCurrentMeshEntity();
1079 oss << " error occurred while meshing entity:\n" <<
1080 " tag=" << e->tag() << "\n" <<
1081 " dimension=" << e->dim() << "\n" <<
1082 " native pointer=" << e->getNativePtr();
1083 //if(e->geomType() != GEntity::CompoundCurve and e->geomType() != GEntity::CompoundSurface)
1085 //SMESH_subMesh *sm = _mesh->GetSubMesh(*((TopoDS_Shape*)e->getNativePtr()));
1086 //SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1087 //SMESH_Comment comment;
1088 //comment << SMESH_Comment(oss.str);
1089 //std::string str = oss.str();
1090 //smError.reset( new SMESH_ComputeError( str ));
1092 // plutot que de faire de la merde ici, on pourait simplement
1093 // remplir une liste pour dire sur quelles entités gmsh se plante
1094 // (puis faire le fillsmesh)
1095 // puis faire une nouvelle routine qui réécrit les messages d'erreur
1096 // probleme : gmsh peut planté en Fatal, dans ce cas pas de fillsmesh
1099 if (level == "Fatal")
1101 CTX::instance()->lock = 0;
1105 printf("%s\n", oss.str().c_str());
1109 //=============================================================================
1111 * Here we are going to use the GMSH mesher
1113 //=============================================================================
1115 bool GMSHPlugin_Mesher::Compute()
1117 MESSAGE("GMSHPlugin_Mesher::Compute");
1121 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1122 SetMaxThreadsGmsh();
1124 //RNV: to avoid modification of PATH and PYTHONPATH
1125 char* argv[] = {"-noenv"};
1126 GmshInitialize(1,argv);
1128 _gModel = new GModel();
1130 GmshSetMessageHandler(&msg);
1131 _gModel->importOCCShape((void*)&_shape);
1132 if (_compounds.size() > 0) CreateGmshCompounds();
1135 //Msg::SetVerbosity(100);
1136 //CTX::instance()->mesh.maxNumThreads1D=1;
1138 HideComputedEntities( _gModel );
1140 _gModel->mesh( /*dim=*/ 1 );
1142 Set1DSubMeshes( _gModel );
1144 //_gModel->writeUNV("/tmp/1D.unv", 1,0,0);
1145 //CTX::instance()->mesh.maxNumThreads2D=1;
1147 _gModel->mesh( /*dim=*/ 2 );
1151 Set2DSubMeshes( _gModel );
1153 //CTX::instance()->mesh.maxNumThreads3D=1;
1155 _gModel->mesh( /*dim=*/ 3 );
1157 RestoreVisibility( _gModel );
1159 #ifdef WITH_SMESH_CANCEL_COMPUTE
1163 catch (std::string& str)
1166 std::cerr << "GMSH: exception caught: " << str << std::endl;
1172 std::cerr << "GMSH: Unknown exception caught: " << std::endl;
1173 MESSAGE("Unrecoverable error during Generation of Gmsh Mesh");
1178 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1179 if (_compounds.size() > 0)
1180 SetCompoundMeshVisibility();
1184 GmshSetMessageHandler(nullptr);
1187 MESSAGE("GMSHPlugin_Mesher::Compute:End");
1191 //================================================================================
1193 * \brief Set 1D sub-meshes to GModel. GMSH 1D mesh is made by now.
1194 * \param [inout] _gModel - GMSH model
1196 //================================================================================
1198 void GMSHPlugin_Mesher::Set1DSubMeshes( GModel* gModel )
1200 SMESHDS_Mesh* meshDS = _mesh->GetMeshDS();
1202 for(GModel::eiter it = gModel->firstEdge(); it != gModel->lastEdge(); ++it)
1206 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1207 if ( !gEdge->haveParametrization())
1209 if ( gEdge->geomType() == GEntity::CompoundCurve )
1213 TopoDS_Edge topoEdge = *((TopoDS_Edge*)gEdge->getNativePtr());
1214 if ( !HasSubMesh( topoEdge ))
1215 continue; // empty sub-mesh
1217 gEdge->deleteMesh();
1219 // get node parameters on topoEdge
1220 StdMeshers_FaceSide side( TopoDS_Face(), topoEdge, _mesh, /*fwd=*/true, /*skpMedium=*/true);
1221 const UVPtStructVec& nodeParam = side.GetUVPtStruct();
1222 if ( nodeParam.empty() )
1223 throw std::string("Pb with StdMeshers_FaceSide::GetUVPtStruct()");
1225 // get GMSH mesh vertices on VERTEX'es
1226 std::vector<MVertex *> mVertices( nodeParam.size(), nullptr );
1227 GVertex * gV0 = gEdge->getBeginVertex(), *gV1 = gEdge->getEndVertex();
1228 mVertices[0] = gV0->mesh_vertices[ 0 ];
1229 mVertices.back() = gV1->mesh_vertices[ 0 ];
1230 TopoDS_Vertex v01 = *((TopoDS_Vertex*) gV0->getNativePtr());
1231 TopoDS_Shape v02 = SMESH_MesherHelper::GetSubShapeByNode( nodeParam[0].node, meshDS );
1232 bool reverse = !v01.IsSame( v02 );
1233 if ( mVertices[0] == mVertices.back() )
1234 reverse = ( nodeParam[0].param > nodeParam.back().param );
1235 const SMDS_MeshNode* n0 = reverse ? nodeParam.back().node : nodeParam[0].node;
1236 const SMDS_MeshNode* n1 = reverse ? nodeParam[0].node : nodeParam.back().node;
1237 _nodeMap.insert({ mVertices[ 0 ], n0 });
1238 _nodeMap.insert({ mVertices.back(), n1 });
1240 // create GMSH mesh vertices on gEdge
1241 for ( size_t i = 1; i < nodeParam.size() - 1; ++i )
1243 size_t iN = reverse ? ( nodeParam.size() - 1 - i ) : i;
1244 SMESH_NodeXYZ xyz = nodeParam[ iN ].node;
1245 double lc = segmentSize( nodeParam, iN );
1246 // SVector3 der = gEdge->firstDer(nodeParam[ iN ].param);
1247 // double lc = norm(der) / segmentSize( nodeParam, i );
1249 mVertices[ i ] = new MEdgeVertex( xyz.X(), xyz.Y(), xyz.Z(),
1250 gEdge, nodeParam[ iN ].param, 0, lc);
1251 gEdge->mesh_vertices.push_back( mVertices[ i ]);
1252 _nodeMap.insert({ mVertices[ i ], nodeParam[ iN ].node });
1254 // create GMSH mesh edges
1255 for ( size_t i = 1; i < mVertices.size(); ++i )
1257 gEdge->lines.push_back( new MLine( mVertices[ i - 1 ],
1261 cout << endl << "EDGE " << gEdge->tag() <<
1262 ( topoEdge.Orientation() == TopAbs_FORWARD ? " F" : " R") << endl;
1263 MVertex* mv = gV0->mesh_vertices[ 0 ];
1264 cout << "V0: " << mv->x() << ", " << mv->y() << ", " << mv->z() << ", "<<endl;
1265 for ( size_t i = 0; i < gEdge->mesh_vertices.size(); ++i )
1267 MEdgeVertex* mv = (MEdgeVertex*) gEdge->mesh_vertices[i];
1268 cout << i << ": " << mv->x() << ", " << mv->y() << ", " << mv->z() << ", ";
1270 mv->getParameter(0, t );
1271 cout << ":\t t = "<< t << " lc = " << mv->getLc() << endl;
1273 mv = gV1->mesh_vertices[ 0 ];
1274 cout << "V1: " << mv->x() << ", " << mv->y() << ", " << mv->z() << ", "<<endl;
1280 //================================================================================
1282 * \brief Set 2D sub-meshes to GModel. GMSH 2D mesh is made by now.
1283 * \param [inout] _gModel - GMSH model
1285 //================================================================================
1287 void GMSHPlugin_Mesher::Set2DSubMeshes( GModel* gModel )
1289 if ( _nodeMap.empty() )
1290 return; // no sub-meshes
1292 SMESH_MesherHelper helper( *_mesh );
1294 std::map< const SMDS_MeshNode* , const MVertex * > nodes2mvertMap;
1295 for ( auto & v2n : _nodeMap )
1296 nodes2mvertMap.insert({ v2n.second, v2n.first });
1298 std::vector<MVertex *> mVertices;
1300 for(GModel::fiter it = gModel->firstFace(); it != gModel->lastFace(); ++it)
1304 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1305 if ( !gFace->haveParametrization())
1307 if ( gFace->geomType() == GEntity::CompoundSurface )
1311 TopoDS_Face topoFace = *((TopoDS_Face*)gFace->getNativePtr());
1312 SMESHDS_SubMesh* sm = HasSubMesh( topoFace );
1315 //_gModel->writeUNV("/tmp/befDEl.unv", 1,0,0);
1317 gFace->deleteMesh();
1319 bool reverse = false;
1320 if ( gFace->getRegion(0) )
1322 //GRegion * gRegion = gFace->getRegion(0);
1323 TopoDS_Shape topoSolid = *((TopoDS_Shape*)gFace->getNativePtr());
1324 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( topoSolid, topoFace );
1325 if ( faceOriInSolid >= 0 )
1327 helper.IsReversedSubMesh( TopoDS::Face( topoFace.Oriented( faceOriInSolid )));
1330 for ( SMDS_ElemIteratorPtr fIt = sm->GetElements(); fIt->more(); )
1332 const SMDS_MeshElement* f = fIt->next();
1334 int nbN = f->NbCornerNodes();
1336 throw std::string("Polygon sub-meshes not supported");
1338 mVertices.resize( nbN );
1339 for ( int i = 0; i < nbN; ++i )
1341 const SMDS_MeshNode* n = f->GetNode( i );
1342 MVertex * mv = nullptr;
1343 auto n2v = nodes2mvertMap.find( n );
1344 if ( n2v != nodes2mvertMap.end() )
1346 mv = const_cast< MVertex*>( n2v->second );
1350 if ( n->GetPosition()->GetDim() < 2 )
1351 throw std::string("Wrong mapping of edge nodes to GMSH nodes");
1352 SMESH_NodeXYZ xyz = n;
1354 gp_XY uv = helper.GetNodeUV( topoFace, n, nullptr, &ok );
1355 mv = new MFaceVertex( xyz.X(), xyz.Y(), xyz.Z(), gFace, uv.X(), uv.Y() );
1356 gFace->mesh_vertices.push_back( mv );
1357 nodes2mvertMap.insert({ n, mv });
1358 _nodeMap.insert ({ mv, n });
1360 mVertices[ i ] = mv;
1362 // create GMSH mesh faces
1366 gFace->triangles.push_back (new MTriangle(mVertices[0], mVertices[2], mVertices[1]));
1368 gFace->triangles.push_back (new MTriangle(mVertices[0], mVertices[1], mVertices[2]));
1372 gFace->quadrangles.push_back (new MQuadrangle(mVertices[0], mVertices[3],
1373 mVertices[2], mVertices[1]));
1375 gFace->quadrangles.push_back (new MQuadrangle(mVertices[0], mVertices[1],
1376 mVertices[2], mVertices[3]));
1381 } // loop on GMSH faces
1386 //================================================================================
1388 * \brief Set visibility 0 to already computed geom entities
1389 * to prevent their meshing
1391 //================================================================================
1393 void GMSHPlugin_Mesher::HideComputedEntities( GModel* gModel )
1395 CTX::instance()->mesh.meshOnlyVisible = true;
1397 for(GModel::eiter it = gModel->firstEdge(); it != gModel->lastEdge(); ++it)
1401 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1402 if ( !gEdge->haveParametrization())
1404 if ( gEdge->geomType() == GEntity::CompoundCurve )
1408 TopoDS_Edge topoEdge = *((TopoDS_Edge*)gEdge->getNativePtr());
1409 if ( HasSubMesh( topoEdge ))
1410 gEdge->setVisibility(0);
1414 for(GModel::fiter it = gModel->firstFace(); it != gModel->lastFace(); ++it)
1418 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1419 if ( !gFace->haveParametrization())
1421 if ( gFace->geomType() == GEntity::CompoundSurface )
1425 TopoDS_Face topoFace = *((TopoDS_Face*)gFace->getNativePtr());
1426 if ( HasSubMesh( topoFace ))
1427 gFace->setVisibility(0);
1431 //================================================================================
1433 * \brief Restore visibility of all geom entities
1435 //================================================================================
1437 void GMSHPlugin_Mesher::RestoreVisibility( GModel* gModel )
1439 for(GModel::eiter it = gModel->firstEdge(); it != gModel->lastEdge(); ++it)
1442 gEdge->setVisibility(1);
1444 for(GModel::fiter it = gModel->firstFace(); it != gModel->lastFace(); ++it)
1447 gFace->setVisibility(1);
1452 void GMSHPlugin_Mesher::toPython( GModel* )
1454 const char* pyFile = "/tmp/gMesh.py";
1455 ofstream outfile( pyFile, ios::out );
1456 if ( !outfile ) return;
1458 outfile << "import salome, SMESH" << std::endl
1459 << "from salome.smesh import smeshBuilder" << std::endl
1460 << "smesh = smeshBuilder.New()" << std::endl
1461 << "mesh = smesh.Mesh()" << std::endl << std::endl;
1463 outfile << "## VERTICES" << endl;
1464 for ( GModel::viter it = _gModel->firstVertex(); it != _gModel->lastVertex(); ++it)
1466 GVertex *gVertex = *it;
1468 for(unsigned int i = 0; i < gVertex->mesh_vertices.size(); i++)
1470 MVertex *v = gVertex->mesh_vertices[i];
1471 if ( v->getIndex() >= 0)
1473 outfile << "n" << v->getNum() << " = mesh.AddNode("
1474 << v->x() << ", " << v->y() << ", " << v->z()<< ")"
1475 << " ## tag = " << gVertex->tag() << endl;
1480 for(GModel::eiter it = _gModel->firstEdge(); it != _gModel->lastEdge(); ++it)
1483 outfile << "## GEdge " << gEdge->tag() << endl;
1485 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1486 if(gEdge->haveParametrization())
1488 if ( gEdge->geomType() != GEntity::CompoundCurve )
1490 for ( size_t i = 0; i < gEdge->mesh_vertices.size(); i++ )
1492 MVertex *v = gEdge->mesh_vertices[i];
1493 if ( v->getIndex() >= 0 )
1495 outfile << "n" << v->getNum() << " = mesh.AddNode("
1496 << v->x() << ", " << v->y() << ", " << v->z()<< ")"
1497 << " ## tag = " << gEdge->tag() << endl;
1502 for ( GModel::fiter it = _gModel->firstFace(); it != _gModel->lastFace(); ++it)
1505 if ( _compounds.size() && gFace->geomType() == GEntity::DiscreteSurface )
1507 outfile << "## GFace " << gFace->tag() << endl;
1509 for ( size_t i = 0; i < gFace->mesh_vertices.size(); i++ )
1511 MVertex *v = gFace->mesh_vertices[i];
1512 if ( v->getIndex() >= 0 )
1514 outfile << "n" << v->getNum() << " = mesh.AddNode("
1515 << v->x() << ", " << v->y() << ", " << v->z()<< ")"
1516 << " ## tag = " << gFace->tag() << endl;
1521 std::vector<MVertex*> verts(3);
1522 for(GModel::eiter it = _gModel->firstEdge(); it != _gModel->lastEdge(); ++it)
1525 outfile << "## GEdge " << gEdge->tag() << endl;
1527 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1528 if(gEdge->haveParametrization())
1530 if ( gEdge->geomType() != GEntity::CompoundCurve )
1532 for ( size_t i = 0; i < gEdge->getNumMeshElements(); i++ )
1534 MElement *e = gEdge->getMeshElement(i);
1536 e->getVertices(verts);
1538 outfile << "e" << e->getNum() << " = mesh.AddEdge(["
1539 << "n" << verts[0]->getNum() << ","
1540 << "n" << verts[1]->getNum();
1541 if ( verts.size() == 3 )
1542 outfile << "n" << verts[2]->getNum();
1543 outfile << "])"<< endl;
1547 for ( GModel::fiter it = _gModel->firstFace(); it != _gModel->lastFace(); ++it)
1550 if ( _compounds.size() && gFace->geomType() == GEntity::DiscreteSurface )
1552 outfile << "## GFace " << gFace->tag() << endl;
1554 for ( size_t i = 0; i < gFace->getNumMeshElements(); i++ )
1556 MElement *e = gFace->getMeshElement(i);
1558 e->getVertices(verts);
1560 outfile << "f" << e->getNum() << " = mesh.AddFace([";
1561 for ( size_t j = 0; j < verts.size(); j++)
1563 outfile << "n" << verts[j]->getNum();
1564 if ( j < verts.size()-1 )
1567 outfile << "])" << endl;
1570 std::cout << "Write " << pyFile << std::endl;