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