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 "SMESH_Gen.hxx"
31 #include "SMESH_Mesh.hxx"
32 #include "SMESH_Comment.hxx"
34 #include "SMDS_MeshElement.hxx"
35 #include "SMDS_MeshNode.hxx"
37 #include <TopExp_Explorer.hxx>
38 #include <OSD_File.hxx>
40 #include "utilities.h"
43 #include <sys/sysinfo.h>
46 //#include <Standard_Stream.hxx>
48 #include <BRepGProp.hxx>
49 #include <BRepBndLib.hxx>
50 #include <BRepClass_FaceClassifier.hxx>
51 #include <BRepClass3d_SolidClassifier.hxx>
53 #include <Bnd_Box.hxx>
54 #include <GProp_GProps.hxx>
55 #include <Precision.hxx>
57 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
76 //=============================================================================
80 //=============================================================================
82 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
83 : SMESH_3D_Algo(hypId, studyId, gen)
85 MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
87 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
88 _onlyUnaryInput = false; // Compute() will be called on a compound of solids
93 //=============================================================================
97 //=============================================================================
99 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
101 MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
104 //=============================================================================
108 //=============================================================================
110 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh& aMesh,
111 const TopoDS_Shape& aShape,
112 SMESH_Hypothesis::Hypothesis_Status& aStatus )
114 // MESSAGE("GHS3DPlugin_GHS3D::CheckHypothesis");
115 aStatus = SMESH_Hypothesis::HYP_OK;
119 //================================================================================
121 * \brief Write faces bounding theShape to file
123 //================================================================================
125 static bool writeFaces (ofstream & theFile,
126 SMESHDS_Mesh * theMesh,
127 const TopoDS_Shape& theShape,
128 vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
132 // NB_ELEMS DUMMY_INT
133 // Loop from 1 to NB_ELEMS
134 // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
136 // get all faces bound to theShape
138 // Solids in the ShapeToMesh() can be meshed by different meshers,
139 // so we take faces only from the given shape
140 //TopoDS_Shape theShape = theMesh->ShapeToMesh();
142 list< const SMDS_MeshElement* > faces;
143 // Use TopTools_IndexedMapOfShape in order not to take twice mesh faces from
144 // a geom face shared by two solids
145 TopTools_IndexedMapOfShape faceMap;
146 TopExp::MapShapes( theShape, TopAbs_FACE, faceMap );
148 SMDS_ElemIteratorPtr eIt;
150 const char* space = " ";
151 const int dummyint = 0;
153 list< const SMDS_MeshElement* >::iterator f;
154 map< const SMDS_MeshNode*,int >::iterator it;
155 SMDS_ElemIteratorPtr nodeIt;
156 const SMDS_MeshElement* elem;
160 for ( int i = 0; i < faceMap.Extent(); ++i ) {
161 sm = theMesh->MeshElements( faceMap( i+1 ) );
163 eIt = sm->GetElements();
164 while ( eIt->more() ) {
165 faces.push_back( eIt->next() );
174 cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
176 // NB_ELEMS DUMMY_INT
177 theFile << space << nbFaces << space << dummyint << endl;
179 // Loop from 1 to NB_ELEMS
181 map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
183 for ( ; f != faces.end(); ++f )
187 nbNodes = elem->NbNodes();
188 theFile << space << nbNodes;
190 // NODE_NB_1 NODE_NB_2 ...
191 nodeIt = elem->nodesIterator();
192 while ( nodeIt->more() )
195 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
196 int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
197 it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
198 theFile << space << it->second;
201 // (NB_NODES + 1) times: DUMMY_INT
202 for ( int i=0; i<=nbNodes; i++)
203 theFile << space << dummyint;
208 // put nodes to theNodeByGhs3dId vector
209 theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
210 map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
211 for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
213 theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
219 //=======================================================================
220 //function : writePoints
222 //=======================================================================
224 static bool writePoints (ofstream & theFile,
225 SMESHDS_Mesh * theMesh,
226 const vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
231 // Loop from 1 to NB_NODES
234 //int nbNodes = theMesh->NbNodes();
235 int nbNodes = theNodeByGhs3dId.size();
239 const char* space = " ";
240 const int dummyint = 0;
242 const SMDS_MeshNode* node;
245 theFile << space << nbNodes << endl;
246 cout << nbNodes << " nodes" << endl;
248 // Loop from 1 to NB_NODES
250 vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
251 vector<const SMDS_MeshNode*>::const_iterator after = theNodeByGhs3dId.end();
252 for ( ; nodeIt != after; ++nodeIt )
258 << space << node->X()
259 << space << node->Y()
260 << space << node->Z()
261 << space << dummyint;
269 //=======================================================================
270 //function : findSolid
272 //=======================================================================
274 static TopoDS_Shape findSolid(const SMDS_MeshNode *aNode[],
276 const TopoDS_Shape shape[],
277 const double box[][6],
280 Standard_Real PX, PY, PZ;
283 PX = ( aNode[0]->X() + aNode[1]->X() + aNode[2]->X() + aNode[3]->X() ) / 4.0;
284 PY = ( aNode[0]->Y() + aNode[1]->Y() + aNode[2]->Y() + aNode[3]->Y() ) / 4.0;
285 PZ = ( aNode[0]->Z() + aNode[1]->Z() + aNode[2]->Z() + aNode[3]->Z() ) / 4.0;
286 gp_Pnt aPnt(PX, PY, PZ);
288 BRepClass3d_SolidClassifier SC (aSolid, aPnt, Precision::Confusion());
289 if ( not(SC.State() == TopAbs_IN) ) {
290 for (iShape = 0; iShape < nShape; iShape++) {
291 aSolid = shape[iShape];
292 if ( not( PX < box[iShape][0] || box[iShape][1] < PX ||
293 PY < box[iShape][2] || box[iShape][3] < PY ||
294 PZ < box[iShape][4] || box[iShape][5] < PZ) ) {
295 BRepClass3d_SolidClassifier SC (aSolid, aPnt, Precision::Confusion());
296 if (SC.State() == TopAbs_IN)
304 //=======================================================================
305 //function : readMapIntLine
307 //=======================================================================
309 static char* readMapIntLine(char* ptr, int tab[]) {
313 for ( int i=0; i<17; i++ ) {
314 intVal = strtol(ptr, &ptr, 10);
321 //=======================================================================
322 //function : readLine
324 //=======================================================================
326 #define GHS3DPlugin_BUFLENGTH 256
327 #define GHS3DPlugin_ReadLine(aPtr,aBuf,aFile,aLineNb) \
328 { aPtr = fgets( aBuf, GHS3DPlugin_BUFLENGTH - 2, aFile ); aLineNb++; DUMP(endl); }
331 //=======================================================================
332 //function : readResultFile
334 //=======================================================================
336 static bool readResultFile(const int fileOpen,
337 SMESHDS_Mesh* theMeshDS,
338 TopoDS_Shape tabShape[],
341 vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
351 int nbElems, nbNodes, nbInputNodes;
352 int nodeId, triangleId;
353 int tab[3]/*, tabID[nShape]*/;
355 int ID, shapeID, ghs3dShapeID;
358 vector< int > tabID (nShape, 0);
361 SMDS_MeshNode * aNewNode;
362 const SMDS_MeshNode * node[4];
363 map <int,const SMDS_MeshNode*>::iterator IdNode;
364 SMDS_MeshElement* aTet;
366 // for (int i=0; i<nShape; i++)
369 // Read the file state
370 fileStat = fstat(fileOpen, &status);
371 length = status.st_size;
373 // Mapping the result file into memory
374 ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
377 ptr = readMapIntLine(ptr, tab);
382 nbInputNodes = tab[2];
384 theNodeByGhs3dId.resize( nbNodes );
386 // Reading the nodeId
387 for (int i=0; i < 4*nbElems; i++)
388 nodeId = strtol(ptr, &ptr, 10);
390 // Reading the nodeCoor and update the nodeMap
391 for (int iNode=0; iNode < nbNodes; iNode++) {
392 for (int iCoor=0; iCoor < 3; iCoor++)
393 coord[ iCoor ] = strtod(ptr, &ptr);
394 if ((iNode+1) > nbInputNodes) {
395 aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
396 theNodeByGhs3dId[ iNode ] = aNewNode;
400 // Reading the triangles
401 nbTriangle = strtol(ptr, &ptr, 10);
403 for (int i=0; i < 3*nbTriangle; i++)
404 triangleId = strtol(ptr, &ptr, 10);
408 // Associating the tetrahedrons to the shapes
409 for (int iElem = 0; iElem < nbElems; iElem++) {
410 for (int iNode = 0; iNode < 4; iNode++) {
411 ID = strtol(tetraPtr, &tetraPtr, 10);
412 node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
414 aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
415 ghs3dShapeID = strtol(shapePtr, &shapePtr, 10);
416 if ( !ghs3dShapeID ) ghs3dShapeID = 1;
417 if ( tabID[ ghs3dShapeID - 1 ] == 0 ) {
419 aSolid = tabShape[0];
420 aSolid = findSolid(node, aSolid, tabShape, tabBox, nShape /*nbTriangle*/);
421 shapeID = theMeshDS->ShapeToIndex( aSolid );
422 tabID[ ghs3dShapeID - 1] = shapeID;
425 shapeID = tabID[ ghs3dShapeID - 1];
427 theMeshDS->SetMeshElementOnShape( aTet, shapeID );
428 // set new nodes on to the shape
429 SMDS_ElemIteratorPtr nodeIt = aTet->nodesIterator();
430 while ( nodeIt->more() ) {
431 const SMDS_MeshNode * n = castToNode( nodeIt->next() );
432 if ( !n->GetPosition()->GetShapeId() )
433 theMeshDS->SetNodeInVolume( n, shapeID );
437 cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
438 munmap(mapPtr, length);
443 //=======================================================================
444 //function : getTmpDir
446 //=======================================================================
448 static TCollection_AsciiString getTmpDir()
450 TCollection_AsciiString aTmpDir;
452 char *Tmp_dir = getenv("SALOME_TMP_DIR");
453 if(Tmp_dir != NULL) {
456 if(aTmpDir.Value(aTmpDir.Length()) != '\\') aTmpDir+='\\';
458 if(aTmpDir.Value(aTmpDir.Length()) != '/') aTmpDir+='/';
463 aTmpDir = TCollection_AsciiString("C:\\");
465 aTmpDir = TCollection_AsciiString("/tmp/");
471 //================================================================================
473 * \brief Look for a line containing a text in a file
474 * \retval bool - true if the line is found
476 //================================================================================
478 static bool findLineContaing(const TCollection_AsciiString& theText,
479 const TCollection_AsciiString& theFile,
480 TCollection_AsciiString & theFoundLine)
483 if ( FILE * aFile = fopen( theFile.ToCString(), "r" ))
486 char aBuffer[ GHS3DPlugin_BUFLENGTH ];
489 GHS3DPlugin_ReadLine( aPtr, aBuffer, aFile, aLineNb );
492 found = theFoundLine.Search( theText ) >= 0;
494 } while ( aPtr && !found );
501 //=============================================================================
503 *Here we are going to use the GHS3D mesher
505 //=============================================================================
507 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh& theMesh,
508 const TopoDS_Shape& theShape)
510 // theShape is a compound of solids as _onlyUnaryInput = false
512 SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
515 /*if (_iShape == 0 && _nbShape == 0)*/ {
517 cout << "Ghs3d execution..." << endl;
520 //TopExp_Explorer exp (meshDS->ShapeToMesh(), TopAbs_SOLID);
521 TopExp_Explorer exp (theShape, TopAbs_SOLID);
522 for (; exp.More(); exp.Next())
528 /*if ( _iShape == _nbShape )*/ {
530 // create bounding box for every shape
533 TopoDS_Shape tabShape[_nbShape];
534 double tabBox[_nbShape][6];
535 Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
537 TopExp_Explorer expBox (theShape, TopAbs_SOLID);
538 for (; expBox.More(); expBox.Next()) {
539 tabShape[iShape] = expBox.Current();
541 BRepBndLib::Add(expBox.Current(), BoundingBox);
542 BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
543 tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
544 tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
545 tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
549 // make a unique working file name
550 // to avoid access to the same files by eg different users
552 TCollection_AsciiString aGenericName, aTmpDir = getTmpDir();
553 aGenericName = aTmpDir + "GHS3D_";
555 aGenericName += GetCurrentProcessId();
557 aGenericName += getpid();
560 aGenericName += meshDS->ShapeToIndex( theShape );
562 TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
563 TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
564 aFacesFileName = aGenericName + ".faces"; // in faces
565 aPointsFileName = aGenericName + ".points"; // in points
566 aResultFileName = aGenericName + ".noboite";// out points and volumes
567 aBadResFileName = aGenericName + ".boite"; // out bad result
568 aBbResFileName = aGenericName + ".bb"; // out vertex stepsize
569 aLogFileName = aGenericName + ".log"; // log
575 ofstream aFacesFile ( aFacesFileName.ToCString() , ios::out);
576 ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
580 aFacesFile->is_open() && aPointsFile->is_open();
582 aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
585 return error(SMESH_Comment("Can't write into ") << aTmpDir);
587 vector <const SMDS_MeshNode*> aNodeByGhs3dId;
589 Ok = ( writeFaces ( aFacesFile, meshDS, theShape, aNodeByGhs3dId ) &&
590 writePoints( aPointsFile, meshDS, aNodeByGhs3dId ));
596 if ( !getenv("GHS3D_KEEP_FILES") ) {
597 OSD_File( aFacesFileName ).Remove();
598 OSD_File( aPointsFileName ).Remove();
600 return error(COMPERR_BAD_INPUT_MESH);
604 // run ghs3d mesher WIN32???
607 // ghs3d need to know amount of memory it may use (MB).
608 // Default memory is defined at ghs3d installation but it may be not enough,
609 // so allow to use about all available memory
611 TCollection_AsciiString memory;
616 int err = sysinfo( &si );
618 int freeMem = si.totalram * si.mem_unit / 1024 / 1024;
620 memory += int( 0.7 * freeMem );
624 OSD_File( aResultFileName ).Remove(); // old file prevents writing a new one
626 TCollection_AsciiString cmd( "ghs3d " ); // command to run
629 " -c 0" // 0 - mesh all components, 1 - keep only the main component
630 " -f " + aGenericName + // file to read
631 " 1>" + aLogFileName; // dump into file
633 MESSAGE("GHS3DPlugin_GHS3D::Compute() " << cmd );
634 system( cmd.ToCString() ); // run
637 cout << "End of Ghs3d execution !" << endl;
643 // Mapping the result file
646 fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
647 if ( fileOpen < 0 ) {
649 cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
654 Ok = readResultFile( fileOpen, meshDS, tabShape, tabBox, _nbShape, aNodeByGhs3dId );
657 // ---------------------
658 // remove working files
659 // ---------------------
663 OSD_File( aLogFileName ).Remove();
665 else if ( OSD_File( aLogFileName ).Size() > 0 )
667 // get problem description from the log file
668 SMESH_Comment comment;
669 TCollection_AsciiString foundLine;
670 if ( findLineContaing( "has expired",aLogFileName,foundLine) &&
671 foundLine.Search("Licence") >= 0)
673 foundLine.LeftAdjust();
674 comment << foundLine;
676 if ( findLineContaing( "%% ERROR",aLogFileName,foundLine))
678 foundLine.LeftAdjust();
679 comment << foundLine;
681 if ( findLineContaing( "%% NO SAVING OPERATION",aLogFileName,foundLine))
683 comment << "Too many elements generated for a trial version.\n";
685 if ( comment.empty() )
686 comment << "See " << aLogFileName << " for problem description";
688 comment << "See " << aLogFileName << " for more information";
689 error(COMPERR_ALGO_FAILED, comment);
693 // the log file is empty
694 OSD_File( aLogFileName ).Remove();
695 error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
698 if ( !getenv("GHS3D_KEEP_FILES") )
700 OSD_File( aFacesFileName ).Remove();
701 OSD_File( aPointsFileName ).Remove();
702 OSD_File( aResultFileName ).Remove();
703 OSD_File( aBadResFileName ).Remove();
704 OSD_File( aBbResFileName ).Remove();
706 /*if ( _iShape == _nbShape )*/ {
707 cout << aResultFileName.ToCString() << " Output file ";
710 cout << "treated !" << endl;