1 // Copyright (C) 2007-2015 CEA/DEN, EDF R&D
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 // Lesser General Public License for more details.
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 // File : HexoticPlugin_Hexotic.cxx
22 // Author : Lioka RAZAFINDRAZAKA (CEA)
25 #include "HexoticPlugin_Hexotic.hxx"
26 #include "HexoticPlugin_Hypothesis.hxx"
28 #include "utilities.h"
31 #include <sys/sysinfo.h>
44 #include <SMESHDS_Mesh.hxx>
45 #include <SMESHDS_GroupBase.hxx>
46 #include <SMESH_ComputeError.hxx>
47 #include <SMESH_File.hxx>
48 #include <SMESH_Gen.hxx>
49 #include <SMESH_HypoFilter.hxx>
50 #include <SMESH_MesherHelper.hxx>
51 #include <SMESH_subMesh.hxx>
57 #include <Standard_ProgramError.hxx>
59 #include <BRep_Tool.hxx>
60 #include <BRepBndLib.hxx>
61 #include <BRepClass3d_SolidClassifier.hxx>
62 #include <BRepBuilderAPI_MakeVertex.hxx>
63 #include <BRepExtrema_DistShapeShape.hxx>
64 #include <OSD_File.hxx>
65 #include <Precision.hxx>
66 #include <TopExp_Explorer.hxx>
67 #include <TopTools_MapOfShape.hxx>
71 #include <GCPnts_UniformAbscissa.hxx>
72 #include <IntCurvesFace_Intersector.hxx>
73 #include <BRepMesh_IncrementalMesh.hxx>
74 #include <Poly_Triangulation.hxx>
76 #include <GeomAdaptor_Curve.hxx>
78 #include <Basics_Utils.hxx>
79 #include "GEOMImpl_Types.hxx"
80 #include "GEOM_wrap.hxx"
82 static void removeFile( const TCollection_AsciiString& fileName )
85 OSD_File( fileName ).Remove();
87 catch ( Standard_ProgramError ) {
88 MESSAGE("Can't remove file: " << fileName.ToCString() << " ; file does not exist or permission denied");
92 //=============================================================================
96 //=============================================================================
98 HexoticPlugin_Hexotic::HexoticPlugin_Hexotic(int hypId, int studyId, SMESH_Gen* gen)
99 : SMESH_3D_Algo(hypId, studyId, gen)
101 MESSAGE("HexoticPlugin_Hexotic::HexoticPlugin_Hexotic");
103 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
104 // _onlyUnaryInput = false;
105 _requireShape = false;
108 _hexoticFilesKept=false;
109 _compatibleHypothesis.push_back( HexoticPlugin_Hypothesis::GetHypType() );
110 #ifdef WITH_BLSURFPLUGIN
113 _compute_canceled = false;
115 // Copy of what is done in BLSURFPLugin TODO : share the code
116 smeshGen_i = SMESH_Gen_i::GetSMESHGen();
117 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
118 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
121 myStudy = aStudyMgr->GetStudyByID(_studyId);
123 MESSAGE("myStudy->StudyId() = " << myStudy->StudyId());
126 //=============================================================================
130 //=============================================================================
132 HexoticPlugin_Hexotic::~HexoticPlugin_Hexotic()
134 MESSAGE("HexoticPlugin_Hexotic::~HexoticPlugin_Hexotic");
138 #ifdef WITH_BLSURFPLUGIN
139 bool HexoticPlugin_Hexotic::CheckBLSURFHypothesis( SMESH_Mesh& aMesh,
140 const TopoDS_Shape& aShape )
142 // MESSAGE("HexoticPlugin_Hexotic::CheckBLSURFHypothesis");
145 std::list<const SMESHDS_Hypothesis*>::const_iterator itl;
146 const SMESHDS_Hypothesis* theHyp;
148 // If a BLSURF hypothesis is applied, get it
149 SMESH_HypoFilter blsurfFilter;
150 blsurfFilter.Init( blsurfFilter.HasName( BLSURFPlugin_Hypothesis::GetHypType() ));
151 std::list<const SMESHDS_Hypothesis *> appliedHyps;
152 aMesh.GetHypotheses( aShape, blsurfFilter, appliedHyps, false );
154 if ( appliedHyps.size() > 0 ) {
155 itl = appliedHyps.begin();
156 theHyp = (*itl); // use only the first hypothesis
157 std::string hypName = theHyp->GetName();
158 if (hypName == BLSURFPlugin_Hypothesis::GetHypType()) {
159 _blsurfHypo = static_cast<const BLSURFPlugin_Hypothesis*> (theHyp);
168 //=============================================================================
172 //=============================================================================
174 bool HexoticPlugin_Hexotic::CheckHypothesis( SMESH_Mesh& aMesh,
175 const TopoDS_Shape& aShape,
176 SMESH_Hypothesis::Hypothesis_Status& aStatus )
178 // MESSAGE("HexoticPlugin_Hexotic::CheckHypothesis");
181 std::list<const SMESHDS_Hypothesis*>::const_iterator itl;
182 const SMESHDS_Hypothesis* theHyp;
184 const std::list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape, false);
185 int nbHyp = hyps.size();
187 aStatus = SMESH_Hypothesis::HYP_OK;
188 return true; // can work with no hypothesis
192 theHyp = (*itl); // use only the first hypothesis
194 std::string hypName = theHyp->GetName();
195 if (hypName == HexoticPlugin_Hypothesis::GetHypType() ) {
196 _hypothesis = static_cast<const HexoticPlugin_Hypothesis*> (theHyp);
198 aStatus = SMESH_Hypothesis::HYP_OK;
201 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
203 #ifdef WITH_BLSURFPLUGIN
204 CheckBLSURFHypothesis(aMesh, aShape);
207 return aStatus == SMESH_Hypothesis::HYP_OK;
210 //=======================================================================
211 //function : findShape
213 //=======================================================================
215 static TopoDS_Shape findShape(SMDS_MeshNode** t_Node,
217 const TopoDS_Shape* t_Shape,
222 int iShape, nbNode = 8;
224 for ( int i=0; i<3; i++ ) {
226 for ( int j=0; j<nbNode; j++ ) {
227 if ( i == 0) pntCoor[i] += t_Node[j]->X();
228 if ( i == 1) pntCoor[i] += t_Node[j]->Y();
229 if ( i == 2) pntCoor[i] += t_Node[j]->Z();
231 pntCoor[i] /= nbNode;
233 gp_Pnt aPnt(pntCoor[0], pntCoor[1], pntCoor[2]);
235 if ( aShape.IsNull() ) aShape = t_Shape[0];
236 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
237 if ( !(SC.State() == TopAbs_IN) ) {
239 for (iShape = 0; iShape < nShape && aShape.IsNull(); iShape++) {
240 if ( !( pntCoor[0] < t_Box[iShape][0] || t_Box[iShape][1] < pntCoor[0] ||
241 pntCoor[1] < t_Box[iShape][2] || t_Box[iShape][3] < pntCoor[1] ||
242 pntCoor[2] < t_Box[iShape][4] || t_Box[iShape][5] < pntCoor[2]) ) {
243 BRepClass3d_SolidClassifier SC (t_Shape[iShape], aPnt, Precision::Confusion());
244 if (SC.State() == TopAbs_IN)
245 aShape = t_Shape[iShape];
252 //=======================================================================
253 //function : findEdge
255 //=======================================================================
257 static int findEdge(const SMDS_MeshNode* aNode,
258 const SMESHDS_Mesh* aMesh,
260 const TopoDS_Shape* t_Edge) {
262 TopoDS_Shape aPntShape, foundEdge;
263 TopoDS_Vertex aVertex;
264 gp_Pnt aPnt( aNode->X(), aNode->Y(), aNode->Z() );
267 double nearest = RealLast(), *t_Dist;
268 double epsilon = Precision::Confusion();
270 t_Dist = new double[ nEdge ];
271 aPntShape = BRepBuilderAPI_MakeVertex( aPnt ).Shape();
272 aVertex = TopoDS::Vertex( aPntShape );
274 for ( ind=0; ind < nEdge; ind++ ) {
275 BRepExtrema_DistShapeShape aDistance ( aVertex, t_Edge[ind] );
276 t_Dist[ind] = aDistance.Value();
277 if ( t_Dist[ind] < nearest ) {
278 nearest = t_Dist[ind];
279 foundEdge = t_Edge[ind];
281 if ( nearest < epsilon )
287 return aMesh->ShapeToIndex( foundEdge );
290 //=======================================================================
291 //function : getNbShape
293 //=======================================================================
295 static int getNbShape(std::string aFile, std::string aString, int defaultValue=0) {
296 int number = defaultValue;
298 std::ifstream file(aFile.c_str());
299 while ( !file.eof() ) {
300 getline( file, aLine);
301 if ( aLine == aString ) {
302 getline( file, aLine);
303 std::istringstream stringFlux( aLine );
304 stringFlux >> number;
305 number = ( number + defaultValue + std::abs(number - defaultValue) ) / 2;
313 //=======================================================================
314 //function : countShape
316 //=======================================================================
318 template < class Mesh, class Shape >
319 static int countShape( Mesh* mesh, Shape shape ) {
320 TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
321 TopTools_MapOfShape mapShape;
323 for ( ; expShape.More(); expShape.Next() ) {
324 if (mapShape.Add(expShape.Current())) {
331 //=======================================================================
332 //function : getShape
334 //=======================================================================
336 template < class Mesh, class Shape, class Tab >
337 void getShape(Mesh* mesh, Shape shape, Tab *t_Shape) {
338 TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
339 TopTools_MapOfShape mapShape;
340 for ( int i=0; expShape.More(); expShape.Next() ) {
341 if (mapShape.Add(expShape.Current())) {
342 t_Shape[i] = expShape.Current();
349 //=======================================================================
350 //function : printWarning
352 //=======================================================================
354 static void printWarning(const int nbExpected, std::string aString, const int nbFound) {
356 cout << "WARNING : " << nbExpected << " " << aString << " expected, MG-Hexa has found " << nbFound << std::endl;
357 cout << "=======" << std::endl;
362 //=======================================================================
363 //function : removeHexoticFiles
365 //=======================================================================
367 static void removeHexoticFiles(TCollection_AsciiString file_In, TCollection_AsciiString file_Out) {
368 removeFile( file_In );
369 removeFile( file_Out );
372 //=======================================================================
373 //function : splitQuads
374 //purpose : splits all quadrangles into triangles
375 //=======================================================================
377 static void splitQuads(SMESH_Mesh& aMesh)
379 SMESH_MeshEditor spliter( &aMesh );
381 TIDSortedElemSet elems;
382 SMDS_ElemIteratorPtr eIt = aMesh.GetMeshDS()->elementsIterator();
384 elems.insert( elems.end(), eIt->next() );
386 spliter.QuadToTri ( elems, /*the13Diag=*/true);
389 //=======================================================================
390 //function : readResult
391 //purpose : Read GMF file in case of a mesh with geometry
392 //=======================================================================
394 static bool readResult(std::string theFile,
395 HexoticPlugin_Hexotic* theAlgo,
396 SMESHDS_Mesh* theMesh,
398 const TopoDS_Shape* tabShape,
401 // ---------------------------------
402 // Optimisation of the plugin ...
403 // Retrieve the correspondance edge --> shape
404 // (which is very costly) only when user
405 // has defined at least one group of edges
406 // which should be rare for a 3d mesh !
407 // ---------------------------------
409 bool retrieve_edges = false;
410 const std::set<SMESHDS_GroupBase*>& aGroups = theMesh->GetGroups();
411 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
412 for (; GrIt != aGroups.end(); GrIt++)
414 SMESHDS_GroupBase* aGrp = *GrIt;
417 if ( aGrp->GetType() == SMDSAbs_Edge )
419 retrieve_edges = true;
424 // ---------------------------------
425 // Read generated elements and nodes
426 // ---------------------------------
429 TopoDS_Vertex aVertex;
431 int EndOfFile = 0, nbElem = 0, nField = 9, nbRef = 0;
432 int aHexoticNodeID = 0, shapeID, hexoticShapeID;
433 const int IdShapeRef = 2;
434 int *tabID, *tabRef, *nodeAssigne;
435 bool *tabDummy, hasDummy = false;
436 double epsilon = Precision::Confusion();
437 std::map <std::string,int> mapField;
438 SMDS_MeshNode** HexoticNode;
439 TopoDS_Shape *tabCorner, *tabEdge;
441 const int nbDomains = countShape( theMesh, TopAbs_SHELL );
442 const int holeID = -1;
444 // tabID = new int[nbShape];
445 tabID = new int[nbDomains];
446 tabRef = new int[nField];
447 tabDummy = new bool[nField];
449 for (int i=0; i<nbDomains; i++)
451 if ( nbDomains == 1 )
452 tabID[0] = theMesh->ShapeToIndex( tabShape[0] );
454 mapField["MeshVersionFormatted"] = 0; tabRef[0] = 0; tabDummy[0] = false;
455 mapField["Dimension"] = 1; tabRef[1] = 0; tabDummy[1] = false;
456 mapField["Vertices"] = 2; tabRef[2] = 3; tabDummy[2] = true;
457 mapField["Corners"] = 3; tabRef[3] = 1; tabDummy[3] = false;
458 mapField["Edges"] = 4; tabRef[4] = 2; tabDummy[4] = true;
459 mapField["Ridges"] = 5; tabRef[5] = 1; tabDummy[5] = false;
460 mapField["Quadrilaterals"] = 6; tabRef[6] = 4; tabDummy[6] = true;
461 mapField["Hexahedra"] = 7; tabRef[7] = 8; tabDummy[7] = true;
462 mapField["End"] = 8; tabRef[8] = 0; tabDummy[0] = false;
464 SMDS_NodeIteratorPtr itOnHexoticInputNode = theMesh->nodesIterator();
465 while ( itOnHexoticInputNode->more() )
466 theMesh->RemoveNode( itOnHexoticInputNode->next() );
468 int nbVertices = getNbShape(theFile, "Vertices");
469 int nbCorners = getNbShape(theFile, "Corners", countShape( theMesh, TopAbs_VERTEX ));
470 int nbShapeEdge = countShape( theMesh, TopAbs_EDGE );
472 tabCorner = new TopoDS_Shape[ nbCorners ];
473 tabEdge = new TopoDS_Shape[ nbShapeEdge ];
474 nodeAssigne = new int[ nbVertices + 1 ];
475 HexoticNode = new SMDS_MeshNode*[ nbVertices + 1 ];
477 getShape(theMesh, TopAbs_VERTEX, tabCorner);
478 getShape(theMesh, TopAbs_EDGE, tabEdge);
480 MESSAGE("Read " << theFile << " file");
481 std::ifstream fileRes(theFile.c_str());
484 while ( EndOfFile == 0 ) {
488 if (mapField.count(token)) {
489 nField = mapField[token];
490 nbRef = tabRef[nField];
491 hasDummy = tabDummy[nField];
499 if ( nField < (mapField.size() - 1) && nField >= 0 )
503 case 0: { // "MeshVersionFormatted"
504 MESSAGE(token << " " << nbElem);
507 case 1: { // "Dimension"
508 MESSAGE("Mesh dimension " << nbElem << "D");
511 case 2: { // "Vertices"
512 MESSAGE("Read " << nbElem << " " << token);
515 SMDS_MeshNode * aHexoticNode;
517 coord = new double[nbRef];
518 for ( int iElem = 0; iElem < nbElem; iElem++ ) {
519 if(theAlgo->computeCanceled())
523 aHexoticID = iElem + 1;
524 for ( int iCoord = 0; iCoord < 3; iCoord++ )
525 fileRes >> coord[ iCoord ];
527 aHexoticNode = theMesh->AddNode(coord[0], coord[1], coord[2]);
528 HexoticNode[ aHexoticID ] = aHexoticNode;
529 nodeAssigne[ aHexoticID ] = 0;
537 case 6: // "Quadrilaterals"
538 case 7: { // "Hexahedra"
539 MESSAGE("Read " << nbElem << " " << token);
540 SMDS_MeshNode** node;
541 int nodeDim, *nodeID;
542 SMDS_MeshElement * aHexoticElement = 0;
544 node = new SMDS_MeshNode*[ nbRef ];
545 nodeID = new int[ nbRef ];
546 for ( int iElem = 0; iElem < nbElem; iElem++ ) {
547 if(theAlgo->computeCanceled())
551 for ( int iRef = 0; iRef < nbRef; iRef++ ) {
552 fileRes >> aHexoticNodeID; // read nbRef aHexoticNodeID
553 node[ iRef ] = HexoticNode[ aHexoticNodeID ];
554 nodeID[ iRef ] = aHexoticNodeID;
559 case 3: { // "Corners"
561 gp_Pnt HexoticPnt ( node[0]->X(), node[0]->Y(), node[0]->Z() );
562 for ( int i=0; i<nbElem; i++ ) {
563 aVertex = TopoDS::Vertex( tabCorner[i] );
564 gp_Pnt aPnt = BRep_Tool::Pnt( aVertex );
565 if ( aPnt.Distance( HexoticPnt ) < epsilon )
572 aHexoticElement = theMesh->AddEdge( node[0], node[1] );
574 if ( nodeAssigne[ nodeID[0] ] == 0 || nodeAssigne[ nodeID[0] ] == 2 )
577 shapeID = findEdge( node[iNode], theMesh, nbShapeEdge, tabEdge );
582 case 5: { // "Ridges"
585 case 6: { // "Quadrilaterals"
587 aHexoticElement = theMesh->AddFace( node[0], node[1], node[2], node[3] );
591 case 7: { // "Hexahedra"
593 if ( nbDomains > 1 ) {
594 hexoticShapeID = dummy - IdShapeRef;
595 if ( tabID[ hexoticShapeID ] == 0 ) {
596 aShape = findShape(node, aShape, tabShape, tabBox, nbShape);
597 shapeID = aShape.IsNull() ? holeID : theMesh->ShapeToIndex( aShape );
598 tabID[ hexoticShapeID ] = shapeID;
601 shapeID = tabID[ hexoticShapeID ];
602 if ( iElem == (nbElem - 1) ) {
603 int shapeAssociated = 0;
604 for ( int i=0; i<nbDomains; i++ ) {
606 shapeAssociated += 1;
608 if ( shapeAssociated != nbShape )
609 printWarning(nbShape, "domains", shapeAssociated);
615 if ( shapeID != holeID )
616 aHexoticElement = theMesh->AddVolume( node[0], node[3], node[2], node[1], node[4], node[7], node[6], node[5] );
621 if ( token != "Ridges" && ( shapeID > 0 || token == "Corners")) {
622 for ( int i=0; i<nbRef; i++ ) {
623 if ( nodeAssigne[ nodeID[i] ] == 0 ) {
624 if ( token == "Corners" ) theMesh->SetNodeOnVertex( node[0], aVertex );
625 else if ( token == "Edges" ) theMesh->SetNodeOnEdge( node[i], shapeID );
626 else if ( token == "Quadrilaterals" ) theMesh->SetNodeOnFace( node[i], shapeID );
627 else if ( token == "Hexahedra" ) theMesh->SetNodeInVolume( node[i], shapeID );
628 nodeAssigne[ nodeID[i] ] = nodeDim;
631 if ( token != "Corners" && aHexoticElement )
632 theMesh->SetMeshElementOnShape( aHexoticElement, shapeID );
641 MESSAGE("End of " << theFile << " file");
645 MESSAGE("Unknown Token: " << token);
651 // remove nodes in holes
654 SMESHDS_SubMesh* subMesh;
655 for ( int i = 1; i <= nbVertices; ++i )
656 if ( HexoticNode[i]->NbInverseElements() == 0 )
658 subMesh = HexoticNode[i]->getshapeId() > 0 ? theMesh->MeshElements(HexoticNode[i]->getshapeId() ) : 0;
659 theMesh->RemoveFreeNode( HexoticNode[i], subMesh, /*fromGroups=*/false );
667 delete [] nodeAssigne;
668 delete [] HexoticNode;
673 //=======================================================================
674 //function : readResult
675 //purpose : Read GMF file in case of a mesh w/o geometry
676 //=======================================================================
678 static bool readResult(std::string theFile,
679 HexoticPlugin_Hexotic* theAlgo,
680 SMESH_MesherHelper* theHelper)
682 SMESHDS_Mesh* theMesh = theHelper->GetMeshDS();
684 // ---------------------------------
685 // Read generated elements and nodes
686 // ---------------------------------
689 const int nbField = 9;
690 int nField, EndOfFile = 0, nbElem = 0, nbRef = 0;
691 int aHexoticNodeID = 0, shapeID;
692 int tabRef[nbField], *nodeAssigne;
693 bool tabDummy[nbField], hasDummy = false;
694 std::map <std::string,int> mapField;
695 SMDS_MeshNode** HexoticNode;
697 mapField["MeshVersionFormatted"] = 0; tabRef[0] = 0; tabDummy[0] = false;
698 mapField["Dimension"] = 1; tabRef[1] = 0; tabDummy[1] = false;
699 mapField["Vertices"] = 2; tabRef[2] = 3; tabDummy[2] = true;
700 mapField["Corners"] = 3; tabRef[3] = 1; tabDummy[3] = false;
701 mapField["Edges"] = 4; tabRef[4] = 2; tabDummy[4] = true;
702 mapField["Ridges"] = 5; tabRef[5] = 1; tabDummy[5] = false;
703 mapField["Quadrilaterals"] = 6; tabRef[6] = 4; tabDummy[6] = true;
704 mapField["Hexahedra"] = 7; tabRef[7] = 8; tabDummy[7] = true;
705 mapField["End"] = 8; tabRef[8] = 0; tabDummy[8] = false;
708 // theMesh->Clear(); -- this does not remove imported mesh
709 SMDS_ElemIteratorPtr eIt = theMesh->elementsIterator();
711 theMesh->RemoveFreeElement( eIt->next(), /*sm=*/0 );
712 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator();
713 while ( nIt->more() )
714 theMesh->RemoveFreeNode( nIt->next(), /*sm=*/0 );
717 int nbVertices = getNbShape(theFile, "Vertices");
718 HexoticNode = new SMDS_MeshNode*[ nbVertices + 1 ];
719 nodeAssigne = new int[ nbVertices + 1 ];
721 MESSAGE("Read " << theFile << " file");
722 std::ifstream fileRes(theFile.c_str());
730 if (mapField.count(token)) {
731 nField = mapField[token];
732 nbRef = tabRef[nField];
733 hasDummy = tabDummy[nField];
741 if ( nField < (mapField.size() - 1) && nField >= 0 )
745 case 0: { // "MeshVersionFormatted"
746 MESSAGE(token << " " << nbElem);
749 case 1: { // "Dimension"
750 MESSAGE("Mesh dimension " << nbElem << "D");
753 case 2: { // "Vertices"
754 MESSAGE("Read " << nbElem << " " << token);
757 SMDS_MeshNode * aHexoticNode;
759 for ( int iElem = 0; iElem < nbElem; iElem++ ) {
760 if(theAlgo->computeCanceled())
764 aHexoticID = iElem + 1;
765 for ( int iCoord = 0; iCoord < 3; iCoord++ )
766 fileRes >> coord[ iCoord ];
768 aHexoticNode = theMesh->AddNode(coord[0], coord[1], coord[2]);
769 HexoticNode[ aHexoticID ] = aHexoticNode;
770 nodeAssigne[ aHexoticID ] = 0;
777 case 6: // "Quadrilaterals"
778 case 7: { // "Hexahedra"
779 MESSAGE("Read " << nbElem << " " << token);
780 std::vector< SMDS_MeshNode* > node( nbRef );
781 std::vector< int > nodeID( nbRef );
783 for ( int iElem = 0; iElem < nbElem; iElem++ )
785 if(theAlgo->computeCanceled())
789 for ( int iRef = 0; iRef < nbRef; iRef++ )
791 fileRes >> aHexoticNodeID; // read nbRef aHexoticNodeID
792 node [ iRef ] = HexoticNode[ aHexoticNodeID ];
793 nodeID[ iRef ] = aHexoticNodeID;
800 theHelper->AddEdge( node[0], node[1] ); break;
801 case 6: // "Quadrilaterals"
802 theMesh->AddFace( node[0], node[1], node[2], node[3] ); break;
803 case 7: // "Hexahedra"
804 theHelper->AddVolume( node[0], node[3], node[2], node[1],
805 node[4], node[7], node[6], node[5] ); break;
809 for ( int iRef = 0; iRef < nbRef; iRef++ )
810 nodeAssigne[ nodeID[ iRef ]] = 1;
816 MESSAGE("End of " << theFile << " file");
820 MESSAGE("Unknown Token: " << token);
826 shapeID = theHelper->GetSubShapeID();
827 for ( int i = 0; i < nbVertices; ++i )
828 if ( !nodeAssigne[ i+1 ])
829 theMesh->SetNodeInVolume( HexoticNode[ i+1 ], shapeID );
831 delete [] HexoticNode;
832 delete [] nodeAssigne;
836 //=============================================================================
838 * Pass parameters to MG-Hexa
840 //=============================================================================
842 void HexoticPlugin_Hexotic::SetParameters(const HexoticPlugin_Hypothesis* hyp) {
844 MESSAGE("HexoticPlugin_Hexotic::SetParameters");
846 _hexesMinLevel = hyp->GetHexesMinLevel();
847 _hexesMaxLevel = hyp->GetHexesMaxLevel();
848 _hexesMinSize = hyp->GetMinSize();
849 _hexesMaxSize = hyp->GetMaxSize();
850 _hexoticIgnoreRidges = hyp->GetHexoticIgnoreRidges();
851 _hexoticInvalidElements = hyp->GetHexoticInvalidElements();
852 _hexoticSharpAngleThreshold = hyp->GetHexoticSharpAngleThreshold();
853 _hexoticNbProc = hyp->GetHexoticNbProc();
854 _hexoticWorkingDirectory = hyp->GetHexoticWorkingDirectory();
855 _hexoticVerbosity = hyp->GetHexoticVerbosity();
856 _hexoticMaxMemory = hyp->GetHexoticMaxMemory();
857 _hexoticSdMode = hyp->GetHexoticSdMode();
858 _sizeMaps = hyp->GetSizeMaps();
859 _nbLayers = hyp->GetNbLayers();
860 _firstLayerSize = hyp->GetFirstLayerSize();
861 _direction = hyp->GetDirection();
862 _growth = hyp->GetGrowth();
863 _facesWithLayers = hyp->GetFacesWithLayers();
864 _imprintedFaces = hyp->GetImprintedFaces();
868 cout << "WARNING : The MG-Hexa default parameters are taken into account" << std::endl;
869 cout << "=======" << std::endl;
870 _hexesMinLevel = hyp->GetDefaultHexesMinLevel();
871 _hexesMaxLevel = hyp->GetDefaultHexesMaxLevel();
872 _hexesMinSize = hyp->GetDefaultMinSize();
873 _hexesMaxSize = hyp->GetDefaultMaxSize();
874 _hexoticIgnoreRidges = hyp->GetDefaultHexoticIgnoreRidges();
875 _hexoticInvalidElements = hyp->GetDefaultHexoticInvalidElements();
876 _hexoticSharpAngleThreshold = hyp->GetDefaultHexoticSharpAngleThreshold();
877 _hexoticNbProc = hyp->GetDefaultHexoticNbProc();
878 _hexoticWorkingDirectory = hyp->GetDefaultHexoticWorkingDirectory();
879 _hexoticVerbosity = hyp->GetDefaultHexoticVerbosity();
880 _hexoticMaxMemory = hyp->GetDefaultHexoticMaxMemory();
881 _hexoticSdMode = hyp->GetDefaultHexoticSdMode();
882 _sizeMaps = hyp->GetDefaultHexoticSizeMaps();
883 _nbLayers = hyp->GetDefaultNbLayers();
884 _firstLayerSize = hyp->GetDefaultFirstLayerSize();
885 _direction = hyp->GetDefaultDirection();
886 _growth = hyp->GetDefaultGrowth();
887 _facesWithLayers = hyp->GetDefaultFacesWithLayers();
888 _imprintedFaces = hyp->GetDefaultImprintedFaces();
892 //=======================================================================
893 //function : getTmpDir
895 //=======================================================================
897 static TCollection_AsciiString getTmpDir()
899 TCollection_AsciiString aTmpDir;
901 char *Tmp_dir = getenv("SALOME_TMP_DIR");
903 if(Tmp_dir == NULL) {
904 Tmp_dir = getenv("TEMP");
906 Tmp_dir = getenv("TMP");
910 if(Tmp_dir != NULL) {
913 if(aTmpDir.Value(aTmpDir.Length()) != '\\') aTmpDir+='\\';
915 if(aTmpDir.Value(aTmpDir.Length()) != '/') aTmpDir+='/';
920 aTmpDir = TCollection_AsciiString("C:\\");
922 aTmpDir = TCollection_AsciiString("/tmp/");
928 //=======================================================================
929 //function : getSuffix
930 //purpose : Returns a suffix that will be unique for the current process
931 //=======================================================================
933 static TCollection_AsciiString getSuffix()
935 TCollection_AsciiString aSuffix = "";
938 aSuffix += getenv("USER");
940 aSuffix += getenv("USERNAME");
943 aSuffix += Kernel_Utils::GetHostname().c_str();
948 aSuffix += _getpid();
954 //================================================================================
956 * \brief Returns a command to run MG-Hexa mesher
958 //================================================================================
960 std::string HexoticPlugin_Hexotic::getHexoticCommand(const TCollection_AsciiString& Hexotic_In,
961 const TCollection_AsciiString& Hexotic_Out,
962 const TCollection_AsciiString& Hexotic_SizeMap_Prefix) const
965 cout << "MG-Hexa execution..." << std::endl;
966 cout << _name << " parameters :" << std::endl;
967 cout << " " << _name << " Verbosity = " << _hexoticVerbosity << std::endl;
968 cout << " " << _name << " Max Memory = " << _hexoticMaxMemory << std::endl;
969 cout << " " << _name << " Segments Min Level = " << _hexesMinLevel << std::endl;
970 cout << " " << _name << " Segments Max Level = " << _hexesMaxLevel << std::endl;
971 cout << " " << _name << " Segments Min Size = " << _hexesMinSize << std::endl;
972 cout << " " << _name << " Segments Max Size = " << _hexesMaxSize << std::endl;
973 cout << " " << "MG-Hexa can ignore ridges : " << (_hexoticIgnoreRidges ? "yes":"no") << std::endl;
974 cout << " " << "MG-Hexa authorize invalide elements : " << ( _hexoticInvalidElements ? "yes":"no") << std::endl;
975 cout << " " << _name << " Sharp angle threshold = " << _hexoticSharpAngleThreshold << " degrees" << std::endl;
976 cout << " " << _name << " Number of threads = " << _hexoticNbProc << std::endl;
977 cout << " " << _name << " Working directory = \"" << _hexoticWorkingDirectory << "\"" << std::endl;
978 cout << " " << _name << " Sub. Dom mode = " << _hexoticSdMode << std::endl;
979 cout << " " << _name << " Number of layers = " << _nbLayers << std::endl;
980 cout << " " << _name << " Size of the first layer = " << _firstLayerSize << std::endl;
981 cout << " " << _name << " Direction of the layers = " << ( _direction ? "Inward" : "Outward" ) << std::endl;
982 cout << " " << _name << " Growth = " << _growth << std::endl;
983 if (!_facesWithLayers.empty()) {
984 cout << " " << _name << " Faces with layers = ";
985 for (int i = 0; i < _facesWithLayers.size(); i++)
987 cout << _facesWithLayers.at(i);
988 if ((i + 1) != _facesWithLayers.size())
993 if (!_imprintedFaces.empty()) {
994 cout << " " << _name << " Imprinted faces = ";
995 for (int i = 0; i < _imprintedFaces.size(); i++)
997 cout << _imprintedFaces.at(i);
998 if ((i + 1) != _imprintedFaces.size())
1004 TCollection_AsciiString run_Hexotic( "mg-hexa.exe" );
1006 TCollection_AsciiString minl = " --min_level ", maxl = " --max_level ", angle = " --ridge_angle ";
1007 TCollection_AsciiString mins = " --min_size ", maxs = " --max_size ";
1008 TCollection_AsciiString in = " --in ", out = " --out ";
1009 TCollection_AsciiString sizeMap = " --read_sizemap ";
1010 TCollection_AsciiString ignoreRidges = " --compute_ridges no ", invalideElements = " --allow_invalid_elements yes ";
1011 TCollection_AsciiString subdom = " --components ";
1013 TCollection_AsciiString proc = " --max_number_of_threads ";
1015 TCollection_AsciiString verb = " --verbose ";
1016 TCollection_AsciiString maxmem = " --max_memory ";
1018 TCollection_AsciiString comNbLayers = " --number_of_boundary_layers ";
1019 TCollection_AsciiString comFirstLayerSize = " --height_of_the_first_layer ";
1020 TCollection_AsciiString comDirection = " --boundary_layers_subdomain_direction ";
1021 TCollection_AsciiString comGrowth = " --boundary_layers_geometric_progression ";
1022 TCollection_AsciiString comFacesWithLayers = " --boundary_layers_surface_ids ";
1023 TCollection_AsciiString comImptintedFaces = " --imprinted_surface_ids ";
1025 TCollection_AsciiString minLevel, maxLevel, minSize, maxSize, sharpAngle, mode, nbproc, verbosity, maxMemory,
1026 nbLayers, firstLayerSize, direction, growth, facesWithLayers, imprintedFaces;
1027 minLevel = _hexesMinLevel;
1028 maxLevel = _hexesMaxLevel;
1029 minSize = _hexesMinSize;
1030 maxSize = _hexesMaxSize;
1031 sharpAngle = _hexoticSharpAngleThreshold;
1032 // Mode translation for mg-tetra 1.1
1033 switch ( _hexoticSdMode )
1036 mode = "outside_skin_only";
1039 mode = "outside_components";
1045 mode = "all --manifold_geometry no";
1048 nbproc = _hexoticNbProc;
1049 verbosity = _hexoticVerbosity;
1050 maxMemory = _hexoticMaxMemory;
1051 nbLayers = _nbLayers;
1052 firstLayerSize = _firstLayerSize;
1053 direction = _direction ? "1" : "-1";
1055 for (int i = 0; i < _facesWithLayers.size(); i++)
1057 facesWithLayers += _facesWithLayers[i];
1058 if ((i + 1) != _facesWithLayers.size())
1059 facesWithLayers += ",";
1061 for (int i = 0; i < _imprintedFaces.size(); i++)
1063 imprintedFaces += _imprintedFaces[i];
1064 if ((i + 1) != _imprintedFaces.size())
1065 imprintedFaces += ",";
1068 if (_hexoticIgnoreRidges)
1069 run_Hexotic += ignoreRidges;
1071 if (_hexoticInvalidElements)
1072 run_Hexotic += invalideElements;
1074 if (_hexesMinSize > 0)
1075 run_Hexotic += mins + minSize;
1077 if (_hexesMaxSize > 0)
1078 run_Hexotic += maxs + maxSize;
1080 if (_hexesMinLevel > 0)
1081 run_Hexotic += minl + minLevel;
1083 if (_hexesMaxLevel > 0)
1084 run_Hexotic += maxl + maxLevel;
1086 if (_hexoticSharpAngleThreshold > 0)
1087 run_Hexotic += angle + sharpAngle;
1089 if (_sizeMaps.begin() != _sizeMaps.end())
1090 run_Hexotic += sizeMap + Hexotic_SizeMap_Prefix;
1093 run_Hexotic += comNbLayers + nbLayers;
1095 if (_firstLayerSize > 0)
1096 run_Hexotic += comFirstLayerSize + firstLayerSize;
1098 run_Hexotic += comDirection + direction;
1101 run_Hexotic += comGrowth + growth;
1103 if (!_facesWithLayers.empty())
1104 run_Hexotic += comFacesWithLayers + facesWithLayers;
1106 if (!_imprintedFaces.empty())
1107 run_Hexotic += comImptintedFaces + imprintedFaces;
1109 run_Hexotic += in + Hexotic_In + out + Hexotic_Out;
1110 run_Hexotic += subdom + mode;
1112 run_Hexotic += proc + nbproc;
1114 run_Hexotic += verb + verbosity;
1115 run_Hexotic += maxmem + maxMemory;
1117 return run_Hexotic.ToCString();
1120 // TODO : this is a duplication of some code found in BLSURFPlugin_BLSURF find a proper
1122 TopoDS_Shape HexoticPlugin_Hexotic::entryToShape(std::string entry)
1124 MESSAGE("HexoticPlugin_Hexotic::entryToShape "<<entry );
1125 GEOM::GEOM_Object_var aGeomObj;
1126 TopoDS_Shape S = TopoDS_Shape();
1127 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
1128 if (!aSObj->_is_nil()) {
1129 CORBA::Object_var obj = aSObj->GetObject();
1130 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
1131 aSObj->UnRegister();
1133 if ( !aGeomObj->_is_nil() )
1134 S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
1138 //================================================================================
1140 * \brief Produces a .mesh file with the size maps informations to give to Hexotic
1142 //================================================================================
1143 std::vector<std::string> HexoticPlugin_Hexotic::writeSizeMapFile( std::string sizeMapPrefix )
1145 HexoticPlugin_Hypothesis::THexoticSizeMaps::iterator it ;
1147 // The GMF driver will be used to write the size map file
1148 DriverGMF_Write aWriter;
1149 aWriter.SetSizeMapPrefix( sizeMapPrefix );
1151 std::vector<Control_Pnt> points;
1152 // Iterate on the size maps
1153 for (it=_sizeMaps.begin(); it!=_sizeMaps.end(); it++)
1155 // Step 1 : Get the GEOM object entry and the size
1156 // from the _sizeMaps infos
1157 std::string anEntry = it->first;
1158 double aLocalSize = it->second;
1159 TopoDS_Shape aShape = entryToShape( anEntry );
1161 // Step 2 : Create the points
1162 createControlPoints( aShape, aLocalSize, points );
1164 // Write the .mesh size map file
1165 aWriter.PerformSizeMap(points);
1166 return aWriter.GetSizeMapFiles();
1169 //================================================================================
1171 * \brief Fills a vector of points from which a size map input file can be written
1173 //================================================================================
1174 void HexoticPlugin_Hexotic::createControlPoints( const TopoDS_Shape& aShape,
1175 const double& theSize,
1176 std::vector<Control_Pnt>& thePoints )
1178 if ( aShape.ShapeType() == TopAbs_VERTEX )
1180 gp_Pnt aPnt = BRep_Tool::Pnt( TopoDS::Vertex(aShape) );
1181 Control_Pnt aControl_Pnt( aPnt, theSize );
1182 thePoints.push_back( aControl_Pnt );
1184 if ( aShape.ShapeType() == TopAbs_EDGE )
1186 createPointsSampleFromEdge( aShape, theSize, thePoints );
1188 else if ( aShape.ShapeType() == TopAbs_WIRE )
1191 for (Ex.Init(aShape,TopAbs_EDGE); Ex.More(); Ex.Next())
1193 createPointsSampleFromEdge( Ex.Current(), theSize, thePoints );
1196 else if ( aShape.ShapeType() == TopAbs_FACE )
1198 createPointsSampleFromFace( aShape, theSize, thePoints );
1200 else if ( aShape.ShapeType() == TopAbs_SOLID )
1202 createPointsSampleFromSolid( aShape, theSize, thePoints );
1204 else if ( aShape.ShapeType() == TopAbs_COMPOUND )
1206 TopoDS_Iterator it( aShape );
1207 for(; it.More(); it.Next())
1209 createControlPoints( it.Value(), theSize, thePoints );
1214 //================================================================================
1216 * \brief Fills a vector of points with point samples approximately
1217 * \brief spaced with a given size
1219 //================================================================================
1220 void HexoticPlugin_Hexotic::createPointsSampleFromEdge( const TopoDS_Shape& aShape,
1221 const double& theSize,
1222 std::vector<Control_Pnt>& thePoints )
1224 double step = theSize;
1226 Handle( Geom_Curve ) aCurve = BRep_Tool::Curve( TopoDS::Edge( aShape ), first, last );
1227 GeomAdaptor_Curve C ( aCurve );
1228 GCPnts_UniformAbscissa DiscretisationAlgo(C, step , first, last, Precision::Confusion());
1229 int nbPoints = DiscretisationAlgo.NbPoints();
1231 for ( int i = 1; i <= nbPoints; i++ )
1233 double param = DiscretisationAlgo.Parameter( i );
1235 aCurve->D0( param, aPnt );
1236 aPnt.SetSize(theSize);
1237 thePoints.push_back( aPnt );
1241 //================================================================================
1243 * \brief Fills a vector of points with point samples approximately
1244 * \brief spaced with a given size
1246 //================================================================================
1247 void HexoticPlugin_Hexotic::createPointsSampleFromFace( const TopoDS_Shape& aShape,
1248 const double& theSize,
1249 std::vector<Control_Pnt>& thePoints )
1251 BRepMesh_IncrementalMesh M(aShape, 0.01, Standard_True);
1252 TopLoc_Location aLocation;
1253 TopoDS_Face aFace = TopoDS::Face(aShape);
1255 // Triangulate the face
1256 Handle(Poly_Triangulation) aTri = BRep_Tool::Triangulation (aFace, aLocation);
1258 // Get the transformation associated to the face location
1259 gp_Trsf aTrsf = aLocation.Transformation();
1262 int nbTriangles = aTri->NbTriangles();
1263 Poly_Array1OfTriangle triangles(1,nbTriangles);
1264 triangles=aTri->Triangles();
1267 int nbNodes = aTri->NbNodes();
1268 TColgp_Array1OfPnt nodes(1,nbNodes);
1269 nodes = aTri->Nodes();
1271 // Iterate on triangles and subdivide them
1272 for(int i=1; i<=nbTriangles; i++)
1274 Poly_Triangle aTriangle = triangles.Value(i);
1275 gp_Pnt p1 = nodes.Value(aTriangle.Value(1));
1276 gp_Pnt p2 = nodes.Value(aTriangle.Value(2));
1277 gp_Pnt p3 = nodes.Value(aTriangle.Value(3));
1279 p1.Transform(aTrsf);
1280 p2.Transform(aTrsf);
1281 p3.Transform(aTrsf);
1283 subdivideTriangle(p1, p2, p3, theSize, thePoints);
1287 //================================================================================
1289 * \brief Fills a vector of points with point samples approximately
1290 * \brief spaced with a given size
1292 //================================================================================
1293 void HexoticPlugin_Hexotic::createPointsSampleFromSolid( const TopoDS_Shape& aShape,
1294 const double& theSize,
1295 std::vector<Control_Pnt>& thePoints )
1297 // Compute the bounding box
1298 double Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1300 BRepBndLib::Add(aShape, B);
1301 B.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1303 // Create the points
1304 double step = theSize;
1306 for ( double x=Xmin; x-Xmax<Precision::Confusion(); x=x+step )
1308 for ( double y=Ymin; y-Ymax<Precision::Confusion(); y=y+step )
1310 // Step1 : generate the Zmin -> Zmax line
1311 gp_Pnt startPnt(x, y, Zmin);
1312 gp_Pnt endPnt(x, y, Zmax);
1313 gp_Vec aVec(startPnt, endPnt);
1314 gp_Lin aLine(startPnt, aVec);
1315 double endParam = Zmax - Zmin;
1317 // Step2 : for each face of aShape:
1318 std::set<double> intersections;
1319 std::set<double>::iterator it = intersections.begin();
1322 for (Ex.Init(aShape,TopAbs_FACE); Ex.More(); Ex.Next())
1324 // check if there is an intersection
1325 IntCurvesFace_Intersector anIntersector(TopoDS::Face(Ex.Current()), Precision::Confusion());
1326 anIntersector.Perform(aLine, 0, endParam);
1328 // get the intersection's parameter and store it
1329 int nbPoints = anIntersector.NbPnt();
1330 for(int i = 0 ; i < nbPoints ; i++ )
1332 it = intersections.insert( it, anIntersector.WParameter(i+1) );
1335 // Step3 : go through the line chunk by chunk
1336 if ( intersections.begin() != intersections.end() )
1338 std::set<double>::iterator intersectionsIterator=intersections.begin();
1339 double first = *intersectionsIterator;
1340 intersectionsIterator++;
1341 bool innerPoints = true;
1342 for ( ; intersectionsIterator!=intersections.end() ; intersectionsIterator++ )
1344 double second = *intersectionsIterator;
1347 // If the last chunk was outside of the shape or this is the first chunk
1348 // add the points in the range [first, second] to the points vector
1349 double localStep = (second -first) / ceil( (second - first) / step );
1350 for ( double z = Zmin + first; z < Zmin + second; z = z + localStep )
1352 thePoints.push_back(Control_Pnt( x, y, z, theSize ));
1354 thePoints.push_back(Control_Pnt( x, y, Zmin + second, theSize ));
1357 innerPoints = !innerPoints;
1364 //================================================================================
1366 * \brief Subdivides a triangle until it reaches a certain size (recursive function)
1368 //================================================================================
1369 void HexoticPlugin_Hexotic::subdivideTriangle( const gp_Pnt& p1,
1372 const double& theSize,
1373 std::vector<Control_Pnt>& thePoints)
1375 // Size threshold to stop subdividing
1376 // This value ensures that two control points are distant no more than 2*theSize
1379 // The greater distance D of the mass center M to each Edge is 1/3 * Median
1380 // and Median < sqrt(3/4) * a where a is the greater side (by using Apollonius' thorem).
1381 // So D < 1/3 * sqrt(3/4) * a and if a < sqrt(3) * S then D < S/2
1382 // and the distance between two mass centers of two neighbouring triangles
1383 // sharing an edge is < 2 * 1/2 * S = S
1384 // If the traingles share a Vertex and no Edge the distance of the mass centers
1385 // to the Vertices is 2*D < S so the mass centers are distant of less than 2*S
1387 double threshold = sqrt( 3. ) * theSize;
1389 if ( (p1.Distance(p2) > threshold ||
1390 p2.Distance(p3) > threshold ||
1391 p3.Distance(p1) > threshold))
1393 std::vector<gp_Pnt> midPoints = computePointsForSplitting(p1, p2, p3);
1395 subdivideTriangle( midPoints[0], midPoints[1], midPoints[2], theSize, thePoints );
1396 subdivideTriangle( midPoints[0], p2, midPoints[1], theSize, thePoints );
1397 subdivideTriangle( midPoints[2], midPoints[1], p3, theSize, thePoints );
1398 subdivideTriangle( p1, midPoints[0], midPoints[2], theSize, thePoints );
1402 double x = (p1.X() + p2.X() + p3.X()) / 3 ;
1403 double y = (p1.Y() + p2.Y() + p3.Y()) / 3 ;
1404 double z = (p1.Z() + p2.Z() + p3.Z()) / 3 ;
1406 Control_Pnt massCenter( x ,y ,z, theSize );
1407 thePoints.push_back( massCenter );
1411 //================================================================================
1413 * \brief Returns the appropriate points for splitting a triangle
1414 * \brief the tangency points of the incircle are used in order to have mostly
1415 * \brief well-shaped sub-triangles
1417 //================================================================================
1418 std::vector<gp_Pnt> HexoticPlugin_Hexotic::computePointsForSplitting( const gp_Pnt& p1,
1422 std::vector<gp_Pnt> midPoints;
1423 //Change coordinates
1424 gp_Trsf Trsf_1; // Identity transformation
1425 gp_Ax3 reference_system(gp::Origin(), gp::DZ(), gp::DX()); // OXY
1428 gp_Vec Vaux(p1, p2);
1431 gp_Dir Dz = Dx.Crossed(Daux);
1432 gp_Ax3 current_system(p1, Dz, Dx);
1434 Trsf_1.SetTransformation( reference_system, current_system );
1436 gp_Pnt A = p1.Transformed(Trsf_1);
1437 gp_Pnt B = p2.Transformed(Trsf_1);
1438 gp_Pnt C = p3.Transformed(Trsf_1);
1440 double a = B.Distance(C) ;
1441 double b = A.Distance(C) ;
1442 double c = B.Distance(A) ;
1444 // Incenter coordinates
1445 // see http://mathworld.wolfram.com/Incenter.html
1446 double Xi = ( b*B.X() + c*C.X() ) / ( a + b + c );
1447 double Yi = ( b*B.Y() ) / ( a + b + c );
1448 gp_Pnt Center(Xi, Yi, 0);
1450 // Calculate the tangency points of the incircle
1451 gp_Pnt T1 = tangencyPoint( A, B, Center);
1452 gp_Pnt T2 = tangencyPoint( B, C, Center);
1453 gp_Pnt T3 = tangencyPoint( C, A, Center);
1455 gp_Pnt p1_2 = T1.Transformed(Trsf_1.Inverted());
1456 gp_Pnt p2_3 = T2.Transformed(Trsf_1.Inverted());
1457 gp_Pnt p3_1 = T3.Transformed(Trsf_1.Inverted());
1459 midPoints.push_back(p1_2);
1460 midPoints.push_back(p2_3);
1461 midPoints.push_back(p3_1);
1466 //================================================================================
1468 * \brief Computes the tangency points of the circle of center Center with
1469 * \brief the straight line (p1 p2)
1471 //================================================================================
1472 gp_Pnt HexoticPlugin_Hexotic::tangencyPoint(const gp_Pnt& p1,
1474 const gp_Pnt& Center)
1479 // The tangency point is the intersection of the straight line (p1 p2)
1480 // and the straight line (Center T) which is orthogonal to (p1 p2)
1481 if ( fabs(p1.X() - p2.X()) <= Precision::Confusion() )
1483 Xt=p1.X(); // T is on (p1 p2)
1484 Yt=Center.Y(); // (Center T) is orthogonal to (p1 p2)
1486 else if ( fabs(p1.Y() - p2.Y()) <= Precision::Confusion() )
1488 Yt=p1.Y(); // T is on (p1 p2)
1489 Xt=Center.X(); // (Center T) is orthogonal to (p1 p2)
1493 // First straight line coefficients (equation y=a*x+b)
1494 double a = (p2.Y() - p1.Y()) / (p2.X() - p1.X()) ;
1495 double b = p1.Y() - a*p1.X(); // p1 is on this straight line
1497 // Second straight line coefficients (equation y=c*x+d)
1498 double c = -1 / a; // The 2 lines are orthogonal
1499 double d = Center.Y() - c*Center.X(); // Center is on this straight line
1501 Xt = (d - b) / (a - c);
1505 return gp_Pnt( Xt, Yt, 0 );
1508 //=============================================================================
1510 * Here we are going to use the MG-Hexa mesher
1512 //=============================================================================
1514 bool HexoticPlugin_Hexotic::Compute(SMESH_Mesh& aMesh,
1515 const TopoDS_Shape& aShape)
1517 _compute_canceled = false;
1519 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1520 TCollection_AsciiString hexahedraMessage;
1522 if (_iShape == 0 && _nbShape == 0) {
1523 _nbShape = countShape( meshDS, TopAbs_SOLID ); // we count the number of shapes
1526 // to prevent from displaying error message after computing,
1527 // SetIsAlwaysComputed( true ) to empty sub-meshes
1528 vector< SMESH_subMesh* > subMeshesAlwaysComp;
1529 for ( int i = 0; i < _nbShape; ++i )
1530 if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( aShape ))
1532 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/true,
1533 /*complexShapeFirst=*/false);
1534 while ( smIt->more() )
1537 if ( !sm->IsMeshComputed() )
1539 sm->SetIsAlwaysComputed( true );
1540 subMeshesAlwaysComp.push_back( sm );
1547 if (_iShape == _nbShape ) {
1549 // create bounding box for each shape of the compound
1552 TopoDS_Shape *tabShape;
1555 tabShape = new TopoDS_Shape[_nbShape];
1556 tabBox = new double*[_nbShape];
1557 for (int i=0; i<_nbShape; i++)
1558 tabBox[i] = new double[6];
1559 double Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1561 TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
1562 for (; expBox.More(); expBox.Next()) {
1563 tabShape[iShape] = expBox.Current();
1564 Bnd_Box BoundingBox;
1565 BRepBndLib::Add(expBox.Current(), BoundingBox);
1566 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1567 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1568 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1569 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1573 SetParameters(_hypothesis);
1575 // TCollection_AsciiString aTmpDir = getTmpDir();
1576 TCollection_AsciiString aTmpDir = _hexoticWorkingDirectory.c_str();
1578 if ( aTmpDir.Value(aTmpDir.Length()) != '\\' ) aTmpDir += '\\';
1580 if ( aTmpDir.Value(aTmpDir.Length()) != '/' ) aTmpDir += '/';
1582 TCollection_AsciiString Hexotic_In(""), Hexotic_Out, Hexotic_SizeMap_Prefix;
1583 TCollection_AsciiString modeFile_In( "chmod 666 " ), modeFile_Out( "chmod 666 " ); TCollection_AsciiString aLogFileName = aTmpDir + "Hexotic"+getSuffix()+".log"; // log
1585 std::map <int,int> aSmdsToHexoticIdMap;
1586 std::map <int,const SMDS_MeshNode*> aHexoticIdToNodeMap;
1588 Hexotic_Out = aTmpDir + "Hexotic"+getSuffix()+"_Out.mesh";
1589 #ifdef WITH_BLSURFPLUGIN
1590 bool defaultInputFile = true;
1591 if (_blsurfHypo && !_blsurfHypo->GetQuadAllowed()) {
1592 Hexotic_In = TCollection_AsciiString(_blsurfHypo->GetGMFFile().c_str());
1593 if (Hexotic_In != "")
1594 defaultInputFile = false;
1596 if (defaultInputFile) {
1598 Hexotic_In = aTmpDir + "Hexotic"+getSuffix()+"_In.mesh";
1599 removeHexoticFiles(Hexotic_In, Hexotic_Out);
1600 splitQuads(aMesh); // quadrangles are no longer acceptable as input
1602 cout << "Creating MG-Hexa input mesh file : " << Hexotic_In << std::endl;
1603 aMesh.ExportGMF(Hexotic_In.ToCString(), meshDS, true);
1604 #ifdef WITH_BLSURFPLUGIN
1607 removeFile( Hexotic_Out );
1611 Hexotic_SizeMap_Prefix = aTmpDir + "Hexotic_SizeMap" + getSuffix();
1612 std::vector<std::string> sizeMapFiles = writeSizeMapFile( Hexotic_SizeMap_Prefix.ToCString() );
1614 std::string run_Hexotic = getHexoticCommand(Hexotic_In, Hexotic_Out, Hexotic_SizeMap_Prefix);
1615 run_Hexotic += std::string(" 1> ") + aLogFileName.ToCString(); // dump into file
1618 cout << "MG-Hexa command : " << run_Hexotic << std::endl;
1621 modeFile_In += Hexotic_In;
1622 system( modeFile_In.ToCString() );
1624 aSmdsToHexoticIdMap.clear();
1625 aHexoticIdToNodeMap.clear();
1627 MESSAGE("HexoticPlugin_Hexotic::Compute");
1629 int status = system( run_Hexotic.data() );
1635 std::ifstream fileRes( Hexotic_Out.ToCString() );
1637 modeFile_Out += Hexotic_Out;
1638 system( modeFile_Out.ToCString() );
1640 if ( ! fileRes.fail() ) {
1641 Ok = readResult( Hexotic_Out.ToCString(),
1643 meshDS, _nbShape, tabShape, tabBox );
1645 /*********************
1646 // TODO: Detect and remove elements in holes in case of sd mode = 4
1647 // Remove previous nodes and elements
1648 SMDS_ElemIteratorPtr itElement = meshDS->elementsIterator();
1649 SMDS_NodeIteratorPtr itNode = meshDS->nodesIterator();
1651 while ( itElement->more() )
1652 meshDS->RemoveElement( itElement->next() );
1653 while ( itNode->more() )
1654 meshDS->RemoveNode( itNode->next() );
1656 SMESH_ComputeErrorPtr myError = aMesh.GMFToMesh(Hexotic_Out.ToCString());
1659 hexahedraMessage = "success";
1661 removeFile(Hexotic_Out);
1662 removeFile(Hexotic_In);
1663 removeFile(aLogFileName);
1664 for( int i=0; i<sizeMapFiles.size(); i++)
1666 removeFile( TCollection_AsciiString( sizeMapFiles[i].c_str() ) );
1671 hexahedraMessage = "failed";
1675 hexahedraMessage = "failed";
1676 cout << "Problem with MG-Hexa output file " << Hexotic_Out.ToCString() << std::endl;
1679 SMESH_File logFile( aLogFileName.ToCString() );
1680 if ( !logFile.eof() )
1682 char msgLic[] = " Dlim ";
1683 const char* fileBeg = logFile.getPos(), *fileEnd = fileBeg + logFile.size();
1684 if ( std::search( fileBeg, fileEnd, msgLic, msgLic+strlen(msgLic)) != fileEnd )
1685 error("Licence problems.");
1688 if ( status > 0 && WEXITSTATUS(status) == 127 )
1689 error("mg-hexa.exe: command not found");
1692 if ( status == 0 && err == ENOENT ) {
1693 error("mg-hexa.exe: command not found");
1697 cout << "Hexahedra meshing " << hexahedraMessage << std::endl;
1700 // restore "always computed" flag of sub-meshes (0022127)
1701 for ( size_t iSM = 0; iSM < subMeshesAlwaysComp.size(); ++iSM )
1702 subMeshesAlwaysComp[ iSM ]->SetIsAlwaysComputed( false );
1705 for (int i=0; i<_nbShape; i++)
1706 delete [] tabBox[i];
1711 if(_compute_canceled)
1712 return error(SMESH_Comment("interruption initiated by user"));
1716 //=============================================================================
1718 * \brief Computes mesh without geometry
1719 * \param aMesh - the mesh
1720 * \param aHelper - helper that must be used for adding elements to \aaMesh
1721 * \retval bool - is a success
1723 * The method is called if ( !aMesh->HasShapeToMesh() )
1725 //=============================================================================
1727 bool HexoticPlugin_Hexotic::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelper)
1729 _compute_canceled = false;
1731 SMESH_ComputeErrorPtr myError = SMESH_ComputeError::New();
1734 TCollection_AsciiString hexahedraMessage;
1736 SetParameters(_hypothesis);
1738 TCollection_AsciiString aTmpDir = _hexoticWorkingDirectory.c_str();//getTmpDir();
1739 TCollection_AsciiString Hexotic_In, Hexotic_Out, Hexotic_SizeMap_Prefix;
1740 TCollection_AsciiString modeFile_In( "chmod 666 " ), modeFile_Out( "chmod 666 " );
1741 TCollection_AsciiString aLogFileName = aTmpDir + "Hexotic"+getSuffix()+".log"; // log
1743 std::map <int,int> aSmdsToHexoticIdMap;
1744 std::map <int,const SMDS_MeshNode*> aHexoticIdToNodeMap;
1746 Hexotic_In = aTmpDir + "Hexotic"+getSuffix()+"_In.mesh";
1747 Hexotic_Out = aTmpDir + "Hexotic"+getSuffix()+"_Out.mesh";
1748 Hexotic_SizeMap_Prefix = aTmpDir + "Hexotic_SizeMap";
1751 std::vector<std::string> sizeMapFiles = writeSizeMapFile( Hexotic_SizeMap_Prefix.ToCString() );
1753 std::string run_Hexotic = getHexoticCommand(Hexotic_In, Hexotic_Out, Hexotic_SizeMap_Prefix);
1754 run_Hexotic += std::string(" 1> ") + aLogFileName.ToCString(); // dump into file
1756 removeHexoticFiles(Hexotic_In, Hexotic_Out);
1758 splitQuads(aMesh); // quadrangles are no longer acceptable as input
1761 cout << "Creating MG-Hexa input mesh file : " << Hexotic_In << std::endl;
1762 aMesh.ExportGMF(Hexotic_In.ToCString(), aHelper->GetMeshDS());
1764 modeFile_In += Hexotic_In;
1765 system( modeFile_In.ToCString() );
1767 aSmdsToHexoticIdMap.clear();
1768 aHexoticIdToNodeMap.clear();
1770 MESSAGE("HexoticPlugin_Hexotic::Compute");
1773 cout << "MG-Hexa command : " << run_Hexotic << std::endl;
1774 system( run_Hexotic.data() );
1780 std::ifstream fileRes( Hexotic_Out.ToCString() );
1781 modeFile_Out += Hexotic_Out;
1782 system( modeFile_Out.ToCString() );
1783 if ( ! fileRes.fail() ) {
1784 Ok = readResult( Hexotic_Out.ToCString(),
1789 // Remove previous nodes and elements
1790 SMDS_ElemIteratorPtr itElement = aHelper->GetMeshDS()->elementsIterator();
1791 SMDS_NodeIteratorPtr itNode = aHelper->GetMeshDS()->nodesIterator();
1793 while ( itElement->more() )
1794 aHelper->GetMeshDS()->RemoveElement( itElement->next() );
1795 while ( itNode->more() )
1796 aHelper->GetMeshDS()->RemoveNode( itNode->next() );
1799 myError = aMesh.GMFToMesh(Hexotic_Out.ToCString());
1801 itElement = aHelper->GetMeshDS()->elementsIterator();
1802 itNode = aHelper->GetMeshDS()->nodesIterator();
1804 // Assign nodes and elements to the pseudo shape
1805 while ( itNode->more() )
1806 aHelper->GetMeshDS()->SetNodeInVolume(itNode->next(), 1);
1807 while ( itElement->more() )
1808 aHelper->GetMeshDS()->SetMeshElementOnShape(itElement->next(), 1);
1812 hexahedraMessage = "success";
1814 hexahedraMessage = "failed";
1818 myError->myName = COMPERR_EXCEPTION;
1820 hexahedraMessage = "failed";
1821 cout << "Problem with MG-Hexa output file " << Hexotic_Out << std::endl;
1823 SMESH_File logFile( aLogFileName.ToCString() );
1824 if ( !logFile.eof() )
1826 char msgLic[] = " Dlim ";
1827 const char* fileBeg = logFile.getPos(), *fileEnd = fileBeg + logFile.size();
1828 if ( std::search( fileBeg, fileEnd, msgLic, msgLic+strlen(msgLic)) != fileEnd )
1829 return error("Licence problems.");
1831 return error(SMESH_Comment("Problem with MG-Hexa output file ")<<Hexotic_Out);
1833 cout << "Hexahedra meshing " << hexahedraMessage << std::endl;
1836 if(_compute_canceled)
1837 return error(SMESH_Comment("interruption initiated by user"));
1838 removeFile(Hexotic_Out);
1839 removeFile(Hexotic_In);
1840 removeFile(aLogFileName);
1841 for( int i=0; i<sizeMapFiles.size(); i++)
1843 removeFile( TCollection_AsciiString(sizeMapFiles[i].c_str()) );
1847 return myError->IsOK();
1851 //=============================================================================
1855 //=============================================================================
1857 bool HexoticPlugin_Hexotic::Evaluate(SMESH_Mesh& aMesh,
1858 const TopoDS_Shape& aShape,
1859 MapShapeNbElems& aResMap)
1861 std::vector<int> aResVec(SMDSEntity_Last);
1862 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
1863 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
1864 aResMap.insert(std::make_pair(sm,aResVec));
1866 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1867 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Evaluation is not implemented",this));
1872 void HexoticPlugin_Hexotic::CancelCompute()
1874 _compute_canceled = true;
1877 TCollection_AsciiString aTmpDir = _hexoticWorkingDirectory.c_str(); //getTmpDir();
1878 TCollection_AsciiString Hexotic_In = aTmpDir + "Hexotic_In.mesh";
1879 TCollection_AsciiString cmd = TCollection_AsciiString("ps ux | grep ") + Hexotic_In;
1880 cmd += TCollection_AsciiString(" | grep -v grep | awk '{print $2}' | xargs kill -9 > /dev/null 2>&1");
1881 system( cmd.ToCString() );