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