1 // Copyright (C) 2012-2015 ALNEOS
2 // Copyright (C) 2016-2023 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>
35 #include <StdMeshers_QuadToTriaAdaptor.hxx>
38 #include <BRep_Tool.hxx>
39 #include <GeomAPI_ProjectPointOnCurve.hxx>
46 #include <TopExp_Explorer.hxx>
50 #include <MTriangle.h>
51 #include <MQuadrangle.h>
52 #if GMSH_MAJOR_VERSION >=4
53 #include <GmshGlobal.h>
54 #include <gmsh/Context.h>
57 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
67 SBoundingBox3d _bounds;
71 //================================================================================
73 * \brief Retrieve ShapeBounds from a compound GEdge
75 //================================================================================
77 bool getBoundsOfShapes( GEdge* gEdge,
78 std::vector< ShapeBounds > & topoEdges )
81 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
82 for ( size_t i = 0; i < gEdge->compound.size(); ++i )
84 GEdge* gE = static_cast< GEdge* >( gEdge->compound[ i ]);
85 topoEdges.push_back( ShapeBounds{ gE->bounds(), *((TopoDS_Edge*)gE->getNativePtr()) });
88 for ( size_t i = 0; i < gEdge->_compound.size(); ++i )
90 GEdge* gE = static_cast< GEdge* >( gEdge->_compound[ i ]);
91 topoEdges.push_back( ShapeBounds{ gE->bounds(), *((TopoDS_Edge*)gE->getNativePtr()) });
94 return topoEdges.size();
97 //================================================================================
99 * \brief Retrieve ShapeBounds from a compound GFace
101 //================================================================================
103 bool getBoundsOfShapes( GFace* gFace,
104 std::vector< ShapeBounds > & topoFaces )
107 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
108 for ( size_t i = 0; i < gFace->compound.size(); ++i )
110 GFace* gF = static_cast< GFace* >( gFace->compound[ i ]);
111 topoFaces.push_back( ShapeBounds{ gF->bounds(), *((TopoDS_Face*)gF->getNativePtr()) });
114 for ( size_t i = 0; i < gFace->_compound.size(); ++i )
116 GFace* gF = static_cast< GFace* >( gFace->_compound[ i ]);
117 topoFaces.push_back( ShapeBounds{ gF->bounds(), *((TopoDS_Face*)gF->getNativePtr()) });
120 return topoFaces.size();
122 //================================================================================
124 * \brief Find a shape whose bounding box includes a given point
126 //================================================================================
128 TopoDS_Shape getShapeAtPoint( const SPoint3& point, const std::vector< ShapeBounds > & shapes )
131 float distmin = std::numeric_limits<float>::max();
132 for ( size_t i = 0; i < shapes.size(); ++i )
134 float dist = GMSHPlugin_Mesher::DistBoundingBox( shapes[i]._bounds, point );
137 shape = shapes[i]._shape;
146 double segmentSize( const UVPtStructVec& nodeParam, size_t i )
148 double l1 = SMESH_NodeXYZ( nodeParam[i].node ).Distance( nodeParam[i-1].node );
149 double l2 = SMESH_NodeXYZ( nodeParam[i].node ).Distance( nodeParam[i+1].node );
150 return 0.5 * ( l1 + l2 );
154 //=============================================================================
158 //=============================================================================
160 GMSHPlugin_Mesher::GMSHPlugin_Mesher (SMESH_Mesh* mesh,
161 const TopoDS_Shape& aShape,
169 // il faudra peut être mettre un truc par defaut si l'utilisateur ne rentre rien en para
170 //defaultParameters();
173 //void GMSHPlugin_Mesher::defaultParameters(){}
175 void GMSHPlugin_Mesher::SetParameters(const GMSHPlugin_Hypothesis* hyp)
179 _algo2d = hyp->Get2DAlgo();
180 _algo3d = hyp->Get3DAlgo();
181 _recomb2DAlgo = hyp->GetRecomb2DAlgo();
182 _recombineAll = hyp->GetRecombineAll();
183 _subdivAlgo = hyp->GetSubdivAlgo();
184 _remeshAlgo = hyp->GetRemeshAlgo();
185 _remeshPara = hyp->GetRemeshPara();
186 _smouthSteps = hyp->GetSmouthSteps();
187 _sizeFactor = hyp->GetSizeFactor();
188 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=10
189 _meshCurvatureSize = hyp->GetMeshCurvatureSize();
191 _minSize = hyp->GetMinSize();
192 _maxSize = hyp->GetMaxSize();
193 _secondOrder = hyp->GetSecondOrder();
194 _useIncomplElem = hyp->GetUseIncomplElem();
195 _compounds = hyp->GetCompoundOnEntries();
202 _recombineAll = false;
208 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=10
209 _meshCurvatureSize = 0;
213 _secondOrder = false;
214 _useIncomplElem = true;
220 //================================================================================
222 * \brief Set maximum number of threads to be used by Gmsh
224 //================================================================================
226 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
227 void GMSHPlugin_Mesher::SetMaxThreadsGmsh()
229 MESSAGE("GMSHPlugin_Mesher::SetMaxThreadsGmsh");
230 // compound meshing (_compounds.size() > 0) and quad meshing (_algo2d >= 5) will
231 // not be multi-threaded
232 if (_compounds.size() > 0 || _algo2d >= 5){
236 _maxThreads = omp_get_max_threads();
240 //================================================================================
242 * \brief Initialize GMSH model
244 //================================================================================
245 void GMSHPlugin_Mesher::FillGMSHMesh()
248 gmsh::model::add("mesh");
250 SMESHDS_Mesh* meshDS = _mesh->GetMeshDS();
255 std::vector<int> aTrinagle(3, 0);
257 const int invalid_ID = -1;
260 typedef map< const SMDS_MeshNode*, int, TIDCompare > TNodeToIDMap;
261 typedef TNodeToIDMap::value_type TN2ID;
262 typedef map<std::pair<int, int>, int> TLineToIDMap;
263 TNodeToIDMap aNodeToID;
264 TLineToIDMap aLineToID;
266 TopAbs_ShapeEnum aMainType = _mesh->GetShapeToMesh().ShapeType();
267 bool aCheckReverse = (aMainType == TopAbs_COMPOUND || aMainType == TopAbs_COMPSOLID);
269 SMESH_MesherHelper aHelper(*_mesh);
270 SMESH_ProxyMesh::Ptr aProxyMesh(new SMESH_ProxyMesh(*_mesh));
271 if (_mesh->NbQuadrangles() > 0)
273 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
274 Adaptor->Compute(*_mesh, _shape, aProxyMesh.get());
275 aProxyMesh.reset(Adaptor);
278 std::map<const SMDS_MeshElement*, bool, TIDCompare> listElements;
279 for (TopExp_Explorer exFa(_shape, TopAbs_FACE); exFa.More(); exFa.Next())
281 const TopoDS_Shape& aShapeFace = exFa.Current();
282 int faceID = meshDS->ShapeToIndex(aShapeFace);
285 if (aCheckReverse && aHelper.NbAncestors(aShapeFace, *_mesh, _shape.ShapeType()) > 1)
286 // IsReversedSubMesh() can work wrong on strongly curved faces,
287 // so we use it as less as possible
288 isRev = aHelper.IsReversedSubMesh(TopoDS::Face(aShapeFace));
290 const SMESHDS_SubMesh* aSubMeshDSFace = aProxyMesh->GetSubMesh(aShapeFace);
291 if (!aSubMeshDSFace) continue;
293 SMDS_ElemIteratorPtr iteratorElem = aSubMeshDSFace->GetElements();
294 if (aHelper.IsQuadraticSubMesh(_shape) &&
295 dynamic_cast<const SMESH_ProxyMesh::SubMesh*>(aSubMeshDSFace))
297 // add medium nodes of proxy triangles to helper
298 while (iteratorElem->more())
299 aHelper.AddTLinks(static_cast<const SMDS_MeshFace*>(iteratorElem->next()));
301 iteratorElem = aSubMeshDSFace->GetElements();
303 while (iteratorElem->more()) // loop on elements on a geom face
306 const SMDS_MeshElement* elem = iteratorElem->next();
309 if (elem->NbCornerNodes() != 3)
311 listElements[elem] = isRev;
315 for (auto const& ls : listElements)
317 // Add nodes of triangles and triangles them-selves to netgen mesh
318 // add three nodes of
319 bool hasDegen = false;
320 for (int iN = 0; iN < 3; ++iN)
322 const SMDS_MeshNode* aNode = ls.first->GetNode(iN);
323 const int shapeID = aNode->getshapeId();
324 if (aNode->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE &&
325 aHelper.IsDegenShape(shapeID))
327 // ignore all nodes on degeneraged edge and use node on its vertex instead
328 TopoDS_Shape vertex = TopoDS_Iterator(meshDS->IndexToShape(shapeID)).Value();
329 aNode = SMESH_Algo::VertexNode(TopoDS::Vertex(vertex), meshDS);
332 int& ngID = aNodeToID.insert(TN2ID(aNode, invalid_ID)).first->second;
333 if (ngID == invalid_ID)
336 gmsh::model::occ::addPoint(aNode->X(), aNode->Y(), aNode->Z(), 1.e-2, ngID);
338 aTrinagle[ls.second ? 2 - iN : iN] = ngID;
341 if (hasDegen && (aTrinagle[0] == aTrinagle[1] ||
342 aTrinagle[0] == aTrinagle[2] ||
343 aTrinagle[2] == aTrinagle[1]))
347 std::vector<int> LinesID(3, 0);
348 for (int anIndex = 0; anIndex < 3; ++anIndex)
350 int aNextIndex = (anIndex + 1) % 3;
351 if (aLineToID.find({ aTrinagle[anIndex], aTrinagle[aNextIndex] }) == aLineToID.end()
352 && aLineToID.find({ aTrinagle[aNextIndex], aTrinagle[anIndex] }) == aLineToID.end())
354 LinesID[anIndex] = aLineToID.insert({ { aTrinagle[aNextIndex], aTrinagle[anIndex] }, ++aNbOfLines }).first->second;
355 gmsh::model::occ::addLine(aTrinagle[anIndex], aTrinagle[aNextIndex], LinesID[anIndex]);
359 LinesID[anIndex] = aLineToID.find({ aTrinagle[anIndex], aTrinagle[aNextIndex] })->second;
360 if (LinesID[anIndex] == 0)
361 LinesID[anIndex] = aLineToID.find({ aTrinagle[aNextIndex], aTrinagle[anIndex] })->second;
365 if (!aProxyMesh->IsTemporary(ls.first))
366 swap(aTrinagle[1], aTrinagle[2]);
368 int aTag = gmsh::model::occ::addCurveLoop(LinesID);
371 // Generate 1D and 2D mesh
372 _gModel->mesh( /*dim=*/ 1);
373 Set1DSubMeshes(_gModel);
374 _gModel->mesh( /*dim=*/ 2);
377 //================================================================================
379 * \brief Set Gmsh Options
381 //================================================================================
382 void GMSHPlugin_Mesher::SetGmshOptions()
384 MESSAGE("GMSHPlugin_Mesher::SetGmshOptions");
386 printf("We chose _algo2d %d \n", _algo2d );
387 printf("We chose _algo3d %d \n", _algo3d );
388 printf("We chose _recomb2DAlgo %d \n", _recomb2DAlgo );
389 printf("We chose _recombineAll %d \n", (_recombineAll)?1:0);
390 printf("We chose _subdivAlgo %d \n", _subdivAlgo );
391 printf("We chose _remeshAlgo %d \n", _remeshAlgo );
392 printf("We chose _remeshPara %d \n", _remeshPara );
393 printf("We chose _smouthSteps %e \n", _smouthSteps );
394 printf("We chose _sizeFactor %e \n", _sizeFactor );
395 printf("We chose _minSize %e \n", _minSize );
396 printf("We chose _maxSize %e \n", _maxSize );
397 printf("We chose _secondOrder %d \n", (_secondOrder)?1:0);
398 printf("We chose _useIncomplElem %d \n", (_useIncomplElem)?1:0);
399 printf("We are in dimension %d \n", (_is2d)?2:3);
402 std::map <int,double> mapAlgo2d;
403 mapAlgo2d[0]=2; // Automatic
404 mapAlgo2d[1]=1; // MeshAdapt
405 mapAlgo2d[2]=5; // Delaunay
406 mapAlgo2d[3]=6; // Frontal-Delaunay
407 mapAlgo2d[4]=8; // DelQuad (Frontal-Delaunay for Quads)
408 mapAlgo2d[5]=9; // Packing of parallelograms
409 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=10
410 mapAlgo2d[6]=11;// Quasistructured quads with cross-fields
413 std::map <int,double> mapAlgo3d;
414 mapAlgo3d[0]=1; // Delaunay
415 mapAlgo3d[1]=4; // Frontal
416 mapAlgo3d[2]=7; // MMG3D
417 mapAlgo3d[3]=9; // R-tree
418 mapAlgo3d[4]=10;// HXT
421 ok = GmshSetOption("Mesh", "Algorithm" , mapAlgo2d[_algo2d]) ;
425 ok = GmshSetOption("Mesh", "Algorithm3D" , mapAlgo3d[_algo3d]) ;
428 ok = GmshSetOption("Mesh", "RecombinationAlgorithm" , (double)_recomb2DAlgo) ;
430 ok = GmshSetOption("Mesh", "RecombineAll" , (_recombineAll)?1.:0.) ;
432 ok = GmshSetOption("Mesh", "SubdivisionAlgorithm" , (double)_subdivAlgo) ;
434 ok = GmshSetOption("Mesh", "RemeshAlgorithm" , (double)_remeshAlgo) ;
436 ok = GmshSetOption("Mesh", "RemeshParametrization" , (double)_remeshPara) ;
438 ok = GmshSetOption("Mesh", "Smoothing" , (double)_smouthSteps) ;
440 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=10
441 ok = GmshSetOption("Mesh", "MeshSizeFromCurvature" , _meshCurvatureSize) ;
444 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
445 ok = GmshSetOption("Mesh", "MeshSizeFactor" , _sizeFactor) ;
447 ok = GmshSetOption("Mesh", "MeshSizeMin" , _minSize) ;
449 ok = GmshSetOption("Mesh", "MeshSizeMax" , _maxSize) ;
452 ok = GmshSetOption("Mesh", "CharacteristicLengthFactor" , _sizeFactor) ;
454 ok = GmshSetOption("Mesh", "CharacteristicLengthMin" , _minSize) ;
456 ok = GmshSetOption("Mesh", "CharacteristicLengthMax" , _maxSize) ;
459 ok = GmshSetOption("Mesh", "ElementOrder" , (_secondOrder)?2.:1.) ;
463 ok = GmshSetOption("Mesh", "SecondOrderIncomplete" ,(_useIncomplElem)?1.:0.);
467 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
468 /*ok = GmshSetOption("Mesh", "MaxNumThreads1D" , 0. ) ; // Coarse-grain algo threads
470 ok = GmshSetOption("Mesh", "MaxNumThreads2D" , 0. ) ; // Coarse-grain algo threads
472 ok = GmshSetOption("Mesh", "MaxNumThreads3D" , 0. ) ; // Fine-grain algo threads HXT
475 ok = GmshSetOption("General", "NumThreads" , _maxThreads ) ; // system default i.e. OMP_NUM_THREADS
478 if ( GMSHPlugin_Hypothesis::Algo3D::hxt == _algo3d ){
479 MESSAGE("GMSHPlugin_Mesher::SetGmshOptions: HXT algorithm is being used. Setting number of threads to 1.");
480 ok = GmshSetOption("Mesh", "MaxNumThreads3D" , 1. );
487 //================================================================================
489 * \brief Create and add Compounds into GModel _gModel.
491 //================================================================================
493 void GMSHPlugin_Mesher::CreateGmshCompounds()
495 MESSAGE("GMSHPlugin_Mesher::CreateGmshCompounds");
497 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
499 OCC_Internals* occgeo = _gModel->getOCCInternals();
500 bool toSynchronize = false;
502 for(std::set<std::string>::const_iterator its = _compounds.begin();its != _compounds.end(); ++its )
504 GEOM::GEOM_Object_var aGeomObj;
505 TopoDS_Shape geomShape = TopoDS_Shape();
506 SALOMEDS::SObject_var aSObj = SMESH_Gen_i::GetSMESHGen()->getStudyServant()->FindObjectID( (*its).c_str() );
507 SALOMEDS::GenericAttribute_var anAttr;
508 if (!aSObj->_is_nil() && aSObj->FindAttribute(anAttr, "AttributeIOR"))
510 SALOMEDS::AttributeIOR_var anIOR = SALOMEDS::AttributeIOR::_narrow(anAttr);
511 CORBA::String_var aVal = anIOR->Value();
512 CORBA::Object_var obj = SMESH_Gen_i::GetSMESHGen()->getStudyServant()->ConvertIORToObject(aVal);
513 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
515 geomShape = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
516 if ( geomShape.IsNull() )
519 TopAbs_ShapeEnum geomType = geomShape.ShapeType();
520 if ( geomType == TopAbs_COMPOUND)// voir s'il ne faut pas mettre une erreur dans le cas contraire
522 MESSAGE("shapeType == TopAbs_COMPOUND");
523 TopoDS_Iterator it(geomShape);
526 TopAbs_ShapeEnum shapeType = it.Value().ShapeType();
527 std::vector< std::pair< int, int > > dimTags;
528 for ( ; it.More(); it.Next())
530 const TopoDS_Shape& topoShape = it.Value();
531 ASSERT(topoShape.ShapeType() == shapeType);
532 if ( _mesh->GetMeshDS()->ShapeToIndex( topoShape ) > 0 )
533 occgeo->importShapes( &topoShape, false, dimTags );
536 TopoDS_Shape face = TopExp_Explorer( _shape, shapeType ).Current();
537 SMESH_subMesh* sm = _mesh->GetSubMesh( face );
538 sm->GetComputeError() =
539 SMESH_ComputeError::New
540 ( COMPERR_WARNING, "Compound shape does not belong to the main geometry. Ignored");
543 std::vector<int> tags;
544 int dim = ( shapeType == TopAbs_EDGE ) ? 1 : 2;
545 for ( size_t i = 0; i < dimTags.size(); ++i )
547 if ( dimTags[i].first == dim )
548 tags.push_back( dimTags[i].second );
552 _gModel->getGEOInternals()->setCompoundMesh( dim, tags );
553 toSynchronize = true;
556 _gModel->getGEOInternals()->synchronize(_gModel);
561 //================================================================================
563 * \brief For a compound mesh set the mesh components to be transmitted to SMESH
565 //================================================================================
567 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
568 void GMSHPlugin_Mesher::SetCompoundMeshVisibility()
571 // Loop over all faces, if the face belongs to a compound entry then
572 // for all (boundary) edges whithin the face visibility is set to 0,
573 // if the face doesn't belong to a compound entry then visibility is
574 // set to 1 for all its (boundary) edges. Later, in FillSMesh() func
575 // getVisibility() (returns either 1 or 0) is used to decide weather
576 // the mesh of edge should be transmitted to SMESH or not.
578 for ( GModel::fiter itF = _gModel->firstFace(); itF != _gModel->lastFace(); ++itF )
580 std::vector< GEdge *> faceEdges = (*itF)->edges();
582 for ( auto itE = faceEdges.begin(); itE != faceEdges.end(); ++itE )
584 if ( ((*itF)->compound.size()) )
585 (*itE)->setVisibility(0);
587 (*itE)->setVisibility(1);
592 // Loop over all edges, if the edge belongs to a compound entry then
593 // for all (boundary) vertices whithin the edge visibility is set to
594 // 0, if the edge doesn't belong to a compound entry then visibility
595 // is set to 1 for all its (boundary) vertices. Later, in FillSMesh()
596 // func getVisibility() (returns either 1 or 0) is used to decide we-
597 // ather the mesh of vertices should be transmitted to SMESH or not.
599 for ( GModel::eiter itE = _gModel->firstEdge(); itE != _gModel->lastEdge(); ++itE )
601 std::vector<GVertex *> bndVerticies = (*itE)->vertices();
603 for( auto itV = bndVerticies.begin(); itV != bndVerticies.end(); ++itV )
605 if(((*itE)->compound.size()))
606 (*itV)->setVisibility(0);
608 (*itV)->setVisibility(1);
615 //================================================================================
617 * \brief Get a node by a GMSH mesh vertex
619 //================================================================================
621 const SMDS_MeshNode* GMSHPlugin_Mesher::Node( const MVertex* v )
623 std::map< const MVertex *, const SMDS_MeshNode* >::iterator v2n = _nodeMap.find( v );
624 if ( v2n != _nodeMap.end() )
630 //================================================================================
632 * \brief Return a corresponding sub-mesh if a shape is meshed
634 //================================================================================
636 SMESHDS_SubMesh* GMSHPlugin_Mesher::HasSubMesh( const TopoDS_Shape& s )
638 if ( SMESHDS_SubMesh* sm = _mesh->GetMeshDS()->MeshElements( s ))
640 if ( s.ShapeType() == TopAbs_VERTEX )
641 return ( sm->NbNodes() > 0 ) ? sm : nullptr;
643 return ( sm->NbElements() > 0 ) ? sm : nullptr;
648 //================================================================================
650 * \brief Write mesh from GModel instance to SMESH instance
652 //================================================================================
654 void GMSHPlugin_Mesher::FillSMesh()
656 SMESHDS_Mesh* meshDS = _mesh->GetMeshDS();
659 for ( GModel::viter it = _gModel->firstVertex(); it != _gModel->lastVertex(); ++it)
661 GVertex *gVertex = *it;
663 // GET topoVertex CORRESPONDING TO gVertex
664 TopoDS_Vertex topoVertex = *((TopoDS_Vertex*)gVertex->getNativePtr());
666 if (gVertex->getVisibility() == 0) // belongs to a compound
668 SMESH_subMesh* sm = _mesh->GetSubMesh(topoVertex);
669 sm->SetIsAlwaysComputed(true); // prevent from displaying errors
673 // FILL SMESH FOR topoVertex
675 for( size_t i = 0; i < gVertex->mesh_vertices.size(); i++)
677 MVertex *v = gVertex->mesh_vertices[i];
678 if(v->getIndex() >= 0)
680 if ( SMESHDS_SubMesh* sm = HasSubMesh( topoVertex ))
682 const SMDS_MeshNode *node = sm->GetNodes()->next();
683 _nodeMap.insert({ v, node });
687 SMDS_MeshNode *node = meshDS->AddNode( v->x(),v->y(),v->z() );
688 meshDS->SetNodeOnVertex( node, topoVertex );
689 _nodeMap.insert({ v, node });
693 // WE DON'T ADD 0D ELEMENTS because it does not follow the salome meshers philosophy
695 // for(unsigned int i = 0; i < gVertex->getNumMeshElements(); i++)
697 // MElement *e = gVertex->getMeshElement(i);
698 // std::vector<MVertex*> verts;
699 // e->getVertices(verts);
700 // ASSERT(verts.size()==1);
701 // SMDS_Mesh0DElement* zeroDElement = 0;
702 // zeroDElement = meshDS->Add0DElementWithID(verts[0]->getNum(),e->getNum());
703 // meshDS->SetMeshElementOnShape(zeroDElement, topoVertex);
708 for(GModel::eiter it = _gModel->firstEdge(); it != _gModel->lastEdge(); ++it)
712 // GET topoEdge CORRESPONDING TO gEdge
713 TopoDS_Edge topoEdge;
714 std::vector< ShapeBounds > topoEdges;
715 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
716 if(gEdge->haveParametrization())
718 if ( gEdge->geomType() != GEntity::CompoundCurve )
721 topoEdge = *((TopoDS_Edge*)gEdge->getNativePtr());
722 if (gEdge->getVisibility() == 0) // belongs to a compound
724 SMESH_subMesh* sm = _mesh->GetSubMesh(topoEdge);
725 sm->SetIsAlwaysComputed(true); // prevent from displaying errors
728 if ( HasSubMesh( topoEdge ))
729 continue; // a meshed sub-mesh
731 bool isCompound = getBoundsOfShapes( gEdge, topoEdges );
733 // FILL SMESH FOR topoEdge
735 for ( size_t i = 0; i < gEdge->mesh_vertices.size(); i++ )
737 MVertex *v = gEdge->mesh_vertices[i];
738 if ( v->getIndex() >= 0 )
740 SMDS_MeshNode *node = meshDS->AddNode( v->x(),v->y(),v->z() );
743 topoEdge = TopoDS::Edge( getShapeAtPoint( v->point(), topoEdges ));
745 // Based on BLSURFPlugin_BLSURF
746 gp_Pnt point3D( v->x(),v->y(),v->z() );
747 Standard_Real p0 = 0.0;
748 Standard_Real p1 = 1.0;
750 Handle(Geom_Curve) curve = BRep_Tool::Curve(topoEdge, loc, p0, p1);
752 if ( !curve.IsNull() )
754 if ( !loc.IsIdentity() )
755 point3D.Transform( loc.Transformation().Inverted() );
757 GeomAPI_ProjectPointOnCurve proj(point3D, curve, p0, p1);
760 if ( proj.NbPoints() > 0 )
761 pa = (double)proj.LowerDistanceParameter();
763 meshDS->SetNodeOnEdge( node, topoEdge, pa );
767 meshDS->SetNodeOnEdge( node, topoEdge );
769 //END on BLSURFPlugin_BLSURF
772 _nodeMap.insert({ v, node });
777 for ( GModel::eiter it = _gModel->firstEdge(); it != _gModel->lastEdge(); ++it )
780 if ( gEdge->getVisibility() == 0) // belongs to a compound
783 TopoDS_Edge topoEdge;
784 std::vector< ShapeBounds > topoEdges;
785 bool isCompound = getBoundsOfShapes( gEdge, topoEdges );
787 topoEdge = *((TopoDS_Edge*)gEdge->getNativePtr());
789 if ( HasSubMesh( topoEdge ))
790 continue; // a meshed sub-mesh
793 std::vector<MVertex*> verts(3);
794 for ( size_t i = 0; i < gEdge->getNumMeshElements(); i++ )
796 MElement *e = gEdge->getMeshElement(i);
798 e->getVertices(verts);
800 // if a node wasn't set, it is assigned here
801 for ( size_t j = 0; j < verts.size(); j++ )
803 if ( verts[j]->onWhat()->getVisibility() == 0 )
805 SMDS_MeshNode *node = meshDS->AddNode(verts[j]->x(),verts[j]->y(),verts[j]->z() );
807 gp_Pnt point3D( verts[j]->x(),verts[j]->y(),verts[j]->z() );
808 Standard_Real p0 = 0.0;
809 Standard_Real p1 = 1.0;
811 Handle(Geom_Curve) curve = BRep_Tool::Curve(topoEdge, loc, p0, p1);
813 if ( !curve.IsNull() )
815 if ( !loc.IsIdentity() )
816 point3D.Transform( loc.Transformation().Inverted() );
818 GeomAPI_ProjectPointOnCurve proj(point3D, curve, p0, p1);
821 if ( proj.NbPoints() > 0 )
822 pa = (double)proj.LowerDistanceParameter();
824 meshDS->SetNodeOnEdge( node, topoEdge, pa );
828 meshDS->SetNodeOnEdge( node, topoEdge );
831 verts[j]->setEntity(gEdge);
832 _nodeMap.insert({ verts[j], node });
836 SMDS_MeshEdge* edge = 0;
837 switch (verts.size())
840 edge = meshDS->AddEdge(Node( verts[0]),
844 edge = meshDS->AddEdge(Node( verts[0]),
853 topoEdge = TopoDS::Edge( getShapeAtPoint( e->barycenter(), topoEdges ));
855 meshDS->SetMeshElementOnShape( edge, topoEdge );
860 for ( GModel::fiter it = _gModel->firstFace(); it != _gModel->lastFace(); ++it)
864 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
865 // Gmsh since version 4.3 is now producing extra surface and mesh when
866 // compounds are involved. Since in Gmsh meshing procedure needs acess
867 // to each of the original topology and the meshed topology. Hence we
868 // bypass the additional mesh in case of compounds. Note, similar cri-
869 // teria also occurs in the following 'for' loop.
870 if ( _compounds.size() && gFace->geomType() == GEntity::DiscreteSurface )
874 // GET topoFace CORRESPONDING TO gFace
875 TopoDS_Face topoFace;
876 std::vector< ShapeBounds > topoFaces;
878 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
879 if(gFace->haveParametrization())
881 if ( gFace->geomType() != GEntity::CompoundSurface )
884 topoFace = *((TopoDS_Face*)gFace->getNativePtr());
885 if (gFace->getVisibility() == 0) // belongs to a compound
887 SMESH_subMesh* sm = _mesh->GetSubMesh(topoFace);
888 sm->SetIsAlwaysComputed(true); // prevent from displaying errors
891 if ( HasSubMesh( topoFace ))
892 continue; // a meshed sub-mesh
894 bool isCompound = getBoundsOfShapes( gFace, topoFaces );
896 // FILL SMESH FOR topoFace
898 for ( size_t i = 0; i < gFace->mesh_vertices.size(); i++ )
900 MVertex *v = gFace->mesh_vertices[i];
901 if ( v->getIndex() >= 0 )
903 SMDS_MeshNode *node = meshDS->AddNode( v->x(),v->y(),v->z() );
906 topoFace = TopoDS::Face( getShapeAtPoint( v->point(), topoFaces ));
908 meshDS->SetNodeOnFace( node, topoFace );
909 _nodeMap.insert({ v, node });
914 for ( GModel::fiter it = _gModel->firstFace(); it != _gModel->lastFace(); ++it)
918 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
919 if ( _compounds.size() && gFace->geomType() == GEntity::DiscreteSurface )
922 bool isCompound = (!gFace->haveParametrization());
924 bool isCompound = ( gFace->geomType() == GEntity::CompoundSurface );
927 if ( !isCompound && gFace->getVisibility() == 0 )
928 continue; // belongs to a compound
930 TopoDS_Face topoFace;
931 std::vector< ShapeBounds > topoFaces;
933 getBoundsOfShapes( gFace, topoFaces );
935 topoFace = *((TopoDS_Face*)gFace->getNativePtr());
937 if ( HasSubMesh( topoFace ))
938 continue; // a meshed sub-mesh
941 std::vector<MVertex*> verts;
942 for ( size_t i = 0; i < gFace->getNumMeshElements(); i++ )
944 MElement *e = gFace->getMeshElement(i);
946 e->getVertices(verts);
947 SMDS_MeshFace* face = 0;
949 // if a node wasn't set, it is assigned here
950 for ( size_t j = 0; j < verts.size(); j++)
952 if(verts[j]->onWhat()->getVisibility() == 0)
954 SMDS_MeshNode *node = meshDS->AddNode(verts[j]->x(),verts[j]->y(),verts[j]->z());
955 meshDS->SetNodeOnFace( node, topoFace );
956 _nodeMap.insert({ verts[j], node });
957 verts[j]->setEntity(gFace);
960 switch (verts.size())
963 face = meshDS->AddFace(Node( verts[0]),
968 face = meshDS->AddFace(Node( verts[0]),
974 face = meshDS->AddFace(Node( verts[0]),
982 face = meshDS->AddFace(Node( verts[0]),
992 face = meshDS->AddFace(Node( verts[0]),
1008 topoFace = TopoDS::Face( getShapeAtPoint( e->barycenter(), topoFaces ));
1010 meshDS->SetMeshElementOnShape(face, topoFace);
1015 for ( GModel::riter it = _gModel->firstRegion(); it != _gModel->lastRegion(); ++it)
1017 GRegion *gRegion = *it;
1018 if (gRegion->getVisibility() == 0)
1021 // GET topoSolid CORRESPONDING TO gRegion
1022 TopoDS_Solid topoSolid = *((TopoDS_Solid*)gRegion->getNativePtr());
1024 // FILL SMESH FOR topoSolid
1027 for( size_t i = 0; i < gRegion->mesh_vertices.size(); i++)
1029 MVertex *v = gRegion->mesh_vertices[i];
1030 if(v->getIndex() >= 0)
1032 SMDS_MeshNode *node = meshDS->AddNode( v->x(),v->y(),v->z() );
1033 meshDS->SetNodeInVolume( node, topoSolid );
1034 _nodeMap.insert({ v, node });
1039 std::vector<MVertex*> verts;
1040 for( size_t i = 0; i < gRegion->getNumMeshElements(); i++)
1042 MElement *e = gRegion->getMeshElement(i);
1044 e->getVertices(verts);
1045 SMDS_MeshVolume* volume = 0;
1046 switch (verts.size()){
1048 volume = meshDS->AddVolume(Node( verts[0]),
1054 volume = meshDS->AddVolume(Node( verts[0]),
1061 volume = meshDS->AddVolume(Node( verts[0]),
1069 volume = meshDS->AddVolume(Node( verts[0]),
1079 volume = meshDS->AddVolume(Node( verts[0]),
1091 volume = meshDS->AddVolume(Node( verts[0] ),
1105 case 14: // same as case 13, because no pyra14 in smesh
1106 volume = meshDS->AddVolume(Node( verts[0] ),
1121 volume = meshDS->AddVolume(Node( verts[0] ),
1137 case 18: // same as case 15, because no penta18 in smesh
1138 volume = meshDS->AddVolume(Node( verts[0] ),
1155 volume = meshDS->AddVolume(Node( verts[0] ),
1177 volume = meshDS->AddVolume(Node( verts[0] ),
1209 meshDS->SetMeshElementOnShape(volume, topoSolid);
1216 //================================================================================
1218 * \brief Find if SPoint point is in SBoundingBox3d bounds
1220 //================================================================================
1222 float GMSHPlugin_Mesher::DistBoundingBox(const SBoundingBox3d& bounds, const SPoint3& point)
1224 SPoint3 min = bounds.min();
1225 SPoint3 max = bounds.max();
1229 if (point.x() < min.x())
1230 x = min.x()-point.x();
1231 else if (point.x() > max.x())
1232 x = point.x()-max.x();
1236 if (point.y() < min.y())
1237 y = min.y()-point.y();
1238 else if (point.y() > max.y())
1239 y = point.y()-max.y();
1243 if (point.z() < min.z())
1244 z = min.z()-point.z();
1245 else if (point.z() > max.z())
1246 z = point.z()-max.z();
1252 //================================================================================
1254 * \brief Reimplemented GmshMessage call. Actions done if errors occurs
1255 * during gmsh meshing. We define here what to display and what to do.
1257 //================================================================================
1258 void GMSHPlugin_Mesher::mymsg::operator()(std::string level, std::string msg)
1260 //MESSAGE("level="<< level.c_str() << ", msg=" << msg.c_str()<< "\n");
1261 printf("level=%s msg=%s\n", level.c_str(), msg.c_str());
1263 if(level == "Fatal" || level == "Error")
1265 std::ostringstream oss;
1266 if (level == "Fatal")
1267 oss << "Fatal error during Generation of Gmsh Mesh\n";
1269 oss << "Error during Generation of Gmsh Mesh\n";
1270 oss << " " << msg.c_str() << "\n";
1271 GEntity *e = _gModel->getCurrentMeshEntity();
1274 oss << " error occurred while meshing entity:\n" <<
1275 " tag=" << e->tag() << "\n" <<
1276 " dimension=" << e->dim() << "\n" <<
1277 " native pointer=" << e->getNativePtr();
1278 //if(e->geomType() != GEntity::CompoundCurve and e->geomType() != GEntity::CompoundSurface)
1280 //SMESH_subMesh *sm = _mesh->GetSubMesh(*((TopoDS_Shape*)e->getNativePtr()));
1281 //SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1282 //SMESH_Comment comment;
1283 //comment << SMESH_Comment(oss.str);
1284 //std::string str = oss.str();
1285 //smError.reset( new SMESH_ComputeError( str ));
1287 // plutot que de faire de la merde ici, on pourait simplement
1288 // remplir une liste pour dire sur quelles entités gmsh se plante
1289 // (puis faire le fillsmesh)
1290 // puis faire une nouvelle routine qui réécrit les messages d'erreur
1291 // probleme : gmsh peut planté en Fatal, dans ce cas pas de fillsmesh
1294 if (level == "Fatal")
1296 CTX::instance()->lock = 0;
1300 printf("%s\n", oss.str().c_str());
1304 //=============================================================================
1306 * Here we are going to use the GMSH mesher
1308 //=============================================================================
1310 bool GMSHPlugin_Mesher::Compute()
1312 MESSAGE("GMSHPlugin_Mesher::Compute");
1316 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1317 SetMaxThreadsGmsh();
1319 //RNV: to avoid modification of PATH and PYTHONPATH
1320 char* argv[] = {"-noenv"};
1321 GmshInitialize(1,argv);
1323 _gModel = new GModel();
1325 GmshSetMessageHandler(&msg);
1326 _gModel->importOCCShape((void*)&_shape);
1327 if (_compounds.size() > 0) CreateGmshCompounds();
1331 HideComputedEntities(_gModel);
1336 Set2DSubMeshes(_gModel);
1338 _gModel->mesh( /*dim=*/ 3);
1342 //Msg::SetVerbosity(100);
1343 //CTX::instance()->mesh.maxNumThreads1D=1;
1345 _gModel->mesh( /*dim=*/ 1);
1347 Set1DSubMeshes(_gModel);
1349 //_gModel->writeUNV("/tmp/1D.unv", 1,0,0);
1350 //CTX::instance()->mesh.maxNumThreads2D=1;
1352 _gModel->mesh( /*dim=*/ 2);
1356 Set2DSubMeshes(_gModel);
1358 //CTX::instance()->mesh.maxNumThreads3D=1;
1360 _gModel->mesh( /*dim=*/ 3);
1362 RestoreVisibility(_gModel);
1364 #ifdef WITH_SMESH_CANCEL_COMPUTE
1368 catch (std::string& str)
1371 std::cerr << "GMSH: exception caught: " << str << std::endl;
1377 std::cerr << "GMSH: Unknown exception caught: " << std::endl;
1378 MESSAGE("Unrecoverable error during Generation of Gmsh Mesh");
1383 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1384 if (_compounds.size() > 0)
1385 SetCompoundMeshVisibility();
1389 GmshSetMessageHandler(nullptr);
1392 MESSAGE("GMSHPlugin_Mesher::Compute:End");
1396 //================================================================================
1398 * \brief Set 1D sub-meshes to GModel. GMSH 1D mesh is made by now.
1399 * \param [inout] _gModel - GMSH model
1401 //================================================================================
1403 void GMSHPlugin_Mesher::Set1DSubMeshes( GModel* gModel )
1405 SMESHDS_Mesh* meshDS = _mesh->GetMeshDS();
1407 for(GModel::eiter it = gModel->firstEdge(); it != gModel->lastEdge(); ++it)
1411 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1412 if ( !gEdge->haveParametrization())
1414 if ( gEdge->geomType() == GEntity::CompoundCurve )
1418 TopoDS_Edge topoEdge = *((TopoDS_Edge*)gEdge->getNativePtr());
1419 if ( !HasSubMesh( topoEdge ))
1420 continue; // empty sub-mesh
1422 gEdge->deleteMesh();
1424 // get node parameters on topoEdge
1425 StdMeshers_FaceSide side( TopoDS_Face(), topoEdge, _mesh, /*fwd=*/true, /*skpMedium=*/true);
1426 const UVPtStructVec& nodeParam = side.GetUVPtStruct();
1427 if ( nodeParam.empty() )
1428 throw std::string("Pb with StdMeshers_FaceSide::GetUVPtStruct()");
1430 // get GMSH mesh vertices on VERTEX'es
1431 std::vector<MVertex *> mVertices( nodeParam.size(), nullptr );
1432 GVertex * gV0 = gEdge->getBeginVertex(), *gV1 = gEdge->getEndVertex();
1433 mVertices[0] = gV0->mesh_vertices[ 0 ];
1434 mVertices.back() = gV1->mesh_vertices[ 0 ];
1435 TopoDS_Vertex v01 = *((TopoDS_Vertex*) gV0->getNativePtr());
1436 TopoDS_Shape v02 = SMESH_MesherHelper::GetSubShapeByNode( nodeParam[0].node, meshDS );
1437 bool reverse = !v01.IsSame( v02 );
1438 if ( mVertices[0] == mVertices.back() )
1439 reverse = ( nodeParam[0].param > nodeParam.back().param );
1440 const SMDS_MeshNode* n0 = reverse ? nodeParam.back().node : nodeParam[0].node;
1441 const SMDS_MeshNode* n1 = reverse ? nodeParam[0].node : nodeParam.back().node;
1442 _nodeMap.insert({ mVertices[ 0 ], n0 });
1443 _nodeMap.insert({ mVertices.back(), n1 });
1445 // create GMSH mesh vertices on gEdge
1446 for ( size_t i = 1; i < nodeParam.size() - 1; ++i )
1448 size_t iN = reverse ? ( nodeParam.size() - 1 - i ) : i;
1449 SMESH_NodeXYZ xyz = nodeParam[ iN ].node;
1450 double lc = segmentSize( nodeParam, iN );
1451 // SVector3 der = gEdge->firstDer(nodeParam[ iN ].param);
1452 // double lc = norm(der) / segmentSize( nodeParam, i );
1454 mVertices[ i ] = new MEdgeVertex( xyz.X(), xyz.Y(), xyz.Z(),
1455 gEdge, nodeParam[ iN ].param, 0, lc);
1456 gEdge->mesh_vertices.push_back( mVertices[ i ]);
1457 _nodeMap.insert({ mVertices[ i ], nodeParam[ iN ].node });
1459 // create GMSH mesh edges
1460 for ( size_t i = 1; i < mVertices.size(); ++i )
1462 gEdge->lines.push_back( new MLine( mVertices[ i - 1 ],
1466 cout << endl << "EDGE " << gEdge->tag() <<
1467 ( topoEdge.Orientation() == TopAbs_FORWARD ? " F" : " R") << endl;
1468 MVertex* mv = gV0->mesh_vertices[ 0 ];
1469 cout << "V0: " << mv->x() << ", " << mv->y() << ", " << mv->z() << ", "<<endl;
1470 for ( size_t i = 0; i < gEdge->mesh_vertices.size(); ++i )
1472 MEdgeVertex* mv = (MEdgeVertex*) gEdge->mesh_vertices[i];
1473 cout << i << ": " << mv->x() << ", " << mv->y() << ", " << mv->z() << ", ";
1475 mv->getParameter(0, t );
1476 cout << ":\t t = "<< t << " lc = " << mv->getLc() << endl;
1478 mv = gV1->mesh_vertices[ 0 ];
1479 cout << "V1: " << mv->x() << ", " << mv->y() << ", " << mv->z() << ", "<<endl;
1485 //================================================================================
1487 * \brief Set 2D sub-meshes to GModel. GMSH 2D mesh is made by now.
1488 * \param [inout] _gModel - GMSH model
1490 //================================================================================
1492 void GMSHPlugin_Mesher::Set2DSubMeshes( GModel* gModel )
1494 if ( _nodeMap.empty() )
1495 return; // no sub-meshes
1497 SMESH_MesherHelper helper( *_mesh );
1499 std::map< const SMDS_MeshNode* , const MVertex * > nodes2mvertMap;
1500 for ( auto & v2n : _nodeMap )
1501 nodes2mvertMap.insert({ v2n.second, v2n.first });
1503 std::vector<MVertex *> mVertices;
1505 for(GModel::fiter it = gModel->firstFace(); it != gModel->lastFace(); ++it)
1509 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1510 if ( !gFace->haveParametrization())
1512 if ( gFace->geomType() == GEntity::CompoundSurface )
1516 TopoDS_Face topoFace = *((TopoDS_Face*)gFace->getNativePtr());
1517 SMESHDS_SubMesh* sm = HasSubMesh( topoFace );
1520 //_gModel->writeUNV("/tmp/befDEl.unv", 1,0,0);
1522 gFace->deleteMesh();
1524 bool reverse = false;
1525 if ( gFace->getRegion(0) )
1527 //GRegion * gRegion = gFace->getRegion(0);
1528 TopoDS_Shape topoSolid = *((TopoDS_Shape*)gFace->getNativePtr());
1529 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( topoSolid, topoFace );
1530 if ( faceOriInSolid >= 0 )
1532 helper.IsReversedSubMesh( TopoDS::Face( topoFace.Oriented( faceOriInSolid )));
1535 for ( SMDS_ElemIteratorPtr fIt = sm->GetElements(); fIt->more(); )
1537 const SMDS_MeshElement* f = fIt->next();
1539 int nbN = f->NbCornerNodes();
1541 throw std::string("Polygon sub-meshes not supported");
1543 mVertices.resize( nbN );
1544 for ( int i = 0; i < nbN; ++i )
1546 const SMDS_MeshNode* n = f->GetNode( i );
1547 MVertex * mv = nullptr;
1548 auto n2v = nodes2mvertMap.find( n );
1549 if ( n2v != nodes2mvertMap.end() )
1551 mv = const_cast< MVertex*>( n2v->second );
1555 if ( n->GetPosition()->GetDim() < 2 )
1556 throw std::string("Wrong mapping of edge nodes to GMSH nodes");
1557 SMESH_NodeXYZ xyz = n;
1559 gp_XY uv = helper.GetNodeUV( topoFace, n, nullptr, &ok );
1560 mv = new MFaceVertex( xyz.X(), xyz.Y(), xyz.Z(), gFace, uv.X(), uv.Y() );
1561 gFace->mesh_vertices.push_back( mv );
1562 nodes2mvertMap.insert({ n, mv });
1563 _nodeMap.insert ({ mv, n });
1565 mVertices[ i ] = mv;
1567 // create GMSH mesh faces
1571 gFace->triangles.push_back (new MTriangle(mVertices[0], mVertices[2], mVertices[1]));
1573 gFace->triangles.push_back (new MTriangle(mVertices[0], mVertices[1], mVertices[2]));
1577 gFace->quadrangles.push_back (new MQuadrangle(mVertices[0], mVertices[3],
1578 mVertices[2], mVertices[1]));
1580 gFace->quadrangles.push_back (new MQuadrangle(mVertices[0], mVertices[1],
1581 mVertices[2], mVertices[3]));
1586 } // loop on GMSH faces
1591 //================================================================================
1593 * \brief Set visibility 0 to already computed geom entities
1594 * to prevent their meshing
1596 //================================================================================
1598 void GMSHPlugin_Mesher::HideComputedEntities( GModel* gModel )
1600 CTX::instance()->mesh.meshOnlyVisible = true;
1602 for(GModel::eiter it = gModel->firstEdge(); it != gModel->lastEdge(); ++it)
1606 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1607 if ( !gEdge->haveParametrization())
1609 if ( gEdge->geomType() == GEntity::CompoundCurve )
1613 TopoDS_Edge topoEdge = *((TopoDS_Edge*)gEdge->getNativePtr());
1614 if ( HasSubMesh( topoEdge ))
1615 gEdge->setVisibility(0);
1619 for(GModel::fiter it = gModel->firstFace(); it != gModel->lastFace(); ++it)
1623 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1624 if ( !gFace->haveParametrization())
1626 if ( gFace->geomType() == GEntity::CompoundSurface )
1630 TopoDS_Face topoFace = *((TopoDS_Face*)gFace->getNativePtr());
1631 if ( HasSubMesh( topoFace ))
1632 gFace->setVisibility(0);
1636 //================================================================================
1638 * \brief Restore visibility of all geom entities
1640 //================================================================================
1642 void GMSHPlugin_Mesher::RestoreVisibility( GModel* gModel )
1644 for(GModel::eiter it = gModel->firstEdge(); it != gModel->lastEdge(); ++it)
1647 gEdge->setVisibility(1);
1649 for(GModel::fiter it = gModel->firstFace(); it != gModel->lastFace(); ++it)
1652 gFace->setVisibility(1);
1657 void GMSHPlugin_Mesher::toPython( GModel* )
1659 const char* pyFile = "/tmp/gMesh.py";
1660 ofstream outfile( pyFile, ios::out );
1661 if ( !outfile ) return;
1663 outfile << "import salome, SMESH" << std::endl
1664 << "from salome.smesh import smeshBuilder" << std::endl
1665 << "smesh = smeshBuilder.New()" << std::endl
1666 << "mesh = smesh.Mesh()" << std::endl << std::endl;
1668 outfile << "## VERTICES" << endl;
1669 for ( GModel::viter it = _gModel->firstVertex(); it != _gModel->lastVertex(); ++it)
1671 GVertex *gVertex = *it;
1673 for(unsigned int i = 0; i < gVertex->mesh_vertices.size(); i++)
1675 MVertex *v = gVertex->mesh_vertices[i];
1676 if ( v->getIndex() >= 0)
1678 outfile << "n" << v->getNum() << " = mesh.AddNode("
1679 << v->x() << ", " << v->y() << ", " << v->z()<< ")"
1680 << " ## tag = " << gVertex->tag() << endl;
1685 for(GModel::eiter it = _gModel->firstEdge(); it != _gModel->lastEdge(); ++it)
1688 outfile << "## GEdge " << gEdge->tag() << endl;
1690 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1691 if(gEdge->haveParametrization())
1693 if ( gEdge->geomType() != GEntity::CompoundCurve )
1695 for ( size_t i = 0; i < gEdge->mesh_vertices.size(); i++ )
1697 MVertex *v = gEdge->mesh_vertices[i];
1698 if ( v->getIndex() >= 0 )
1700 outfile << "n" << v->getNum() << " = mesh.AddNode("
1701 << v->x() << ", " << v->y() << ", " << v->z()<< ")"
1702 << " ## tag = " << gEdge->tag() << endl;
1707 for ( GModel::fiter it = _gModel->firstFace(); it != _gModel->lastFace(); ++it)
1710 if ( _compounds.size() && gFace->geomType() == GEntity::DiscreteSurface )
1712 outfile << "## GFace " << gFace->tag() << endl;
1714 for ( size_t i = 0; i < gFace->mesh_vertices.size(); i++ )
1716 MVertex *v = gFace->mesh_vertices[i];
1717 if ( v->getIndex() >= 0 )
1719 outfile << "n" << v->getNum() << " = mesh.AddNode("
1720 << v->x() << ", " << v->y() << ", " << v->z()<< ")"
1721 << " ## tag = " << gFace->tag() << endl;
1726 std::vector<MVertex*> verts(3);
1727 for(GModel::eiter it = _gModel->firstEdge(); it != _gModel->lastEdge(); ++it)
1730 outfile << "## GEdge " << gEdge->tag() << endl;
1732 #if GMSH_MAJOR_VERSION >=4 && GMSH_MINOR_VERSION >=8
1733 if(gEdge->haveParametrization())
1735 if ( gEdge->geomType() != GEntity::CompoundCurve )
1737 for ( size_t i = 0; i < gEdge->getNumMeshElements(); i++ )
1739 MElement *e = gEdge->getMeshElement(i);
1741 e->getVertices(verts);
1743 outfile << "e" << e->getNum() << " = mesh.AddEdge(["
1744 << "n" << verts[0]->getNum() << ","
1745 << "n" << verts[1]->getNum();
1746 if ( verts.size() == 3 )
1747 outfile << "n" << verts[2]->getNum();
1748 outfile << "])"<< endl;
1752 for ( GModel::fiter it = _gModel->firstFace(); it != _gModel->lastFace(); ++it)
1755 if ( _compounds.size() && gFace->geomType() == GEntity::DiscreteSurface )
1757 outfile << "## GFace " << gFace->tag() << endl;
1759 for ( size_t i = 0; i < gFace->getNumMeshElements(); i++ )
1761 MElement *e = gFace->getMeshElement(i);
1763 e->getVertices(verts);
1765 outfile << "f" << e->getNum() << " = mesh.AddFace([";
1766 for ( size_t j = 0; j < verts.size(); j++)
1768 outfile << "n" << verts[j]->getNum();
1769 if ( j < verts.size()-1 )
1772 outfile << "])" << endl;
1775 std::cout << "Write " << pyFile << std::endl;