1 // Copyright (C) 2005 CEA/DEN, EDF R&D, OPEN CASCADE, PRINCIPIA 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.
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
19 //=============================================================================
20 // File : GHS3DPlugin_GHS3D.cxx
22 // Author : Edward AGAPOV, modified by Lioka RAZAFINDRAZAKA (CEA) 09/02/2007
24 // Copyright : CEA 2003
26 //=============================================================================
29 #include "GHS3DPlugin_GHS3D.hxx"
30 #include "GHS3DPlugin_Hypothesis.hxx"
32 #include "SMESH_Gen.hxx"
33 #include "SMESH_Mesh.hxx"
34 #include "SMESH_Comment.hxx"
35 #include "SMESH_MesherHelper.hxx"
36 #include "SMESH_MeshEditor.hxx"
38 #include "SMDS_MeshElement.hxx"
39 #include "SMDS_MeshNode.hxx"
41 #include <BRepAdaptor_Surface.hxx>
42 #include <BRepBndLib.hxx>
43 #include <BRepClass3d_SolidClassifier.hxx>
44 #include <BRepTools.hxx>
45 #include <BRep_Tool.hxx>
46 #include <Bnd_Box.hxx>
47 #include <GeomAPI_ProjectPointOnSurf.hxx>
48 #include <OSD_File.hxx>
49 #include <Precision.hxx>
50 #include <Quantity_Parameter.hxx>
51 #include <Standard_ErrorHandler.hxx>
52 #include <Standard_Failure.hxx>
54 #include <TopExp_Explorer.hxx>
55 #include <TopTools_IndexedMapOfShape.hxx>
56 #include <TopTools_ListIteratorOfListOfShape.hxx>
58 //#include <BRepClass_FaceClassifier.hxx>
59 //#include <BRepGProp.hxx>
60 //#include <GProp_GProps.hxx>
62 #include "utilities.h"
65 #include <sys/sysinfo.h>
68 //#include <Standard_Stream.hxx>
71 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
92 //=============================================================================
96 //=============================================================================
98 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
99 : SMESH_3D_Algo(hypId, studyId, gen)
101 MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
103 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
104 _onlyUnaryInput = false; // Compute() will be called on a compound of solids
107 _compatibleHypothesis.push_back("GHS3D_Parameters");
110 //=============================================================================
114 //=============================================================================
116 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
118 MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
121 //=============================================================================
125 //=============================================================================
127 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh& aMesh,
128 const TopoDS_Shape& aShape,
129 Hypothesis_Status& aStatus )
131 aStatus = SMESH_Hypothesis::HYP_OK;
133 // there is only one compatible Hypothesis so far
136 const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
138 _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
140 _keepFiles = _hyp->GetKeepFiles();
145 //=======================================================================
146 //function : findShape
148 //=======================================================================
150 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
152 const TopoDS_Shape shape[],
155 TopAbs_State * state = 0)
158 int j, iShape, nbNode = 4;
160 for ( j=0; j<nbNode; j++ )
161 aPnt += gp_XYZ( aNode[j]->X(), aNode[j]->Y(), aNode[j]->Z() );
164 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
165 if (state) *state = SC.State();
166 if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
167 for (iShape = 0; iShape < nShape; iShape++) {
168 aShape = shape[iShape];
169 if ( !( aPnt.X() < box[iShape][0] || box[iShape][1] < aPnt.X() ||
170 aPnt.Y() < box[iShape][2] || box[iShape][3] < aPnt.Y() ||
171 aPnt.Z() < box[iShape][4] || box[iShape][5] < aPnt.Z()) ) {
172 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
173 if (state) *state = SC.State();
174 if (SC.State() == TopAbs_IN)
182 //=======================================================================
183 //function : readMapIntLine
185 //=======================================================================
187 static char* readMapIntLine(char* ptr, int tab[]) {
191 for ( int i=0; i<17; i++ ) {
192 intVal = strtol(ptr, &ptr, 10);
199 //=======================================================================
200 //function : readLine
202 //=======================================================================
204 #define GHS3DPlugin_BUFLENGTH 256
205 #define GHS3DPlugin_ReadLine(aPtr,aBuf,aFile,aLineNb) \
206 { aPtr = fgets( aBuf, GHS3DPlugin_BUFLENGTH - 2, aFile ); aLineNb++; DUMP(endl); }
208 //=======================================================================
209 //function : countShape
211 //=======================================================================
213 template < class Mesh, class Shape >
214 static int countShape( Mesh* mesh, Shape shape ) {
215 TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
217 for ( ; expShape.More(); expShape.Next() ) {
223 //=======================================================================
224 //function : writeFaces
226 //=======================================================================
228 static bool writeFaces (ofstream & theFile,
229 SMESHDS_Mesh * theMesh,
230 const map <int,int> & theSmdsToGhs3dIdMap)
234 // NB_ELEMS DUMMY_INT
235 // Loop from 1 to NB_ELEMS
236 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
238 int nbShape = countShape( theMesh, TopAbs_FACE );
240 int *tabID; tabID = new int[nbShape];
241 TopoDS_Shape *tabShape; tabShape = new TopoDS_Shape[nbShape];
243 SMESHDS_SubMesh* theSubMesh;
244 const SMDS_MeshElement* aFace;
245 const char* space = " ";
246 const int dummyint = 0;
247 map<int,int>::const_iterator itOnMap;
248 SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
249 int shapeID, nbNodes, aSmdsID;
252 cout << " " << theMesh->NbFaces() << " shapes of 2D dimension" << endl;
255 theFile << space << theMesh->NbFaces() << space << dummyint << endl;
257 TopExp_Explorer expface( theMesh->ShapeToMesh(), TopAbs_FACE );
258 for ( int i = 0; expface.More(); expface.Next(), i++ ) {
260 aShape = expface.Current();
261 shapeID = theMesh->ShapeToIndex( aShape );
263 for ( int j=0; j<=i; j++) {
264 if ( shapeID == tabID[j] ) {
271 tabShape[i] = aShape;
274 for ( int i =0; i < nbShape; i++ ) {
275 if ( tabID[i] != 0 ) {
276 aShape = tabShape[i];
278 theSubMesh = theMesh->MeshElements( aShape );
279 if ( !theSubMesh ) continue;
280 itOnSubMesh = theSubMesh->GetElements();
281 while ( itOnSubMesh->more() ) {
282 aFace = itOnSubMesh->next();
283 nbNodes = aFace->NbNodes();
285 theFile << space << nbNodes;
287 itOnSubFace = aFace->nodesIterator();
288 while ( itOnSubFace->more() ) {
290 aSmdsID = itOnSubFace->next()->GetID();
291 itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
292 ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
294 theFile << space << (*itOnMap).second;
297 // (NB_NODES + 1) times: DUMMY_INT
298 for ( int j=0; j<=nbNodes; j++)
299 theFile << space << dummyint;
312 //=======================================================================
313 //function : writeFaces
314 //purpose : Write Faces in case if generate 3D mesh w/o geometry
315 //=======================================================================
317 static bool writeFaces (ofstream & theFile,
318 SMESHDS_Mesh * theMesh,
319 vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
323 // NB_ELEMS DUMMY_INT
324 // Loop from 1 to NB_ELEMS
325 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
329 list< const SMDS_MeshElement* > faces;
330 list< const SMDS_MeshElement* >::iterator f;
331 map< const SMDS_MeshNode*,int >::iterator it;
332 SMDS_ElemIteratorPtr nodeIt;
333 const SMDS_MeshElement* elem;
336 const char* space = " ";
337 const int dummyint = 0;
339 //get all faces from mesh
340 SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
341 while ( eIt->more() ) {
342 const SMDS_MeshElement* elem = eIt->next();
345 faces.push_back( elem );
352 cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
354 // NB_ELEMS DUMMY_INT
355 theFile << space << nbFaces << space << dummyint << endl;
357 // Loop from 1 to NB_ELEMS
359 map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
361 for ( ; f != faces.end(); ++f )
365 nbNodes = elem->NbNodes();
366 theFile << space << nbNodes;
368 // NODE_NB_1 NODE_NB_2 ...
369 nodeIt = elem->nodesIterator();
370 while ( nodeIt->more() )
373 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
374 int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
375 it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
376 theFile << space << it->second;
379 // (NB_NODES + 1) times: DUMMY_INT
380 for ( int i=0; i<=nbNodes; i++)
381 theFile << space << dummyint;
386 // put nodes to theNodeByGhs3dId vector
387 theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
388 map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
389 for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
391 theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
397 //=======================================================================
398 //function : writePoints
400 //=======================================================================
402 static bool writePoints (ofstream & theFile,
403 SMESHDS_Mesh * theMesh,
404 map <int,int> & theSmdsToGhs3dIdMap,
405 map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap)
410 // Loop from 1 to NB_NODES
413 int nbNodes = theMesh->NbNodes();
417 const char* space = " ";
418 const int dummyint = 0;
421 SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
422 const SMDS_MeshNode* node;
425 theFile << space << nbNodes << endl;
427 cout << "The initial 2D mesh contains :" << endl;
428 cout << " " << nbNodes << " vertices" << endl;
430 // Loop from 1 to NB_NODES
435 theSmdsToGhs3dIdMap.insert( map <int,int>::value_type( node->GetID(), aGhs3dID ));
436 theGhs3dIdToNodeMap.insert (map <int,const SMDS_MeshNode*>::value_type( aGhs3dID, node ));
441 << space << node->X()
442 << space << node->Y()
443 << space << node->Z()
444 << space << dummyint;
452 //=======================================================================
453 //function : writePoints
455 //=======================================================================
457 static bool writePoints (ofstream & theFile,
458 SMESHDS_Mesh * theMesh,
459 const vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
464 // Loop from 1 to NB_NODES
467 //int nbNodes = theMesh->NbNodes();
468 int nbNodes = theNodeByGhs3dId.size();
472 const char* space = " ";
473 const int dummyint = 0;
475 const SMDS_MeshNode* node;
478 theFile << space << nbNodes << endl;
479 cout << nbNodes << " nodes" << endl;
481 // Loop from 1 to NB_NODES
483 vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
484 vector<const SMDS_MeshNode*>::const_iterator after = theNodeByGhs3dId.end();
485 for ( ; nodeIt != after; ++nodeIt )
491 << space << node->X()
492 << space << node->Y()
493 << space << node->Z()
494 << space << dummyint;
502 //=======================================================================
503 //function : findShapeID
504 //purpose : find the solid corresponding to GHS3D sub-domain following
505 // the technique proposed in GHS3D manual in chapter
506 // "B.4 Subdomain (sub-region) assignment"
507 //=======================================================================
509 static int findShapeID(SMESH_Mesh& mesh,
510 const SMDS_MeshNode* node1,
511 const SMDS_MeshNode* node2,
512 const SMDS_MeshNode* node3,
513 const bool toMeshHoles)
515 const int invalidID = 0;
516 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
518 // face th enodes belong to
519 const SMDS_MeshElement * face = meshDS->FindFace(node1,node2,node3);
523 // geom face the face assigned to
524 SMESH_MeshEditor editor(&mesh);
525 int geomFaceID = editor.FindShape( face );
528 TopoDS_Shape shape = meshDS->IndexToShape( geomFaceID );
529 if ( shape.IsNull() || shape.ShapeType() != TopAbs_FACE )
531 TopoDS_Face geomFace = TopoDS::Face( shape );
533 // solids bounded by geom face
534 TopTools_IndexedMapOfShape solids, shells;
535 TopTools_ListIteratorOfListOfShape ansIt = mesh.GetAncestors(geomFace);
536 for ( ; ansIt.More(); ansIt.Next() ) {
537 switch ( ansIt.Value().ShapeType() ) {
539 solids.Add( ansIt.Value() ); break;
541 shells.Add( ansIt.Value() ); break;
545 // analyse found solids
546 if ( solids.Extent() == 0 || shells.Extent() == 0)
549 const TopoDS_Solid& solid1 = TopoDS::Solid( solids(1) );
550 if ( solids.Extent() == 1 )
553 return meshDS->ShapeToIndex( solid1 );
555 // - are we at a hole boundary face?
556 if ( shells(1).IsSame( BRepTools::OuterShell( solid1 )) )
557 return meshDS->ShapeToIndex( solid1 ); // - no
560 // find orientation of geom face within the first solid
561 TopExp_Explorer fExp( solid1, TopAbs_FACE );
562 for ( ; fExp.More(); fExp.Next() )
563 if ( geomFace.IsSame( fExp.Current() )) {
564 geomFace = TopoDS::Face( fExp.Current() );
568 return invalidID; // face not found
570 // find UV of node1 on geomFace
571 SMESH_MesherHelper helper( mesh );
572 gp_XY uv = helper.GetNodeUV( geomFace, node1 );
574 // check that uv is correct
575 gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
576 double tol = BRep_Tool::Tolerance( geomFace );
577 BRepAdaptor_Surface surface( geomFace );
578 if ( node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol ) {
579 // project node1 onto geomFace to get right UV
580 GeomAPI_ProjectPointOnSurf projector( node1Pnt, surface.Surface().Surface() );
581 if ( !projector.IsDone() || projector.NbPoints() < 1 )
583 Quantity_Parameter U,V;
584 projector.LowerDistanceParameters(U,V);
587 // normale to face at node1
588 gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
589 gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
590 gp_Vec vec12( node1Pnt, node2Pnt );
591 gp_Vec vec13( node1Pnt, node3Pnt );
592 gp_Vec meshNormal = vec12 ^ vec13;
593 if ( meshNormal.SquareMagnitude() < DBL_MIN )
596 // normale to geomFace at UV
598 surface.D1( uv.X(), uv.Y(), node1Pnt, du, dv );
599 gp_Vec geomNormal = du ^ dv;
600 if ( geomNormal.SquareMagnitude() < DBL_MIN )
601 return findShapeID( mesh, node2, node3, node1, toMeshHoles );
602 if ( geomFace.Orientation() == TopAbs_REVERSED )
603 geomNormal.Reverse();
606 bool isReverse = ( meshNormal * geomNormal ) < 0;
608 return meshDS->ShapeToIndex( solid1 );
610 if ( solids.Extent() == 1 )
611 return HOLE_ID; // we are inside a hole
613 return meshDS->ShapeToIndex( solids(2) );
616 //=======================================================================
617 //function : readResultFile
619 //=======================================================================
621 static bool readResultFile(const int fileOpen,
623 TopoDS_Shape tabShape[],
626 map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
636 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
639 int nbElems, nbNodes, nbInputNodes;
640 int nodeId/*, triangleId*/;
642 int ID, shapeID, ghs3dShapeID;
645 nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
647 int *tab, *tabID, *nodeID, *nodeAssigne;
649 const SMDS_MeshNode **node;
652 //tabID = new int[nbShape];
654 coord = new double[3];
655 node = new const SMDS_MeshNode*[4];
658 SMDS_MeshNode * aNewNode;
659 map <int,const SMDS_MeshNode*>::iterator itOnNode;
660 SMDS_MeshElement* aTet;
665 // Read the file state
666 fileStat = fstat(fileOpen, &status);
667 length = status.st_size;
669 // Mapping the result file into memory
670 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
673 ptr = readMapIntLine(ptr, tab);
678 nbInputNodes = tab[2];
680 nodeAssigne = new int[ nbNodes+1 ];
683 aSolid = tabShape[0];
685 // Reading the nodeId
686 for (int i=0; i < 4*nbElems; i++)
687 nodeId = strtol(ptr, &ptr, 10);
689 // Reading the nodeCoor and update the nodeMap
690 for (int iNode=1; iNode <= nbNodes; iNode++) {
691 for (int iCoor=0; iCoor < 3; iCoor++)
692 coord[ iCoor ] = strtod(ptr, &ptr);
693 nodeAssigne[ iNode ] = 1;
694 if ( iNode > nbInputNodes ) {
695 nodeAssigne[ iNode ] = 0;
696 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
697 theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
701 // Reading the number of triangles which corresponds to the number of sub-domains
702 nbTriangle = strtol(ptr, &ptr, 10);
704 tabID = new int[nbTriangle];
705 for (int i=0; i < nbTriangle; i++) {
707 // find the solid corresponding to GHS3D sub-domain following
708 // the technique proposed in GHS3D manual in chapter
709 // "B.4 Subdomain (sub-region) assignment"
710 int nodeId1 = strtol(ptr, &ptr, 10);
711 int nodeId2 = strtol(ptr, &ptr, 10);
712 int nodeId3 = strtol(ptr, &ptr, 10);
713 if ( nbTriangle > 1 ) {
714 const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
715 const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
716 const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
719 tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles );
721 cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << endl;
723 } catch ( Standard_Failure ) {
730 if ( nbTriangle <= nbShape ) // no holes
731 toMeshHoles = true; // not avoid creating tetras in holes
733 // Associating the tetrahedrons to the shapes
734 shapeID = compoundID;
735 for (int iElem = 0; iElem < nbElems; iElem++) {
736 for (int iNode = 0; iNode < 4; iNode++) {
737 ID = strtol(tetraPtr, &tetraPtr, 10);
738 itOnNode = theGhs3dIdToNodeMap.find(ID);
739 node[ iNode ] = itOnNode->second;
740 nodeID[ iNode ] = ID;
742 // We always run GHS3D with "to mesh holes'==TRUE but we must not create
743 // tetras within holes depending on hypo option,
744 // so we first check if aTet is inside a hole and then create it
745 //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
746 if ( nbTriangle > 1 ) {
747 shapeID = HOLE_ID; // negative shapeID means not to create tetras if !toMeshHoles
748 ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
749 if ( tabID[ ghs3dShapeID ] == 0 ) {
751 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
752 if ( toMeshHoles || state == TopAbs_IN )
753 shapeID = theMeshDS->ShapeToIndex( aSolid );
754 tabID[ ghs3dShapeID ] = shapeID;
757 shapeID = tabID[ ghs3dShapeID ];
759 else if ( nbShape > 1 ) {
760 // Case where nbTriangle == 1 while nbShape == 2 encountered
761 // with compound of 2 boxes and "To mesh holes"==False,
762 // so there are no subdomains specified for each tetrahedron.
763 // Try to guess a solid by a node already bound to shape
765 for ( int i=0; i<4 && shapeID==0; i++ ) {
766 if ( nodeAssigne[ nodeID[i] ] == 1 &&
767 node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
768 node[i]->GetPosition()->GetShapeId() > 1 )
770 shapeID = node[i]->GetPosition()->GetShapeId();
774 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
775 shapeID = theMeshDS->ShapeToIndex( aSolid );
778 // set new nodes and tetrahedron onto the shape
779 for ( int i=0; i<4; i++ ) {
780 if ( nodeAssigne[ nodeID[i] ] == 0 ) {
781 if ( shapeID != HOLE_ID )
782 theMeshDS->SetNodeInVolume( node[i], shapeID );
783 nodeAssigne[ nodeID[i] ] = shapeID;
786 if ( toMeshHoles || shapeID != HOLE_ID ) {
787 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
788 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
791 shapeIDs.insert( shapeID );
794 // Remove nodes of tetras inside holes if !toMeshHoles
795 if ( !toMeshHoles ) {
796 itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
797 for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
798 ID = itOnNode->first;
799 if ( nodeAssigne[ ID ] == HOLE_ID )
800 theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
805 cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
806 munmap(mapPtr, length);
814 delete [] nodeAssigne;
817 if ( shapeIDs.size() != nbShape ) {
818 cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << endl;
819 for (int i=0; i<nbShape; i++) {
820 shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
821 if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
822 cout << " Solid #" << shapeID << " not found" << endl;
830 //=======================================================================
831 //function : readResultFile
833 //=======================================================================
835 static bool readResultFile(const int fileOpen,
836 SMESHDS_Mesh* theMeshDS,
838 vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
848 int nbElems, nbNodes, nbInputNodes;
849 int nodeId, triangleId;
855 const SMDS_MeshNode **node;
858 coord = new double[3];
859 node = new const SMDS_MeshNode*[4];
861 SMDS_MeshNode * aNewNode;
862 map <int,const SMDS_MeshNode*>::iterator IdNode;
863 SMDS_MeshElement* aTet;
865 // Read the file state
866 fileStat = fstat(fileOpen, &status);
867 length = status.st_size;
869 // Mapping the result file into memory
870 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
873 ptr = readMapIntLine(ptr, tab);
878 nbInputNodes = tab[2];
880 theNodeByGhs3dId.resize( nbNodes );
882 // Reading the nodeId
883 for (int i=0; i < 4*nbElems; i++)
884 nodeId = strtol(ptr, &ptr, 10);
886 // Reading the nodeCoor and update the nodeMap
887 shapeID = theMeshDS->ShapeToIndex( aSolid );
888 for (int iNode=0; iNode < nbNodes; iNode++) {
889 for (int iCoor=0; iCoor < 3; iCoor++)
890 coord[ iCoor ] = strtod(ptr, &ptr);
891 if ((iNode+1) > nbInputNodes) {
892 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
893 theMeshDS->SetNodeInVolume( aNewNode, shapeID );
894 theNodeByGhs3dId[ iNode ] = aNewNode;
898 // Reading the triangles
899 nbTriangle = strtol(ptr, &ptr, 10);
901 for (int i=0; i < 3*nbTriangle; i++)
902 triangleId = strtol(ptr, &ptr, 10);
906 // Associating the tetrahedrons to the shapes
907 for (int iElem = 0; iElem < nbElems; iElem++) {
908 for (int iNode = 0; iNode < 4; iNode++) {
909 ID = strtol(tetraPtr, &tetraPtr, 10);
910 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
912 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
913 shapeID = theMeshDS->ShapeToIndex( aSolid );
914 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
917 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
918 munmap(mapPtr, length);
928 //================================================================================
930 * \brief Look for a line containing a text in a file
931 * \retval bool - true if the line is found
933 //================================================================================
935 static bool findLineContaing(const TCollection_AsciiString& theText,
936 const TCollection_AsciiString& theFile,
937 TCollection_AsciiString & theFoundLine)
940 if ( FILE * aFile = fopen( theFile.ToCString(), "r" ))
943 char aBuffer[ GHS3DPlugin_BUFLENGTH ];
946 GHS3DPlugin_ReadLine( aPtr, aBuffer, aFile, aLineNb );
949 found = theFoundLine.Search( theText ) >= 0;
951 } while ( aPtr && !found );
958 //================================================================================
960 * \brief Provide human readable text by error code reported by ghs3d
962 //================================================================================
964 static TCollection_AsciiString translateError(const TCollection_AsciiString& errorLine)
966 int beg = errorLine.Location("ERR ", 1, errorLine.Length());
967 if ( !beg ) return errorLine;
969 TCollection_AsciiString errCodeStr = errorLine.SubString( beg + 4 , errorLine.Length());
970 if ( !errCodeStr.IsIntegerValue() )
972 Standard_Integer errCode = errCodeStr.IntegerValue();
974 int codeEnd = beg + 7;
975 while ( codeEnd < errorLine.Length() && isdigit( errorLine.Value(codeEnd+1) ))
977 TCollection_AsciiString error = errorLine.SubString( beg, codeEnd ) + ": ";
981 return "The surface mesh includes a face of type type other than edge, "
982 "triangle or quadrilateral. This face type is not supported.";
984 return error + "Not enough memory for the face table";
986 return error + "Not enough memory";
988 return error + "Not enough memory";
990 return error + "Face is ignored";
992 return error + errorLine.SubString( beg + 5, errorLine.Length()) +
993 ": End of file. Some data are missing in the file.";
995 return error + errorLine.SubString( beg + 5, errorLine.Length()) +
996 ": Read error on the file. There are wrong data in the file.";
998 return error + "the metric file is inadequate (dimension other than 3).";
1000 return error + "the metric file is inadequate (values not per vertices).";
1002 return error + "the metric file contains more than one field.";
1004 return error + "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1005 "value of number of mesh vertices in the \".noboite\" file.";
1007 return error + "Too many sub-domains";
1009 return error + "the number of vertices is negative or null.";
1011 return error + "the number of faces is negative or null.";
1013 return error + "A facehas a null vertex.";
1015 return error + "incompatible data";
1017 return error + "the number of vertices is negative or null.";
1019 return error + "the number of vertices is negative or null (in the \".mesh\" file).";
1021 return error + "the number of faces is negative or null.";
1023 return error + "A face appears more than once in the input surface mesh.";
1025 return error + "An edge appears more than once in the input surface mesh.";
1027 return error + "A face has a vertex negative or null.";
1029 return error + "NOT ENOUGH MEMORY";
1031 return error + "Not enough available memory.";
1033 return error + "Some initial points cannot be inserted. The surface mesh is probably very bad "
1034 "in terms of quality or the input list of points is wrong";
1036 return error + "Some vertices are too close to one another or coincident.";
1038 return error + "Some vertices are too close to one another or coincident.";
1040 return error + "A vertex cannot be inserted.";
1042 return error + "There are at least two points considered as coincident";
1044 return error + "Some vertices are too close to one another or coincident";
1046 return error + "The surface mesh regeneration step has failed";
1048 return error + "Constrained edge cannot be enforced";
1050 return error + "Constrained face cannot be enforced";
1052 return error + "Missing faces";
1054 return error + "No guess to start the definition of the connected component(s)";
1056 return error + "The surface mesh includes at least one hole. The domain is not well defined";
1058 return error + "Impossible to define a component";
1060 return error + "The surface edge intersects another surface edge";
1062 return error + "The surface edge intersects the surface face";
1064 return error + "One boundary point lies within a surface face";
1066 return error + "One surface edge intersects a surface face";
1068 return error + "One boundary point lies within a surface edge";
1070 return error + "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1071 "to too many swaps";
1073 return error + "Edge is unique (i.e., bounds a hole in the surface)";
1075 return error + "Presumably, the surface mesh is not compatible with the domain being processed";
1077 return error + "Too many components, too many sub-domain";
1079 return error + "The surface mesh includes at least one hole. "
1080 "Therefore there is no domain properly defined";
1082 return error + "Statistics.";
1084 return error + "Statistics.";
1086 return error + "Warning, it is dramatically tedious to enforce the boundary items";
1088 return error + "Not enough memory at this time, nevertheless, the program continues. "
1089 "The expected mesh will be correct but not really as large as required";
1091 return error + "see above error code, resulting quality may be poor.";
1093 return error + "Not enough memory at this time, nevertheless, the program continues (warning)";
1095 return error + "Unknown face type.";
1098 return error + errorLine.SubString( beg + 5, errorLine.Length()) +
1099 ": End of file. Some data are missing in the file.";
1101 return error + "A too small volume element is detected";
1103 return error + "There exists at least a null or negative volume element";
1105 return error + "There exist null or negative volume elements";
1107 return error + "A too small volume element is detected. A face is considered being degenerated";
1109 return error + "Some element is suspected to be very bad shaped or wrong";
1111 return error + "A too bad quality face is detected. This face is considered degenerated";
1113 return error + "A too bad quality face is detected. This face is degenerated";
1115 return error + "Presumably, the surface mesh is not compatible with the domain being processed";
1117 return error + "Abnormal error occured, contact hotline";
1119 return error + "Not enough memory for the face table";
1121 return error + "The algorithm cannot run further. "
1122 "The surface mesh is probably very bad in terms of quality.";
1124 return error + "Bad vertex number";
1129 //=============================================================================
1131 *Here we are going to use the GHS3D mesher
1133 //=============================================================================
1135 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1136 const TopoDS_Shape& theShape)
1139 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1141 _nbShape = countShape( meshDS, TopAbs_SOLID ); // we count the number of shapes
1143 // create bounding box for every shape inside the compound
1146 TopoDS_Shape* tabShape;
1148 tabShape = new TopoDS_Shape[_nbShape];
1149 tabBox = new double*[_nbShape];
1150 for (int i=0; i<_nbShape; i++)
1151 tabBox[i] = new double[6];
1152 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1154 TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
1155 for (; expBox.More(); expBox.Next()) {
1156 tabShape[iShape] = expBox.Current();
1157 Bnd_Box BoundingBox;
1158 BRepBndLib::Add(expBox.Current(), BoundingBox);
1159 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1160 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1161 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1162 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1166 // a unique working file name
1167 // to avoid access to the same files by eg different users
1168 TCollection_AsciiString aGenericName
1169 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1171 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1172 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1173 aFacesFileName = aGenericName + ".faces"; // in faces
1174 aPointsFileName = aGenericName + ".points"; // in points
1175 aResultFileName = aGenericName + ".noboite";// out points and volumes
1176 aBadResFileName = aGenericName + ".boite"; // out bad result
1177 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1178 aLogFileName = aGenericName + ".log"; // log
1180 // -----------------
1182 // -----------------
1184 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1185 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1189 aFacesFile->is_open() && aPointsFile->is_open();
1191 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1194 INFOS( "Can't write into " << aFacesFileName);
1195 return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1197 map <int,int> aSmdsToGhs3dIdMap;
1198 map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1200 Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap ) &&
1201 writeFaces ( aFacesFile, meshDS, aSmdsToGhs3dIdMap );
1204 aPointsFile.close();
1207 if ( !_keepFiles ) {
1208 OSD_File( aFacesFileName ).Remove();
1209 OSD_File( aPointsFileName ).Remove();
1211 return error(COMPERR_BAD_INPUT_MESH);
1213 OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1215 // -----------------
1217 // -----------------
1219 TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1220 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1221 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1224 cout << "Ghs3d execution..." << endl;
1225 cout << cmd << endl;
1227 system( cmd.ToCString() ); // run
1230 cout << "End of Ghs3d execution !" << endl;
1236 // Mapping the result file
1239 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1240 if ( fileOpen < 0 ) {
1242 cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << endl;
1247 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1248 Ok = readResultFile( fileOpen, theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1252 // ---------------------
1253 // remove working files
1254 // ---------------------
1259 OSD_File( aLogFileName ).Remove();
1261 else if ( OSD_File( aLogFileName ).Size() > 0 )
1263 INFOS( "GHS3D Error, see the " << aLogFileName.ToCString() << " file" );
1265 // get problem description from the log file
1266 SMESH_Comment comment;
1267 TCollection_AsciiString foundLine;
1268 if ( findLineContaing( "has expired",aLogFileName,foundLine) &&
1269 foundLine.Search("Licence") >= 0)
1271 foundLine.LeftAdjust();
1272 comment << foundLine;
1274 if ( findLineContaing( "%% ERROR",aLogFileName,foundLine) ||
1275 findLineContaing( " ERR ",aLogFileName,foundLine))
1277 foundLine.LeftAdjust();
1278 comment << translateError( foundLine );
1280 else if ( findLineContaing( "%% NO SAVING OPERATION",aLogFileName,foundLine))
1282 comment << "Too many elements generated for a trial version.\n";
1284 if ( comment.empty() )
1285 comment << "See " << aLogFileName << " for problem description";
1287 comment << "\nSee " << aLogFileName << " for more information";
1288 error(COMPERR_ALGO_FAILED, comment);
1292 // the log file is empty
1293 OSD_File( aLogFileName ).Remove();
1294 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1295 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1298 if ( !_keepFiles ) {
1299 OSD_File( aFacesFileName ).Remove();
1300 OSD_File( aPointsFileName ).Remove();
1301 OSD_File( aResultFileName ).Remove();
1302 OSD_File( aBadResFileName ).Remove();
1303 OSD_File( aBbResFileName ).Remove();
1305 cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1308 cout << "treated !" << endl;
1311 _nbShape = 0; // re-initializing _nbShape for the next Compute() method call
1318 //=============================================================================
1320 *Here we are going to use the GHS3D mesher w/o geometry
1322 //=============================================================================
1323 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1324 SMESH_MesherHelper* aHelper)
1326 MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1328 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1329 TopoDS_Shape theShape = aHelper->GetSubShape();
1331 // a unique working file name
1332 // to avoid access to the same files by eg different users
1333 TCollection_AsciiString aGenericName
1334 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1336 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1337 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1338 aFacesFileName = aGenericName + ".faces"; // in faces
1339 aPointsFileName = aGenericName + ".points"; // in points
1340 aResultFileName = aGenericName + ".noboite";// out points and volumes
1341 aBadResFileName = aGenericName + ".boite"; // out bad result
1342 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1343 aLogFileName = aGenericName + ".log"; // log
1345 // -----------------
1347 // -----------------
1349 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1350 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1353 aFacesFile->is_open() && aPointsFile->is_open();
1355 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1359 return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1361 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1363 Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1364 writePoints( aPointsFile, meshDS, aNodeByGhs3dId));
1367 aPointsFile.close();
1370 if ( !_keepFiles ) {
1371 OSD_File( aFacesFileName ).Remove();
1372 OSD_File( aPointsFileName ).Remove();
1374 return error(COMPERR_BAD_INPUT_MESH);
1376 OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1378 // -----------------
1380 // -----------------
1382 TCollection_AsciiString cmd =
1383 (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1384 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1385 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1387 system( cmd.ToCString() ); // run
1393 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1394 if ( fileOpen < 0 ) {
1396 cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
1401 Ok = readResultFile( fileOpen, meshDS, theShape ,aNodeByGhs3dId );
1404 // ---------------------
1405 // remove working files
1406 // ---------------------
1411 OSD_File( aLogFileName ).Remove();
1413 else if ( OSD_File( aLogFileName ).Size() > 0 )
1415 INFOS( "GHS3D Error, see the " << aLogFileName.ToCString() << " file" );
1417 // get problem description from the log file
1418 SMESH_Comment comment;
1419 TCollection_AsciiString foundLine;
1420 if ( findLineContaing( "has expired",aLogFileName,foundLine) &&
1421 foundLine.Search("Licence") >= 0)
1423 foundLine.LeftAdjust();
1424 comment << foundLine;
1426 if ( findLineContaing( "%% ERROR",aLogFileName,foundLine) ||
1427 findLineContaing( " ERR ",aLogFileName,foundLine))
1429 foundLine.LeftAdjust();
1430 comment << translateError( foundLine );
1432 else if ( findLineContaing( "%% NO SAVING OPERATION",aLogFileName,foundLine))
1434 comment << "Too many elements generated for a trial version.\n";
1436 if ( comment.empty() )
1437 comment << "See " << aLogFileName << " for problem description";
1439 comment << "\nSee " << aLogFileName << " for more information";
1440 error(COMPERR_ALGO_FAILED, comment);
1443 // the log file is empty
1444 OSD_File( aLogFileName ).Remove();
1445 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1446 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1451 OSD_File( aFacesFileName ).Remove();
1452 OSD_File( aPointsFileName ).Remove();
1453 OSD_File( aResultFileName ).Remove();
1454 OSD_File( aBadResFileName ).Remove();
1455 OSD_File( aBbResFileName ).Remove();