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 <BRep_Tool.hxx>
45 #include <Bnd_Box.hxx>
46 #include <GeomAPI_ProjectPointOnSurf.hxx>
47 #include <OSD_File.hxx>
48 #include <Precision.hxx>
49 #include <Quantity_Parameter.hxx>
50 #include <Standard_ErrorHandler.hxx>
51 #include <Standard_Failure.hxx>
53 #include <TopExp_Explorer.hxx>
54 #include <TopTools_IndexedMapOfShape.hxx>
55 #include <TopTools_ListIteratorOfListOfShape.hxx>
57 //#include <BRepClass_FaceClassifier.hxx>
58 //#include <BRepGProp.hxx>
59 //#include <GProp_GProps.hxx>
61 #include "utilities.h"
64 #include <sys/sysinfo.h>
67 //#include <Standard_Stream.hxx>
70 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
89 //=============================================================================
93 //=============================================================================
95 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
96 : SMESH_3D_Algo(hypId, studyId, gen)
98 MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
100 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
101 _onlyUnaryInput = false; // Compute() will be called on a compound of solids
104 _compatibleHypothesis.push_back("GHS3D_Parameters");
107 //=============================================================================
111 //=============================================================================
113 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
115 MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
118 //=============================================================================
122 //=============================================================================
124 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh& aMesh,
125 const TopoDS_Shape& aShape,
126 Hypothesis_Status& aStatus )
128 aStatus = SMESH_Hypothesis::HYP_OK;
130 // there is only one compatible Hypothesis so far
133 const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
135 _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
137 _keepFiles = _hyp->GetKeepFiles();
142 //=======================================================================
143 //function : findShape
145 //=======================================================================
147 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
149 const TopoDS_Shape shape[],
152 TopAbs_State * state = 0)
155 int j, iShape, nbNode = 4;
157 for ( j=0; j<nbNode; j++ )
158 aPnt += gp_XYZ( aNode[j]->X(), aNode[j]->Y(), aNode[j]->Z() );
161 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
162 if (state) *state = SC.State();
163 if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
164 for (iShape = 0; iShape < nShape; iShape++) {
165 aShape = shape[iShape];
166 if ( !( aPnt.X() < box[iShape][0] || box[iShape][1] < aPnt.X() ||
167 aPnt.Y() < box[iShape][2] || box[iShape][3] < aPnt.Y() ||
168 aPnt.Z() < box[iShape][4] || box[iShape][5] < aPnt.Z()) ) {
169 BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
170 if (state) *state = SC.State();
171 if (SC.State() == TopAbs_IN)
179 //=======================================================================
180 //function : readMapIntLine
182 //=======================================================================
184 static char* readMapIntLine(char* ptr, int tab[]) {
188 for ( int i=0; i<17; i++ ) {
189 intVal = strtol(ptr, &ptr, 10);
196 //=======================================================================
197 //function : readLine
199 //=======================================================================
201 #define GHS3DPlugin_BUFLENGTH 256
202 #define GHS3DPlugin_ReadLine(aPtr,aBuf,aFile,aLineNb) \
203 { aPtr = fgets( aBuf, GHS3DPlugin_BUFLENGTH - 2, aFile ); aLineNb++; DUMP(endl); }
205 //=======================================================================
206 //function : countShape
208 //=======================================================================
210 template < class Mesh, class Shape >
211 static int countShape( Mesh* mesh, Shape shape ) {
212 TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
214 for ( ; expShape.More(); expShape.Next() ) {
220 //=======================================================================
221 //function : writeFaces
223 //=======================================================================
225 static bool writeFaces (ofstream & theFile,
226 SMESHDS_Mesh * theMesh,
227 const map <int,int> & theSmdsToGhs3dIdMap)
231 // NB_ELEMS DUMMY_INT
232 // Loop from 1 to NB_ELEMS
233 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
235 int nbShape = countShape( theMesh, TopAbs_FACE );
237 int *tabID; tabID = new int[nbShape];
238 TopoDS_Shape *tabShape; tabShape = new TopoDS_Shape[nbShape];
240 SMESHDS_SubMesh* theSubMesh;
241 const SMDS_MeshElement* aFace;
242 const char* space = " ";
243 const int dummyint = 0;
244 map<int,int>::const_iterator itOnMap;
245 SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
246 int shapeID, nbNodes, aSmdsID;
249 cout << " " << theMesh->NbFaces() << " shapes of 2D dimension" << endl;
252 theFile << space << theMesh->NbFaces() << space << dummyint << endl;
254 TopExp_Explorer expface( theMesh->ShapeToMesh(), TopAbs_FACE );
255 for ( int i = 0; expface.More(); expface.Next(), i++ ) {
257 aShape = expface.Current();
258 shapeID = theMesh->ShapeToIndex( aShape );
260 for ( int j=0; j<=i; j++) {
261 if ( shapeID == tabID[j] ) {
268 tabShape[i] = aShape;
271 for ( int i =0; i < nbShape; i++ ) {
272 if ( tabID[i] != 0 ) {
273 aShape = tabShape[i];
275 theSubMesh = theMesh->MeshElements( aShape );
276 if ( !theSubMesh ) continue;
277 itOnSubMesh = theSubMesh->GetElements();
278 while ( itOnSubMesh->more() ) {
279 aFace = itOnSubMesh->next();
280 nbNodes = aFace->NbNodes();
282 theFile << space << nbNodes;
284 itOnSubFace = aFace->nodesIterator();
285 while ( itOnSubFace->more() ) {
287 aSmdsID = itOnSubFace->next()->GetID();
288 itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
289 ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
291 theFile << space << (*itOnMap).second;
294 // (NB_NODES + 1) times: DUMMY_INT
295 for ( int j=0; j<=nbNodes; j++)
296 theFile << space << dummyint;
309 //=======================================================================
310 //function : writeFaces
311 //purpose : Write Faces in case if generate 3D mesh w/o geometry
312 //=======================================================================
314 static bool writeFaces (ofstream & theFile,
315 SMESHDS_Mesh * theMesh,
316 vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
320 // NB_ELEMS DUMMY_INT
321 // Loop from 1 to NB_ELEMS
322 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
326 list< const SMDS_MeshElement* > faces;
327 list< const SMDS_MeshElement* >::iterator f;
328 map< const SMDS_MeshNode*,int >::iterator it;
329 SMDS_ElemIteratorPtr nodeIt;
330 const SMDS_MeshElement* elem;
333 const char* space = " ";
334 const int dummyint = 0;
336 //get all faces from mesh
337 SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
338 while ( eIt->more() ) {
339 const SMDS_MeshElement* elem = eIt->next();
342 faces.push_back( elem );
349 cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
351 // NB_ELEMS DUMMY_INT
352 theFile << space << nbFaces << space << dummyint << endl;
354 // Loop from 1 to NB_ELEMS
356 map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
358 for ( ; f != faces.end(); ++f )
362 nbNodes = elem->NbNodes();
363 theFile << space << nbNodes;
365 // NODE_NB_1 NODE_NB_2 ...
366 nodeIt = elem->nodesIterator();
367 while ( nodeIt->more() )
370 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
371 int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
372 it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
373 theFile << space << it->second;
376 // (NB_NODES + 1) times: DUMMY_INT
377 for ( int i=0; i<=nbNodes; i++)
378 theFile << space << dummyint;
383 // put nodes to theNodeByGhs3dId vector
384 theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
385 map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
386 for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
388 theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
394 //=======================================================================
395 //function : writePoints
397 //=======================================================================
399 static bool writePoints (ofstream & theFile,
400 SMESHDS_Mesh * theMesh,
401 map <int,int> & theSmdsToGhs3dIdMap,
402 map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap)
407 // Loop from 1 to NB_NODES
410 int nbNodes = theMesh->NbNodes();
414 const char* space = " ";
415 const int dummyint = 0;
418 SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
419 const SMDS_MeshNode* node;
422 theFile << space << nbNodes << endl;
424 cout << "The initial 2D mesh contains :" << endl;
425 cout << " " << nbNodes << " vertices" << endl;
427 // Loop from 1 to NB_NODES
432 theSmdsToGhs3dIdMap.insert( map <int,int>::value_type( node->GetID(), aGhs3dID ));
433 theGhs3dIdToNodeMap.insert (map <int,const SMDS_MeshNode*>::value_type( aGhs3dID, node ));
438 << space << node->X()
439 << space << node->Y()
440 << space << node->Z()
441 << space << dummyint;
449 //=======================================================================
450 //function : writePoints
452 //=======================================================================
454 static bool writePoints (ofstream & theFile,
455 SMESHDS_Mesh * theMesh,
456 const vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
461 // Loop from 1 to NB_NODES
464 //int nbNodes = theMesh->NbNodes();
465 int nbNodes = theNodeByGhs3dId.size();
469 const char* space = " ";
470 const int dummyint = 0;
472 const SMDS_MeshNode* node;
475 theFile << space << nbNodes << endl;
476 cout << nbNodes << " nodes" << endl;
478 // Loop from 1 to NB_NODES
480 vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
481 vector<const SMDS_MeshNode*>::const_iterator after = theNodeByGhs3dId.end();
482 for ( ; nodeIt != after; ++nodeIt )
488 << space << node->X()
489 << space << node->Y()
490 << space << node->Z()
491 << space << dummyint;
499 //=======================================================================
500 //function : findShapeID
501 //purpose : find the solid corresponding to GHS3D sub-domain following
502 // the technique proposed in GHS3D manual in chapter
503 // "B.4 Subdomain (sub-region) assignment"
504 //=======================================================================
506 static int findShapeID(SMESH_Mesh& mesh,
507 const SMDS_MeshNode* node1,
508 const SMDS_MeshNode* node2,
509 const SMDS_MeshNode* node3)
511 const int invalidID = 0;
512 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
514 // face th enodes belong to
515 const SMDS_MeshElement * face = meshDS->FindFace(node1,node2,node3);
519 // geom face the face assigned to
520 SMESH_MeshEditor editor(&mesh);
521 int geomFaceID = editor.FindShape( face );
524 TopoDS_Shape shape = meshDS->IndexToShape( geomFaceID );
525 if ( shape.IsNull() || shape.ShapeType() != TopAbs_FACE )
527 TopoDS_Face geomFace = TopoDS::Face( shape );
529 // solids bounded by geom face
530 TopTools_IndexedMapOfShape solids;
531 TopTools_ListIteratorOfListOfShape ansIt = mesh.GetAncestors(geomFace);
532 for ( ; ansIt.More(); ansIt.Next() ) {
533 if ( ansIt.Value().ShapeType() == TopAbs_SOLID )
534 solids.Add( ansIt.Value() );
536 // analyse found solids
537 if ( solids.Extent() == 0 )
539 if ( solids.Extent() == 1 )
540 return meshDS->ShapeToIndex( solids(1) );
542 // find orientation of geom face within the first solid
543 TopoDS_Shape solid1 = solids(1);
544 TopExp_Explorer fExp( solid1, TopAbs_FACE );
545 for ( ; fExp.More(); fExp.Next() )
546 if ( geomFace.IsSame( fExp.Current() )) {
547 geomFace = TopoDS::Face( fExp.Current() );
551 return invalidID; // face not found
553 // find UV of node1 on geomFace
554 SMESH_MesherHelper helper( mesh );
555 gp_XY uv = helper.GetNodeUV( geomFace, node1 );
557 // check that uv is correct
558 gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
559 double tol = BRep_Tool::Tolerance( geomFace );
560 BRepAdaptor_Surface surface( geomFace );
561 if ( node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol ) {
562 // project node1 onto geomFace to get right UV
563 GeomAPI_ProjectPointOnSurf projector( node1Pnt, surface.Surface().Surface() );
564 if ( !projector.IsDone() || projector.NbPoints() < 1 )
566 Quantity_Parameter U,V;
567 projector.LowerDistanceParameters(U,V);
570 // normale to face at node1
571 gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
572 gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
573 gp_Vec vec12( node1Pnt, node2Pnt );
574 gp_Vec vec13( node1Pnt, node3Pnt );
575 gp_Vec meshNormal = vec12 ^ vec13;
576 if ( meshNormal.SquareMagnitude() < DBL_MIN )
579 // normale to geomFace at UV
581 surface.D1( uv.X(), uv.Y(), node1Pnt, du, dv );
582 gp_Vec geomNormal = du ^ dv;
583 if ( geomNormal.SquareMagnitude() < DBL_MIN )
584 return findShapeID( mesh, node2, node3, node1 );
585 if ( geomFace.Orientation() == TopAbs_REVERSED )
586 geomNormal.Reverse();
589 bool isReverse = ( meshNormal * geomNormal ) < 0;
591 return meshDS->ShapeToIndex( solids(2) );
593 return meshDS->ShapeToIndex( solid1 );
596 //=======================================================================
597 //function : readResultFile
599 //=======================================================================
601 static bool readResultFile(const int fileOpen,
603 TopoDS_Shape tabShape[],
606 map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
616 SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
619 int nbElems, nbNodes, nbInputNodes;
620 int nodeId/*, triangleId*/;
622 int ID, shapeID, ghs3dShapeID;
625 nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
627 int *tab, *tabID, *nodeID, *nodeAssigne;
629 const SMDS_MeshNode **node;
632 //tabID = new int[nbShape];
634 coord = new double[3];
635 node = new const SMDS_MeshNode*[4];
638 SMDS_MeshNode * aNewNode;
639 map <int,const SMDS_MeshNode*>::iterator itOnNode;
640 SMDS_MeshElement* aTet;
645 // Read the file state
646 fileStat = fstat(fileOpen, &status);
647 length = status.st_size;
649 // Mapping the result file into memory
650 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
653 ptr = readMapIntLine(ptr, tab);
658 nbInputNodes = tab[2];
660 nodeAssigne = new int[ nbNodes+1 ];
663 aSolid = tabShape[0];
665 // Reading the nodeId
666 for (int i=0; i < 4*nbElems; i++)
667 nodeId = strtol(ptr, &ptr, 10);
669 // Reading the nodeCoor and update the nodeMap
670 for (int iNode=1; iNode <= nbNodes; iNode++) {
671 for (int iCoor=0; iCoor < 3; iCoor++)
672 coord[ iCoor ] = strtod(ptr, &ptr);
673 nodeAssigne[ iNode ] = 1;
674 if ( iNode > nbInputNodes ) {
675 nodeAssigne[ iNode ] = 0;
676 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
677 theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
681 // Reading the number of triangles which corresponds to the number of sub-domains
682 nbTriangle = strtol(ptr, &ptr, 10);
684 tabID = new int[nbTriangle];
685 for (int i=0; i < nbTriangle; i++) {
687 // find the solid corresponding to GHS3D sub-domain following
688 // the technique proposed in GHS3D manual in chapter
689 // "B.4 Subdomain (sub-region) assignment"
690 int nodeId1 = strtol(ptr, &ptr, 10);
691 int nodeId2 = strtol(ptr, &ptr, 10);
692 int nodeId3 = strtol(ptr, &ptr, 10);
693 const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
694 const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
695 const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
698 tabID[i] = findShapeID( theMesh, n1, n2, n3 );
700 cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << endl;
702 } catch ( Standard_Failure ) {
708 if ( nbTriangle <= nbShape ) // no holes
709 toMeshHoles = true; // not avoid creating tetras in holes
711 // Associating the tetrahedrons to the shapes
712 shapeID = compoundID;
713 for (int iElem = 0; iElem < nbElems; iElem++) {
714 for (int iNode = 0; iNode < 4; iNode++) {
715 ID = strtol(tetraPtr, &tetraPtr, 10);
716 itOnNode = theGhs3dIdToNodeMap.find(ID);
717 node[ iNode ] = itOnNode->second;
718 nodeID[ iNode ] = ID;
720 // We always run GHS3D with "to mesh holes'==TRUE but we must not create
721 // tetras within holes depending on hypo option,
722 // so we first check if aTet is inside a hole and then create it
723 //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
724 if ( nbTriangle > 1 ) {
725 shapeID = -1; // negative shapeID means not to create tetras if !toMeshHoles
726 ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
727 if ( tabID[ ghs3dShapeID ] == 0 ) {
729 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
730 if ( toMeshHoles || state == TopAbs_IN )
731 shapeID = theMeshDS->ShapeToIndex( aSolid );
732 tabID[ ghs3dShapeID ] = shapeID;
735 shapeID = tabID[ ghs3dShapeID ];
737 else if ( nbShape > 1 ) {
738 // Case where nbTriangle == 1 while nbShape == 2 encountered
739 // with compound of 2 boxes and "To mesh holes"==False,
740 // so there are no subdomains specified for each tetrahedron.
741 // Try to guess a solid by a node already bound to shape
743 for ( int i=0; i<4 && shapeID==0; i++ ) {
744 if ( nodeAssigne[ nodeID[i] ] == 1 &&
745 node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
746 node[i]->GetPosition()->GetShapeId() > 1 )
748 shapeID = node[i]->GetPosition()->GetShapeId();
752 aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
753 shapeID = theMeshDS->ShapeToIndex( aSolid );
756 // set new nodes and tetrahedron onto the shape
757 for ( int i=0; i<4; i++ ) {
758 if ( nodeAssigne[ nodeID[i] ] == 0 ) {
760 theMeshDS->SetNodeInVolume( node[i], shapeID );
761 nodeAssigne[ nodeID[i] ] = shapeID;
764 if ( toMeshHoles || shapeID > 0 ) {
765 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
766 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
769 shapeIDs.insert( shapeID );
772 // Remove nodes of tetras inside holes if !toMeshHoles
773 if ( !toMeshHoles ) {
774 itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
775 for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
776 ID = itOnNode->first;
777 if ( nodeAssigne[ ID ] < 0 )
778 theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
783 cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
784 munmap(mapPtr, length);
792 delete [] nodeAssigne;
795 if ( shapeIDs.size() != nbShape ) {
796 cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << endl;
797 for (int i=0; i<nbShape; i++) {
798 shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
799 if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
800 cout << " Solid #" << shapeID << " not found" << endl;
808 //=======================================================================
809 //function : readResultFile
811 //=======================================================================
813 static bool readResultFile(const int fileOpen,
814 SMESHDS_Mesh* theMeshDS,
816 vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
826 int nbElems, nbNodes, nbInputNodes;
827 int nodeId, triangleId;
833 const SMDS_MeshNode **node;
836 coord = new double[3];
837 node = new const SMDS_MeshNode*[4];
839 SMDS_MeshNode * aNewNode;
840 map <int,const SMDS_MeshNode*>::iterator IdNode;
841 SMDS_MeshElement* aTet;
843 // Read the file state
844 fileStat = fstat(fileOpen, &status);
845 length = status.st_size;
847 // Mapping the result file into memory
848 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
851 ptr = readMapIntLine(ptr, tab);
856 nbInputNodes = tab[2];
858 theNodeByGhs3dId.resize( nbNodes );
860 // Reading the nodeId
861 for (int i=0; i < 4*nbElems; i++)
862 nodeId = strtol(ptr, &ptr, 10);
864 // Reading the nodeCoor and update the nodeMap
865 shapeID = theMeshDS->ShapeToIndex( aSolid );
866 for (int iNode=0; iNode < nbNodes; iNode++) {
867 for (int iCoor=0; iCoor < 3; iCoor++)
868 coord[ iCoor ] = strtod(ptr, &ptr);
869 if ((iNode+1) > nbInputNodes) {
870 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
871 theMeshDS->SetNodeInVolume( aNewNode, shapeID );
872 theNodeByGhs3dId[ iNode ] = aNewNode;
876 // Reading the triangles
877 nbTriangle = strtol(ptr, &ptr, 10);
879 for (int i=0; i < 3*nbTriangle; i++)
880 triangleId = strtol(ptr, &ptr, 10);
884 // Associating the tetrahedrons to the shapes
885 for (int iElem = 0; iElem < nbElems; iElem++) {
886 for (int iNode = 0; iNode < 4; iNode++) {
887 ID = strtol(tetraPtr, &tetraPtr, 10);
888 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
890 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
891 shapeID = theMeshDS->ShapeToIndex( aSolid );
892 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
895 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
896 munmap(mapPtr, length);
906 //================================================================================
908 * \brief Look for a line containing a text in a file
909 * \retval bool - true if the line is found
911 //================================================================================
913 static bool findLineContaing(const TCollection_AsciiString& theText,
914 const TCollection_AsciiString& theFile,
915 TCollection_AsciiString & theFoundLine)
918 if ( FILE * aFile = fopen( theFile.ToCString(), "r" ))
921 char aBuffer[ GHS3DPlugin_BUFLENGTH ];
924 GHS3DPlugin_ReadLine( aPtr, aBuffer, aFile, aLineNb );
927 found = theFoundLine.Search( theText ) >= 0;
929 } while ( aPtr && !found );
936 //================================================================================
938 * \brief Provide human readable text by error code reported by ghs3d
940 //================================================================================
942 static TCollection_AsciiString translateError(const TCollection_AsciiString& errorLine)
944 int beg = errorLine.Location("ERR ", 1, errorLine.Length());
945 if ( !beg ) return errorLine;
947 TCollection_AsciiString errCodeStr = errorLine.SubString( beg + 4 , errorLine.Length());
948 if ( !errCodeStr.IsIntegerValue() )
950 Standard_Integer errCode = errCodeStr.IntegerValue();
952 int codeEnd = beg + 7;
953 while ( codeEnd < errorLine.Length() && isdigit( errorLine.Value(codeEnd+1) ))
955 TCollection_AsciiString error = errorLine.SubString( beg, codeEnd ) + ": ";
959 return "The surface mesh includes a face of type type other than edge, "
960 "triangle or quadrilateral. This face type is not supported.";
962 return error + "Not enough memory for the face table";
964 return error + "Not enough memory";
966 return error + "Not enough memory";
968 return error + "Face is ignored";
970 return error + errorLine.SubString( beg + 5, errorLine.Length()) +
971 ": End of file. Some data are missing in the file.";
973 return error + errorLine.SubString( beg + 5, errorLine.Length()) +
974 ": Read error on the file. There are wrong data in the file.";
976 return error + "the metric file is inadequate (dimension other than 3).";
978 return error + "the metric file is inadequate (values not per vertices).";
980 return error + "the metric file contains more than one field.";
982 return error + "the number of values in the \".bb\" (metric file) is incompatible with the expected"
983 "value of number of mesh vertices in the \".noboite\" file.";
985 return error + "Too many sub-domains";
987 return error + "the number of vertices is negative or null.";
989 return error + "the number of faces is negative or null.";
991 return error + "A facehas a null vertex.";
993 return error + "incompatible data";
995 return error + "the number of vertices is negative or null.";
997 return error + "the number of vertices is negative or null (in the \".mesh\" file).";
999 return error + "the number of faces is negative or null.";
1001 return error + "A face appears more than once in the input surface mesh.";
1003 return error + "An edge appears more than once in the input surface mesh.";
1005 return error + "A face has a vertex negative or null.";
1007 return error + "NOT ENOUGH MEMORY";
1009 return error + "Not enough available memory.";
1011 return error + "Some initial points cannot be inserted. The surface mesh is probably very bad "
1012 "in terms of quality or the input list of points is wrong";
1014 return error + "Some vertices are too close to one another or coincident.";
1016 return error + "Some vertices are too close to one another or coincident.";
1018 return error + "A vertex cannot be inserted.";
1020 return error + "There are at least two points considered as coincident";
1022 return error + "Some vertices are too close to one another or coincident";
1024 return error + "The surface mesh regeneration step has failed";
1026 return error + "Constrained edge cannot be enforced";
1028 return error + "Constrained face cannot be enforced";
1030 return error + "Missing faces";
1032 return error + "No guess to start the definition of the connected component(s)";
1034 return error + "The surface mesh includes at least one hole. The domain is not well defined";
1036 return error + "Impossible to define a component";
1038 return error + "The surface edge intersects another surface edge";
1040 return error + "The surface edge intersects the surface face";
1042 return error + "One boundary point lies within a surface face";
1044 return error + "One surface edge intersects a surface face";
1046 return error + "One boundary point lies within a surface edge";
1048 return error + "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1049 "to too many swaps";
1051 return error + "Edge is unique (i.e., bounds a hole in the surface)";
1053 return error + "Presumably, the surface mesh is not compatible with the domain being processed";
1055 return error + "Too many components, too many sub-domain";
1057 return error + "The surface mesh includes at least one hole. "
1058 "Therefore there is no domain properly defined";
1060 return error + "Statistics.";
1062 return error + "Statistics.";
1064 return error + "Warning, it is dramatically tedious to enforce the boundary items";
1066 return error + "Not enough memory at this time, nevertheless, the program continues. "
1067 "The expected mesh will be correct but not really as large as required";
1069 return error + "see above error code, resulting quality may be poor.";
1071 return error + "Not enough memory at this time, nevertheless, the program continues (warning)";
1073 return error + "Unknown face type.";
1076 return error + errorLine.SubString( beg + 5, errorLine.Length()) +
1077 ": End of file. Some data are missing in the file.";
1079 return error + "A too small volume element is detected";
1081 return error + "There exists at least a null or negative volume element";
1083 return error + "There exist null or negative volume elements";
1085 return error + "A too small volume element is detected. A face is considered being degenerated";
1087 return error + "Some element is suspected to be very bad shaped or wrong";
1089 return error + "A too bad quality face is detected. This face is considered degenerated";
1091 return error + "A too bad quality face is detected. This face is degenerated";
1093 return error + "Presumably, the surface mesh is not compatible with the domain being processed";
1095 return error + "Abnormal error occured, contact hotline";
1097 return error + "Not enough memory for the face table";
1099 return error + "The algorithm cannot run further. "
1100 "The surface mesh is probably very bad in terms of quality.";
1102 return error + "Bad vertex number";
1107 //=============================================================================
1109 *Here we are going to use the GHS3D mesher
1111 //=============================================================================
1113 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1114 const TopoDS_Shape& theShape)
1117 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1119 _nbShape = countShape( meshDS, TopAbs_SOLID ); // we count the number of shapes
1121 // create bounding box for every shape inside the compound
1124 TopoDS_Shape* tabShape;
1126 tabShape = new TopoDS_Shape[_nbShape];
1127 tabBox = new double*[_nbShape];
1128 for (int i=0; i<_nbShape; i++)
1129 tabBox[i] = new double[6];
1130 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1132 TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
1133 for (; expBox.More(); expBox.Next()) {
1134 tabShape[iShape] = expBox.Current();
1135 Bnd_Box BoundingBox;
1136 BRepBndLib::Add(expBox.Current(), BoundingBox);
1137 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1138 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1139 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1140 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1144 // a unique working file name
1145 // to avoid access to the same files by eg different users
1146 TCollection_AsciiString aGenericName
1147 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1149 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1150 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1151 aFacesFileName = aGenericName + ".faces"; // in faces
1152 aPointsFileName = aGenericName + ".points"; // in points
1153 aResultFileName = aGenericName + ".noboite";// out points and volumes
1154 aBadResFileName = aGenericName + ".boite"; // out bad result
1155 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1156 aLogFileName = aGenericName + ".log"; // log
1158 // -----------------
1160 // -----------------
1162 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1163 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1167 aFacesFile->is_open() && aPointsFile->is_open();
1169 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1172 INFOS( "Can't write into " << aFacesFileName);
1173 return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1175 map <int,int> aSmdsToGhs3dIdMap;
1176 map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1178 Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap ) &&
1179 writeFaces ( aFacesFile, meshDS, aSmdsToGhs3dIdMap );
1182 aPointsFile.close();
1185 if ( !_keepFiles ) {
1186 OSD_File( aFacesFileName ).Remove();
1187 OSD_File( aPointsFileName ).Remove();
1189 return error(COMPERR_BAD_INPUT_MESH);
1191 OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1193 // -----------------
1195 // -----------------
1197 TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1198 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1199 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1202 cout << "Ghs3d execution..." << endl;
1203 cout << cmd << endl;
1205 system( cmd.ToCString() ); // run
1208 cout << "End of Ghs3d execution !" << endl;
1214 // Mapping the result file
1217 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1218 if ( fileOpen < 0 ) {
1220 cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << endl;
1225 _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1226 Ok = readResultFile( fileOpen, theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1230 // ---------------------
1231 // remove working files
1232 // ---------------------
1237 OSD_File( aLogFileName ).Remove();
1239 else if ( OSD_File( aLogFileName ).Size() > 0 )
1241 INFOS( "GHS3D Error, see the " << aLogFileName.ToCString() << " file" );
1243 // get problem description from the log file
1244 SMESH_Comment comment;
1245 TCollection_AsciiString foundLine;
1246 if ( findLineContaing( "has expired",aLogFileName,foundLine) &&
1247 foundLine.Search("Licence") >= 0)
1249 foundLine.LeftAdjust();
1250 comment << foundLine;
1252 if ( findLineContaing( "%% ERROR",aLogFileName,foundLine) ||
1253 findLineContaing( " ERR ",aLogFileName,foundLine))
1255 foundLine.LeftAdjust();
1256 comment << translateError( foundLine );
1258 else if ( findLineContaing( "%% NO SAVING OPERATION",aLogFileName,foundLine))
1260 comment << "Too many elements generated for a trial version.\n";
1262 if ( comment.empty() )
1263 comment << "See " << aLogFileName << " for problem description";
1265 comment << "\nSee " << aLogFileName << " for more information";
1266 error(COMPERR_ALGO_FAILED, comment);
1270 // the log file is empty
1271 OSD_File( aLogFileName ).Remove();
1272 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1273 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1276 if ( !_keepFiles ) {
1277 OSD_File( aFacesFileName ).Remove();
1278 OSD_File( aPointsFileName ).Remove();
1279 OSD_File( aResultFileName ).Remove();
1280 OSD_File( aBadResFileName ).Remove();
1281 OSD_File( aBbResFileName ).Remove();
1283 cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1286 cout << "treated !" << endl;
1289 _nbShape = 0; // re-initializing _nbShape for the next Compute() method call
1296 //=============================================================================
1298 *Here we are going to use the GHS3D mesher w/o geometry
1300 //=============================================================================
1301 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
1302 SMESH_MesherHelper* aHelper)
1304 MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1306 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1307 TopoDS_Shape theShape = aHelper->GetSubShape();
1309 // a unique working file name
1310 // to avoid access to the same files by eg different users
1311 TCollection_AsciiString aGenericName
1312 = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1314 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1315 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1316 aFacesFileName = aGenericName + ".faces"; // in faces
1317 aPointsFileName = aGenericName + ".points"; // in points
1318 aResultFileName = aGenericName + ".noboite";// out points and volumes
1319 aBadResFileName = aGenericName + ".boite"; // out bad result
1320 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
1321 aLogFileName = aGenericName + ".log"; // log
1323 // -----------------
1325 // -----------------
1327 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
1328 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1331 aFacesFile->is_open() && aPointsFile->is_open();
1333 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1337 return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1339 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1341 Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1342 writePoints( aPointsFile, meshDS, aNodeByGhs3dId));
1345 aPointsFile.close();
1348 if ( !_keepFiles ) {
1349 OSD_File( aFacesFileName ).Remove();
1350 OSD_File( aPointsFileName ).Remove();
1352 return error(COMPERR_BAD_INPUT_MESH);
1354 OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1356 // -----------------
1358 // -----------------
1360 TCollection_AsciiString cmd =
1361 (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1362 cmd += TCollection_AsciiString(" -f ") + aGenericName; // file to read
1363 cmd += TCollection_AsciiString(" 1>" ) + aLogFileName; // dump into file
1365 system( cmd.ToCString() ); // run
1371 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1372 if ( fileOpen < 0 ) {
1374 cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
1379 Ok = readResultFile( fileOpen, meshDS, theShape ,aNodeByGhs3dId );
1382 // ---------------------
1383 // remove working files
1384 // ---------------------
1389 OSD_File( aLogFileName ).Remove();
1391 else if ( OSD_File( aLogFileName ).Size() > 0 )
1393 INFOS( "GHS3D Error, see the " << aLogFileName.ToCString() << " file" );
1395 // get problem description from the log file
1396 SMESH_Comment comment;
1397 TCollection_AsciiString foundLine;
1398 if ( findLineContaing( "has expired",aLogFileName,foundLine) &&
1399 foundLine.Search("Licence") >= 0)
1401 foundLine.LeftAdjust();
1402 comment << foundLine;
1404 if ( findLineContaing( "%% ERROR",aLogFileName,foundLine) ||
1405 findLineContaing( " ERR ",aLogFileName,foundLine))
1407 foundLine.LeftAdjust();
1408 comment << translateError( foundLine );
1410 else if ( findLineContaing( "%% NO SAVING OPERATION",aLogFileName,foundLine))
1412 comment << "Too many elements generated for a trial version.\n";
1414 if ( comment.empty() )
1415 comment << "See " << aLogFileName << " for problem description";
1417 comment << "\nSee " << aLogFileName << " for more information";
1418 error(COMPERR_ALGO_FAILED, comment);
1421 // the log file is empty
1422 OSD_File( aLogFileName ).Remove();
1423 INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1424 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1429 OSD_File( aFacesFileName ).Remove();
1430 OSD_File( aPointsFileName ).Remove();
1431 OSD_File( aResultFileName ).Remove();
1432 OSD_File( aBadResFileName ).Remove();
1433 OSD_File( aBbResFileName ).Remove();