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