]> SALOME platform Git repositories - plugins/ghs3dplugin.git/blob - src/GHS3DPlugin_GHS3D.cxx
Salome HOME
0020042: EDF 864 SMESH: Mesh of holes (GHS3D/BLSurf)
[plugins/ghs3dplugin.git] / src / GHS3DPlugin_GHS3D.cxx
1 //  Copyright (C) 2004-2008  CEA/DEN, EDF 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 // $Header$
25 //=============================================================================
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 #include "SMESH_MeshEditor.hxx"
37
38 #include "SMDS_MeshElement.hxx"
39 #include "SMDS_MeshNode.hxx"
40 #include "SMDS_FaceOfNodes.hxx"
41 #include "SMDS_VolumeOfNodes.hxx"
42
43 #include <BRepAdaptor_Surface.hxx>
44 #include <BRepBndLib.hxx>
45 #include <BRepClass3d_SolidClassifier.hxx>
46 #include <BRepTools.hxx>
47 #include <BRep_Tool.hxx>
48 #include <Bnd_Box.hxx>
49 #include <GeomAPI_ProjectPointOnSurf.hxx>
50 #include <OSD_File.hxx>
51 #include <Precision.hxx>
52 #include <Quantity_Parameter.hxx>
53 #include <Standard_ErrorHandler.hxx>
54 #include <Standard_Failure.hxx>
55 #include <TopExp.hxx>
56 #include <TopExp_Explorer.hxx>
57 #include <TopTools_IndexedMapOfShape.hxx>
58 #include <TopTools_ListIteratorOfListOfShape.hxx>
59 #include <TopoDS.hxx>
60 //#include <BRepClass_FaceClassifier.hxx>
61 //#include <BRepGProp.hxx>
62 //#include <GProp_GProps.hxx>
63
64 #include "utilities.h"
65
66 #ifndef WIN32
67 #include <sys/sysinfo.h>
68 #endif
69
70 //#include <Standard_Stream.hxx>
71
72
73 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
74
75 #ifdef _DEBUG_
76 #define DUMP(txt) \
77 //  cout << txt
78 #else
79 #define DUMP(txt)
80 #endif
81
82 extern "C"
83 {
84 #ifndef WNT
85 #include <unistd.h>
86 #include <sys/mman.h>
87 #endif
88 #include <sys/stat.h>
89 #include <fcntl.h>
90 }
91
92 #define HOLE_ID -1
93
94 //=============================================================================
95 /*!
96  *  
97  */
98 //=============================================================================
99
100 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
101   : SMESH_3D_Algo(hypId, studyId, gen)
102 {
103   MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
104   _name = "GHS3D_3D";
105   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
106   _onlyUnaryInput = false; // Compute() will be called on a compound of solids
107   _iShape=0;
108   _nbShape=0;
109   _compatibleHypothesis.push_back("GHS3D_Parameters");
110   _requireShape = false; // can work without shape
111 }
112
113 //=============================================================================
114 /*!
115  *  
116  */
117 //=============================================================================
118
119 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
120 {
121   MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
122 }
123
124 //=============================================================================
125 /*!
126  *  
127  */
128 //=============================================================================
129
130 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh&         aMesh,
131                                           const TopoDS_Shape& aShape,
132                                           Hypothesis_Status&  aStatus )
133 {
134   aStatus = SMESH_Hypothesis::HYP_OK;
135
136   // there is only one compatible Hypothesis so far
137   _hyp = 0;
138   _keepFiles = false;
139   const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
140   if ( !hyps.empty() )
141     _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
142   if ( _hyp )
143     _keepFiles = _hyp->GetKeepFiles();
144
145   return true;
146 }
147
148 //=======================================================================
149 //function : findShape
150 //purpose  : 
151 //=======================================================================
152
153 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
154                               TopoDS_Shape        aShape,
155                               const TopoDS_Shape  shape[],
156                               double**            box,
157                               const int           nShape,
158                               TopAbs_State *      state = 0)
159 {
160   gp_XYZ aPnt(0,0,0);
161   int j, iShape, nbNode = 4;
162
163   for ( j=0; j<nbNode; j++ )
164     aPnt += gp_XYZ( aNode[j]->X(), aNode[j]->Y(), aNode[j]->Z() );
165   aPnt /= nbNode;
166
167   BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
168   if (state) *state = SC.State();
169   if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
170     for (iShape = 0; iShape < nShape; iShape++) {
171       aShape = shape[iShape];
172       if ( !( aPnt.X() < box[iShape][0] || box[iShape][1] < aPnt.X() ||
173               aPnt.Y() < box[iShape][2] || box[iShape][3] < aPnt.Y() ||
174               aPnt.Z() < box[iShape][4] || box[iShape][5] < aPnt.Z()) ) {
175         BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
176         if (state) *state = SC.State();
177         if (SC.State() == TopAbs_IN)
178           break;
179       }
180     }
181   }
182   return aShape;
183 }
184
185 //=======================================================================
186 //function : readMapIntLine
187 //purpose  : 
188 //=======================================================================
189
190 static char* readMapIntLine(char* ptr, int tab[]) {
191   long int intVal;
192   cout << endl;
193
194   for ( int i=0; i<17; i++ ) {
195     intVal = strtol(ptr, &ptr, 10);
196     if ( i < 3 )
197       tab[i] = intVal;
198   }
199   return ptr;
200 }
201
202 //=======================================================================
203 //function : countShape
204 //purpose  :
205 //=======================================================================
206
207 template < class Mesh, class Shape >
208 static int countShape( Mesh* mesh, Shape shape ) {
209   TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
210   int nbShape = 0;
211   for ( ; expShape.More(); expShape.Next() ) {
212       nbShape++;
213   }
214   return nbShape;
215 }
216
217 //=======================================================================
218 //function : writeFaces
219 //purpose  : 
220 //=======================================================================
221
222 static bool writeFaces (ofstream &            theFile,
223                         SMESHDS_Mesh *        theMesh,
224                         const map <int,int> & theSmdsToGhs3dIdMap)
225 {
226   // record structure:
227   //
228   // NB_ELEMS DUMMY_INT
229   // Loop from 1 to NB_ELEMS
230   // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
231
232   int nbShape = countShape( theMesh, TopAbs_FACE );
233
234   int *tabID;             tabID    = new int[nbShape];
235   TopoDS_Shape *tabShape; tabShape = new TopoDS_Shape[nbShape];
236   TopoDS_Shape aShape;
237   SMESHDS_SubMesh* theSubMesh;
238   const SMDS_MeshElement* aFace;
239   const char* space    = "  ";
240   const int   dummyint = 0;
241   map<int,int>::const_iterator itOnMap;
242   SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
243   int shapeID, nbNodes, aSmdsID;
244   bool idFound;
245
246   cout << "    " << theMesh->NbFaces() << " shapes of 2D dimension" << endl;
247   cout << endl;
248
249   theFile << space << theMesh->NbFaces() << space << dummyint << endl;
250
251   TopExp_Explorer expface( theMesh->ShapeToMesh(), TopAbs_FACE );
252   for ( int i = 0; expface.More(); expface.Next(), i++ ) {
253     tabID[i] = 0;
254     aShape   = expface.Current();
255     shapeID  = theMesh->ShapeToIndex( aShape );
256     idFound  = false;
257     for ( int j=0; j<=i; j++) {
258       if ( shapeID == tabID[j] ) {
259         idFound = true;
260         break;
261       }
262     }
263     if ( ! idFound ) {
264       tabID[i]    = shapeID;
265       tabShape[i] = aShape;
266     }
267   }
268   for ( int i =0; i < nbShape; i++ ) {
269     if ( tabID[i] != 0 ) {
270       aShape      = tabShape[i];
271       shapeID     = tabID[i];
272       theSubMesh  = theMesh->MeshElements( aShape );
273       if ( !theSubMesh ) continue;
274       itOnSubMesh = theSubMesh->GetElements();
275       while ( itOnSubMesh->more() ) {
276         aFace   = itOnSubMesh->next();
277         nbNodes = aFace->NbNodes();
278
279         theFile << space << nbNodes;
280
281         itOnSubFace = aFace->nodesIterator();
282         while ( itOnSubFace->more() ) {
283           // find GHS3D ID
284           aSmdsID = itOnSubFace->next()->GetID();
285           itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
286           ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
287
288           theFile << space << (*itOnMap).second;
289         }
290
291         // (NB_NODES + 1) times: DUMMY_INT
292         for ( int j=0; j<=nbNodes; j++)
293           theFile << space << dummyint;
294
295         theFile << endl;
296       }
297     }
298   }
299
300   delete [] tabID;
301   delete [] tabShape;
302
303   return true;
304 }
305
306 //=======================================================================
307 //function : writeFaces
308 //purpose  : Write Faces in case if generate 3D mesh w/o geometry
309 //=======================================================================
310
311 static bool writeFaces (ofstream &            theFile,
312                         SMESHDS_Mesh *        theMesh,
313                         vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
314 {
315   // record structure:
316   //
317   // NB_ELEMS DUMMY_INT
318   // Loop from 1 to NB_ELEMS
319   //   NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
320
321
322   int nbFaces = 0;
323   list< const SMDS_MeshElement* > faces;
324   list< const SMDS_MeshElement* >::iterator f;
325   map< const SMDS_MeshNode*,int >::iterator it;
326   SMDS_ElemIteratorPtr nodeIt;
327   const SMDS_MeshElement* elem;
328   int nbNodes;
329
330   const char* space    = "  ";
331   const int   dummyint = 0;
332
333   //get all faces from mesh
334   SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
335   while ( eIt->more() ) {
336     const SMDS_MeshElement* elem = eIt->next();
337     if ( !elem )
338       return false;
339     faces.push_back( elem );
340     nbFaces++;
341   }
342
343   if ( nbFaces == 0 )
344     return false;
345   
346   cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
347
348   // NB_ELEMS DUMMY_INT
349   theFile << space << nbFaces << space << dummyint << endl;
350
351   // Loop from 1 to NB_ELEMS
352
353   map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
354   f = faces.begin();
355   for ( ; f != faces.end(); ++f )
356   {
357     // NB_NODES PER FACE
358     elem = *f;
359     nbNodes = elem->NbNodes();
360     theFile << space << nbNodes;
361
362     // NODE_NB_1 NODE_NB_2 ...
363     nodeIt = elem->nodesIterator();
364     while ( nodeIt->more() )
365     {
366       // find GHS3D ID
367       const SMDS_MeshNode* node = castToNode( nodeIt->next() );
368       int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
369       it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
370       theFile << space << it->second;
371     }
372
373     // (NB_NODES + 1) times: DUMMY_INT
374     for ( int i=0; i<=nbNodes; i++)
375       theFile << space << dummyint;
376
377     theFile << endl;
378   }
379
380   // put nodes to theNodeByGhs3dId vector
381   theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
382   map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
383   for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
384   {
385     theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
386   }
387
388   return true;
389 }
390
391 //=======================================================================
392 //function : writePoints
393 //purpose  : 
394 //=======================================================================
395
396 static bool writePoints (ofstream &                       theFile,
397                          SMESHDS_Mesh *                   theMesh,
398                          map <int,int> &                  theSmdsToGhs3dIdMap,
399                          map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap)
400 {
401   // record structure:
402   //
403   // NB_NODES
404   // Loop from 1 to NB_NODES
405   //   X Y Z DUMMY_INT
406
407   int nbNodes = theMesh->NbNodes();
408   if ( nbNodes == 0 )
409     return false;
410
411   const char* space    = "  ";
412   const int   dummyint = 0;
413
414   int aGhs3dID = 1;
415   SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
416   const SMDS_MeshNode* node;
417
418   // NB_NODES
419   theFile << space << nbNodes << endl;
420   cout << endl;
421   cout << "The initial 2D mesh contains :" << endl;
422   cout << "    " << nbNodes << " vertices" << endl;
423
424   // Loop from 1 to NB_NODES
425
426   while ( it->more() )
427   {
428     node = it->next();
429     theSmdsToGhs3dIdMap.insert( make_pair( node->GetID(), aGhs3dID ));
430     theGhs3dIdToNodeMap.insert( make_pair( aGhs3dID, node ));
431     aGhs3dID++;
432
433     // X Y Z DUMMY_INT
434     theFile
435       << space << node->X()
436       << space << node->Y()
437       << space << node->Z()
438       << space << dummyint;
439
440     theFile << endl;
441   }
442
443   return true;
444 }
445
446 //=======================================================================
447 //function : writePoints
448 //purpose  : 
449 //=======================================================================
450
451 static bool writePoints (ofstream &                            theFile,
452                          SMESHDS_Mesh *                        theMesh,
453                          const vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
454 {
455   // record structure:
456   //
457   // NB_NODES
458   // Loop from 1 to NB_NODES
459   //   X Y Z DUMMY_INT
460
461   //int nbNodes = theMesh->NbNodes();
462   int nbNodes = theNodeByGhs3dId.size();
463   if ( nbNodes == 0 )
464     return false;
465
466   const char* space    = "  ";
467   const int   dummyint = 0;
468
469   const SMDS_MeshNode* node;
470
471   // NB_NODES
472   theFile << space << nbNodes << endl;
473   cout << nbNodes << " nodes" << endl;
474
475   // Loop from 1 to NB_NODES
476
477   vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
478   vector<const SMDS_MeshNode*>::const_iterator after  = theNodeByGhs3dId.end();
479   for ( ; nodeIt != after; ++nodeIt )
480   {
481     node = *nodeIt;
482
483     // X Y Z DUMMY_INT
484     theFile
485       << space << node->X()
486       << space << node->Y()
487       << space << node->Z()
488       << space << dummyint;
489
490     theFile << endl;
491   }
492
493   return true;
494 }
495
496 //=======================================================================
497 //function : findShapeID
498 //purpose  : find the solid corresponding to GHS3D sub-domain following
499 //           the technique proposed in GHS3D manual in chapter
500 //           "B.4 Subdomain (sub-region) assignment"
501 //=======================================================================
502
503 static int findShapeID(SMESH_Mesh&          mesh,
504                        const SMDS_MeshNode* node1,
505                        const SMDS_MeshNode* node2,
506                        const SMDS_MeshNode* node3,
507                        const bool           toMeshHoles)
508 {
509   const int invalidID = 0;
510   SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
511
512   // face th enodes belong to
513   const SMDS_MeshElement * face = meshDS->FindFace(node1,node2,node3);
514   if ( !face )
515     return invalidID;
516
517   // geom face the face assigned to
518   SMESH_MeshEditor editor(&mesh);
519   int geomFaceID = editor.FindShape( face );
520   if ( !geomFaceID )
521     return invalidID;
522   TopoDS_Shape shape = meshDS->IndexToShape( geomFaceID );
523   if ( shape.IsNull() || shape.ShapeType() != TopAbs_FACE )
524     return invalidID;
525   TopoDS_Face geomFace = TopoDS::Face( shape );
526
527   // solids bounded by geom face
528   TopTools_IndexedMapOfShape solids, shells;
529   TopTools_ListIteratorOfListOfShape ansIt = mesh.GetAncestors(geomFace);
530   for ( ; ansIt.More(); ansIt.Next() ) {
531     switch ( ansIt.Value().ShapeType() ) {
532     case TopAbs_SOLID:
533       solids.Add( ansIt.Value() ); break;
534     case TopAbs_SHELL:
535       shells.Add( ansIt.Value() ); break;
536     default:;
537     }
538   }
539   // analyse found solids
540   if ( solids.Extent() == 0 || shells.Extent() == 0)
541     return invalidID;
542
543   const TopoDS_Solid& solid1 = TopoDS::Solid( solids(1) );
544   if ( solids.Extent() == 1 )
545   {
546     if ( toMeshHoles )
547       return meshDS->ShapeToIndex( solid1 );
548
549     // - Are we at a hole boundary face?
550     if ( shells(1).IsSame( BRepTools::OuterShell( solid1 )) )
551     { // - No, but maybe a hole is bound by two shapes? Does shells(1) touches another shell?
552       bool touch = false;
553       TopExp_Explorer eExp( shells(1), TopAbs_EDGE );
554       // check if any edge of shells(1) belongs to another shell
555       for ( ; eExp.More() && !touch; eExp.Next() ) {
556         ansIt = mesh.GetAncestors( eExp.Current() );
557         for ( ; ansIt.More() && !touch; ansIt.Next() ) {
558           if ( ansIt.Value().ShapeType() == TopAbs_SHELL )
559             touch = ( !ansIt.Value().IsSame( shells(1) ));
560         }
561       }
562       if (!touch)
563         return meshDS->ShapeToIndex( solid1 );
564     }
565   }
566   // find orientation of geom face within the first solid
567   TopExp_Explorer fExp( solid1, TopAbs_FACE );
568   for ( ; fExp.More(); fExp.Next() )
569     if ( geomFace.IsSame( fExp.Current() )) {
570       geomFace = TopoDS::Face( fExp.Current() );
571       break;
572     }
573   if ( !fExp.More() )
574     return invalidID; // face not found
575
576   // find UV of node1 on geomFace
577   SMESH_MesherHelper helper( mesh );
578   gp_XY uv = helper.GetNodeUV( geomFace, node1 );
579
580   // check that uv is correct
581   gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
582   double tol = BRep_Tool::Tolerance( geomFace );
583   BRepAdaptor_Surface surface( geomFace );
584   if ( node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol ) {
585     // project node1 onto geomFace to get right UV
586     GeomAPI_ProjectPointOnSurf projector( node1Pnt, surface.Surface().Surface() );
587     if ( !projector.IsDone() || projector.NbPoints() < 1 )
588       return invalidID;
589     Quantity_Parameter U,V;
590     projector.LowerDistanceParameters(U,V);
591     uv = gp_XY( U,V );
592   }
593   // normale to face at node1
594   gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
595   gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
596   gp_Vec vec12( node1Pnt, node2Pnt );
597   gp_Vec vec13( node1Pnt, node3Pnt );
598   gp_Vec meshNormal = vec12 ^ vec13;
599   if ( meshNormal.SquareMagnitude() < DBL_MIN )
600     return invalidID;
601   
602   // normale to geomFace at UV
603   gp_Vec du, dv;
604   surface.D1( uv.X(), uv.Y(), node1Pnt, du, dv );
605   gp_Vec geomNormal = du ^ dv;
606   if ( geomNormal.SquareMagnitude() < DBL_MIN )
607     return findShapeID( mesh, node2, node3, node1, toMeshHoles );
608   if ( geomFace.Orientation() == TopAbs_REVERSED )
609     geomNormal.Reverse();
610
611   // compare normals
612   bool isReverse = ( meshNormal * geomNormal ) < 0;
613   if ( !isReverse )
614     return meshDS->ShapeToIndex( solid1 );
615
616   if ( solids.Extent() == 1 )
617     return HOLE_ID; // we are inside a hole
618   else
619     return meshDS->ShapeToIndex( solids(2) );
620 }
621                        
622 //=======================================================================
623 //function : readResultFile
624 //purpose  : 
625 //=======================================================================
626
627 static bool readResultFile(const int                       fileOpen,
628                            SMESH_Mesh&                     theMesh,
629                            TopoDS_Shape                    tabShape[],
630                            double**                        tabBox,
631                            const int                       nbShape,
632                            map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
633                            bool                            toMeshHoles)
634 {
635   struct stat status;
636   size_t      length;
637
638   char *ptr, *mapPtr;
639   char *tetraPtr;
640   char *shapePtr;
641
642   SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
643
644   int fileStat;
645   int nbElems, nbNodes, nbInputNodes;
646   int nodeId/*, triangleId*/;
647   int nbTriangle;
648   int ID, shapeID, ghs3dShapeID;
649   int IdShapeRef = 1;
650   int compoundID =
651     nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
652
653   int *tab, *tabID, *nodeID, *nodeAssigne;
654   double *coord;
655   const SMDS_MeshNode **node;
656
657   tab    = new int[3];
658   //tabID  = new int[nbShape];
659   nodeID = new int[4];
660   coord  = new double[3];
661   node   = new const SMDS_MeshNode*[4];
662
663   TopoDS_Shape aSolid;
664   SMDS_MeshNode * aNewNode;
665   map <int,const SMDS_MeshNode*>::iterator itOnNode;
666   SMDS_MeshElement* aTet;
667 #ifdef _DEBUG_
668   set<int> shapeIDs;
669 #endif
670
671   // Read the file state
672   fileStat = fstat(fileOpen, &status);
673   length   = status.st_size;
674
675   // Mapping the result file into memory
676   ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
677   mapPtr = ptr;
678
679   ptr      = readMapIntLine(ptr, tab);
680   tetraPtr = ptr;
681
682   nbElems      = tab[0];
683   nbNodes      = tab[1];
684   nbInputNodes = tab[2];
685
686   nodeAssigne = new int[ nbNodes+1 ];
687
688   if (nbShape > 0)
689     aSolid = tabShape[0];
690
691   // Reading the nodeId
692   for (int i=0; i < 4*nbElems; i++)
693     nodeId = strtol(ptr, &ptr, 10);
694
695   // Reading the nodeCoor and update the nodeMap
696   for (int iNode=1; iNode <= nbNodes; iNode++) {
697     for (int iCoor=0; iCoor < 3; iCoor++)
698       coord[ iCoor ] = strtod(ptr, &ptr);
699     nodeAssigne[ iNode ] = 1;
700     if ( iNode > nbInputNodes ) {
701       nodeAssigne[ iNode ] = 0;
702       aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
703       theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
704     }
705   }
706
707   // Reading the number of triangles which corresponds to the number of sub-domains
708   nbTriangle = strtol(ptr, &ptr, 10);
709
710   tabID = new int[nbTriangle];
711   for (int i=0; i < nbTriangle; i++) {
712     tabID[i] = 0;
713     // find the solid corresponding to GHS3D sub-domain following
714     // the technique proposed in GHS3D manual in chapter
715     // "B.4 Subdomain (sub-region) assignment"
716     int nodeId1 = strtol(ptr, &ptr, 10);
717     int nodeId2 = strtol(ptr, &ptr, 10);
718     int nodeId3 = strtol(ptr, &ptr, 10);
719     if ( nbTriangle > 1 ) {
720       const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
721       const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
722       const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
723       try {
724         OCC_CATCH_SIGNALS;
725         tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles );
726 #ifdef _DEBUG_
727         cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << endl;
728 #endif
729       } catch ( Standard_Failure ) {
730       } catch (...) {}
731     }
732   }
733
734   shapePtr = ptr;
735
736   if ( nbTriangle <= nbShape ) // no holes
737     toMeshHoles = true; // not avoid creating tetras in holes
738
739   // Associating the tetrahedrons to the shapes
740   shapeID = compoundID;
741   for (int iElem = 0; iElem < nbElems; iElem++) {
742     for (int iNode = 0; iNode < 4; iNode++) {
743       ID = strtol(tetraPtr, &tetraPtr, 10);
744       itOnNode = theGhs3dIdToNodeMap.find(ID);
745       node[ iNode ] = itOnNode->second;
746       nodeID[ iNode ] = ID;
747     }
748     // We always run GHS3D with "to mesh holes'==TRUE but we must not create
749     // tetras within holes depending on hypo option,
750     // so we first check if aTet is inside a hole and then create it 
751     //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
752     if ( nbTriangle > 1 ) {
753       shapeID = HOLE_ID; // negative shapeID means not to create tetras if !toMeshHoles
754       ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
755       if ( tabID[ ghs3dShapeID ] == 0 ) {
756         TopAbs_State state;
757         aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
758         if ( toMeshHoles || state == TopAbs_IN )
759           shapeID = theMeshDS->ShapeToIndex( aSolid );
760         tabID[ ghs3dShapeID ] = shapeID;
761       }
762       else
763         shapeID = tabID[ ghs3dShapeID ];
764     }
765     else if ( nbShape > 1 ) {
766       // Case where nbTriangle == 1 while nbShape == 2 encountered
767       // with compound of 2 boxes and "To mesh holes"==False,
768       // so there are no subdomains specified for each tetrahedron.
769       // Try to guess a solid by a node already bound to shape
770       shapeID = 0;
771       for ( int i=0; i<4 && shapeID==0; i++ ) {
772         if ( nodeAssigne[ nodeID[i] ] == 1 &&
773              node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
774              node[i]->GetPosition()->GetShapeId() > 1 )
775         {
776           shapeID = node[i]->GetPosition()->GetShapeId();
777         }
778       }
779       if ( shapeID==0 ) {
780         aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
781         shapeID = theMeshDS->ShapeToIndex( aSolid );
782       }
783     }
784     // set new nodes and tetrahedron onto the shape
785     for ( int i=0; i<4; i++ ) {
786       if ( nodeAssigne[ nodeID[i] ] == 0 ) {
787         if ( shapeID != HOLE_ID )
788           theMeshDS->SetNodeInVolume( node[i], shapeID );
789         nodeAssigne[ nodeID[i] ] = shapeID;
790       }
791     }
792     if ( toMeshHoles || shapeID != HOLE_ID ) {
793       aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
794       theMeshDS->SetMeshElementOnShape( aTet, shapeID );
795     }
796 #ifdef _DEBUG_
797     shapeIDs.insert( shapeID );
798 #endif
799   }
800   // Remove nodes of tetras inside holes if !toMeshHoles
801   if ( !toMeshHoles ) {
802     itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
803     for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
804       ID = itOnNode->first;
805       if ( nodeAssigne[ ID ] == HOLE_ID )
806         theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
807     }
808   }
809
810   if ( nbElems )
811     cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
812   munmap(mapPtr, length);
813   close(fileOpen);
814
815   delete [] tab;
816   delete [] tabID;
817   delete [] nodeID;
818   delete [] coord;
819   delete [] node;
820   delete [] nodeAssigne;
821
822 #ifdef _DEBUG_
823   if ( shapeIDs.size() != nbShape ) {
824     cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << endl;
825     for (int i=0; i<nbShape; i++) {
826       shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
827       if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
828         cout << "  Solid #" << shapeID << " not found" << endl;
829     }
830   }
831 #endif
832
833   return true;
834 }
835
836 //=======================================================================
837 //function : readResultFile
838 //purpose  : 
839 //=======================================================================
840
841 static bool readResultFile(const int                      fileOpen,
842                            SMESHDS_Mesh*                  theMeshDS,
843                            TopoDS_Shape                   aSolid,
844                            vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
845
846   struct stat  status;
847   size_t       length;
848
849   char *ptr, *mapPtr;
850   char *tetraPtr;
851   char *shapePtr;
852
853   int fileStat;
854   int nbElems, nbNodes, nbInputNodes;
855   int nodeId, triangleId;
856   int nbTriangle;
857   int ID, shapeID;
858
859   int *tab;
860   double *coord;
861   const SMDS_MeshNode **node;
862
863   tab   = new int[3];
864   coord = new double[3];
865   node  = new const SMDS_MeshNode*[4];
866
867   SMDS_MeshNode * aNewNode;
868   map <int,const SMDS_MeshNode*>::iterator IdNode;
869   SMDS_MeshElement* aTet;
870
871   // Read the file state
872   fileStat = fstat(fileOpen, &status);
873   length   = status.st_size;
874
875   // Mapping the result file into memory
876   ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
877   mapPtr = ptr;
878
879   ptr      = readMapIntLine(ptr, tab);
880   tetraPtr = ptr;
881
882   nbElems      = tab[0];
883   nbNodes      = tab[1];
884   nbInputNodes = tab[2];
885
886   theNodeByGhs3dId.resize( nbNodes );
887
888   // Reading the nodeId
889   for (int i=0; i < 4*nbElems; i++)
890     nodeId = strtol(ptr, &ptr, 10);
891
892   // Reading the nodeCoor and update the nodeMap
893   shapeID = theMeshDS->ShapeToIndex( aSolid );
894   for (int iNode=0; iNode < nbNodes; iNode++) {
895     for (int iCoor=0; iCoor < 3; iCoor++)
896       coord[ iCoor ] = strtod(ptr, &ptr);
897     if ((iNode+1) > nbInputNodes) {
898       aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
899       theMeshDS->SetNodeInVolume( aNewNode, shapeID );
900       theNodeByGhs3dId[ iNode ] = aNewNode;
901     }
902   }
903
904   // Reading the triangles
905   nbTriangle = strtol(ptr, &ptr, 10);
906
907   for (int i=0; i < 3*nbTriangle; i++)
908     triangleId = strtol(ptr, &ptr, 10);
909
910   shapePtr = ptr;
911
912   // Associating the tetrahedrons to the shapes
913   for (int iElem = 0; iElem < nbElems; iElem++) {
914     for (int iNode = 0; iNode < 4; iNode++) {
915       ID = strtol(tetraPtr, &tetraPtr, 10);
916       node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
917     }
918     aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
919     shapeID = theMeshDS->ShapeToIndex( aSolid );
920     theMeshDS->SetMeshElementOnShape( aTet, shapeID );
921   }
922   if ( nbElems )
923     cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
924   munmap(mapPtr, length);
925   close(fileOpen);
926
927   delete [] tab;
928   delete [] coord;
929   delete [] node;
930
931   return true;
932 }
933
934 //=============================================================================
935 /*!
936  *Here we are going to use the GHS3D mesher
937  */
938 //=============================================================================
939
940 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh&         theMesh,
941                                 const TopoDS_Shape& theShape)
942 {
943   bool Ok(false);
944   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
945
946   _nbShape = countShape( meshDS, TopAbs_SOLID ); // we count the number of shapes
947
948   // create bounding box for every shape inside the compound
949
950   int iShape = 0;
951   TopoDS_Shape* tabShape;
952   double**      tabBox;
953   tabShape = new TopoDS_Shape[_nbShape];
954   tabBox   = new double*[_nbShape];
955   for (int i=0; i<_nbShape; i++)
956     tabBox[i] = new double[6];
957   Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
958
959   TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
960   for (; expBox.More(); expBox.Next()) {
961     tabShape[iShape] = expBox.Current();
962     Bnd_Box BoundingBox;
963     BRepBndLib::Add(expBox.Current(), BoundingBox);
964     BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
965     tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
966     tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
967     tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
968     iShape++;
969   }
970
971   // a unique working file name
972   // to avoid access to the same files by eg different users
973   TCollection_AsciiString aGenericName
974     = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
975
976   TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
977   TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
978   aFacesFileName  = aGenericName + ".faces";  // in faces
979   aPointsFileName = aGenericName + ".points"; // in points
980   aResultFileName = aGenericName + ".noboite";// out points and volumes
981   aBadResFileName = aGenericName + ".boite";  // out bad result
982   aBbResFileName  = aGenericName + ".bb";     // out vertex stepsize
983   aLogFileName    = aGenericName + ".log";    // log
984
985   // -----------------
986   // make input files
987   // -----------------
988
989   ofstream aFacesFile  ( aFacesFileName.ToCString()  , ios::out);
990   ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
991
992   Ok =
993 #ifdef WIN32
994     aFacesFile->is_open() && aPointsFile->is_open();
995 #else
996     aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
997 #endif
998   if (!Ok) {
999     INFOS( "Can't write into " << aFacesFileName);
1000     return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1001   }
1002   map <int,int> aSmdsToGhs3dIdMap;
1003   map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1004
1005   Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap ) &&
1006        writeFaces ( aFacesFile,  meshDS, aSmdsToGhs3dIdMap );
1007
1008   aFacesFile.close();
1009   aPointsFile.close();
1010
1011   if ( ! Ok ) {
1012     if ( !_keepFiles ) {
1013       OSD_File( aFacesFileName ).Remove();
1014       OSD_File( aPointsFileName ).Remove();
1015     }
1016     return error(COMPERR_BAD_INPUT_MESH);
1017   }
1018   OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1019
1020   // -----------------
1021   // run ghs3d mesher
1022   // -----------------
1023
1024   TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1025   cmd += TCollection_AsciiString(" -f ") + aGenericName;  // file to read
1026   cmd += TCollection_AsciiString(" 1>" ) + aLogFileName;  // dump into file
1027
1028   cout << endl;
1029   cout << "Ghs3d execution..." << endl;
1030   cout << cmd << endl;
1031
1032   system( cmd.ToCString() ); // run
1033
1034   cout << endl;
1035   cout << "End of Ghs3d execution !" << endl;
1036
1037   // --------------
1038   // read a result
1039   // --------------
1040
1041   // Mapping the result file
1042
1043   int fileOpen;
1044   fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1045   if ( fileOpen < 0 ) {
1046     cout << endl;
1047     cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << endl;
1048     cout << "Log: " << aLogFileName << endl;
1049     Ok = false;
1050   }
1051   else {
1052     bool toMeshHoles =
1053       _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1054     Ok = readResultFile( fileOpen, theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1055                          toMeshHoles );
1056   }
1057
1058   // ---------------------
1059   // remove working files
1060   // ---------------------
1061
1062   if ( Ok )
1063   {
1064     if ( !_keepFiles )
1065       OSD_File( aLogFileName ).Remove();
1066   }
1067   else if ( OSD_File( aLogFileName ).Size() > 0 )
1068   {
1069     // get problem description from the log file
1070     _Ghs2smdsConvertor conv( aGhs3dIdToNodeMap );
1071     storeErrorDescription( aLogFileName, conv );
1072   }
1073   else
1074   {
1075     // the log file is empty
1076     OSD_File( aLogFileName ).Remove();
1077     INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1078     error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1079   }
1080
1081   if ( !_keepFiles ) {
1082     OSD_File( aFacesFileName ).Remove();
1083     OSD_File( aPointsFileName ).Remove();
1084     OSD_File( aResultFileName ).Remove();
1085     OSD_File( aBadResFileName ).Remove();
1086     OSD_File( aBbResFileName ).Remove();
1087   }
1088   cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1089   if ( !Ok )
1090     cout << "not ";
1091   cout << "treated !" << endl;
1092   cout << endl;
1093
1094   _nbShape = 0;    // re-initializing _nbShape for the next Compute() method call
1095   delete [] tabShape;
1096   delete [] tabBox;
1097
1098   return Ok;
1099 }
1100
1101 //=============================================================================
1102 /*!
1103  *Here we are going to use the GHS3D mesher w/o geometry
1104  */
1105 //=============================================================================
1106 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh&         theMesh,
1107                                 SMESH_MesherHelper* aHelper)
1108 {
1109   MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1110
1111   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1112   TopoDS_Shape theShape = aHelper->GetSubShape();
1113
1114   // a unique working file name
1115   // to avoid access to the same files by eg different users
1116   TCollection_AsciiString aGenericName
1117     = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1118
1119   TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1120   TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1121   aFacesFileName  = aGenericName + ".faces";  // in faces
1122   aPointsFileName = aGenericName + ".points"; // in points
1123   aResultFileName = aGenericName + ".noboite";// out points and volumes
1124   aBadResFileName = aGenericName + ".boite";  // out bad result
1125   aBbResFileName  = aGenericName + ".bb";     // out vertex stepsize
1126   aLogFileName    = aGenericName + ".log";    // log
1127
1128   // -----------------
1129   // make input files
1130   // -----------------
1131
1132   ofstream aFacesFile  ( aFacesFileName.ToCString()  , ios::out);
1133   ofstream aPointsFile  ( aPointsFileName.ToCString()  , ios::out);
1134   bool Ok =
1135 #ifdef WIN32
1136     aFacesFile->is_open() && aPointsFile->is_open();
1137 #else
1138     aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1139 #endif
1140
1141   if (!Ok)
1142     return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1143
1144   vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1145
1146   Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1147         writePoints( aPointsFile, meshDS, aNodeByGhs3dId));
1148   
1149   aFacesFile.close();
1150   aPointsFile.close();
1151   
1152   if ( ! Ok ) {
1153     if ( !_keepFiles ) {
1154       OSD_File( aFacesFileName ).Remove();
1155       OSD_File( aPointsFileName ).Remove();
1156     }
1157     return error(COMPERR_BAD_INPUT_MESH);
1158   }
1159   OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1160
1161   // -----------------
1162   // run ghs3d mesher
1163   // -----------------
1164
1165   TCollection_AsciiString cmd =
1166     (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1167   cmd += TCollection_AsciiString(" -f ") + aGenericName;  // file to read
1168   cmd += TCollection_AsciiString(" 1>" ) + aLogFileName;  // dump into file
1169   
1170   system( cmd.ToCString() ); // run
1171
1172   // --------------
1173   // read a result
1174   // --------------
1175   int fileOpen;
1176   fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1177   if ( fileOpen < 0 ) {
1178     cout << endl;
1179     cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
1180     cout << "Log: " << aLogFileName << endl;
1181     cout << endl;
1182     Ok = false;
1183   }
1184   else {
1185     Ok = readResultFile( fileOpen, meshDS, theShape ,aNodeByGhs3dId );
1186   }
1187   
1188   // ---------------------
1189   // remove working files
1190   // ---------------------
1191
1192   if ( Ok )
1193   {
1194     if ( !_keepFiles )
1195       OSD_File( aLogFileName ).Remove();
1196   }
1197   else if ( OSD_File( aLogFileName ).Size() > 0 )
1198   {
1199     // get problem description from the log file
1200     _Ghs2smdsConvertor conv( aNodeByGhs3dId );
1201     storeErrorDescription( aLogFileName, conv );
1202   }
1203   else {
1204     // the log file is empty
1205     OSD_File( aLogFileName ).Remove();
1206     INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1207     error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1208   }
1209
1210   if ( !_keepFiles )
1211   {
1212     OSD_File( aFacesFileName ).Remove();
1213     OSD_File( aPointsFileName ).Remove();
1214     OSD_File( aResultFileName ).Remove();
1215     OSD_File( aBadResFileName ).Remove();
1216     OSD_File( aBbResFileName ).Remove();
1217   }
1218   
1219   return Ok;
1220 }
1221
1222 //================================================================================
1223 /*!
1224  * \brief Provide human readable text by error code reported by ghs3d
1225  */
1226 //================================================================================
1227
1228 static string translateError(const int errNum)
1229 {
1230   switch ( errNum ) {
1231   case 0:
1232     return "The surface mesh includes a face of type other than edge, "
1233       "triangle or quadrilateral. This face type is not supported.";
1234   case 1:
1235     return "Not enough memory for the face table.";
1236   case 2:
1237     return "Not enough memory.";
1238   case 3:
1239     return "Not enough memory.";
1240   case 4:
1241     return "Face is ignored.";
1242   case 5:
1243     return "End of file. Some data are missing in the file.";
1244   case 6:
1245     return "Read error on the file. There are wrong data in the file.";
1246   case 7:
1247     return "the metric file is inadequate (dimension other than 3).";
1248   case 8:
1249     return "the metric file is inadequate (values not per vertices).";
1250   case 9:
1251     return "the metric file contains more than one field.";
1252   case 10:
1253     return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1254       "value of number of mesh vertices in the \".noboite\" file.";
1255   case 12:
1256     return "Too many sub-domains.";
1257   case 13:
1258     return "the number of vertices is negative or null.";
1259   case 14:
1260     return "the number of faces is negative or null.";
1261   case 15:
1262     return "A face has a null vertex.";
1263   case 22:
1264     return "incompatible data.";
1265   case 131:
1266     return "the number of vertices is negative or null.";
1267   case 132:
1268     return "the number of vertices is negative or null (in the \".mesh\" file).";
1269   case 133:
1270     return "the number of faces is negative or null.";
1271   case 1000:
1272     return "A face appears more than once in the input surface mesh.";
1273   case 1001:
1274     return "An edge appears more than once in the input surface mesh.";
1275   case 1002:
1276     return "A face has a vertex negative or null.";
1277   case 1003:
1278     return "NOT ENOUGH MEMORY.";
1279   case 2000:
1280     return "Not enough available memory.";
1281   case 2002:
1282     return "Some initial points cannot be inserted. The surface mesh is probably very bad "
1283       "in terms of quality or the input list of points is wrong.";
1284   case 2003:
1285     return "Some vertices are too close to one another or coincident.";
1286   case 2004:
1287     return "Some vertices are too close to one another or coincident.";
1288   case 2012:
1289     return "A vertex cannot be inserted.";
1290   case 2014:
1291     return "There are at least two points considered as coincident.";
1292   case 2103:
1293     return "Some vertices are too close to one another or coincident.";
1294   case 3000:
1295     return "The surface mesh regeneration step has failed.";
1296   case 3009:
1297     return "Constrained edge cannot be enforced.";
1298   case 3019:
1299     return "Constrained face cannot be enforced.";
1300   case 3029:
1301     return "Missing faces.";
1302   case 3100:
1303     return "No guess to start the definition of the connected component(s).";
1304   case 3101:
1305     return "The surface mesh includes at least one hole. The domain is not well defined.";
1306   case 3102:
1307     return "Impossible to define a component.";
1308   case 3103:
1309     return "The surface edge intersects another surface edge.";
1310   case 3104:
1311     return "The surface edge intersects the surface face.";
1312   case 3105:
1313     return "One boundary point lies within a surface face.";
1314   case 3106:
1315     return "One surface edge intersects a surface face.";
1316   case 3107:
1317     return "One boundary point lies within a surface edge.";
1318   case 3108:
1319     return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1320       "to too many swaps.";
1321   case 3109:
1322     return "Edge is unique (i.e., bounds a hole in the surface).";
1323   case 3122:
1324     return "Presumably, the surface mesh is not compatible with the domain being processed.";
1325   case 3123:
1326     return "Too many components, too many sub-domain.";
1327   case 3209:
1328     return "The surface mesh includes at least one hole. "
1329       "Therefore there is no domain properly defined.";
1330   case 3300:
1331     return "Statistics.";
1332   case 3400:
1333     return "Statistics.";
1334   case 3500:
1335     return "Warning, it is dramatically tedious to enforce the boundary items.";
1336   case 4000:
1337     return "Not enough memory at this time, nevertheless, the program continues. "
1338       "The expected mesh will be correct but not really as large as required.";
1339   case 4002:
1340     return "see above error code, resulting quality may be poor.";
1341   case 4003:
1342     return "Not enough memory at this time, nevertheless, the program continues (warning).";
1343   case 8000:
1344     return "Unknown face type.";
1345   case 8005:
1346   case 8006:
1347     return "End of file. Some data are missing in the file.";
1348   case 9000:
1349     return "A too small volume element is detected.";
1350   case 9001:
1351     return "There exists at least a null or negative volume element.";
1352   case 9002:
1353     return "There exist null or negative volume elements.";
1354   case 9003:
1355     return "A too small volume element is detected. A face is considered being degenerated.";
1356   case 9100:
1357     return "Some element is suspected to be very bad shaped or wrong.";
1358   case 9102:
1359     return "A too bad quality face is detected. This face is considered degenerated.";
1360   case 9112:
1361     return "A too bad quality face is detected. This face is degenerated.";
1362   case 9122:
1363     return "Presumably, the surface mesh is not compatible with the domain being processed.";
1364   case 9999:
1365     return "Abnormal error occured, contact hotline.";
1366   case 23600:
1367     return "Not enough memory for the face table.";
1368   case 23601:
1369     return "The algorithm cannot run further. "
1370       "The surface mesh is probably very bad in terms of quality.";
1371   case 23602:
1372     return "Bad vertex number.";
1373   }
1374   return "";
1375 }
1376
1377 //================================================================================
1378 /*!
1379  * \brief Retrieve from a string given number of integers
1380  */
1381 //================================================================================
1382
1383 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
1384 {
1385   ids.clear();
1386   ids.reserve( nbIds );
1387   while ( nbIds )
1388   {
1389     while ( !isdigit( *ptr )) ++ptr;
1390     if ( ptr[-1] == '-' ) --ptr;
1391     ids.push_back( strtol( ptr, &ptr, 10 ));
1392     --nbIds;
1393   }
1394   return ptr;
1395 }
1396
1397 //================================================================================
1398 /*!
1399  * \brief Retrieve problem description form a log file
1400  *  \retval bool - always false
1401  */
1402 //================================================================================
1403
1404 bool GHS3DPlugin_GHS3D::storeErrorDescription(const TCollection_AsciiString& logFile,
1405                                               const _Ghs2smdsConvertor &     toSmdsConvertor )
1406 {
1407   // open file
1408 #ifdef WNT
1409   int file = ::_open (logFile.ToCString(), _O_RDONLY|_O_BINARY);
1410 #else
1411   int file = ::open (logFile.ToCString(), O_RDONLY);
1412 #endif
1413   if ( file < 0 )
1414     return error( SMESH_Comment("See ") << logFile << " for problem description");
1415
1416   // get file size
1417 //   struct stat status;
1418 //   fstat(file, &status);
1419 //   size_t length = status.st_size;
1420   off_t length = lseek( file, 0, SEEK_END);
1421   lseek( file, 0, SEEK_SET);
1422
1423   // read file
1424   vector< char > buf( length );
1425   int nBytesRead = ::read (file, & buf[0], length);
1426   ::close (file);
1427   char* ptr = & buf[0];
1428   char* bufEnd = ptr + nBytesRead;
1429
1430   SMESH_Comment errDescription;
1431
1432   enum { NODE = 1, EDGE, TRIA, VOL, ID = 1 };
1433
1434   // look for errors "ERR #"
1435
1436   set<string> foundErrorStr; // to avoid reporting same error several times
1437   set<int>    elemErrorNums; // not to report different types of errors with bad elements
1438   while ( ++ptr < bufEnd )
1439   {
1440     if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1441       continue;
1442
1443     list<const SMDS_MeshElement*> badElems;
1444     vector<int> nodeIds;
1445
1446     ptr += 4;
1447     char* errBeg = ptr;
1448     int   errNum = strtol(ptr, &ptr, 10);
1449     switch ( errNum ) { // we treat errors enumerated in [SALOME platform 0019316] issue
1450     case 0015:
1451       // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1452       ptr = getIds(ptr, NODE, nodeIds);
1453       ptr = getIds(ptr, TRIA, nodeIds);
1454       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1455       break;
1456     case 1000: // ERR  1000 :  1 3 2
1457       // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1458       ptr = getIds(ptr, TRIA, nodeIds);
1459       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1460       break;
1461     case 1001:
1462       // Edge (e1, e2) appears more than once in the input surface mesh
1463       ptr = getIds(ptr, EDGE, nodeIds);
1464       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1465       break;
1466     case 1002:
1467       // Face (f 1, f 2, f 3) has a vertex negative or null
1468       ptr = getIds(ptr, TRIA, nodeIds);
1469       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1470       break;
1471     case 2004:
1472       // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1473       ptr = getIds(ptr, NODE, nodeIds);
1474       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1475       ptr = getIds(ptr, NODE, nodeIds);
1476       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1477       break;
1478     case 2012:
1479       // Vertex v1 cannot be inserted (warning).
1480       ptr = getIds(ptr, NODE, nodeIds);
1481       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1482       break;
1483     case 2014:
1484       // There are at least two points whose distance is dist, i.e., considered as coincident
1485     case 2103: // ERR  2103 :  16 WITH  3
1486       // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1487       ptr = getIds(ptr, NODE, nodeIds);
1488       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1489       ptr = getIds(ptr, NODE, nodeIds);
1490       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1491       break;
1492     case 3009:
1493       // Constrained edge (e1, e2) cannot be enforced (warning).
1494       ptr = getIds(ptr, EDGE, nodeIds);
1495       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1496       break;
1497     case 3019:
1498       // Constrained face (f 1, f 2, f 3) cannot be enforced
1499       ptr = getIds(ptr, TRIA, nodeIds);
1500       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1501       break;
1502     case 3103: // ERR  3103 :  1 2 WITH  7 3
1503       // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1504       ptr = getIds(ptr, EDGE, nodeIds);
1505       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1506       ptr = getIds(ptr, EDGE, nodeIds);
1507       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1508       break;
1509     case 3104: // ERR  3104 :  9 10 WITH  1 2 3
1510       // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1511       ptr = getIds(ptr, EDGE, nodeIds);
1512       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1513       ptr = getIds(ptr, TRIA, nodeIds);
1514       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1515       break;
1516     case 3105: // ERR  3105 :  8 IN  2 3 5
1517       // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1518       ptr = getIds(ptr, NODE, nodeIds);
1519       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1520       ptr = getIds(ptr, TRIA, nodeIds);
1521       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1522       break;
1523     case 3106:
1524       // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1525       ptr = getIds(ptr, EDGE, nodeIds);
1526       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1527       ptr = getIds(ptr, TRIA, nodeIds);
1528       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1529       break;
1530     case 3107: // ERR  3107 :  2 IN  4 1
1531       // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1532       ptr = getIds(ptr, NODE, nodeIds);
1533       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1534       ptr = getIds(ptr, EDGE, nodeIds);
1535       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1536       break;
1537     case 3109: // ERR  3109 :  EDGE  5 6 UNIQUE
1538       // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1539       ptr = getIds(ptr, EDGE, nodeIds);
1540       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1541       break;
1542     case 9000: // ERR  9000 
1543       //  ELEMENT  261 WITH VERTICES :  7 396 -8 242 
1544       //  VOLUME   : -1.11325045E+11 W.R.T. EPSILON   0.
1545       // A too small volume element is detected. Are reported the index of the element,
1546       // its four vertex indices, its volume and the tolerance threshold value
1547       ptr = getIds(ptr, ID, nodeIds);
1548       ptr = getIds(ptr, VOL, nodeIds);
1549       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1550       // even if all nodes found, volume it most probably invisible,
1551       // add its faces to demenstrate it anyhow
1552       {
1553         vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1554         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1555         faceNodes[2] = nodeIds[3]; // 013
1556         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1557         faceNodes[1] = nodeIds[2]; // 023
1558         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1559         faceNodes[0] = nodeIds[1]; // 123
1560         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1561       }
1562       break;
1563     case 9001: // ERR  9001
1564       //  %% NUMBER OF NEGATIVE VOLUME TETS  :  1 
1565       //  %% THE LARGEST NEGATIVE TET        :   1.75376581E+11 
1566       //  %%  NUMBER OF NULL VOLUME TETS     :  0
1567       // There exists at least a null or negative volume element
1568       break;
1569     case 9002:
1570       // There exist n null or negative volume elements
1571       break;
1572     case 9003:
1573       // A too small volume element is detected
1574       break;
1575     case 9102:
1576       // A too bad quality face is detected. This face is considered degenerated,
1577       // its index, its three vertex indices together with its quality value are reported
1578       break; // same as next
1579     case 9112: // ERR  9112 
1580       //  FACE   2 WITH VERTICES :  4 2 5 
1581       //  SMALL INRADIUS :   0.
1582       // A too bad quality face is detected. This face is degenerated,
1583       // its index, its three vertex indices together with its inradius are reported
1584       ptr = getIds(ptr, ID, nodeIds);
1585       ptr = getIds(ptr, TRIA, nodeIds);
1586       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1587       // add triangle edges as it most probably has zero area and hence invisible
1588       {
1589         vector<int> edgeNodes(2);
1590         edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
1591         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1592         edgeNodes[1] = nodeIds[2]; // 0-2
1593         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1594         edgeNodes[0] = nodeIds[1]; // 1-2
1595         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1596       }
1597       break;
1598     }
1599
1600     bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
1601     if ( !isNewError )
1602       continue; // not to report same error several times
1603
1604 //     const SMDS_MeshElement* nullElem = 0;
1605 //     bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
1606
1607 //     if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
1608 //       bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
1609 //       if ( oneMoreErrorType )
1610 //         continue; // not to report different types of errors with bad elements
1611 //     }
1612
1613     // store bad elements
1614     //if ( allElemsOk ) {
1615       list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
1616       for ( ; elem != badElems.end(); ++elem )
1617         addBadInputElement( *elem );
1618       //}
1619
1620     // make error text
1621     string text = translateError( errNum );
1622     if ( errDescription.find( text ) == text.npos ) {
1623       if ( !errDescription.empty() )
1624         errDescription << "\n";
1625       errDescription << text;
1626     }
1627
1628   } // end while
1629
1630   if ( errDescription.empty() ) { // no errors found
1631     char msg[] = "connection to server failed";
1632     if ( search( &buf[0], bufEnd, msg, msg + strlen(msg)) != bufEnd )
1633       errDescription << "Licence problems.";
1634   }
1635
1636   if ( errDescription.empty() )
1637     errDescription << "See " << logFile << " for problem description";
1638   else
1639     errDescription << "\nSee " << logFile << " for more information";
1640
1641   return error( errDescription );
1642 }
1643
1644 //================================================================================
1645 /*!
1646  * \brief Creates _Ghs2smdsConvertor
1647  */
1648 //================================================================================
1649
1650 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap)
1651   :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
1652 {
1653 }
1654
1655 //================================================================================
1656 /*!
1657  * \brief Creates _Ghs2smdsConvertor
1658  */
1659 //================================================================================
1660
1661 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> &  nodeByGhsId)
1662   : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
1663 {
1664 }
1665
1666 //================================================================================
1667 /*!
1668  * \brief Return SMDS element by ids of GHS3D nodes
1669  */
1670 //================================================================================
1671
1672 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
1673 {
1674   size_t nbNodes = ghsNodes.size();
1675   vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
1676   for ( size_t i = 0; i < nbNodes; ++i ) {
1677     int ghsNode = ghsNodes[ i ];
1678     if ( _ghs2NodeMap ) {
1679       map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
1680       if ( in == _ghs2NodeMap->end() )
1681         return 0;
1682       nodes[ i ] = in->second;
1683     }
1684     else {
1685       if ( ghsNode < 1 || ghsNode > _nodeByGhsId->size() )
1686         return 0;
1687       nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
1688     }
1689   }
1690   if ( nbNodes == 1 )
1691     return nodes[0];
1692
1693   if ( nbNodes == 2 ) {
1694     const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
1695     if ( !edge )
1696       edge = new SMDS_MeshEdge( nodes[0], nodes[1] );
1697     return edge;
1698   }
1699   if ( nbNodes == 3 ) {
1700     const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
1701     if ( !face )
1702       face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
1703     return face;
1704   }
1705   if ( nbNodes == 4 )
1706     return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
1707
1708   return 0;
1709 }
1710