]> SALOME platform Git repositories - plugins/ghs3dplugin.git/blob - src/GHS3DPlugin_GHS3D.cxx
Salome HOME
PAL19680 CEA 4.1.2: Meshers: BLSURF, GHS3D and holed shapes
[plugins/ghs3dplugin.git] / src / GHS3DPlugin_GHS3D.cxx
1 // Copyright (C) 2005  CEA/DEN, EDF R&D, OPEN CASCADE, PRINCIPIA R&D
2 //
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.
7 //
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.
12 //
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
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19 //=============================================================================
20 // File      : GHS3DPlugin_GHS3D.cxx
21 // Created   : 
22 // Author    : Edward AGAPOV, modified by Lioka RAZAFINDRAZAKA (CEA) 09/02/2007
23 // Project   : SALOME
24 // Copyright : CEA 2003
25 // $Header$
26 //=============================================================================
27 using namespace std;
28
29 #include "GHS3DPlugin_GHS3D.hxx"
30 #include "GHS3DPlugin_Hypothesis.hxx"
31
32 #include "SMESH_Gen.hxx"
33 #include "SMESH_Mesh.hxx"
34 #include "SMESH_Comment.hxx"
35 #include "SMESH_MesherHelper.hxx"
36
37 #include "SMDS_MeshElement.hxx"
38 #include "SMDS_MeshNode.hxx"
39
40 #include <TopExp.hxx>
41 #include <TopExp_Explorer.hxx>
42 #include <OSD_File.hxx>
43
44 #include "utilities.h"
45
46 #ifndef WIN32
47 #include <sys/sysinfo.h>
48 #endif
49
50 //#include <Standard_Stream.hxx>
51
52 #include <BRepGProp.hxx>
53 #include <BRepBndLib.hxx>
54 #include <BRepClass_FaceClassifier.hxx>
55 #include <BRepClass3d_SolidClassifier.hxx>
56 #include <TopAbs.hxx>
57 #include <Bnd_Box.hxx>
58 #include <GProp_GProps.hxx>
59 #include <Precision.hxx>
60
61 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
62
63 #ifdef _DEBUG_
64 #define DUMP(txt) \
65 //  cout << txt
66 #else
67 #define DUMP(txt)
68 #endif
69
70 extern "C"
71 {
72 #ifndef WNT
73 #include <unistd.h>
74 #include <sys/mman.h>
75 #endif
76 #include <sys/stat.h>
77 #include <fcntl.h>
78 }
79
80 //=============================================================================
81 /*!
82  *  
83  */
84 //=============================================================================
85
86 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
87   : SMESH_3D_Algo(hypId, studyId, gen)
88 {
89   MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
90   _name = "GHS3D_3D";
91   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
92   _onlyUnaryInput = false; // Compute() will be called on a compound of solids
93   _iShape=0;
94   _nbShape=0;
95   _compatibleHypothesis.push_back("GHS3D_Parameters");
96 }
97
98 //=============================================================================
99 /*!
100  *  
101  */
102 //=============================================================================
103
104 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
105 {
106   MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
107 }
108
109 //=============================================================================
110 /*!
111  *  
112  */
113 //=============================================================================
114
115 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh&         aMesh,
116                                           const TopoDS_Shape& aShape,
117                                           Hypothesis_Status&  aStatus )
118 {
119   aStatus = SMESH_Hypothesis::HYP_OK;
120
121   // there is only one compatible Hypothesis so far
122   _hyp = 0;
123   _keepFiles = false;
124   const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
125   if ( !hyps.empty() )
126     _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
127   if ( _hyp )
128     _keepFiles = _hyp->GetKeepFiles();
129
130   return true;
131 }
132
133 //=======================================================================
134 //function : findShape
135 //purpose  : 
136 //=======================================================================
137
138 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
139                               TopoDS_Shape        aShape,
140                               const TopoDS_Shape  shape[],
141                               double**            box,
142                               const int           nShape) {
143   double *pntCoor;
144   int iShape, nbNode = 4;
145
146   pntCoor = new double[3];
147   for ( int i=0; i<3; i++ ) {
148     pntCoor[i] = 0;
149     for ( int j=0; j<nbNode; j++ ) {
150       if ( i == 0) pntCoor[i] += aNode[j]->X();
151       if ( i == 1) pntCoor[i] += aNode[j]->Y();
152       if ( i == 2) pntCoor[i] += aNode[j]->Z();
153     }
154     pntCoor[i] /= nbNode;
155   }
156
157   gp_Pnt aPnt(pntCoor[0], pntCoor[1], pntCoor[2]);
158   BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
159   if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
160     for (iShape = 0; iShape < nShape; iShape++) {
161       aShape = shape[iShape];
162       if ( !( pntCoor[0] < box[iShape][0] || box[iShape][1] < pntCoor[0] ||
163               pntCoor[1] < box[iShape][2] || box[iShape][3] < pntCoor[1] ||
164               pntCoor[2] < box[iShape][4] || box[iShape][5] < pntCoor[2]) ) {
165         BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
166         if (SC.State() == TopAbs_IN)
167           break;
168       }
169     }
170   }
171   return aShape;
172 }
173
174 //=======================================================================
175 //function : readMapIntLine
176 //purpose  : 
177 //=======================================================================
178
179 static char* readMapIntLine(char* ptr, int tab[]) {
180   long int intVal;
181   cout << endl;
182
183   for ( int i=0; i<17; i++ ) {
184     intVal = strtol(ptr, &ptr, 10);
185     if ( i < 3 )
186       tab[i] = intVal;
187   }
188   return ptr;
189 }
190
191 //=======================================================================
192 //function : readLine
193 //purpose  : 
194 //=======================================================================
195
196 #define GHS3DPlugin_BUFLENGTH 256
197 #define GHS3DPlugin_ReadLine(aPtr,aBuf,aFile,aLineNb) \
198 {  aPtr = fgets( aBuf, GHS3DPlugin_BUFLENGTH - 2, aFile ); aLineNb++; DUMP(endl); }
199
200 //=======================================================================
201 //function : countShape
202 //purpose  :
203 //=======================================================================
204
205 template < class Mesh, class Shape >
206 static int countShape( Mesh* mesh, Shape shape ) {
207   TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
208   int nbShape = 0;
209   for ( ; expShape.More(); expShape.Next() ) {
210       nbShape++;
211   }
212   return nbShape;
213 }
214
215 //=======================================================================
216 //function : writeFaces
217 //purpose  : 
218 //=======================================================================
219
220 static bool writeFaces (ofstream &            theFile,
221                         SMESHDS_Mesh *        theMesh,
222                         const map <int,int> & theSmdsToGhs3dIdMap)
223 {
224   // record structure:
225   //
226   // NB_ELEMS DUMMY_INT
227   // Loop from 1 to NB_ELEMS
228   // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
229
230   int nbShape = countShape( theMesh, TopAbs_FACE );
231
232   int *tabID;             tabID    = new int[nbShape];
233   TopoDS_Shape *tabShape; tabShape = new TopoDS_Shape[nbShape];
234   TopoDS_Shape aShape;
235   SMESHDS_SubMesh* theSubMesh;
236   const SMDS_MeshElement* aFace;
237   const char* space    = "  ";
238   const int   dummyint = 0;
239   map<int,int>::const_iterator itOnMap;
240   SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
241   int shapeID, nbNodes, aSmdsID;
242   bool idFound;
243
244   cout << "    " << theMesh->NbFaces() << " shapes of 2D dimension" << endl;
245   cout << endl;
246
247   theFile << space << theMesh->NbFaces() << space << dummyint << endl;
248
249   TopExp_Explorer expface( theMesh->ShapeToMesh(), TopAbs_FACE );
250   for ( int i = 0; expface.More(); expface.Next(), i++ ) {
251     tabID[i] = 0;
252     aShape   = expface.Current();
253     shapeID  = theMesh->ShapeToIndex( aShape );
254     idFound  = false;
255     for ( int j=0; j<=i; j++) {
256       if ( shapeID == tabID[j] ) {
257         idFound = true;
258         break;
259       }
260     }
261     if ( ! idFound ) {
262       tabID[i]    = shapeID;
263       tabShape[i] = aShape;
264     }
265   }
266   for ( int i =0; i < nbShape; i++ ) {
267     if ( tabID[i] != 0 ) {
268       aShape      = tabShape[i];
269       shapeID     = tabID[i];
270       theSubMesh  = theMesh->MeshElements( aShape );
271       if ( !theSubMesh ) continue;
272       itOnSubMesh = theSubMesh->GetElements();
273       while ( itOnSubMesh->more() ) {
274         aFace   = itOnSubMesh->next();
275         nbNodes = aFace->NbNodes();
276
277         theFile << space << nbNodes;
278
279         itOnSubFace = aFace->nodesIterator();
280         while ( itOnSubFace->more() ) {
281           // find GHS3D ID
282           aSmdsID = itOnSubFace->next()->GetID();
283           itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
284           ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
285
286           theFile << space << (*itOnMap).second;
287         }
288
289         // (NB_NODES + 1) times: DUMMY_INT
290         for ( int j=0; j<=nbNodes; j++)
291           theFile << space << dummyint;
292
293         theFile << endl;
294       }
295     }
296   }
297
298   delete [] tabID;
299   delete [] tabShape;
300
301   return true;
302 }
303
304 //=======================================================================
305 //function : writeFaces
306 //purpose  : Write Faces in case if generate 3D mesh w/o geometry
307 //=======================================================================
308
309 static bool writeFaces (ofstream &            theFile,
310                         SMESHDS_Mesh *        theMesh,
311                         vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
312 {
313   // record structure:
314   //
315   // NB_ELEMS DUMMY_INT
316   // Loop from 1 to NB_ELEMS
317   //   NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
318
319
320   int nbFaces = 0;
321   list< const SMDS_MeshElement* > faces;
322   list< const SMDS_MeshElement* >::iterator f;
323   map< const SMDS_MeshNode*,int >::iterator it;
324   SMDS_ElemIteratorPtr nodeIt;
325   const SMDS_MeshElement* elem;
326   int nbNodes;
327
328   const char* space    = "  ";
329   const int   dummyint = 0;
330
331   //get all faces from mesh
332   SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
333   while ( eIt->more() ) {
334     const SMDS_MeshElement* elem = eIt->next();
335     if ( !elem )
336       return false;
337     faces.push_back( elem );
338     nbFaces++;
339   }
340
341   if ( nbFaces == 0 )
342     return false;
343   
344   cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
345
346   // NB_ELEMS DUMMY_INT
347   theFile << space << nbFaces << space << dummyint << endl;
348
349   // Loop from 1 to NB_ELEMS
350
351   map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
352   f = faces.begin();
353   for ( ; f != faces.end(); ++f )
354   {
355     // NB_NODES PER FACE
356     elem = *f;
357     nbNodes = elem->NbNodes();
358     theFile << space << nbNodes;
359
360     // NODE_NB_1 NODE_NB_2 ...
361     nodeIt = elem->nodesIterator();
362     while ( nodeIt->more() )
363     {
364       // find GHS3D ID
365       const SMDS_MeshNode* node = castToNode( nodeIt->next() );
366       int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
367       it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
368       theFile << space << it->second;
369     }
370
371     // (NB_NODES + 1) times: DUMMY_INT
372     for ( int i=0; i<=nbNodes; i++)
373       theFile << space << dummyint;
374
375     theFile << endl;
376   }
377
378   // put nodes to theNodeByGhs3dId vector
379   theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
380   map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
381   for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
382   {
383     theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
384   }
385
386   return true;
387 }
388
389 //=======================================================================
390 //function : writePoints
391 //purpose  : 
392 //=======================================================================
393
394 static bool writePoints (ofstream &                       theFile,
395                          SMESHDS_Mesh *                   theMesh,
396                          map <int,int> &                  theSmdsToGhs3dIdMap,
397                          map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap)
398 {
399   // record structure:
400   //
401   // NB_NODES
402   // Loop from 1 to NB_NODES
403   //   X Y Z DUMMY_INT
404
405   int nbNodes = theMesh->NbNodes();
406   if ( nbNodes == 0 )
407     return false;
408
409   const char* space    = "  ";
410   const int   dummyint = 0;
411
412   int aGhs3dID = 1;
413   SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
414   const SMDS_MeshNode* node;
415
416   // NB_NODES
417   theFile << space << nbNodes << endl;
418   cout << endl;
419   cout << "The initial 2D mesh contains :" << endl;
420   cout << "    " << nbNodes << " vertices" << endl;
421
422   // Loop from 1 to NB_NODES
423
424   while ( it->more() )
425   {
426     node = it->next();
427     theSmdsToGhs3dIdMap.insert( map <int,int>::value_type( node->GetID(), aGhs3dID ));
428     theGhs3dIdToNodeMap.insert (map <int,const SMDS_MeshNode*>::value_type( aGhs3dID, node ));
429     aGhs3dID++;
430
431     // X Y Z DUMMY_INT
432     theFile
433       << space << node->X()
434       << space << node->Y()
435       << space << node->Z()
436       << space << dummyint;
437
438     theFile << endl;
439   }
440
441   return true;
442 }
443
444 //=======================================================================
445 //function : writePoints
446 //purpose  : 
447 //=======================================================================
448
449 static bool writePoints (ofstream &                            theFile,
450                          SMESHDS_Mesh *                        theMesh,
451                          const vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
452 {
453   // record structure:
454   //
455   // NB_NODES
456   // Loop from 1 to NB_NODES
457   //   X Y Z DUMMY_INT
458
459   //int nbNodes = theMesh->NbNodes();
460   int nbNodes = theNodeByGhs3dId.size();
461   if ( nbNodes == 0 )
462     return false;
463
464   const char* space    = "  ";
465   const int   dummyint = 0;
466
467   const SMDS_MeshNode* node;
468
469   // NB_NODES
470   theFile << space << nbNodes << endl;
471   cout << nbNodes << " nodes" << endl;
472
473   // Loop from 1 to NB_NODES
474
475   vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
476   vector<const SMDS_MeshNode*>::const_iterator after  = theNodeByGhs3dId.end();
477   for ( ; nodeIt != after; ++nodeIt )
478   {
479     node = *nodeIt;
480
481     // X Y Z DUMMY_INT
482     theFile
483       << space << node->X()
484       << space << node->Y()
485       << space << node->Z()
486       << space << dummyint;
487
488     theFile << endl;
489   }
490
491   return true;
492 }
493
494 //=======================================================================
495 //function : readResultFile
496 //purpose  : 
497 //=======================================================================
498
499 static bool readResultFile(const int                       fileOpen,
500                            SMESHDS_Mesh*                   theMeshDS,
501                            TopoDS_Shape                    tabShape[],
502                            double**                        tabBox,
503                            const int                       nbShape,
504                            map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap) {
505
506   struct stat status;
507   size_t      length;
508
509   char *ptr, *mapPtr;
510   char *tetraPtr;
511   char *shapePtr;
512
513   int fileStat;
514   int nbElems, nbNodes, nbInputNodes;
515   int nodeId, triangleId;
516   int nbTriangle;
517   int ID, shapeID, ghs3dShapeID;
518   int IdShapeRef = 1;
519   int compoundID =
520     nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
521
522   int *tab, *tabID, *nodeID, *nodeAssigne;
523   double *coord;
524   const SMDS_MeshNode **node;
525
526   tab    = new int[3];
527   tabID  = new int[nbShape];
528   nodeID = new int[4];
529   coord  = new double[3];
530   node   = new const SMDS_MeshNode*[4];
531
532   TopoDS_Shape aSolid;
533   SMDS_MeshNode * aNewNode;
534   map <int,const SMDS_MeshNode*>::iterator itOnNode;
535   SMDS_MeshElement* aTet;
536 #ifdef _DEBUG_
537   set<int> shapeIDs;
538 #endif
539
540   for (int i=0; i<nbShape; i++)
541     tabID[i] = 0;
542
543   // Read the file state
544   fileStat = fstat(fileOpen, &status);
545   length   = status.st_size;
546
547   // Mapping the result file into memory
548   ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
549   mapPtr = ptr;
550
551   ptr      = readMapIntLine(ptr, tab);
552   tetraPtr = ptr;
553
554   nbElems      = tab[0];
555   nbNodes      = tab[1];
556   nbInputNodes = tab[2];
557
558   nodeAssigne = new int[ nbNodes+1 ];
559
560   // Reading the nodeId
561   for (int i=0; i < 4*nbElems; i++)
562     nodeId = strtol(ptr, &ptr, 10);
563
564   // Reading the nodeCoor and update the nodeMap
565   for (int iNode=0; iNode < nbNodes; iNode++) {
566     for (int iCoor=0; iCoor < 3; iCoor++)
567       coord[ iCoor ] = strtod(ptr, &ptr);
568     nodeAssigne[ iNode+1 ] = 1;
569     if ( (iNode+1) > nbInputNodes ) {
570       nodeAssigne[ iNode+1 ] = 0;
571       aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
572       theGhs3dIdToNodeMap.insert(make_pair( (iNode+1), aNewNode ));
573     }
574   }
575
576   // Reading the number of triangles which corresponds to the number of shapes
577   nbTriangle = strtol(ptr, &ptr, 10);
578
579   for (int i=0; i < 3*nbTriangle; i++)
580     triangleId = strtol(ptr, &ptr, 10);
581
582   shapePtr = ptr;
583
584   // Associating the tetrahedrons to the shapes
585   shapeID = compoundID;
586   for (int iElem = 0; iElem < nbElems; iElem++) {
587     for (int iNode = 0; iNode < 4; iNode++) {
588       ID = strtol(tetraPtr, &tetraPtr, 10);
589       itOnNode = theGhs3dIdToNodeMap.find(ID);
590       node[ iNode ] = itOnNode->second;
591       nodeID[ iNode ] = ID;
592     }
593     aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
594     if ( nbShape > 1 ) {
595       if ( nbTriangle > 1 ) {
596         ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
597         if ( tabID[ ghs3dShapeID ] == 0 ) {
598           if (iElem == 0)
599             aSolid = tabShape[0];
600           aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
601           shapeID = theMeshDS->ShapeToIndex( aSolid );
602           tabID[ ghs3dShapeID ] = shapeID;
603         }
604         else
605           shapeID = tabID[ ghs3dShapeID ];
606       }
607       else {
608         // Case where nbTriangle == 1 while nbShape == 2 encountered
609         // with compound of 2 boxes and "To mesh holes"==False,
610         // so there are no subdomains specified for each tetrahedron.
611         // Try to guess a solid by a node already bound to shape
612         shapeID = 0;
613         for ( int i=0; i<4 && shapeID==0; i++ ) {
614           if ( nodeAssigne[ nodeID[i] ] == 1 &&
615                node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
616                node[i]->GetPosition()->GetShapeId() > 1 )
617           {
618             shapeID = node[i]->GetPosition()->GetShapeId();
619           }
620         }
621         if ( shapeID==0 ) {
622           aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
623           shapeID = theMeshDS->ShapeToIndex( aSolid );
624         }
625       }
626     }
627     // set new nodes and tetrahedron on to the shape
628     for ( int i=0; i<4; i++ ) {
629       if ( nodeAssigne[ nodeID[i] ] == 0 ) {
630         theMeshDS->SetNodeInVolume( node[i], shapeID );
631         nodeAssigne[ nodeID[i] ] = 1;
632       }
633     }
634     theMeshDS->SetMeshElementOnShape( aTet, shapeID );
635 #ifdef _DEBUG_
636     shapeIDs.insert( shapeID );
637 #endif
638   }
639   if ( nbElems )
640     cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
641   munmap(mapPtr, length);
642   close(fileOpen);
643
644   delete [] tab;
645   delete [] tabID;
646   delete [] nodeID;
647   delete [] coord;
648   delete [] node;
649   delete [] nodeAssigne;
650
651 #ifdef _DEBUG_
652   if ( shapeIDs.size() != nbShape ) {
653     cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << endl;
654     for (int i=0; i<nbShape; i++) {
655       shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
656       if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
657         cout << "  Solid #" << shapeID << " not found" << endl;
658     }
659   }
660 #endif
661
662   return true;
663 }
664
665 //=======================================================================
666 //function : readResultFile
667 //purpose  : 
668 //=======================================================================
669
670 static bool readResultFile(const int                      fileOpen,
671                            SMESHDS_Mesh*                  theMeshDS,
672                            TopoDS_Shape                   aSolid,
673                            vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
674
675   struct stat  status;
676   size_t       length;
677
678   char *ptr, *mapPtr;
679   char *tetraPtr;
680   char *shapePtr;
681
682   int fileStat;
683   int nbElems, nbNodes, nbInputNodes;
684   int nodeId, triangleId;
685   int nbTriangle;
686   int ID, shapeID;
687
688   int *tab;
689   double *coord;
690   const SMDS_MeshNode **node;
691
692   tab   = new int[3];
693   coord = new double[3];
694   node  = new const SMDS_MeshNode*[4];
695
696   SMDS_MeshNode * aNewNode;
697   map <int,const SMDS_MeshNode*>::iterator IdNode;
698   SMDS_MeshElement* aTet;
699
700   // Read the file state
701   fileStat = fstat(fileOpen, &status);
702   length   = status.st_size;
703
704   // Mapping the result file into memory
705   ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
706   mapPtr = ptr;
707
708   ptr      = readMapIntLine(ptr, tab);
709   tetraPtr = ptr;
710
711   nbElems      = tab[0];
712   nbNodes      = tab[1];
713   nbInputNodes = tab[2];
714
715   theNodeByGhs3dId.resize( nbNodes );
716
717   // Reading the nodeId
718   for (int i=0; i < 4*nbElems; i++)
719     nodeId = strtol(ptr, &ptr, 10);
720
721   // Reading the nodeCoor and update the nodeMap
722   shapeID = theMeshDS->ShapeToIndex( aSolid );
723   for (int iNode=0; iNode < nbNodes; iNode++) {
724     for (int iCoor=0; iCoor < 3; iCoor++)
725       coord[ iCoor ] = strtod(ptr, &ptr);
726     if ((iNode+1) > nbInputNodes) {
727       aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
728       theMeshDS->SetNodeInVolume( aNewNode, shapeID );
729       theNodeByGhs3dId[ iNode ] = aNewNode;
730     }
731   }
732
733   // Reading the triangles
734   nbTriangle = strtol(ptr, &ptr, 10);
735
736   for (int i=0; i < 3*nbTriangle; i++)
737     triangleId = strtol(ptr, &ptr, 10);
738
739   shapePtr = ptr;
740
741   // Associating the tetrahedrons to the shapes
742   for (int iElem = 0; iElem < nbElems; iElem++) {
743     for (int iNode = 0; iNode < 4; iNode++) {
744       ID = strtol(tetraPtr, &tetraPtr, 10);
745       node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
746     }
747     aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
748     shapeID = theMeshDS->ShapeToIndex( aSolid );
749     theMeshDS->SetMeshElementOnShape( aTet, shapeID );
750   }
751   if ( nbElems )
752     cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
753   munmap(mapPtr, length);
754   close(fileOpen);
755
756   delete [] tab;
757   delete [] coord;
758   delete [] node;
759
760   return true;
761 }
762
763 //================================================================================
764 /*!
765  * \brief Look for a line containing a text in a file
766   * \retval bool - true if the line is found
767  */
768 //================================================================================
769
770 static bool findLineContaing(const TCollection_AsciiString& theText,
771                              const TCollection_AsciiString& theFile,
772                              TCollection_AsciiString &      theFoundLine)
773 {
774   bool found = false;
775   if ( FILE * aFile = fopen( theFile.ToCString(), "r" ))
776   {
777     char * aPtr;
778     char aBuffer[ GHS3DPlugin_BUFLENGTH ];
779     int aLineNb = 0;
780     do {
781       GHS3DPlugin_ReadLine( aPtr, aBuffer, aFile, aLineNb );
782       if ( aPtr ) {
783         theFoundLine = aPtr;
784         found = theFoundLine.Search( theText ) >= 0;
785       }
786     } while ( aPtr && !found );
787
788     fclose( aFile );
789   }
790   return found;
791 }
792
793 //================================================================================
794 /*!
795  * \brief Provide human readable text by error code reported by ghs3d
796  */
797 //================================================================================
798
799 static TCollection_AsciiString translateError(const TCollection_AsciiString& errorLine)
800 {
801   int beg = errorLine.Location("ERR ", 1, errorLine.Length());
802   if ( !beg ) return errorLine;
803
804   TCollection_AsciiString errCodeStr = errorLine.SubString( beg + 4 , errorLine.Length());
805   if ( !errCodeStr.IsIntegerValue() )
806     return errorLine;
807   Standard_Integer errCode = errCodeStr.IntegerValue();
808
809   int codeEnd = beg + 7;
810   while ( codeEnd < errorLine.Length() && isdigit( errorLine.Value(codeEnd+1) ))
811     codeEnd++;
812   TCollection_AsciiString error = errorLine.SubString( beg, codeEnd ) + ": ";
813
814   switch ( errCode ) {
815   case 0:
816     return "The surface mesh includes a face of type type other than edge, "
817       "triangle or quadrilateral. This face type is not supported.";
818   case 1:
819     return error + "Not enough memory for the face table";
820   case 2:
821     return error + "Not enough memory";
822   case 3:
823     return error + "Not enough memory";
824   case 4:
825     return error + "Face is ignored";
826   case 5:
827     return error + errorLine.SubString( beg + 5, errorLine.Length()) +
828       ": End of file. Some data are missing in the file.";
829   case 6:
830     return error + errorLine.SubString( beg + 5, errorLine.Length()) +
831       ": Read error on the file. There are wrong data in the file.";
832   case 7:
833     return error + "the metric file is inadequate (dimension other than 3).";
834   case 8:
835     return error + "the metric file is inadequate (values not per vertices).";
836   case 9:
837       return error + "the metric file contains more than one field.";
838   case 10:
839     return error + "the number of values in the \".bb\" (metric file) is incompatible with the expected"
840       "value of number of mesh vertices in the \".noboite\" file.";
841   case 12:
842     return error + "Too many sub-domains";
843   case 13:
844     return error + "the number of vertices is negative or null.";
845   case 14:
846     return error + "the number of faces is negative or null.";
847   case 15:
848     return error + "A facehas a null vertex.";
849   case 22:
850     return error + "incompatible data";
851   case 131:
852     return error + "the number of vertices is negative or null.";
853   case 132:
854     return error + "the number of vertices is negative or null (in the \".mesh\" file).";
855   case 133:
856     return error + "the number of faces is negative or null.";
857   case 1000:
858     return error + "A face appears more than once in the input surface mesh.";
859   case 1001:
860     return error + "An edge appears more than once in the input surface mesh.";
861   case 1002:
862     return error + "A face has a vertex negative or null.";
863   case 1003:
864     return error + "NOT ENOUGH MEMORY";
865   case 2000:
866     return error + "Not enough available memory.";
867   case 2002:
868     return error + "Some initial points cannot be inserted. The surface mesh is probably very bad "
869       "in terms of quality or the input list of points is wrong";
870   case 2003:
871     return error + "Some vertices are too close to one another or coincident.";
872   case 2004:
873     return error + "Some vertices are too close to one another or coincident.";
874   case 2012:
875     return error + "A vertex cannot be inserted.";
876   case 2014:
877     return error + "There are at least two points considered as coincident";
878   case 2103:
879     return error + "Some vertices are too close to one another or coincident";
880   case 3000:
881     return error + "The surface mesh regeneration step has failed";
882   case 3009:
883     return error + "Constrained edge cannot be enforced";
884   case 3019:
885     return error + "Constrained face cannot be enforced";
886   case 3029:
887     return error + "Missing faces";
888   case 3100:
889     return error + "No guess to start the definition of the connected component(s)";
890   case 3101:
891     return error + "The surface mesh includes at least one hole. The domain is not well defined";
892   case 3102:
893     return error + "Impossible to define a component";
894   case 3103:
895     return error + "The surface edge intersects another surface edge";
896   case 3104:
897     return error + "The surface edge intersects the surface face";
898   case 3105:
899     return error + "One boundary point lies within a surface face";
900   case 3106:
901     return error + "One surface edge intersects a surface face";
902   case 3107:
903     return error + "One boundary point lies within a surface edge";
904   case 3108:
905     return error + "Insufficient memory ressources detected due to a bad quality surface mesh leading "
906       "to too many swaps";
907   case 3109:
908     return error + "Edge is unique (i.e., bounds a hole in the surface)";
909   case 3122:
910     return error + "Presumably, the surface mesh is not compatible with the domain being processed";
911   case 3123:
912     return error + "Too many components, too many sub-domain";
913   case 3209:
914     return error + "The surface mesh includes at least one hole. "
915       "Therefore there is no domain properly defined";
916   case 3300:
917     return error + "Statistics.";
918   case 3400:
919     return error + "Statistics.";
920   case 3500:
921     return error + "Warning, it is dramatically tedious to enforce the boundary items";
922   case 4000:
923     return error + "Not enough memory at this time, nevertheless, the program continues. "
924       "The expected mesh will be correct but not really as large as required";
925   case 4002:
926     return error + "see above error code, resulting quality may be poor.";
927   case 4003:
928     return error + "Not enough memory at this time, nevertheless, the program continues (warning)";
929   case 8000:
930     return error + "Unknown face type.";
931   case 8005:
932   case 8006:
933     return error + errorLine.SubString( beg + 5, errorLine.Length()) +
934       ": End of file. Some data are missing in the file.";
935   case 9000:
936     return error + "A too small volume element is detected";
937   case 9001:
938     return error + "There exists at least a null or negative volume element";
939   case 9002:
940     return error + "There exist null or negative volume elements";
941   case 9003:
942     return error + "A too small volume element is detected. A face is considered being degenerated";
943   case 9100:
944     return error + "Some element is suspected to be very bad shaped or wrong";
945   case 9102:
946     return error + "A too bad quality face is detected. This face is considered degenerated";
947   case 9112:
948     return error + "A too bad quality face is detected. This face is degenerated";
949   case 9122:
950     return error + "Presumably, the surface mesh is not compatible with the domain being processed";
951   case 9999:
952     return error + "Abnormal error occured, contact hotline";
953   case 23600:
954     return error + "Not enough memory for the face table";
955   case 23601:
956     return error + "The algorithm cannot run further. "
957       "The surface mesh is probably very bad in terms of quality.";
958   case 23602:
959     return error + "Bad vertex number";
960   }
961   return errorLine;
962 }
963
964 //=============================================================================
965 /*!
966  *Here we are going to use the GHS3D mesher
967  */
968 //=============================================================================
969
970 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh&         theMesh,
971                                 const TopoDS_Shape& theShape)
972 {
973   bool Ok(false);
974   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
975
976   _nbShape = countShape( meshDS, TopAbs_SOLID ); // we count the number of shapes
977
978   // create bounding box for every shape inside the compound
979
980   int iShape = 0;
981   TopoDS_Shape* tabShape;
982   double**      tabBox;
983   tabShape = new TopoDS_Shape[_nbShape];
984   tabBox   = new double*[_nbShape];
985   for (int i=0; i<_nbShape; i++)
986     tabBox[i] = new double[6];
987   Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
988
989   TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
990   for (; expBox.More(); expBox.Next()) {
991     tabShape[iShape] = expBox.Current();
992     Bnd_Box BoundingBox;
993     BRepBndLib::Add(expBox.Current(), BoundingBox);
994     BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
995     tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
996     tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
997     tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
998     iShape++;
999   }
1000
1001   // a unique working file name
1002   // to avoid access to the same files by eg different users
1003   TCollection_AsciiString aGenericName
1004     = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1005
1006   TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1007   TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1008   aFacesFileName  = aGenericName + ".faces";  // in faces
1009   aPointsFileName = aGenericName + ".points"; // in points
1010   aResultFileName = aGenericName + ".noboite";// out points and volumes
1011   aBadResFileName = aGenericName + ".boite";  // out bad result
1012   aBbResFileName  = aGenericName + ".bb";     // out vertex stepsize
1013   aLogFileName    = aGenericName + ".log";    // log
1014
1015   // -----------------
1016   // make input files
1017   // -----------------
1018
1019   ofstream aFacesFile  ( aFacesFileName.ToCString()  , ios::out);
1020   ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1021
1022   Ok =
1023 #ifdef WIN32
1024     aFacesFile->is_open() && aPointsFile->is_open();
1025 #else
1026     aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1027 #endif
1028   if (!Ok) {
1029     INFOS( "Can't write into " << aFacesFileName);
1030     return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1031   }
1032   map <int,int> aSmdsToGhs3dIdMap;
1033   map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1034
1035   Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap ) &&
1036        writeFaces ( aFacesFile,  meshDS, aSmdsToGhs3dIdMap );
1037
1038   aFacesFile.close();
1039   aPointsFile.close();
1040
1041   if ( ! Ok ) {
1042     if ( !_keepFiles ) {
1043       OSD_File( aFacesFileName ).Remove();
1044       OSD_File( aPointsFileName ).Remove();
1045     }
1046     return error(COMPERR_BAD_INPUT_MESH);
1047   }
1048   OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1049
1050   // -----------------
1051   // run ghs3d mesher
1052   // -----------------
1053
1054   TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1055   cmd += TCollection_AsciiString(" -f ") + aGenericName;  // file to read
1056   cmd += TCollection_AsciiString(" 1>" ) + aLogFileName;  // dump into file
1057
1058   cout << endl;
1059   cout << "Ghs3d execution..." << endl;
1060   cout << cmd << endl;
1061
1062   system( cmd.ToCString() ); // run
1063
1064   cout << endl;
1065   cout << "End of Ghs3d execution !" << endl;
1066
1067   // --------------
1068   // read a result
1069   // --------------
1070
1071   // Mapping the result file
1072
1073   int fileOpen;
1074   fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1075   if ( fileOpen < 0 ) {
1076     cout << endl;
1077     cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << endl;
1078     Ok = false;
1079   }
1080   else
1081     Ok = readResultFile( fileOpen, meshDS, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap );
1082
1083   // ---------------------
1084   // remove working files
1085   // ---------------------
1086
1087   if ( Ok )
1088   {
1089     if ( !_keepFiles )
1090       OSD_File( aLogFileName ).Remove();
1091   }
1092   else if ( OSD_File( aLogFileName ).Size() > 0 )
1093   {
1094     INFOS( "GHS3D Error, see the " << aLogFileName.ToCString() << " file" );
1095
1096     // get problem description from the log file
1097     SMESH_Comment comment;
1098     TCollection_AsciiString foundLine;
1099     if ( findLineContaing( "has expired",aLogFileName,foundLine) &&
1100          foundLine.Search("Licence") >= 0)
1101     {
1102       foundLine.LeftAdjust();
1103       comment << foundLine;
1104     }
1105     if ( findLineContaing( "%% ERROR",aLogFileName,foundLine) ||
1106          findLineContaing( " ERR ",aLogFileName,foundLine))
1107     {
1108       foundLine.LeftAdjust();
1109       comment << translateError( foundLine );
1110     }
1111     else if ( findLineContaing( "%% NO SAVING OPERATION",aLogFileName,foundLine))
1112     {
1113       comment << "Too many elements generated for a trial version.\n";
1114     }
1115     if ( comment.empty() )
1116       comment << "See " << aLogFileName << " for problem description";
1117     else
1118       comment << "\nSee " << aLogFileName << " for more information";
1119     error(COMPERR_ALGO_FAILED, comment);
1120   }
1121   else
1122   {
1123     // the log file is empty
1124     OSD_File( aLogFileName ).Remove();
1125     INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1126     error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1127   }
1128
1129   if ( !_keepFiles ) {
1130     OSD_File( aFacesFileName ).Remove();
1131     OSD_File( aPointsFileName ).Remove();
1132     OSD_File( aResultFileName ).Remove();
1133     OSD_File( aBadResFileName ).Remove();
1134     OSD_File( aBbResFileName ).Remove();
1135   }
1136   cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1137   if ( !Ok )
1138     cout << "not ";
1139   cout << "treated !" << endl;
1140   cout << endl;
1141
1142   _nbShape = 0;    // re-initializing _nbShape for the next Compute() method call
1143   delete [] tabShape;
1144   delete [] tabBox;
1145
1146   return Ok;
1147 }
1148
1149 //=============================================================================
1150 /*!
1151  *Here we are going to use the GHS3D mesher w/o geometry
1152  */
1153 //=============================================================================
1154 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh&         theMesh,
1155                                 SMESH_MesherHelper* aHelper)
1156 {
1157   MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1158
1159   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1160   TopoDS_Shape theShape = aHelper->GetSubShape();
1161
1162   // a unique working file name
1163   // to avoid access to the same files by eg different users
1164   TCollection_AsciiString aGenericName
1165     = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1166
1167   TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1168   TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1169   aFacesFileName  = aGenericName + ".faces";  // in faces
1170   aPointsFileName = aGenericName + ".points"; // in points
1171   aResultFileName = aGenericName + ".noboite";// out points and volumes
1172   aBadResFileName = aGenericName + ".boite";  // out bad result
1173   aBbResFileName  = aGenericName + ".bb";     // out vertex stepsize
1174   aLogFileName    = aGenericName + ".log";    // log
1175
1176   // -----------------
1177   // make input files
1178   // -----------------
1179
1180   ofstream aFacesFile  ( aFacesFileName.ToCString()  , ios::out);
1181   ofstream aPointsFile  ( aPointsFileName.ToCString()  , ios::out);
1182   bool Ok =
1183 #ifdef WIN32
1184     aFacesFile->is_open() && aPointsFile->is_open();
1185 #else
1186     aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1187 #endif
1188
1189   if (!Ok)
1190     return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1191
1192   vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1193
1194   Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1195         writePoints( aPointsFile, meshDS, aNodeByGhs3dId));
1196   
1197   aFacesFile.close();
1198   aPointsFile.close();
1199   
1200   if ( ! Ok ) {
1201     if ( !_keepFiles ) {
1202       OSD_File( aFacesFileName ).Remove();
1203       OSD_File( aPointsFileName ).Remove();
1204     }
1205     return error(COMPERR_BAD_INPUT_MESH);
1206   }
1207   OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1208
1209   // -----------------
1210   // run ghs3d mesher
1211   // -----------------
1212
1213   TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1214   cmd += TCollection_AsciiString(" -f ") + aGenericName;  // file to read
1215   cmd += TCollection_AsciiString(" 1>" ) + aLogFileName;  // dump into file
1216   
1217   system( cmd.ToCString() ); // run
1218
1219   // --------------
1220   // read a result
1221   // --------------
1222   int fileOpen;
1223   fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1224   if ( fileOpen < 0 ) {
1225     cout << endl;
1226     cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
1227     cout << endl;
1228     Ok = false;
1229   }
1230   else {
1231     Ok = readResultFile( fileOpen, meshDS, theShape ,aNodeByGhs3dId );
1232   }
1233   
1234   // ---------------------
1235   // remove working files
1236   // ---------------------
1237
1238   if ( Ok )
1239   {
1240     if ( !_keepFiles )
1241       OSD_File( aLogFileName ).Remove();
1242   }
1243   else if ( OSD_File( aLogFileName ).Size() > 0 )
1244   {
1245     INFOS( "GHS3D Error, see the " << aLogFileName.ToCString() << " file" );
1246
1247     // get problem description from the log file
1248     SMESH_Comment comment;
1249     TCollection_AsciiString foundLine;
1250     if ( findLineContaing( "has expired",aLogFileName,foundLine) &&
1251          foundLine.Search("Licence") >= 0)
1252     {
1253       foundLine.LeftAdjust();
1254       comment << foundLine;
1255     }
1256     if ( findLineContaing( "%% ERROR",aLogFileName,foundLine) ||
1257          findLineContaing( " ERR ",aLogFileName,foundLine))
1258     {
1259       foundLine.LeftAdjust();
1260       comment << translateError( foundLine );
1261     }
1262     else if ( findLineContaing( "%% NO SAVING OPERATION",aLogFileName,foundLine))
1263     {
1264       comment << "Too many elements generated for a trial version.\n";
1265     }
1266     if ( comment.empty() )
1267       comment << "See " << aLogFileName << " for problem description";
1268     else
1269       comment << "\nSee " << aLogFileName << " for more information";
1270     error(COMPERR_ALGO_FAILED, comment);
1271   }
1272   else {
1273     // the log file is empty
1274     OSD_File( aLogFileName ).Remove();
1275     INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1276     error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1277   }
1278
1279   if ( !_keepFiles )
1280   {
1281     OSD_File( aFacesFileName ).Remove();
1282     OSD_File( aPointsFileName ).Remove();
1283     OSD_File( aResultFileName ).Remove();
1284     OSD_File( aBadResFileName ).Remove();
1285     OSD_File( aBbResFileName ).Remove();
1286   }
1287   
1288   return Ok;
1289 }