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