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