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