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