]> SALOME platform Git repositories - plugins/ghs3dplugin.git/blob - src/GHS3DPlugin_GHS3D.cxx
Salome HOME
Merge from BR_V5_DEV 16Feb09
[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       return meshDS->ShapeToIndex( solid1 ); // - no
552   }
553
554   // find orientation of geom face within the first solid
555   TopExp_Explorer fExp( solid1, TopAbs_FACE );
556   for ( ; fExp.More(); fExp.Next() )
557     if ( geomFace.IsSame( fExp.Current() )) {
558       geomFace = TopoDS::Face( fExp.Current() );
559       break;
560     }
561   if ( !fExp.More() )
562     return invalidID; // face not found
563
564   // find UV of node1 on geomFace
565   SMESH_MesherHelper helper( mesh );
566   gp_XY uv = helper.GetNodeUV( geomFace, node1 );
567
568   // check that uv is correct
569   gp_Pnt node1Pnt ( node1->X(), node1->Y(), node1->Z() );
570   double tol = BRep_Tool::Tolerance( geomFace );
571   BRepAdaptor_Surface surface( geomFace );
572   if ( node1Pnt.Distance( surface.Value( uv.X(), uv.Y() )) > 2 * tol ) {
573     // project node1 onto geomFace to get right UV
574     GeomAPI_ProjectPointOnSurf projector( node1Pnt, surface.Surface().Surface() );
575     if ( !projector.IsDone() || projector.NbPoints() < 1 )
576       return invalidID;
577     Quantity_Parameter U,V;
578     projector.LowerDistanceParameters(U,V);
579     uv = gp_XY( U,V );
580   }
581   // normale to face at node1
582   gp_Pnt node2Pnt ( node2->X(), node2->Y(), node2->Z() );
583   gp_Pnt node3Pnt ( node3->X(), node3->Y(), node3->Z() );
584   gp_Vec vec12( node1Pnt, node2Pnt );
585   gp_Vec vec13( node1Pnt, node3Pnt );
586   gp_Vec meshNormal = vec12 ^ vec13;
587   if ( meshNormal.SquareMagnitude() < DBL_MIN )
588     return invalidID;
589   
590   // normale to geomFace at UV
591   gp_Vec du, dv;
592   surface.D1( uv.X(), uv.Y(), node1Pnt, du, dv );
593   gp_Vec geomNormal = du ^ dv;
594   if ( geomNormal.SquareMagnitude() < DBL_MIN )
595     return findShapeID( mesh, node2, node3, node1, toMeshHoles );
596   if ( geomFace.Orientation() == TopAbs_REVERSED )
597     geomNormal.Reverse();
598
599   // compare normals
600   bool isReverse = ( meshNormal * geomNormal ) < 0;
601   if ( !isReverse )
602     return meshDS->ShapeToIndex( solid1 );
603
604   if ( solids.Extent() == 1 )
605     return HOLE_ID; // we are inside a hole
606   else
607     return meshDS->ShapeToIndex( solids(2) );
608 }
609                        
610 //=======================================================================
611 //function : readResultFile
612 //purpose  : 
613 //=======================================================================
614
615 static bool readResultFile(const int                       fileOpen,
616                            SMESH_Mesh&                     theMesh,
617                            TopoDS_Shape                    tabShape[],
618                            double**                        tabBox,
619                            const int                       nbShape,
620                            map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
621                            bool                            toMeshHoles)
622 {
623   struct stat status;
624   size_t      length;
625
626   char *ptr, *mapPtr;
627   char *tetraPtr;
628   char *shapePtr;
629
630   SMESHDS_Mesh* theMeshDS = theMesh.GetMeshDS();
631
632   int fileStat;
633   int nbElems, nbNodes, nbInputNodes;
634   int nodeId/*, triangleId*/;
635   int nbTriangle;
636   int ID, shapeID, ghs3dShapeID;
637   int IdShapeRef = 1;
638   int compoundID =
639     nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
640
641   int *tab, *tabID, *nodeID, *nodeAssigne;
642   double *coord;
643   const SMDS_MeshNode **node;
644
645   tab    = new int[3];
646   //tabID  = new int[nbShape];
647   nodeID = new int[4];
648   coord  = new double[3];
649   node   = new const SMDS_MeshNode*[4];
650
651   TopoDS_Shape aSolid;
652   SMDS_MeshNode * aNewNode;
653   map <int,const SMDS_MeshNode*>::iterator itOnNode;
654   SMDS_MeshElement* aTet;
655 #ifdef _DEBUG_
656   set<int> shapeIDs;
657 #endif
658
659   // Read the file state
660   fileStat = fstat(fileOpen, &status);
661   length   = status.st_size;
662
663   // Mapping the result file into memory
664   ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
665   mapPtr = ptr;
666
667   ptr      = readMapIntLine(ptr, tab);
668   tetraPtr = ptr;
669
670   nbElems      = tab[0];
671   nbNodes      = tab[1];
672   nbInputNodes = tab[2];
673
674   nodeAssigne = new int[ nbNodes+1 ];
675
676   if (nbShape > 0)
677     aSolid = tabShape[0];
678
679   // Reading the nodeId
680   for (int i=0; i < 4*nbElems; i++)
681     nodeId = strtol(ptr, &ptr, 10);
682
683   // Reading the nodeCoor and update the nodeMap
684   for (int iNode=1; iNode <= nbNodes; iNode++) {
685     for (int iCoor=0; iCoor < 3; iCoor++)
686       coord[ iCoor ] = strtod(ptr, &ptr);
687     nodeAssigne[ iNode ] = 1;
688     if ( iNode > nbInputNodes ) {
689       nodeAssigne[ iNode ] = 0;
690       aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
691       theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
692     }
693   }
694
695   // Reading the number of triangles which corresponds to the number of sub-domains
696   nbTriangle = strtol(ptr, &ptr, 10);
697
698   tabID = new int[nbTriangle];
699   for (int i=0; i < nbTriangle; i++) {
700     tabID[i] = 0;
701     // find the solid corresponding to GHS3D sub-domain following
702     // the technique proposed in GHS3D manual in chapter
703     // "B.4 Subdomain (sub-region) assignment"
704     int nodeId1 = strtol(ptr, &ptr, 10);
705     int nodeId2 = strtol(ptr, &ptr, 10);
706     int nodeId3 = strtol(ptr, &ptr, 10);
707     if ( nbTriangle > 1 ) {
708       const SMDS_MeshNode* n1 = theGhs3dIdToNodeMap[ nodeId1 ];
709       const SMDS_MeshNode* n2 = theGhs3dIdToNodeMap[ nodeId2 ];
710       const SMDS_MeshNode* n3 = theGhs3dIdToNodeMap[ nodeId3 ];
711       try {
712         OCC_CATCH_SIGNALS;
713         tabID[i] = findShapeID( theMesh, n1, n2, n3, toMeshHoles );
714 #ifdef _DEBUG_
715         cout << i+1 << " subdomain: findShapeID() returns " << tabID[i] << endl;
716 #endif
717       } catch ( Standard_Failure ) {
718       } catch (...) {}
719     }
720   }
721
722   shapePtr = ptr;
723
724   if ( nbTriangle <= nbShape ) // no holes
725     toMeshHoles = true; // not avoid creating tetras in holes
726
727   // Associating the tetrahedrons to the shapes
728   shapeID = compoundID;
729   for (int iElem = 0; iElem < nbElems; iElem++) {
730     for (int iNode = 0; iNode < 4; iNode++) {
731       ID = strtol(tetraPtr, &tetraPtr, 10);
732       itOnNode = theGhs3dIdToNodeMap.find(ID);
733       node[ iNode ] = itOnNode->second;
734       nodeID[ iNode ] = ID;
735     }
736     // We always run GHS3D with "to mesh holes'==TRUE but we must not create
737     // tetras within holes depending on hypo option,
738     // so we first check if aTet is inside a hole and then create it 
739     //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
740     if ( nbTriangle > 1 ) {
741       shapeID = HOLE_ID; // negative shapeID means not to create tetras if !toMeshHoles
742       ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
743       if ( tabID[ ghs3dShapeID ] == 0 ) {
744         TopAbs_State state;
745         aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
746         if ( toMeshHoles || state == TopAbs_IN )
747           shapeID = theMeshDS->ShapeToIndex( aSolid );
748         tabID[ ghs3dShapeID ] = shapeID;
749       }
750       else
751         shapeID = tabID[ ghs3dShapeID ];
752     }
753     else if ( nbShape > 1 ) {
754       // Case where nbTriangle == 1 while nbShape == 2 encountered
755       // with compound of 2 boxes and "To mesh holes"==False,
756       // so there are no subdomains specified for each tetrahedron.
757       // Try to guess a solid by a node already bound to shape
758       shapeID = 0;
759       for ( int i=0; i<4 && shapeID==0; i++ ) {
760         if ( nodeAssigne[ nodeID[i] ] == 1 &&
761              node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
762              node[i]->GetPosition()->GetShapeId() > 1 )
763         {
764           shapeID = node[i]->GetPosition()->GetShapeId();
765         }
766       }
767       if ( shapeID==0 ) {
768         aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
769         shapeID = theMeshDS->ShapeToIndex( aSolid );
770       }
771     }
772     // set new nodes and tetrahedron onto the shape
773     for ( int i=0; i<4; i++ ) {
774       if ( nodeAssigne[ nodeID[i] ] == 0 ) {
775         if ( shapeID != HOLE_ID )
776           theMeshDS->SetNodeInVolume( node[i], shapeID );
777         nodeAssigne[ nodeID[i] ] = shapeID;
778       }
779     }
780     if ( toMeshHoles || shapeID != HOLE_ID ) {
781       aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
782       theMeshDS->SetMeshElementOnShape( aTet, shapeID );
783     }
784 #ifdef _DEBUG_
785     shapeIDs.insert( shapeID );
786 #endif
787   }
788   // Remove nodes of tetras inside holes if !toMeshHoles
789   if ( !toMeshHoles ) {
790     itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
791     for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
792       ID = itOnNode->first;
793       if ( nodeAssigne[ ID ] == HOLE_ID )
794         theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
795     }
796   }
797
798   if ( nbElems )
799     cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
800   munmap(mapPtr, length);
801   close(fileOpen);
802
803   delete [] tab;
804   delete [] tabID;
805   delete [] nodeID;
806   delete [] coord;
807   delete [] node;
808   delete [] nodeAssigne;
809
810 #ifdef _DEBUG_
811   if ( shapeIDs.size() != nbShape ) {
812     cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << endl;
813     for (int i=0; i<nbShape; i++) {
814       shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
815       if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
816         cout << "  Solid #" << shapeID << " not found" << endl;
817     }
818   }
819 #endif
820
821   return true;
822 }
823
824 //=======================================================================
825 //function : readResultFile
826 //purpose  : 
827 //=======================================================================
828
829 static bool readResultFile(const int                      fileOpen,
830                            SMESHDS_Mesh*                  theMeshDS,
831                            TopoDS_Shape                   aSolid,
832                            vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
833
834   struct stat  status;
835   size_t       length;
836
837   char *ptr, *mapPtr;
838   char *tetraPtr;
839   char *shapePtr;
840
841   int fileStat;
842   int nbElems, nbNodes, nbInputNodes;
843   int nodeId, triangleId;
844   int nbTriangle;
845   int ID, shapeID;
846
847   int *tab;
848   double *coord;
849   const SMDS_MeshNode **node;
850
851   tab   = new int[3];
852   coord = new double[3];
853   node  = new const SMDS_MeshNode*[4];
854
855   SMDS_MeshNode * aNewNode;
856   map <int,const SMDS_MeshNode*>::iterator IdNode;
857   SMDS_MeshElement* aTet;
858
859   // Read the file state
860   fileStat = fstat(fileOpen, &status);
861   length   = status.st_size;
862
863   // Mapping the result file into memory
864   ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
865   mapPtr = ptr;
866
867   ptr      = readMapIntLine(ptr, tab);
868   tetraPtr = ptr;
869
870   nbElems      = tab[0];
871   nbNodes      = tab[1];
872   nbInputNodes = tab[2];
873
874   theNodeByGhs3dId.resize( nbNodes );
875
876   // Reading the nodeId
877   for (int i=0; i < 4*nbElems; i++)
878     nodeId = strtol(ptr, &ptr, 10);
879
880   // Reading the nodeCoor and update the nodeMap
881   shapeID = theMeshDS->ShapeToIndex( aSolid );
882   for (int iNode=0; iNode < nbNodes; iNode++) {
883     for (int iCoor=0; iCoor < 3; iCoor++)
884       coord[ iCoor ] = strtod(ptr, &ptr);
885     if ((iNode+1) > nbInputNodes) {
886       aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
887       theMeshDS->SetNodeInVolume( aNewNode, shapeID );
888       theNodeByGhs3dId[ iNode ] = aNewNode;
889     }
890   }
891
892   // Reading the triangles
893   nbTriangle = strtol(ptr, &ptr, 10);
894
895   for (int i=0; i < 3*nbTriangle; i++)
896     triangleId = strtol(ptr, &ptr, 10);
897
898   shapePtr = ptr;
899
900   // Associating the tetrahedrons to the shapes
901   for (int iElem = 0; iElem < nbElems; iElem++) {
902     for (int iNode = 0; iNode < 4; iNode++) {
903       ID = strtol(tetraPtr, &tetraPtr, 10);
904       node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
905     }
906     aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
907     shapeID = theMeshDS->ShapeToIndex( aSolid );
908     theMeshDS->SetMeshElementOnShape( aTet, shapeID );
909   }
910   if ( nbElems )
911     cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
912   munmap(mapPtr, length);
913   close(fileOpen);
914
915   delete [] tab;
916   delete [] coord;
917   delete [] node;
918
919   return true;
920 }
921
922 //=============================================================================
923 /*!
924  *Here we are going to use the GHS3D mesher
925  */
926 //=============================================================================
927
928 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh&         theMesh,
929                                 const TopoDS_Shape& theShape)
930 {
931   bool Ok(false);
932   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
933
934   _nbShape = countShape( meshDS, TopAbs_SOLID ); // we count the number of shapes
935
936   // create bounding box for every shape inside the compound
937
938   int iShape = 0;
939   TopoDS_Shape* tabShape;
940   double**      tabBox;
941   tabShape = new TopoDS_Shape[_nbShape];
942   tabBox   = new double*[_nbShape];
943   for (int i=0; i<_nbShape; i++)
944     tabBox[i] = new double[6];
945   Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
946
947   TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
948   for (; expBox.More(); expBox.Next()) {
949     tabShape[iShape] = expBox.Current();
950     Bnd_Box BoundingBox;
951     BRepBndLib::Add(expBox.Current(), BoundingBox);
952     BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
953     tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
954     tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
955     tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
956     iShape++;
957   }
958
959   // a unique working file name
960   // to avoid access to the same files by eg different users
961   TCollection_AsciiString aGenericName
962     = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
963
964   TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
965   TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
966   aFacesFileName  = aGenericName + ".faces";  // in faces
967   aPointsFileName = aGenericName + ".points"; // in points
968   aResultFileName = aGenericName + ".noboite";// out points and volumes
969   aBadResFileName = aGenericName + ".boite";  // out bad result
970   aBbResFileName  = aGenericName + ".bb";     // out vertex stepsize
971   aLogFileName    = aGenericName + ".log";    // log
972
973   // -----------------
974   // make input files
975   // -----------------
976
977   ofstream aFacesFile  ( aFacesFileName.ToCString()  , ios::out);
978   ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
979
980   Ok =
981 #ifdef WIN32
982     aFacesFile->is_open() && aPointsFile->is_open();
983 #else
984     aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
985 #endif
986   if (!Ok) {
987     INFOS( "Can't write into " << aFacesFileName);
988     return error(SMESH_Comment("Can't write into ") << aFacesFileName);
989   }
990   map <int,int> aSmdsToGhs3dIdMap;
991   map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
992
993   Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap ) &&
994        writeFaces ( aFacesFile,  meshDS, aSmdsToGhs3dIdMap );
995
996   aFacesFile.close();
997   aPointsFile.close();
998
999   if ( ! Ok ) {
1000     if ( !_keepFiles ) {
1001       OSD_File( aFacesFileName ).Remove();
1002       OSD_File( aPointsFileName ).Remove();
1003     }
1004     return error(COMPERR_BAD_INPUT_MESH);
1005   }
1006   OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1007
1008   // -----------------
1009   // run ghs3d mesher
1010   // -----------------
1011
1012   TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1013   cmd += TCollection_AsciiString(" -f ") + aGenericName;  // file to read
1014   cmd += TCollection_AsciiString(" 1>" ) + aLogFileName;  // dump into file
1015
1016   cout << endl;
1017   cout << "Ghs3d execution..." << endl;
1018   cout << cmd << endl;
1019
1020   system( cmd.ToCString() ); // run
1021
1022   cout << endl;
1023   cout << "End of Ghs3d execution !" << endl;
1024
1025   // --------------
1026   // read a result
1027   // --------------
1028
1029   // Mapping the result file
1030
1031   int fileOpen;
1032   fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1033   if ( fileOpen < 0 ) {
1034     cout << endl;
1035     cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << endl;
1036     cout << "Log: " << aLogFileName << endl;
1037     Ok = false;
1038   }
1039   else {
1040     bool toMeshHoles =
1041       _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1042     Ok = readResultFile( fileOpen, theMesh, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1043                          toMeshHoles );
1044   }
1045
1046   // ---------------------
1047   // remove working files
1048   // ---------------------
1049
1050   if ( Ok )
1051   {
1052     if ( !_keepFiles )
1053       OSD_File( aLogFileName ).Remove();
1054   }
1055   else if ( OSD_File( aLogFileName ).Size() > 0 )
1056   {
1057     // get problem description from the log file
1058     _Ghs2smdsConvertor conv( aGhs3dIdToNodeMap );
1059     storeErrorDescription( aLogFileName, conv );
1060   }
1061   else
1062   {
1063     // the log file is empty
1064     OSD_File( aLogFileName ).Remove();
1065     INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1066     error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1067   }
1068
1069   if ( !_keepFiles ) {
1070     OSD_File( aFacesFileName ).Remove();
1071     OSD_File( aPointsFileName ).Remove();
1072     OSD_File( aResultFileName ).Remove();
1073     OSD_File( aBadResFileName ).Remove();
1074     OSD_File( aBbResFileName ).Remove();
1075   }
1076   cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1077   if ( !Ok )
1078     cout << "not ";
1079   cout << "treated !" << endl;
1080   cout << endl;
1081
1082   _nbShape = 0;    // re-initializing _nbShape for the next Compute() method call
1083   delete [] tabShape;
1084   delete [] tabBox;
1085
1086   return Ok;
1087 }
1088
1089 //=============================================================================
1090 /*!
1091  *Here we are going to use the GHS3D mesher w/o geometry
1092  */
1093 //=============================================================================
1094 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh&         theMesh,
1095                                 SMESH_MesherHelper* aHelper)
1096 {
1097   MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1098
1099   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1100   TopoDS_Shape theShape = aHelper->GetSubShape();
1101
1102   // a unique working file name
1103   // to avoid access to the same files by eg different users
1104   TCollection_AsciiString aGenericName
1105     = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1106
1107   TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1108   TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1109   aFacesFileName  = aGenericName + ".faces";  // in faces
1110   aPointsFileName = aGenericName + ".points"; // in points
1111   aResultFileName = aGenericName + ".noboite";// out points and volumes
1112   aBadResFileName = aGenericName + ".boite";  // out bad result
1113   aBbResFileName  = aGenericName + ".bb";     // out vertex stepsize
1114   aLogFileName    = aGenericName + ".log";    // log
1115
1116   // -----------------
1117   // make input files
1118   // -----------------
1119
1120   ofstream aFacesFile  ( aFacesFileName.ToCString()  , ios::out);
1121   ofstream aPointsFile  ( aPointsFileName.ToCString()  , ios::out);
1122   bool Ok =
1123 #ifdef WIN32
1124     aFacesFile->is_open() && aPointsFile->is_open();
1125 #else
1126     aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1127 #endif
1128
1129   if (!Ok)
1130     return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1131
1132   vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1133
1134   Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1135         writePoints( aPointsFile, meshDS, aNodeByGhs3dId));
1136   
1137   aFacesFile.close();
1138   aPointsFile.close();
1139   
1140   if ( ! Ok ) {
1141     if ( !_keepFiles ) {
1142       OSD_File( aFacesFileName ).Remove();
1143       OSD_File( aPointsFileName ).Remove();
1144     }
1145     return error(COMPERR_BAD_INPUT_MESH);
1146   }
1147   OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1148
1149   // -----------------
1150   // run ghs3d mesher
1151   // -----------------
1152
1153   TCollection_AsciiString cmd =
1154     (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1155   cmd += TCollection_AsciiString(" -f ") + aGenericName;  // file to read
1156   cmd += TCollection_AsciiString(" 1>" ) + aLogFileName;  // dump into file
1157   
1158   system( cmd.ToCString() ); // run
1159
1160   // --------------
1161   // read a result
1162   // --------------
1163   int fileOpen;
1164   fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1165   if ( fileOpen < 0 ) {
1166     cout << endl;
1167     cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
1168     cout << "Log: " << aLogFileName << endl;
1169     cout << endl;
1170     Ok = false;
1171   }
1172   else {
1173     Ok = readResultFile( fileOpen, meshDS, theShape ,aNodeByGhs3dId );
1174   }
1175   
1176   // ---------------------
1177   // remove working files
1178   // ---------------------
1179
1180   if ( Ok )
1181   {
1182     if ( !_keepFiles )
1183       OSD_File( aLogFileName ).Remove();
1184   }
1185   else if ( OSD_File( aLogFileName ).Size() > 0 )
1186   {
1187     // get problem description from the log file
1188     _Ghs2smdsConvertor conv( aNodeByGhs3dId );
1189     storeErrorDescription( aLogFileName, conv );
1190   }
1191   else {
1192     // the log file is empty
1193     OSD_File( aLogFileName ).Remove();
1194     INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1195     error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1196   }
1197
1198   if ( !_keepFiles )
1199   {
1200     OSD_File( aFacesFileName ).Remove();
1201     OSD_File( aPointsFileName ).Remove();
1202     OSD_File( aResultFileName ).Remove();
1203     OSD_File( aBadResFileName ).Remove();
1204     OSD_File( aBbResFileName ).Remove();
1205   }
1206   
1207   return Ok;
1208 }
1209
1210 //================================================================================
1211 /*!
1212  * \brief Provide human readable text by error code reported by ghs3d
1213  */
1214 //================================================================================
1215
1216 static string translateError(const int errNum)
1217 {
1218   switch ( errNum ) {
1219   case 0:
1220     return "The surface mesh includes a face of type other than edge, "
1221       "triangle or quadrilateral. This face type is not supported.";
1222   case 1:
1223     return "Not enough memory for the face table.";
1224   case 2:
1225     return "Not enough memory.";
1226   case 3:
1227     return "Not enough memory.";
1228   case 4:
1229     return "Face is ignored.";
1230   case 5:
1231     return "End of file. Some data are missing in the file.";
1232   case 6:
1233     return "Read error on the file. There are wrong data in the file.";
1234   case 7:
1235     return "the metric file is inadequate (dimension other than 3).";
1236   case 8:
1237     return "the metric file is inadequate (values not per vertices).";
1238   case 9:
1239     return "the metric file contains more than one field.";
1240   case 10:
1241     return "the number of values in the \".bb\" (metric file) is incompatible with the expected"
1242       "value of number of mesh vertices in the \".noboite\" file.";
1243   case 12:
1244     return "Too many sub-domains.";
1245   case 13:
1246     return "the number of vertices is negative or null.";
1247   case 14:
1248     return "the number of faces is negative or null.";
1249   case 15:
1250     return "A face has a null vertex.";
1251   case 22:
1252     return "incompatible data.";
1253   case 131:
1254     return "the number of vertices is negative or null.";
1255   case 132:
1256     return "the number of vertices is negative or null (in the \".mesh\" file).";
1257   case 133:
1258     return "the number of faces is negative or null.";
1259   case 1000:
1260     return "A face appears more than once in the input surface mesh.";
1261   case 1001:
1262     return "An edge appears more than once in the input surface mesh.";
1263   case 1002:
1264     return "A face has a vertex negative or null.";
1265   case 1003:
1266     return "NOT ENOUGH MEMORY.";
1267   case 2000:
1268     return "Not enough available memory.";
1269   case 2002:
1270     return "Some initial points cannot be inserted. The surface mesh is probably very bad "
1271       "in terms of quality or the input list of points is wrong.";
1272   case 2003:
1273     return "Some vertices are too close to one another or coincident.";
1274   case 2004:
1275     return "Some vertices are too close to one another or coincident.";
1276   case 2012:
1277     return "A vertex cannot be inserted.";
1278   case 2014:
1279     return "There are at least two points considered as coincident.";
1280   case 2103:
1281     return "Some vertices are too close to one another or coincident.";
1282   case 3000:
1283     return "The surface mesh regeneration step has failed.";
1284   case 3009:
1285     return "Constrained edge cannot be enforced.";
1286   case 3019:
1287     return "Constrained face cannot be enforced.";
1288   case 3029:
1289     return "Missing faces.";
1290   case 3100:
1291     return "No guess to start the definition of the connected component(s).";
1292   case 3101:
1293     return "The surface mesh includes at least one hole. The domain is not well defined.";
1294   case 3102:
1295     return "Impossible to define a component.";
1296   case 3103:
1297     return "The surface edge intersects another surface edge.";
1298   case 3104:
1299     return "The surface edge intersects the surface face.";
1300   case 3105:
1301     return "One boundary point lies within a surface face.";
1302   case 3106:
1303     return "One surface edge intersects a surface face.";
1304   case 3107:
1305     return "One boundary point lies within a surface edge.";
1306   case 3108:
1307     return "Insufficient memory ressources detected due to a bad quality surface mesh leading "
1308       "to too many swaps.";
1309   case 3109:
1310     return "Edge is unique (i.e., bounds a hole in the surface).";
1311   case 3122:
1312     return "Presumably, the surface mesh is not compatible with the domain being processed.";
1313   case 3123:
1314     return "Too many components, too many sub-domain.";
1315   case 3209:
1316     return "The surface mesh includes at least one hole. "
1317       "Therefore there is no domain properly defined.";
1318   case 3300:
1319     return "Statistics.";
1320   case 3400:
1321     return "Statistics.";
1322   case 3500:
1323     return "Warning, it is dramatically tedious to enforce the boundary items.";
1324   case 4000:
1325     return "Not enough memory at this time, nevertheless, the program continues. "
1326       "The expected mesh will be correct but not really as large as required.";
1327   case 4002:
1328     return "see above error code, resulting quality may be poor.";
1329   case 4003:
1330     return "Not enough memory at this time, nevertheless, the program continues (warning).";
1331   case 8000:
1332     return "Unknown face type.";
1333   case 8005:
1334   case 8006:
1335     return "End of file. Some data are missing in the file.";
1336   case 9000:
1337     return "A too small volume element is detected.";
1338   case 9001:
1339     return "There exists at least a null or negative volume element.";
1340   case 9002:
1341     return "There exist null or negative volume elements.";
1342   case 9003:
1343     return "A too small volume element is detected. A face is considered being degenerated.";
1344   case 9100:
1345     return "Some element is suspected to be very bad shaped or wrong.";
1346   case 9102:
1347     return "A too bad quality face is detected. This face is considered degenerated.";
1348   case 9112:
1349     return "A too bad quality face is detected. This face is degenerated.";
1350   case 9122:
1351     return "Presumably, the surface mesh is not compatible with the domain being processed.";
1352   case 9999:
1353     return "Abnormal error occured, contact hotline.";
1354   case 23600:
1355     return "Not enough memory for the face table.";
1356   case 23601:
1357     return "The algorithm cannot run further. "
1358       "The surface mesh is probably very bad in terms of quality.";
1359   case 23602:
1360     return "Bad vertex number.";
1361   }
1362   return "";
1363 }
1364
1365 //================================================================================
1366 /*!
1367  * \brief Retrieve from a string given number of integers
1368  */
1369 //================================================================================
1370
1371 static char* getIds( char* ptr, int nbIds, vector<int>& ids )
1372 {
1373   ids.clear();
1374   ids.reserve( nbIds );
1375   while ( nbIds )
1376   {
1377     while ( !isdigit( *ptr )) ++ptr;
1378     if ( ptr[-1] == '-' ) --ptr;
1379     ids.push_back( strtol( ptr, &ptr, 10 ));
1380     --nbIds;
1381   }
1382   return ptr;
1383 }
1384
1385 //================================================================================
1386 /*!
1387  * \brief Retrieve problem description form a log file
1388  *  \retval bool - always false
1389  */
1390 //================================================================================
1391
1392 bool GHS3DPlugin_GHS3D::storeErrorDescription(const TCollection_AsciiString& logFile,
1393                                               const _Ghs2smdsConvertor &     toSmdsConvertor )
1394 {
1395   // open file
1396 #ifdef WNT
1397   int file = ::_open (logFile.ToCString(), _O_RDONLY|_O_BINARY);
1398 #else
1399   int file = ::open (logFile.ToCString(), O_RDONLY);
1400 #endif
1401   if ( file < 0 )
1402     return error( SMESH_Comment("See ") << logFile << " for problem description");
1403
1404   // get file size
1405 //   struct stat status;
1406 //   fstat(file, &status);
1407 //   size_t length = status.st_size;
1408   off_t length = lseek( file, 0, SEEK_END);
1409   lseek( file, 0, SEEK_SET);
1410
1411   // read file
1412   vector< char > buf( length );
1413   int nBytesRead = ::read (file, & buf[0], length);
1414   ::close (file);
1415   char* ptr = & buf[0];
1416   char* bufEnd = ptr + nBytesRead;
1417
1418   SMESH_Comment errDescription;
1419
1420   enum { NODE = 1, EDGE, TRIA, VOL, ID = 1 };
1421
1422   // look for errors "ERR #"
1423
1424   set<string> foundErrorStr; // to avoid reporting same error several times
1425   set<int>    elemErrorNums; // not to report different types of errors with bad elements
1426   while ( ++ptr < bufEnd )
1427   {
1428     if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1429       continue;
1430
1431     list<const SMDS_MeshElement*> badElems;
1432     vector<int> nodeIds;
1433
1434     ptr += 4;
1435     char* errBeg = ptr;
1436     int   errNum = strtol(ptr, &ptr, 10);
1437     switch ( errNum ) { // we treat errors enumerated in [SALOME platform 0019316] issue
1438     case 0015:
1439       // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1440       ptr = getIds(ptr, NODE, nodeIds);
1441       ptr = getIds(ptr, TRIA, nodeIds);
1442       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1443       break;
1444     case 1000: // ERR  1000 :  1 3 2
1445       // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1446       ptr = getIds(ptr, TRIA, nodeIds);
1447       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1448       break;
1449     case 1001:
1450       // Edge (e1, e2) appears more than once in the input surface mesh
1451       ptr = getIds(ptr, EDGE, nodeIds);
1452       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1453       break;
1454     case 1002:
1455       // Face (f 1, f 2, f 3) has a vertex negative or null
1456       ptr = getIds(ptr, TRIA, nodeIds);
1457       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1458       break;
1459     case 2004:
1460       // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1461       ptr = getIds(ptr, NODE, nodeIds);
1462       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1463       ptr = getIds(ptr, NODE, nodeIds);
1464       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1465       break;
1466     case 2012:
1467       // Vertex v1 cannot be inserted (warning).
1468       ptr = getIds(ptr, NODE, nodeIds);
1469       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1470       break;
1471     case 2014:
1472       // There are at least two points whose distance is dist, i.e., considered as coincident
1473     case 2103: // ERR  2103 :  16 WITH  3
1474       // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1475       ptr = getIds(ptr, NODE, nodeIds);
1476       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1477       ptr = getIds(ptr, NODE, nodeIds);
1478       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1479       break;
1480     case 3009:
1481       // Constrained edge (e1, e2) cannot be enforced (warning).
1482       ptr = getIds(ptr, EDGE, nodeIds);
1483       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1484       break;
1485     case 3019:
1486       // Constrained face (f 1, f 2, f 3) cannot be enforced
1487       ptr = getIds(ptr, TRIA, nodeIds);
1488       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1489       break;
1490     case 3103: // ERR  3103 :  1 2 WITH  7 3
1491       // The surface edge (e1, e2) intersects another surface edge (e3, e4)
1492       ptr = getIds(ptr, EDGE, nodeIds);
1493       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1494       ptr = getIds(ptr, EDGE, nodeIds);
1495       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1496       break;
1497     case 3104: // ERR  3104 :  9 10 WITH  1 2 3
1498       // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
1499       ptr = getIds(ptr, EDGE, nodeIds);
1500       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1501       ptr = getIds(ptr, TRIA, nodeIds);
1502       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1503       break;
1504     case 3105: // ERR  3105 :  8 IN  2 3 5
1505       // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
1506       ptr = getIds(ptr, NODE, nodeIds);
1507       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1508       ptr = getIds(ptr, TRIA, nodeIds);
1509       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1510       break;
1511     case 3106:
1512       // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
1513       ptr = getIds(ptr, EDGE, nodeIds);
1514       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1515       ptr = getIds(ptr, TRIA, nodeIds);
1516       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1517       break;
1518     case 3107: // ERR  3107 :  2 IN  4 1
1519       // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
1520       ptr = getIds(ptr, NODE, nodeIds);
1521       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1522       ptr = getIds(ptr, EDGE, nodeIds);
1523       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1524       break;
1525     case 3109: // ERR  3109 :  EDGE  5 6 UNIQUE
1526       // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1527       ptr = getIds(ptr, EDGE, nodeIds);
1528       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1529       break;
1530     case 9000: // ERR  9000 
1531       //  ELEMENT  261 WITH VERTICES :  7 396 -8 242 
1532       //  VOLUME   : -1.11325045E+11 W.R.T. EPSILON   0.
1533       // A too small volume element is detected. Are reported the index of the element,
1534       // its four vertex indices, its volume and the tolerance threshold value
1535       ptr = getIds(ptr, ID, nodeIds);
1536       ptr = getIds(ptr, VOL, nodeIds);
1537       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1538       // even if all nodes found, volume it most probably invisible,
1539       // add its faces to demenstrate it anyhow
1540       {
1541         vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
1542         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1543         faceNodes[2] = nodeIds[3]; // 013
1544         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1545         faceNodes[1] = nodeIds[2]; // 023
1546         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1547         faceNodes[0] = nodeIds[1]; // 123
1548         badElems.push_back( toSmdsConvertor.getElement(faceNodes));
1549       }
1550       break;
1551     case 9001: // ERR  9001
1552       //  %% NUMBER OF NEGATIVE VOLUME TETS  :  1 
1553       //  %% THE LARGEST NEGATIVE TET        :   1.75376581E+11 
1554       //  %%  NUMBER OF NULL VOLUME TETS     :  0
1555       // There exists at least a null or negative volume element
1556       break;
1557     case 9002:
1558       // There exist n null or negative volume elements
1559       break;
1560     case 9003:
1561       // A too small volume element is detected
1562       break;
1563     case 9102:
1564       // A too bad quality face is detected. This face is considered degenerated,
1565       // its index, its three vertex indices together with its quality value are reported
1566       break; // same as next
1567     case 9112: // ERR  9112 
1568       //  FACE   2 WITH VERTICES :  4 2 5 
1569       //  SMALL INRADIUS :   0.
1570       // A too bad quality face is detected. This face is degenerated,
1571       // its index, its three vertex indices together with its inradius are reported
1572       ptr = getIds(ptr, ID, nodeIds);
1573       ptr = getIds(ptr, TRIA, nodeIds);
1574       badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1575       // add triangle edges as it most probably has zero area and hence invisible
1576       {
1577         vector<int> edgeNodes(2);
1578         edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
1579         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1580         edgeNodes[1] = nodeIds[2]; // 0-2
1581         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1582         edgeNodes[0] = nodeIds[1]; // 1-2
1583         badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1584       }
1585       break;
1586     }
1587
1588     bool isNewError = foundErrorStr.insert( string( errBeg, ptr )).second;
1589     if ( !isNewError )
1590       continue; // not to report same error several times
1591
1592 //     const SMDS_MeshElement* nullElem = 0;
1593 //     bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
1594
1595 //     if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
1596 //       bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
1597 //       if ( oneMoreErrorType )
1598 //         continue; // not to report different types of errors with bad elements
1599 //     }
1600
1601     // store bad elements
1602     //if ( allElemsOk ) {
1603       list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
1604       for ( ; elem != badElems.end(); ++elem )
1605         addBadInputElement( *elem );
1606       //}
1607
1608     // make error text
1609     string text = translateError( errNum );
1610     if ( errDescription.find( text ) == text.npos ) {
1611       if ( !errDescription.empty() )
1612         errDescription << "\n";
1613       errDescription << text;
1614     }
1615
1616   } // end while
1617
1618   if ( errDescription.empty() ) { // no errors found
1619     char msg[] = "connection to server failed";
1620     if ( search( &buf[0], bufEnd, msg, msg + strlen(msg)) != bufEnd )
1621       errDescription << "Licence problems.";
1622   }
1623
1624   if ( errDescription.empty() )
1625     errDescription << "See " << logFile << " for problem description";
1626   else
1627     errDescription << "\nSee " << logFile << " for more information";
1628
1629   return error( errDescription );
1630 }
1631
1632 //================================================================================
1633 /*!
1634  * \brief Creates _Ghs2smdsConvertor
1635  */
1636 //================================================================================
1637
1638 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const map <int,const SMDS_MeshNode*> & ghs2NodeMap)
1639   :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
1640 {
1641 }
1642
1643 //================================================================================
1644 /*!
1645  * \brief Creates _Ghs2smdsConvertor
1646  */
1647 //================================================================================
1648
1649 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const vector <const SMDS_MeshNode*> &  nodeByGhsId)
1650   : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
1651 {
1652 }
1653
1654 //================================================================================
1655 /*!
1656  * \brief Return SMDS element by ids of GHS3D nodes
1657  */
1658 //================================================================================
1659
1660 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const vector<int>& ghsNodes) const
1661 {
1662   size_t nbNodes = ghsNodes.size();
1663   vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
1664   for ( size_t i = 0; i < nbNodes; ++i ) {
1665     int ghsNode = ghsNodes[ i ];
1666     if ( _ghs2NodeMap ) {
1667       map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
1668       if ( in == _ghs2NodeMap->end() )
1669         return 0;
1670       nodes[ i ] = in->second;
1671     }
1672     else {
1673       if ( ghsNode < 1 || ghsNode > _nodeByGhsId->size() )
1674         return 0;
1675       nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
1676     }
1677   }
1678   if ( nbNodes == 1 )
1679     return nodes[0];
1680
1681   if ( nbNodes == 2 ) {
1682     const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
1683     if ( !edge )
1684       edge = new SMDS_MeshEdge( nodes[0], nodes[1] );
1685     return edge;
1686   }
1687   if ( nbNodes == 3 ) {
1688     const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
1689     if ( !face )
1690       face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
1691     return face;
1692   }
1693   if ( nbNodes == 4 )
1694     return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
1695
1696   return 0;
1697 }
1698