Salome HOME
Merge from BR_V5_DEV 17Feb09
[plugins/hexoticplugin.git] / src / HexoticPlugin / HexoticPlugin_Hexotic.cxx
1 //  Copyright (C) 2007-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   : HexoticPlugin_Hexotic.cxx
21 // Author : Lioka RAZAFINDRAZAKA (CEA)
22 // ---
23 //
24 #include "HexoticPlugin_Hexotic.hxx"
25 #include "HexoticPlugin_Hypothesis.hxx"
26
27 #include "SMDS_MeshElement.hxx"
28 #include "SMDS_MeshNode.hxx"
29
30 #include <TopoDS.hxx>
31 #include <TopExp_Explorer.hxx>
32 #include <OSD_File.hxx>
33
34 #include "utilities.h"
35
36 #ifndef WIN32
37 #include <sys/sysinfo.h>
38 #endif
39
40 #ifdef _DEBUG_
41 #define DUMP(txt) \
42 //  cout << txt
43 #else
44 #define DUMP(txt)
45 #endif
46
47 #include <SMESH_Gen.hxx>
48 #include <SMESHDS_Mesh.hxx>
49 #include <SMESH_ControlsDef.hxx>
50
51 #include <list>
52 #include <cstdlib>
53 #include <iostream>
54
55 #include <BRepClass3d_SolidClassifier.hxx>
56 #include <Bnd_Box.hxx>
57 #include <BRepBndLib.hxx>
58 #include <Precision.hxx>
59
60 #include <BRepBuilderAPI_MakeVertex.hxx>
61 #include <BRep_Tool.hxx>
62 #include <TopTools_MapOfShape.hxx>
63 #include <TopTools_Array1OfShape.hxx>
64 #include <BRepExtrema_DistShapeShape.hxx>
65 #include <TColStd_Array1OfReal.hxx>
66 #include <BRepExtrema_DistShapeShape.hxx>
67
68 //=============================================================================
69 /*!
70  *  
71  */
72 //=============================================================================
73
74 HexoticPlugin_Hexotic::HexoticPlugin_Hexotic(int hypId, int studyId, SMESH_Gen* gen)
75   : SMESH_3D_Algo(hypId, studyId, gen)
76 {
77   MESSAGE("HexoticPlugin_Hexotic::HexoticPlugin_Hexotic");
78   _name = "Hexotic_3D";
79   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
80 //   _onlyUnaryInput = false;
81   _iShape=0;
82   _nbShape=0;
83   _hexoticFilesKept=false;
84   _compatibleHypothesis.push_back("Hexotic_Parameters");
85 }
86
87 //=============================================================================
88 /*!
89  *  
90  */
91 //=============================================================================
92
93 HexoticPlugin_Hexotic::~HexoticPlugin_Hexotic()
94 {
95   MESSAGE("HexoticPlugin_Hexotic::~HexoticPlugin_Hexotic");
96 }
97
98 //=============================================================================
99 /*!
100  *  
101  */
102 //=============================================================================
103
104 bool HexoticPlugin_Hexotic::CheckHypothesis( SMESH_Mesh&                          aMesh,
105                                              const TopoDS_Shape&                  aShape,
106                                              SMESH_Hypothesis::Hypothesis_Status& aStatus )
107 {
108   // MESSAGE("HexoticPlugin_Hexotic::CheckHypothesis");
109   _hypothesis = NULL;
110
111   std::list<const SMESHDS_Hypothesis*>::const_iterator itl;
112   const SMESHDS_Hypothesis* theHyp;
113
114   const std::list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
115   int nbHyp = hyps.size();
116   if (!nbHyp) {
117     aStatus = SMESH_Hypothesis::HYP_OK;
118     return true;  // can work with no hypothesis
119   }
120
121   itl = hyps.begin();
122   theHyp = (*itl); // use only the first hypothesis
123
124   std::string hypName = theHyp->GetName();
125   if (hypName == "Hexotic_Parameters") {
126     _hypothesis = static_cast<const HexoticPlugin_Hypothesis*> (theHyp);
127     ASSERT(_hypothesis);
128     aStatus = SMESH_Hypothesis::HYP_OK;
129   }
130   else
131     aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
132
133   return aStatus == SMESH_Hypothesis::HYP_OK;
134 }
135
136 //=======================================================================
137 //function : findShape
138 //purpose  :
139 //=======================================================================
140
141 static TopoDS_Shape findShape(SMDS_MeshNode**     t_Node,
142                               TopoDS_Shape        aShape,
143                               const TopoDS_Shape* t_Shape,
144                               double**            t_Box,
145                               const int           nShape) {
146   double *pntCoor;
147   int iShape, nbNode = 8;
148
149   pntCoor = new double[3];
150   for ( int i=0; i<3; i++ ) {
151     pntCoor[i] = 0;
152     for ( int j=0; j<nbNode; j++ ) {
153       if ( i == 0) pntCoor[i] += t_Node[j]->X();
154       if ( i == 1) pntCoor[i] += t_Node[j]->Y();
155       if ( i == 2) pntCoor[i] += t_Node[j]->Z();
156     }
157     pntCoor[i] /= nbNode;
158   }
159
160   gp_Pnt aPnt(pntCoor[0], pntCoor[1], pntCoor[2]);
161   BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
162   if ( not(SC.State() == TopAbs_IN) ) {
163     for (iShape = 0; iShape < nShape; iShape++) {
164       aShape = t_Shape[iShape];
165       if ( not( pntCoor[0] < t_Box[iShape][0] || t_Box[iShape][1] < pntCoor[0] ||
166                 pntCoor[1] < t_Box[iShape][2] || t_Box[iShape][3] < pntCoor[1] ||
167                 pntCoor[2] < t_Box[iShape][4] || t_Box[iShape][5] < pntCoor[2]) ) {
168         BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
169         if (SC.State() == TopAbs_IN)
170           break;
171       }
172     }
173   }
174   delete [] pntCoor;
175   return aShape;
176 }
177
178 //=======================================================================
179 //function : findEdge
180 //purpose  :
181 //=======================================================================
182
183 static int findEdge(const SMDS_MeshNode* aNode,
184                     const SMESHDS_Mesh*  theMesh,
185                     const int            nEdge,
186                     const TopoDS_Shape*  t_Edge) {
187
188   TopoDS_Shape aPntShape, foundEdge;
189   TopoDS_Vertex aVertex;
190   gp_Pnt aPnt( aNode->X(), aNode->Y(), aNode->Z() );
191
192   int foundInd, ind;
193   double nearest = RealLast(), *t_Dist;
194   double epsilon = Precision::Confusion();
195
196   t_Dist = new double[ nEdge ];
197   aPntShape = BRepBuilderAPI_MakeVertex( aPnt ).Shape();
198   aVertex   = TopoDS::Vertex( aPntShape );
199
200   for ( ind=0; ind < nEdge; ind++ ) {
201     BRepExtrema_DistShapeShape aDistance ( aVertex, t_Edge[ind] );
202     t_Dist[ind] = aDistance.Value();
203     if ( t_Dist[ind] < nearest ) {
204       nearest   = t_Dist[ind];
205       foundEdge = t_Edge[ind];
206       foundInd  = ind;
207       if ( nearest < epsilon )
208         ind = nEdge;
209     }
210   }
211
212   delete [] t_Dist;
213   return theMesh->ShapeToIndex( foundEdge );
214 }
215
216 //=======================================================================
217 //function : getNbShape
218 //purpose  :
219 //=======================================================================
220
221 static int getNbShape(std::string aFile, std::string aString) {
222   int number;
223   std::string aLine;
224   std::ifstream file(aFile.c_str());
225   while ( !file.eof() ) {
226     getline( file, aLine);
227     if ( aLine == aString ) {
228       getline( file, aLine);
229       std::istringstream stringFlux( aLine );
230       stringFlux >> number;
231       break;
232     }
233   }
234   file.close();
235   return number;
236 }
237
238 //=======================================================================
239 //function : countShape
240 //purpose  :
241 //=======================================================================
242
243 template < class Mesh, class Shape >
244 static int countShape( Mesh* mesh, Shape shape ) {
245   TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
246   TopTools_MapOfShape mapShape;
247   int nbShape = 0;
248   for ( ; expShape.More(); expShape.Next() ) {
249     if (mapShape.Add(expShape.Current())) {
250       nbShape++;
251     }
252   }
253   return nbShape;
254 }
255
256 //=======================================================================
257 //function : getShape
258 //purpose  :
259 //=======================================================================
260
261 template < class Mesh, class Shape, class Tab >
262 void getShape(Mesh* mesh, Shape shape, Tab *t_Shape) {
263   TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
264   TopTools_MapOfShape mapShape;
265   for ( int i=0; expShape.More(); expShape.Next() ) {
266     if (mapShape.Add(expShape.Current())) {
267       t_Shape[i] = expShape.Current();
268       i++;
269     }
270   }
271   return;
272 }
273
274 //=======================================================================
275 //function : printWarning
276 //purpose  :
277 //=======================================================================
278
279 static void printWarning(const int nbExpected, std::string aString, const int nbFound) {
280   cout << std::endl;
281   cout << "WARNING : " << nbExpected << " " << aString << " expected, Hexotic has found " << nbFound << std::endl;
282   cout << "=======" << std::endl;
283   cout << std::endl;
284   return;
285 }
286
287 //=======================================================================
288 //function : removeHexoticFiles
289 //purpose  :
290 //=======================================================================
291
292 static void removeHexoticFiles(TCollection_AsciiString file_In, TCollection_AsciiString file_Out) {
293   OSD_File( file_In  ).Remove();
294   OSD_File( file_Out ).Remove();
295 }
296
297 //=======================================================================
298 //function : writeHexoticFile
299 //purpose  : 
300 //=======================================================================
301
302 static bool writeHexoticFile (std::ofstream&                       theFile,
303                               const SMESHDS_Mesh*                  theMesh,
304                               std::map <int,int>&                  theSmdsToHexoticIdMap,
305                               std::map <int,const SMDS_MeshNode*>& theHexoticIdToNodeMap,
306                               const TCollection_AsciiString&       Hexotic_In) {
307   cout << std::endl;
308   cout << "Creating Hexotic processed mesh file : " << Hexotic_In << std::endl;
309
310   int nbShape = 0;
311
312   TopExp_Explorer expface(theMesh->ShapeToMesh(), TopAbs_FACE);
313   for ( ; expface.More(); expface.Next() )
314     nbShape++;
315
316   int *tabID;
317   int *tabNodeId;
318   TopoDS_Shape *tabShape, aShape;
319
320   int shapeID;
321   bool idFound;
322
323   int nbVertices    = 0;
324   int nbTriangles   = 0;
325   const char* space = "  ";
326   int dummy_1D      = 0;
327   int dummy_2D;
328
329   int aSmdsNodeID = 1;
330   const SMDS_MeshNode* aNode;
331   SMDS_NodeIteratorPtr itOnNode;
332
333   std::list< const SMDS_MeshElement* > faces;
334   std::list< const SMDS_MeshElement* >::iterator itListFace;
335   const SMDS_MeshElement* aFace;
336   SMESHDS_SubMesh* theSubMesh;
337   std::map<int,int>::const_iterator itOnSmdsNode;
338   SMDS_ElemIteratorPtr itOnSubNode, itOnSubFace;
339
340 // Writing SMESH points into Hexotic File
341
342   nbVertices = theMesh->NbNodes();
343
344   theFile << "MeshVersionFormatted 1" << std::endl;
345   theFile << std::endl;
346   theFile << "Dimension" << std::endl;
347   theFile << 3 << std::endl;
348   theFile << "# Set of mesh vertices" << std::endl;
349   theFile << "Vertices" << std::endl;
350   theFile << nbVertices << std::endl;
351
352   tabID     = new int[nbShape];
353   tabNodeId = new int[ nbVertices ];
354   tabShape  = new TopoDS_Shape[nbShape];
355
356   itOnNode = theMesh->nodesIterator();
357   while ( itOnNode->more() ) {
358       aNode = itOnNode->next();
359       dummy_1D = aNode->GetPosition()->GetShapeId();
360       tabNodeId[ aSmdsNodeID - 1 ] = 0;
361       idFound  = false;
362       for ( int j=0; j< aSmdsNodeID; j++ ) {
363         if ( dummy_1D == tabNodeId[j] ) {
364           idFound = true;
365           break;
366         }
367       }
368       if ( not idFound )
369         tabNodeId[ aSmdsNodeID - 1 ] = dummy_1D;
370       theSmdsToHexoticIdMap.insert(std::map <int,int>::value_type( aNode->GetID(), aSmdsNodeID ));
371       theHexoticIdToNodeMap.insert(std::map <int,const SMDS_MeshNode*>::value_type( aSmdsNodeID, aNode ));
372       aSmdsNodeID++;
373       theFile << aNode->X() << space << aNode->Y() << space << aNode->Z() << space << dummy_1D << std::endl;
374       }
375
376 // Writing SMESH faces into Hexotic File
377
378   nbTriangles = theMesh->NbFaces();
379
380   theFile << std::endl;
381   theFile << "# Set of mesh triangles (v1,v2,v3,tag)" << std::endl;
382   theFile << "Triangles" << std::endl;
383   theFile << nbTriangles << std::endl;
384
385   expface.ReInit();
386   for ( int i = 0; expface.More(); expface.Next(), i++ ) {
387     tabID[i] = 0;
388     aShape   = expface.Current();
389     shapeID  = theMesh->ShapeToIndex( aShape );
390     idFound  = false;
391     for ( int j=0; j<=i; j++) {
392       if ( shapeID == tabID[j] ) {
393         idFound = true;
394         break;
395       }
396     }
397     if ( not idFound ) {
398       tabID[i]    = shapeID;
399       tabShape[i] = aShape;
400     }
401   }
402   for ( int i=0; i<nbShape; i++ ) {
403     if ( not (tabID[i] == 0) ) {
404       aShape      = tabShape[i];
405       shapeID     = tabID[i];
406       theSubMesh  = theMesh->MeshElements( aShape );
407       itOnSubFace = theSubMesh->GetElements();
408       while ( itOnSubFace->more() ) {
409         aFace    = itOnSubFace->next();
410         dummy_2D = shapeID;
411         itOnSubNode = aFace->nodesIterator();
412         while ( itOnSubNode->more() ) {
413           aSmdsNodeID  = itOnSubNode->next()->GetID();
414           itOnSmdsNode = theSmdsToHexoticIdMap.find( aSmdsNodeID );
415           ASSERT( itOnSmdsNode != theSmdsToHexoticIdMap.end() );
416           theFile << (*itOnSmdsNode).second << space;
417         }
418         theFile << dummy_2D << std::endl;
419       }
420     }
421   }
422
423   theFile << std::endl;
424   theFile << "End" << std::endl;
425
426   cout << "Processed mesh file created, it contains :" << std::endl;
427   cout << "    " << nbVertices  << " vertices"  << std::endl;
428   cout << "    " << nbTriangles << " triangles" << std::endl;
429   cout << std::endl;
430
431   delete [] tabID;
432   delete [] tabNodeId;
433   delete [] tabShape;
434
435   return true;
436 }
437
438 //=======================================================================
439 //function : readResult
440 //purpose  : 
441 //=======================================================================
442
443 static bool readResult(std::string         theFile,
444                        SMESHDS_Mesh*       theMesh,
445                        const int           nbShape,
446                        const TopoDS_Shape* tabShape,
447                        double**            tabBox)
448 {
449   // ---------------------------------
450   // Read generated elements and nodes
451   // ---------------------------------
452
453   TopoDS_Shape aShape;
454   TopoDS_Vertex aVertex;
455   std::string token;
456   int EndOfFile = 0, nbElem = 0, nField = 9, nbRef = 0;
457   int aHexoticNodeID = 0, shapeID, hexoticShapeID;
458   int IdShapeRef = 2;
459   int *tabID, *tabRef, *nodeAssigne;
460   bool *tabDummy, hasDummy = false;
461   double epsilon = Precision::Confusion();
462   std::map <std::string,int> mapField;
463   SMDS_MeshNode** HexoticNode;
464   TopoDS_Shape *tabCorner, *tabEdge;
465
466   tabID    = new int[nbShape];
467   tabRef   = new int[nField];
468   tabDummy = new bool[nField];
469
470   for (int i=0; i<nbShape; i++)
471     tabID[i] = 0;
472   if ( nbShape == 1 )
473     tabID[0] = theMesh->ShapeToIndex( tabShape[0] );
474
475   mapField["MeshVersionFormatted"] = 0; tabRef[0] = 0; tabDummy[0] = false;
476   mapField["Dimension"]            = 1; tabRef[1] = 0; tabDummy[1] = false;
477   mapField["Vertices"]             = 2; tabRef[2] = 3; tabDummy[2] = true;
478   mapField["Corners"]              = 3; tabRef[3] = 1; tabDummy[3] = false;
479   mapField["Edges"]                = 4; tabRef[4] = 2; tabDummy[4] = true;
480   mapField["Ridges"]               = 5; tabRef[5] = 1; tabDummy[5] = false;
481   mapField["Quadrilaterals"]       = 6; tabRef[6] = 4; tabDummy[6] = true;
482   mapField["Hexahedra"]            = 7; tabRef[7] = 8; tabDummy[7] = true;
483   mapField["End"]                  = 8; tabRef[8] = 0; tabDummy[0] = false;
484
485   SMDS_NodeIteratorPtr itOnHexoticInputNode = theMesh->nodesIterator();
486   while ( itOnHexoticInputNode->more() )
487     theMesh->RemoveNode( itOnHexoticInputNode->next() );
488
489   int nbVertices   = getNbShape(theFile, "Vertices");
490   int nbHexCorners = getNbShape(theFile, "Corners");
491   int nbCorners    = countShape( theMesh, TopAbs_VERTEX );
492   int nbShapeEdge  = countShape( theMesh, TopAbs_EDGE );
493
494   if ( nbHexCorners != nbCorners ) {
495     printWarning(nbCorners, "corners", nbHexCorners);
496     if ( nbHexCorners > nbCorners )
497       nbCorners = nbHexCorners;
498   }
499
500   tabCorner   = new TopoDS_Shape[ nbCorners ];
501   tabEdge     = new TopoDS_Shape[ nbShapeEdge ];
502   nodeAssigne = new int[ nbVertices + 1 ];
503   HexoticNode = new SMDS_MeshNode*[ nbVertices + 1 ];
504
505   getShape(theMesh, TopAbs_VERTEX, tabCorner);
506   getShape(theMesh, TopAbs_EDGE,   tabEdge);
507
508   MESSAGE("Read " << theFile << " file");
509   std::ifstream fileRes(theFile.c_str());
510   ASSERT(fileRes);
511
512   while ( EndOfFile == 0  ) {
513     int dummy;
514     fileRes >> token;
515
516     if (mapField.count(token)) {
517       nField   = mapField[token];
518       nbRef    = tabRef[nField];
519       hasDummy = tabDummy[nField];
520     }
521     else {
522       nField = -1;
523       nbRef = 0;
524     }
525
526     nbElem = 0;
527     if ( nField < (mapField.size() - 1) && nField >= 0 )
528       fileRes >> nbElem;
529
530     switch (nField) {
531       case 0: { // "MeshVersionFormatted"
532         MESSAGE(token << " " << nbElem);
533         break;
534       }
535       case 1: { // "Dimension"
536         MESSAGE("Mesh dimension " << nbElem << "D");
537         break;
538       }
539       case 2: { // "Vertices"
540         MESSAGE("Read " << nbElem << " " << token);
541         int aHexoticID;
542         double *coord;
543         SMDS_MeshNode * aHexoticNode;
544
545         coord = new double[nbRef];
546         for ( int iElem = 0; iElem < nbElem; iElem++ ) {
547           aHexoticID = iElem + 1;
548           for ( int iCoord = 0; iCoord < 3; iCoord++ )
549             fileRes >> coord[ iCoord ];
550           fileRes >> dummy;
551           aHexoticNode = theMesh->AddNode(coord[0], coord[1], coord[2]);
552           HexoticNode[ aHexoticID ] = aHexoticNode;
553           nodeAssigne[ aHexoticID ] = 0;
554         }
555         delete [] coord;
556         break;
557       }
558       case 3: // "Corners"
559       case 4: // "Edges"
560       case 5: // "Ridges"
561       case 6: // "Quadrilaterals"
562       case 7: { // "Hexahedra"
563         MESSAGE("Read " << nbElem << " " << token);
564         SMDS_MeshNode** node;
565         int nodeDim, *nodeID;
566         SMDS_MeshElement * aHexoticElement;
567
568         node   = new SMDS_MeshNode*[ nbRef ];
569         nodeID = new int[ nbRef ];
570         for ( int iElem = 0; iElem < nbElem; iElem++ ) {
571           for ( int iRef = 0; iRef < nbRef; iRef++ ) {
572             fileRes >> aHexoticNodeID;                          // read nbRef aHexoticNodeID
573             node[ iRef ]   = HexoticNode[ aHexoticNodeID ];
574             nodeID[ iRef ] = aHexoticNodeID;
575           }
576           if ( hasDummy )
577             fileRes >> dummy;
578           switch (nField) {
579             case 3: { // "Corners"
580               nodeDim = 1;
581               gp_Pnt HexoticPnt ( node[0]->X(), node[0]->Y(), node[0]->Z() );
582               for ( int i=0; i<nbElem; i++ ) {
583                 aVertex = TopoDS::Vertex( tabCorner[i] );
584                 gp_Pnt aPnt = BRep_Tool::Pnt( aVertex );
585                 if ( aPnt.Distance( HexoticPnt ) < epsilon )
586                   break;
587               }
588               break;
589             }
590             case 4: { // "Edges"
591               nodeDim = 2;
592               aHexoticElement = theMesh->AddEdge( node[0], node[1] );
593               int iNode = 1;
594               if ( nodeAssigne[ nodeID[0] ] == 0 || nodeAssigne[ nodeID[0] ] == 2 )
595                 iNode = 0;
596               shapeID = findEdge( node[iNode], theMesh, nbShapeEdge, tabEdge );
597               break;
598             }
599             case 5: { // "Ridges"
600               break;
601             }
602             case 6: { // "Quadrilaterals"
603               nodeDim = 3;
604               aHexoticElement = theMesh->AddFace( node[0], node[1], node[2], node[3] );
605               shapeID = dummy;
606               break;
607             }
608             case 7: { // "Hexahedra"
609               nodeDim = 4;
610               aHexoticElement = theMesh->AddVolume( node[0], node[3], node[2], node[1], node[4], node[7], node[6], node[5] );
611               if ( nbShape > 1 ) {
612                 hexoticShapeID = dummy - IdShapeRef;
613                 if ( tabID[ hexoticShapeID ] == 0 ) {
614                   if (iElem == 0)
615                     aShape = tabShape[0];
616                   aShape = findShape(node, aShape, tabShape, tabBox, nbShape);
617                   shapeID = theMesh->ShapeToIndex( aShape );
618                   tabID[ hexoticShapeID ] = shapeID;
619                 }
620                 else
621                   shapeID = tabID[ hexoticShapeID ];
622                 if ( iElem == (nbElem - 1) ) {
623                   int shapeAssociated = 0;
624                   for ( int i=0; i<nbShape; i++ ) {
625                     if (tabID[i] != 0 )
626                       shapeAssociated += 1;
627                   }
628                   if ( shapeAssociated != nbShape )
629                     printWarning(nbShape, "domains", shapeAssociated);
630                 }
631               }
632               else {
633                 shapeID = tabID[0];
634               }
635               break;
636             }
637           }
638           if ( token != "Ridges" ) {
639             for ( int i=0; i<nbRef; i++ ) {
640               if ( nodeAssigne[ nodeID[i] ] == 0 ) {
641                 if      ( token == "Corners" )        theMesh->SetNodeOnVertex( node[0], aVertex );
642                 else if ( token == "Edges" )          theMesh->SetNodeOnEdge( node[i], shapeID );
643                 else if ( token == "Quadrilaterals" ) theMesh->SetNodeOnFace( node[i], shapeID );
644                 else if ( token == "Hexahedra" )      theMesh->SetNodeInVolume( node[i], shapeID );
645                 nodeAssigne[ nodeID[i] ] = nodeDim;
646               }
647             }
648             if ( token != "Corners" )
649               theMesh->SetMeshElementOnShape( aHexoticElement, shapeID );
650           }
651         }
652         delete [] node;
653         delete [] nodeID;
654         break;
655       }
656       case 8: { // "End"
657         EndOfFile = 1;
658         MESSAGE("End of " << theFile << " file");
659         break;
660       }
661       default: {
662         MESSAGE("Unknown Token: " << token);
663       }
664     }
665   }
666   cout << std::endl;
667   delete [] tabID;
668   delete [] tabRef;
669   delete [] tabDummy;
670   delete [] tabCorner;
671   delete [] tabEdge;
672   delete [] nodeAssigne;
673   delete [] HexoticNode;
674   return true;
675 }
676
677 //=============================================================================
678 /*!
679  * Pass parameters to Hexotic
680  */
681 //=============================================================================
682
683 void HexoticPlugin_Hexotic::SetParameters(const HexoticPlugin_Hypothesis* hyp) {
684
685   MESSAGE("HexoticPlugin_Hexotic::SetParameters");
686   if (hyp) {
687     _hexesMinLevel = hyp->GetHexesMinLevel();
688     _hexesMaxLevel = hyp->GetHexesMaxLevel();
689     _hexoticQuadrangles = hyp->GetHexoticQuadrangles();
690     _hexoticIgnoreRidges = hyp->GetHexoticIgnoreRidges();
691     _hexoticInvalidElements = hyp->GetHexoticInvalidElements();
692     _hexoticSharpAngleThreshold = hyp->GetHexoticSharpAngleThreshold();
693   }
694   else {
695     cout << std::endl;
696     cout << "WARNING : The Hexotic default parameters are taken into account" << std::endl;
697     cout << "=======" << std::endl;
698     _hexesMinLevel = hyp->GetDefaultHexesMinLevel();
699     _hexesMaxLevel = hyp->GetDefaultHexesMaxLevel();
700     _hexoticQuadrangles = hyp->GetDefaultHexoticQuadrangles();
701     _hexoticIgnoreRidges = hyp->GetDefaultHexoticIgnoreRidges();
702     _hexoticInvalidElements = hyp->GetDefaultHexoticInvalidElements();
703     _hexoticSharpAngleThreshold = hyp->GetDefaultHexoticSharpAngleThreshold();
704   }
705 }
706
707 //=======================================================================
708 //function : getTmpDir
709 //purpose  :
710 //=======================================================================
711
712 static TCollection_AsciiString getTmpDir()
713 {
714   TCollection_AsciiString aTmpDir;
715
716   char *Tmp_dir = getenv("SALOME_TMP_DIR");
717   if(Tmp_dir != NULL) {
718     aTmpDir = Tmp_dir;
719     #ifdef WIN32
720       if(aTmpDir.Value(aTmpDir.Length()) != '\\') aTmpDir+='\\';
721     #else
722       if(aTmpDir.Value(aTmpDir.Length()) != '/') aTmpDir+='/';
723     #endif      
724   }
725   else {
726     #ifdef WIN32
727       aTmpDir = TCollection_AsciiString("C:\\");
728     #else
729       aTmpDir = TCollection_AsciiString("/tmp/");
730     #endif
731   }
732   return aTmpDir;
733 }
734
735 //=============================================================================
736 /*!
737  * Here we are going to use the Hexotic mesher
738  */
739 //=============================================================================
740
741 bool HexoticPlugin_Hexotic::Compute(SMESH_Mesh&          theMesh,
742                                      const TopoDS_Shape& theShape)
743 {
744   bool Ok = true;
745   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
746   TCollection_AsciiString hexahedraMessage;
747
748   if (_iShape == 0 && _nbShape == 0) {
749     _nbShape = countShape( meshDS, TopAbs_SOLID );  // we count the number of shapes
750     _tabNode = new SMDS_MeshNode*[_nbShape];        // we declare the size of the node array
751   }
752
753   // to prevent from displaying error message after computing,
754   // we need to create one node for each shape theShape.
755
756   _tabNode[_iShape] = meshDS->AddNode(0, 0, 0);
757   meshDS->SetNodeInVolume( _tabNode[_iShape], meshDS->ShapeToIndex(theShape) );
758
759   _iShape++;
760
761   if (_iShape == _nbShape ) {
762
763     for (int i=0; i<_nbShape; i++)        // we destroy the (_nbShape - 1) nodes created and used
764       meshDS->RemoveNode( _tabNode[i] );  // to simulate successful mesh computing.
765     delete [] _tabNode;
766
767     // create bounding box for each shape of the compound
768
769     int iShape = 0;
770     TopoDS_Shape *tabShape;
771     double **tabBox;
772
773     tabShape = new TopoDS_Shape[_nbShape];
774     tabBox   = new double*[_nbShape];
775     for (int i=0; i<_nbShape; i++)
776       tabBox[i] = new double[6];
777     double Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
778
779     TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
780     for (; expBox.More(); expBox.Next()) {
781       tabShape[iShape] = expBox.Current();
782       Bnd_Box BoundingBox;
783       BRepBndLib::Add(expBox.Current(), BoundingBox);
784       BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
785       tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
786       tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
787       tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
788       iShape++;
789     }
790
791     SetParameters(_hypothesis);
792
793     cout << std::endl;
794     cout << "Hexotic execution..." << std::endl;
795     cout << _name << " parameters :" << std::endl;
796     cout << "    " << _name << " Segments Min Level = " << _hexesMinLevel << std::endl;
797     cout << "    " << _name << " Segments Max Level = " << _hexesMaxLevel << std::endl;
798     cout << "    " << "Salome Quadrangles : " << (_hexoticQuadrangles ? "yes":"no") << std::endl;
799     cout << "    " << "Hexotic can ignore ridges : " << (_hexoticIgnoreRidges ? "yes":"no") << std::endl;
800     cout << "    " << "Hexotic authorize invalide elements : " << ( _hexoticInvalidElements ? "yes":"no") << std::endl;
801     cout << "    " << _name << " Sharp angle threshold = " << _hexoticSharpAngleThreshold << " degrees" << std::endl;
802
803     TCollection_AsciiString aTmpDir = getTmpDir();
804     TCollection_AsciiString Hexotic_In, Hexotic_Out;
805     TCollection_AsciiString run_Hexotic( "hexotic" );
806
807     TCollection_AsciiString minl = " -minl ", maxl = " -maxl ", angle = " -ra ";
808     TCollection_AsciiString in   = " -in ",   out  = " -out ";
809     TCollection_AsciiString ignoreRidges = " -nr ", invalideElements = " -inv ";
810     TCollection_AsciiString subdom = " -sd ";
811
812     TCollection_AsciiString minLevel, maxLevel, sharpAngle, mode;
813     minLevel = _hexesMinLevel;
814     maxLevel = _hexesMaxLevel;
815     sharpAngle = _hexoticSharpAngleThreshold;
816     mode = 4;
817
818     std::map <int,int> aSmdsToHexoticIdMap;
819     std::map <int,const SMDS_MeshNode*> aHexoticIdToNodeMap;
820
821     Hexotic_In  = aTmpDir + "Hexotic_In.mesh";
822     Hexotic_Out = aTmpDir + "Hexotic_Out.mesh";
823
824     if (_hexoticIgnoreRidges)
825       run_Hexotic +=  ignoreRidges;
826
827     if (_hexoticInvalidElements)
828       run_Hexotic +=  invalideElements;
829
830     run_Hexotic += angle + sharpAngle + minl + minLevel + maxl + maxLevel + in + Hexotic_In + out + Hexotic_Out;
831     run_Hexotic += subdom + mode;
832
833     cout << std::endl;
834     cout << "Hexotic command : " << run_Hexotic << std::endl;
835
836     removeHexoticFiles(Hexotic_In, Hexotic_Out);
837
838     std::ofstream HexoticFile (Hexotic_In.ToCString(), std::ios::out);
839
840     Ok = ( writeHexoticFile(HexoticFile, meshDS, aSmdsToHexoticIdMap, aHexoticIdToNodeMap, Hexotic_In) );
841
842     HexoticFile.close();
843     aSmdsToHexoticIdMap.clear();
844     aHexoticIdToNodeMap.clear();
845
846     MESSAGE("HexoticPlugin_Hexotic::Compute");
847
848     system( run_Hexotic.ToCString() );
849
850     // --------------
851     // read a result
852     // --------------
853
854     std::ifstream fileRes( Hexotic_Out.ToCString() );
855     if ( ! fileRes.fail() ) {
856       Ok = readResult( Hexotic_Out.ToCString(), meshDS, _nbShape, tabShape, tabBox );
857       hexahedraMessage = "success";
858     }
859     else {
860       hexahedraMessage = "failed";
861       cout << "Problem with Hexotic output file " << Hexotic_Out.ToCString() << std::endl;
862       Ok = false;
863     }
864     cout << "Hexahedra meshing " << hexahedraMessage << std::endl;
865     cout << std::endl;
866
867     delete [] tabShape;
868     for (int i=0; i<_nbShape; i++)
869       delete [] tabBox[i];
870     delete [] tabBox;
871     _nbShape = 0;
872     _iShape  = 0;
873   }
874   return Ok;
875 }
876
877 //=============================================================================
878 /*!
879  *  
880  */
881 //=============================================================================
882
883 std::ostream& HexoticPlugin_Hexotic::SaveTo(std::ostream& save)
884 {
885   return save;
886 }
887
888 //=============================================================================
889 /*!
890  *  
891  */
892 //=============================================================================
893
894 std::istream& HexoticPlugin_Hexotic::LoadFrom(std::istream& load)
895 {
896   return load;
897 }
898
899 //=============================================================================
900 /*!
901  *  
902  */
903 //=============================================================================
904
905 std::ostream& operator << (std::ostream& save, HexoticPlugin_Hexotic& hyp)
906 {
907   return hyp.SaveTo( save );
908 }
909
910 //=============================================================================
911 /*!
912  *  
913  */
914 //=============================================================================
915
916 std::istream& operator >> (std::istream& load, HexoticPlugin_Hexotic& hyp)
917 {
918   return hyp.LoadFrom( load );
919 }