Salome HOME
Merge from BR_V5_DEV 16Feb09
[modules/smesh.git] / src / SMESH / SMESH_Block.cxx
1 //  Copyright (C) 2007-2008  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 //  Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 //  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 //  This library is free software; you can redistribute it and/or
7 //  modify it under the terms of the GNU Lesser General Public
8 //  License as published by the Free Software Foundation; either
9 //  version 2.1 of the License.
10 //
11 //  This library is distributed in the hope that it will be useful,
12 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
13 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 //  Lesser General Public License for more details.
15 //
16 //  You should have received a copy of the GNU Lesser General Public
17 //  License along with this library; if not, write to the Free Software
18 //  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 //  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 // File      : SMESH_Pattern.hxx
23 // Created   : Mon Aug  2 10:30:00 2004
24 // Author    : Edward AGAPOV (eap)
25 //
26 #include "SMESH_Block.hxx"
27
28 #include <BRepAdaptor_Curve.hxx>
29 #include <BRepAdaptor_Curve2d.hxx>
30 #include <BRepAdaptor_Surface.hxx>
31 #include <BRepTools.hxx>
32 #include <BRepTools_WireExplorer.hxx>
33 #include <BRep_Builder.hxx>
34 #include <BRep_Tool.hxx>
35 #include <Bnd_Box.hxx>
36 #include <Extrema_ExtPC.hxx>
37 #include <Extrema_POnCurv.hxx>
38 #include <Geom2d_Curve.hxx>
39 #include <TopExp_Explorer.hxx>
40 #include <TopTools_IndexedDataMapOfShapeListOfShape.hxx>
41 #include <TopTools_ListIteratorOfListOfShape.hxx>
42 #include <TopTools_ListOfShape.hxx>
43 #include <TopoDS.hxx>
44 #include <TopoDS_Compound.hxx>
45 #include <TopoDS_Face.hxx>
46 #include <TopoDS_Iterator.hxx>
47 #include <TopoDS_Wire.hxx>
48 #include <gp_Trsf.hxx>
49 #include <gp_Vec.hxx>
50 #include <math_FunctionSetRoot.hxx>
51 #include <math_Matrix.hxx>
52 #include <math_Vector.hxx>
53
54 #include "SMDS_MeshNode.hxx"
55 #include "SMDS_MeshVolume.hxx"
56 #include "SMDS_VolumeTool.hxx"
57 #include "utilities.h"
58
59 #include <list>
60
61 using namespace std;
62
63 //#define DEBUG_PARAM_COMPUTE
64
65 //================================================================================
66 /*!
67  * \brief Set edge data
68   * \param edgeID - block subshape ID
69   * \param curve - edge geometry
70   * \param isForward - is curve orientation coincides with edge orientation in the block
71  */
72 //================================================================================
73
74 void SMESH_Block::TEdge::Set( const int edgeID, Adaptor3d_Curve* curve, const bool isForward )
75 {
76   myCoordInd = SMESH_Block::GetCoordIndOnEdge( edgeID );
77   if ( myC3d ) delete myC3d;
78   myC3d = curve;
79   myFirst = curve->FirstParameter();
80   myLast = curve->LastParameter();
81   if ( !isForward )
82     std::swap( myFirst, myLast );
83 }
84
85 //================================================================================
86 /*!
87  * \brief Set coordinates of nodes at edge ends to work with mesh block
88   * \param edgeID - block subshape ID
89   * \param node1 - coordinates of node with lower ID
90   * \param node2 - coordinates of node with upper ID
91  */
92 //================================================================================
93
94 void SMESH_Block::TEdge::Set( const int edgeID, const gp_XYZ& node1, const gp_XYZ& node2 )
95 {
96   myCoordInd = SMESH_Block::GetCoordIndOnEdge( edgeID );
97   myNodes[ 0 ] = node1;
98   myNodes[ 1 ] = node2;
99
100   if ( myC3d ) delete myC3d;
101   myC3d = 0;
102 }
103
104 //=======================================================================
105 //function : SMESH_Block::TEdge::GetU
106 //purpose  : 
107 //=======================================================================
108
109 double SMESH_Block::TEdge::GetU( const gp_XYZ& theParams ) const
110 {
111   double u = theParams.Coord( myCoordInd );
112   if ( !myC3d ) // if mesh block
113     return u;
114   return ( 1 - u ) * myFirst + u * myLast;
115 }
116
117 //=======================================================================
118 //function : SMESH_Block::TEdge::Point
119 //purpose  : 
120 //=======================================================================
121
122 gp_XYZ SMESH_Block::TEdge::Point( const gp_XYZ& theParams ) const
123 {
124   double u = GetU( theParams );
125   if ( myC3d ) return myC3d->Value( u ).XYZ();
126   // mesh block
127   return myNodes[0] * ( 1 - u ) + myNodes[1] * u;
128 }
129
130 //================================================================================
131 /*!
132  * \brief Destructor
133  */
134 //================================================================================
135
136 SMESH_Block::TEdge::~TEdge()
137 {
138   if ( myC3d ) delete myC3d;
139 }
140
141 //================================================================================
142 /*!
143  * \brief Set face data
144   * \param faceID - block subshape ID
145   * \param S - face surface geometry
146   * \param c2d - 4 pcurves in the order as returned by GetFaceEdgesIDs(faceID)
147   * \param isForward - orientation of pcurves comparing with block edge direction
148  */
149 //================================================================================
150
151 void SMESH_Block::TFace::Set( const int          faceID,
152                               Adaptor3d_Surface* S,
153                               Adaptor2d_Curve2d* c2D[4],
154                               const bool         isForward[4] )
155 {
156   if ( myS ) delete myS;
157   myS = S;
158   // pcurves
159   vector< int > edgeIdVec;
160   GetFaceEdgesIDs( faceID, edgeIdVec );
161   for ( int iE = 0; iE < edgeIdVec.size(); iE++ ) // loop on 4 edges
162   {
163     myCoordInd[ iE ] = GetCoordIndOnEdge( edgeIdVec[ iE ] );
164     if ( myC2d[ iE ]) delete myC2d[ iE ];
165     myC2d[ iE ] = c2D[ iE ];
166     myFirst[ iE ] = myC2d[ iE ]->FirstParameter();
167     myLast [ iE ] = myC2d[ iE ]->LastParameter();
168     if ( !isForward[ iE ])
169       std::swap( myFirst[ iE ], myLast[ iE ] );
170   }
171   // 2d corners
172   myCorner[ 0 ] = myC2d[ 0 ]->Value( myFirst[0] ).XY();
173   myCorner[ 1 ] = myC2d[ 0 ]->Value( myLast[0] ).XY();
174   myCorner[ 2 ] = myC2d[ 1 ]->Value( myLast[1] ).XY();
175   myCorner[ 3 ] = myC2d[ 1 ]->Value( myFirst[1] ).XY();
176 }
177
178 //================================================================================
179 /*!
180  * \brief Set face data to work with mesh block
181   * \param faceID - block subshape ID
182   * \param edgeU0 - filled data of edge u0 = GetFaceEdgesIDs(faceID)[ 0 ]
183   * \param edgeU1 - filled data of edge u1 = GetFaceEdgesIDs(faceID)[ 1 ]
184  */
185 //================================================================================
186
187 void SMESH_Block::TFace::Set( const int faceID, const TEdge& edgeU0, const TEdge& edgeU1 )
188 {
189   vector< int > edgeIdVec;
190   GetFaceEdgesIDs( faceID, edgeIdVec );
191   myNodes[ 0 ] = edgeU0.NodeXYZ( 1 );
192   myNodes[ 1 ] = edgeU0.NodeXYZ( 0 );
193   myNodes[ 2 ] = edgeU1.NodeXYZ( 0 );
194   myNodes[ 3 ] = edgeU1.NodeXYZ( 1 );
195   myCoordInd[ 0 ] = GetCoordIndOnEdge( edgeIdVec[ 0 ] );
196   myCoordInd[ 1 ] = GetCoordIndOnEdge( edgeIdVec[ 1 ] );
197   myCoordInd[ 2 ] = GetCoordIndOnEdge( edgeIdVec[ 2 ] );
198   myCoordInd[ 3 ] = GetCoordIndOnEdge( edgeIdVec[ 3 ] );
199   if ( myS ) delete myS;
200   myS = 0;
201 }
202
203 //================================================================================
204 /*!
205  * \brief Destructor
206  */
207 //================================================================================
208
209 SMESH_Block::TFace::~TFace()
210 {
211   if ( myS ) delete myS;
212   for ( int i = 0 ; i < 4; ++i )
213     if ( myC2d[ i ]) delete myC2d[ i ];
214 }
215
216 //=======================================================================
217 //function : SMESH_Block::TFace::GetCoefs
218 //purpose  : return coefficients for addition of [0-3]-th edge and vertex
219 //=======================================================================
220
221 void SMESH_Block::TFace::GetCoefs(int           iE,
222                                   const gp_XYZ& theParams,
223                                   double&       Ecoef,
224                                   double&       Vcoef ) const
225 {
226   double dU = theParams.Coord( GetUInd() );
227   double dV = theParams.Coord( GetVInd() );
228   switch ( iE ) {
229   case 0:
230     Ecoef = ( 1 - dV ); // u0
231     Vcoef = ( 1 - dU ) * ( 1 - dV ); break; // 00
232   case 1:
233     Ecoef = dV; // u1
234     Vcoef = dU * ( 1 - dV ); break; // 10
235   case 2:
236     Ecoef = ( 1 - dU ); // 0v
237     Vcoef = dU * dV  ; break; // 11
238   case 3:
239     Ecoef = dU  ; // 1v
240     Vcoef = ( 1 - dU ) * dV  ; break; // 01
241   default: ASSERT(0);
242   }
243 }
244
245 //=======================================================================
246 //function : SMESH_Block::TFace::GetUV
247 //purpose  : 
248 //=======================================================================
249
250 gp_XY SMESH_Block::TFace::GetUV( const gp_XYZ& theParams ) const
251 {
252   gp_XY uv(0.,0.);
253   for ( int iE = 0; iE < 4; iE++ ) // loop on 4 edges
254   {
255     double Ecoef = 0, Vcoef = 0;
256     GetCoefs( iE, theParams, Ecoef, Vcoef );
257     // edge addition
258     double u = theParams.Coord( myCoordInd[ iE ] );
259     u = ( 1 - u ) * myFirst[ iE ] + u * myLast[ iE ];
260     uv += Ecoef * myC2d[ iE ]->Value( u ).XY();
261     // corner addition
262     uv -= Vcoef * myCorner[ iE ];
263   }
264   return uv;
265 }
266
267 //=======================================================================
268 //function : SMESH_Block::TFace::Point
269 //purpose  : 
270 //=======================================================================
271
272 gp_XYZ SMESH_Block::TFace::Point( const gp_XYZ& theParams ) const
273 {
274   gp_XYZ p(0.,0.,0.);
275   if ( !myS ) // if mesh block
276   {
277     for ( int iE = 0; iE < 4; iE++ ) // loop on 4 edges
278     {
279       double Ecoef = 0, Vcoef = 0;
280       GetCoefs( iE, theParams, Ecoef, Vcoef );
281       // edge addition
282       double u = theParams.Coord( myCoordInd[ iE ] );
283       int i1 = 0, i2 = 1;
284       switch ( iE ) {
285       case 1: i1 = 3; i2 = 2; break;
286       case 2: i1 = 1; i2 = 2; break;
287       case 3: i1 = 0; i2 = 3; break;
288       }
289       p += Ecoef * ( myNodes[ i1 ] * ( 1 - u ) + myNodes[ i2 ] * u );
290       // corner addition
291       p -= Vcoef * myNodes[ iE ];
292     }
293     
294   }
295   else // shape block
296   {
297     gp_XY uv = GetUV( theParams );
298     p = myS->Value( uv.X(), uv.Y() ).XYZ();
299   }
300   return p;
301 }
302
303 //=======================================================================
304 //function : GetShapeCoef
305 //purpose  : 
306 //=======================================================================
307
308 double* SMESH_Block::GetShapeCoef (const int theShapeID)
309 {
310   static double shapeCoef[][3] = {
311     //    V000,        V100,        V010,         V110
312     { -1,-1,-1 }, {  1,-1,-1 }, { -1, 1,-1 }, {  1, 1,-1 },
313     //    V001,        V101,        V011,         V111,
314     { -1,-1, 1 }, {  1,-1, 1 }, { -1, 1, 1 }, {  1, 1, 1 },
315     //    Ex00,        Ex10,        Ex01,         Ex11,
316     {  0,-1,-1 }, {  0, 1,-1 }, {  0,-1, 1 }, {  0, 1, 1 },
317     //    E0y0,        E1y0,        E0y1,         E1y1,
318     { -1, 0,-1 }, {  1, 0,-1 }, { -1, 0, 1 }, {  1, 0, 1 },
319     //    E00z,        E10z,        E01z,         E11z,
320     { -1,-1, 0 }, {  1,-1, 0 }, { -1, 1, 0 }, {  1, 1, 0 },
321     //    Fxy0,        Fxy1,        Fx0z,         Fx1z,         F0yz,           F1yz,
322     {  0, 0,-1 }, {  0, 0, 1 }, {  0,-1, 0 }, {  0, 1, 0 }, { -1, 0, 0 }, {  1, 0, 0 },
323     // ID_Shell
324     {  0, 0, 0 }
325   };
326   if ( theShapeID < ID_V000 || theShapeID > ID_F1yz )
327     return shapeCoef[ ID_Shell - 1 ];
328
329   return shapeCoef[ theShapeID - 1 ];
330 }
331
332 //=======================================================================
333 //function : ShellPoint
334 //purpose  : return coordinates of a point in shell
335 //=======================================================================
336
337 bool SMESH_Block::ShellPoint( const gp_XYZ& theParams, gp_XYZ& thePoint ) const
338 {
339   thePoint.SetCoord( 0., 0., 0. );
340   for ( int shapeID = ID_V000; shapeID < ID_Shell; shapeID++ )
341   {
342     // coef
343     double* coefs = GetShapeCoef( shapeID );
344     double k = 1;
345     for ( int iCoef = 0; iCoef < 3; iCoef++ ) {
346       if ( coefs[ iCoef ] != 0 ) {
347         if ( coefs[ iCoef ] < 0 )
348           k *= ( 1. - theParams.Coord( iCoef + 1 ));
349         else
350           k *= theParams.Coord( iCoef + 1 );
351       }
352     }
353     // add point on a shape
354     if ( fabs( k ) > DBL_MIN )
355     {
356       gp_XYZ Ps;
357       if ( shapeID < ID_Ex00 ) // vertex
358         VertexPoint( shapeID, Ps );
359       else if ( shapeID < ID_Fxy0 ) { // edge
360         EdgePoint( shapeID, theParams, Ps );
361         k = -k;
362       } else // face
363         FacePoint( shapeID, theParams, Ps );
364
365       thePoint += k * Ps;
366     }
367   }
368   return true;
369 }
370
371 //=======================================================================
372 //function : ShellPoint
373 //purpose  : computes coordinates of a point in shell by points on sub-shapes;
374 //           thePointOnShape[ subShapeID ] must be a point on a subShape
375 //=======================================================================
376
377 bool SMESH_Block::ShellPoint(const gp_XYZ&         theParams,
378                              const vector<gp_XYZ>& thePointOnShape,
379                              gp_XYZ&               thePoint )
380 {
381   if ( thePointOnShape.size() < ID_F1yz )
382     return false;
383
384   double x = theParams.X(), y = theParams.Y(), z = theParams.Z();
385   double x1 = 1. - x,       y1 = 1. - y,       z1 = 1. - z;
386   const vector<gp_XYZ>& p = thePointOnShape;
387
388   thePoint = 
389     x1 * p[ID_F0yz] + x * p[ID_F1yz] +
390     y1 * p[ID_Fx0z] + y * p[ID_Fx1z] +
391     z1 * p[ID_Fxy0] + z * p[ID_Fxy1] +
392     x1 * (y1 * (z1 * p[ID_V000] + z * p[ID_V001])  +
393           y  * (z1 * p[ID_V010] + z * p[ID_V011])) +
394     x  * (y1 * (z1 * p[ID_V100] + z * p[ID_V101])  +
395           y  * (z1 * p[ID_V110] + z * p[ID_V111]));
396   thePoint -=
397     x1 * (y1 * p[ID_E00z] + y * p[ID_E01z]) + 
398     x  * (y1 * p[ID_E10z] + y * p[ID_E11z]) + 
399     y1 * (z1 * p[ID_Ex00] + z * p[ID_Ex01]) + 
400     y  * (z1 * p[ID_Ex10] + z * p[ID_Ex11]) + 
401     z1 * (x1 * p[ID_E0y0] + x * p[ID_E1y0]) + 
402     z  * (x1 * p[ID_E0y1] + x * p[ID_E1y1]);
403
404   return true;
405 }
406
407 //=======================================================================
408 //function : Constructor
409 //purpose  : 
410 //=======================================================================
411
412 SMESH_Block::SMESH_Block():
413   myNbIterations(0),
414   mySumDist(0.),
415   myTolerance(-1.) // to be re-initialized
416 {
417 }
418
419
420 //=======================================================================
421 //function : NbVariables
422 //purpose  : 
423 //=======================================================================
424
425 Standard_Integer SMESH_Block::NbVariables() const
426 {
427   return 3;
428 }
429
430 //=======================================================================
431 //function : NbEquations
432 //purpose  : 
433 //=======================================================================
434
435 Standard_Integer SMESH_Block::NbEquations() const
436 {
437   return 1;
438 }
439
440 //=======================================================================
441 //function : Value
442 //purpose  : 
443 //=======================================================================
444
445 Standard_Boolean SMESH_Block::Value(const math_Vector& theXYZ, math_Vector& theFxyz) 
446 {
447   gp_XYZ P, params( theXYZ(1), theXYZ(2), theXYZ(3) );
448   if ( params.IsEqual( myParam, DBL_MIN )) { // same param
449     theFxyz( 1 ) = funcValue( myValues[ SQUARE_DIST ]);
450   }
451   else {
452     ShellPoint( params, P );
453     gp_Vec dP( P - myPoint );
454     theFxyz(1) = funcValue( dP.SquareMagnitude() );
455   }
456   return true;
457 }
458
459 //=======================================================================
460 //function : Derivatives
461 //purpose  : 
462 //=======================================================================
463
464 Standard_Boolean SMESH_Block::Derivatives(const math_Vector& XYZ,math_Matrix& Df) 
465 {
466   math_Vector F(1,3);
467   return Values(XYZ,F,Df);
468 }
469
470 //=======================================================================
471 //function : GetStateNumber
472 //purpose  : 
473 //=======================================================================
474
475 Standard_Integer SMESH_Block::GetStateNumber ()
476 {
477   return 0; //myValues[0] < 1e-1;
478 }
479
480 //=======================================================================
481 //function : Values
482 //purpose  : 
483 //=======================================================================
484
485 Standard_Boolean SMESH_Block::Values(const math_Vector& theXYZ,
486                                      math_Vector&       theFxyz,
487                                      math_Matrix&       theDf) 
488 {
489   gp_XYZ P, params( theXYZ(1), theXYZ(2), theXYZ(3) );
490   if ( params.IsEqual( myParam, DBL_MIN )) { // same param
491     theFxyz( 1 )      = funcValue( myValues[ SQUARE_DIST ] );
492     theDf( 1, DRV_1 ) = myValues[ DRV_1 ];
493     theDf( 1, DRV_2 ) = myValues[ DRV_2 ];
494     theDf( 1, DRV_3 ) = myValues[ DRV_3 ];
495     return true;
496   }
497 #ifdef DEBUG_PARAM_COMPUTE
498   MESSAGE ( "PARAM GUESS: " << params.X() << " "<< params.Y() << " "<< params.X() );
499   myNbIterations++; // how many times call ShellPoint()
500 #endif
501   ShellPoint( params, P );
502
503   gp_Vec dP( myPoint, P );
504   double sqDist = dP.SquareMagnitude();
505   theFxyz(1) = funcValue( sqDist );
506
507   if ( sqDist < myTolerance * myTolerance ) { // a solution found
508     myParam = params;
509     myValues[ SQUARE_DIST ] = sqDist;
510     theFxyz(1)  = theDf( 1,1 ) = theDf( 1,2 ) = theDf( 1,3 ) = 0;
511     return true;
512   }
513
514   if ( sqDist < myValues[ SQUARE_DIST ] ) // a better guess
515   {
516     // 3 partial derivatives
517     gp_Vec drv[ 3 ]; // where we move with a small step in each direction
518     for ( int iP = 1; iP <= 3; iP++ ) {
519       if ( iP == myFaceIndex ) {
520         drv[ iP - 1 ] = gp_Vec(0,0,0);
521         continue;
522       }
523       gp_XYZ Pi;
524       bool onEdge = ( theXYZ( iP ) + 0.001 > 1. );
525       if ( onEdge )
526         params.SetCoord( iP, theXYZ( iP ) - 0.001 );
527       else
528         params.SetCoord( iP, theXYZ( iP ) + 0.001 );
529       ShellPoint( params, Pi );
530       params.SetCoord( iP, theXYZ( iP ) ); // restore params
531       gp_Vec dPi ( P, Pi );
532       if ( onEdge ) dPi *= -1.;
533       double mag = dPi.Magnitude();
534       if ( mag > DBL_MIN )
535         dPi /= mag;
536       drv[ iP - 1 ] = dPi;
537     }
538     for ( int iP = 0; iP < 3; iP++ ) {
539 #if 1
540       theDf( 1, iP + 1 ) = dP * drv[iP];
541 #else
542       // Distance from P to plane passing through myPoint and defined
543       // by the 2 other derivative directions:
544       // like IntAna_IntConicQuad::Perform (const gp_Lin& L, const gp_Pln& P)
545       // where L is (P -> myPoint), P is defined by the 2 other derivative direction
546       int iPrev = ( iP ? iP - 1 : 2 );
547       int iNext = ( iP == 2 ? 0 : iP + 1 );
548       gp_Vec plnNorm = drv[ iPrev ].Crossed( drv [ iNext ] );
549       double Direc = plnNorm * drv[ iP ];
550       if ( Abs(Direc) <= DBL_MIN )
551         theDf( 1, iP + 1 ) = dP * drv[ iP ];
552       else {
553         double Dis = plnNorm * P - plnNorm * myPoint;
554         theDf( 1, iP + 1 ) = Dis/Direc;
555       }
556 #endif
557     }
558 #ifdef DEBUG_PARAM_COMPUTE
559     MESSAGE ( "F = " << theFxyz(1) << " DRV: " << theDf(1,1) << " " << theDf(1,2) << " " << theDf(1,3) );
560     myNbIterations +=3; // how many times call ShellPoint()
561 #endif
562
563     // store better values
564     myParam              = params;
565     myValues[SQUARE_DIST]= sqDist;
566     myValues[DRV_1]      = theDf(1,DRV_1);
567     myValues[DRV_2]      = theDf(1,DRV_2);
568     myValues[DRV_3]      = theDf(1,DRV_3);
569   }
570
571   return true;
572 }
573
574 //============================================================================
575 //function : computeParameters
576 //purpose  : compute point parameters in the block using math_FunctionSetRoot
577 //============================================================================
578
579 bool SMESH_Block::computeParameters(const gp_Pnt& thePoint,
580                                     gp_XYZ&       theParams,
581                                     const gp_XYZ& theParamsHint)
582 {
583   myPoint = thePoint.XYZ();
584
585   myParam.SetCoord( -1,-1,-1 );
586   myValues[ SQUARE_DIST ] = 1e100;
587
588   math_Vector low  ( 1, 3, 0.0 );
589   math_Vector up   ( 1, 3, 1.0 );
590   math_Vector tol  ( 1, 3, 1e-4 );
591   math_Vector start( 1, 3, 0.0 );
592   start( 1 ) = theParamsHint.X();
593   start( 2 ) = theParamsHint.Y();
594   start( 3 ) = theParamsHint.Z();
595
596   math_FunctionSetRoot paramSearch( *this, tol );
597
598   mySquareFunc = 0; // large approaching steps
599   //if ( hasHint ) mySquareFunc = 1; // small approaching steps
600
601   double loopTol = 10 * myTolerance;
602   int nbLoops = 0;
603   while ( distance() > loopTol && nbLoops <= 3 )
604   {
605     paramSearch.Perform ( *static_cast<math_FunctionSetWithDerivatives*>(this),
606                           start, low, up );
607     start( 1 ) = myParam.X();
608     start( 2 ) = myParam.Y();
609     start( 3 ) = myParam.Z();
610     mySquareFunc = !mySquareFunc;
611     nbLoops++;
612   }
613 #ifdef DEBUG_PARAM_COMPUTE
614   mySumDist += distance();
615   MESSAGE ( " ------ SOLUTION: ( "<< myParam.X() <<" "<< myParam.Y() <<" "<< myParam.Z() <<" )"<<endl
616          << " ------ DIST : " << distance() << "\t Tol=" << myTolerance << "\t Nb LOOPS=" << nbLoops << endl
617          << " ------ NB IT: " << myNbIterations << ",  SUM DIST: " << mySumDist );
618 #endif
619
620   theParams = myParam;
621
622   if ( myFaceIndex > 0 )
623     theParams.SetCoord( myFaceIndex, myFaceParam );
624
625   return true;
626 }
627
628 //=======================================================================
629 //function : ComputeParameters
630 //purpose  : compute point parameters in the block
631 //=======================================================================
632
633 bool SMESH_Block::ComputeParameters(const gp_Pnt& thePoint,
634                                     gp_XYZ&       theParams,
635                                     const int     theShapeID,
636                                     const gp_XYZ& theParamsHint)
637 {
638   if ( VertexParameters( theShapeID, theParams ))
639     return true;
640
641   if ( IsEdgeID( theShapeID )) {
642     TEdge& e = myEdge[ theShapeID - ID_FirstE ];
643     Adaptor3d_Curve* curve = e.GetCurve();
644     Extrema_ExtPC anExtPC( thePoint, *curve, curve->FirstParameter(), curve->LastParameter() );
645     int i, nb = anExtPC.IsDone() ? anExtPC.NbExt() : 0;
646     for ( i = 1; i <= nb; i++ ) {
647       if ( anExtPC.IsMin( i ))
648         return EdgeParameters( theShapeID, anExtPC.Point( i ).Parameter(), theParams );
649     }
650     return false;
651   }
652
653   const bool isOnFace = IsFaceID( theShapeID );
654   double * coef = GetShapeCoef( theShapeID );
655
656   // Find the first guess paremeters
657
658   gp_XYZ start(0, 0, 0);
659
660   bool hasHint = ( 0 <= theParamsHint.X() && theParamsHint.X() <= 1 &&
661                    0 <= theParamsHint.Y() && theParamsHint.Y() <= 1 &&
662                    0 <= theParamsHint.Y() && theParamsHint.Y() <= 1 );
663   if ( !hasHint && !myGridComputed )
664   {
665     // define the first guess by thePoint projection on lines
666     // connecting vertices
667     bool needGrid = false;
668     gp_XYZ par000( 0, 0, 0 ), par111( 1, 1, 1 );
669     double zero = DBL_MIN * DBL_MIN;
670     for ( int iEdge = 0, iParam = 1; iParam <= 3 && !needGrid; iParam++ )
671     {
672       if ( isOnFace && coef[ iParam - 1 ] != 0 ) {
673         iEdge += 4;
674         continue;
675       }
676       double sumParam = 0;
677       for ( int iE = 0; iE < 4; iE++, iEdge++ ) { // loop on 4 parallel edges
678         gp_Pnt p0 = myEdge[ iEdge ].Point( par000 );
679         gp_Pnt p1 = myEdge[ iEdge ].Point( par111 );
680         gp_Vec v01( p0, p1 ), v0P( p0, thePoint );
681         double len2 = v01.SquareMagnitude();
682         double par = 0;
683         if ( len2 > zero ) {
684           par = v0P.Dot( v01 ) / len2;
685           if ( par < 0 || par > 1 ) { // projection falls out of line ends => needGrid
686             needGrid = true;
687             break;
688           }
689         }
690         sumParam += par;
691       }
692       start.SetCoord( iParam, sumParam / 4.);
693     }
694     if ( needGrid ) {
695       // compute nodes of 3 x 3 x 3 grid
696       int iNode = 0;
697       Bnd_Box box;
698       for ( double x = 0.25; x < 0.9; x += 0.25 )
699         for ( double y = 0.25; y < 0.9; y += 0.25 )
700           for ( double z = 0.25; z < 0.9; z += 0.25 ) {
701             TxyzPair & prmPtn = my3x3x3GridNodes[ iNode++ ];
702             prmPtn.first.SetCoord( x, y, z );
703             ShellPoint( prmPtn.first, prmPtn.second );
704             box.Add( gp_Pnt( prmPtn.second ));
705           }
706       myGridComputed = true;
707       myTolerance = sqrt( box.SquareExtent() ) * 1e-5;
708     }
709   }
710
711   if ( hasHint )
712   {
713     start = theParamsHint;
714   }
715   else if ( myGridComputed )
716   {
717     double minDist = DBL_MAX;
718     gp_XYZ* bestParam = 0;
719     for ( int iNode = 0; iNode < 27; iNode++ ) {
720       TxyzPair & prmPtn = my3x3x3GridNodes[ iNode ];
721       double dist = ( thePoint.XYZ() - prmPtn.second ).SquareModulus();
722       if ( dist < minDist ) {
723         minDist = dist;
724         bestParam = & prmPtn.first;
725       }
726     }
727     start = *bestParam;
728   }
729
730   myFaceIndex = -1;
731   myFaceParam = 0.;
732   if ( isOnFace ) {
733     // put a point on the face
734     for ( int iCoord = 0; iCoord < 3; iCoord++ )
735       if ( coef[ iCoord ] ) {
736         myFaceIndex = iCoord + 1;
737         myFaceParam = ( coef[ iCoord ] < 0.5 ) ? 0.0 : 1.0;
738         start.SetCoord( myFaceIndex, myFaceParam );
739       }
740   }
741
742 #ifdef DEBUG_PARAM_COMPUTE
743   MESSAGE ( " #### POINT " <<thePoint.X()<<" "<<thePoint.Y()<<" "<<thePoint.Z()<<" ####" );
744 #endif
745
746   if ( myTolerance < 0 ) myTolerance = 1e-6;
747
748   const double parDelta = 1e-4;
749   const double sqTolerance = myTolerance * myTolerance;
750
751   gp_XYZ solution = start, params = start;
752   double sqDistance = 1e100; 
753   int nbLoops = 0, nbGetWorst = 0;
754
755   while ( nbLoops <= 100 )
756   {
757     gp_XYZ P, Pi;
758     ShellPoint( params, P );
759
760     gp_Vec dP( thePoint, P );
761     double sqDist = dP.SquareMagnitude();
762
763     if ( sqDist > sqDistance ) { // solution get worse
764       if ( ++nbGetWorst > 2 )
765         return computeParameters( thePoint, theParams, solution );
766     }
767 #ifdef DEBUG_PARAM_COMPUTE
768     MESSAGE ( "PARAMS: ( " << params.X() <<" "<< params.Y() <<" "<< params.Z() <<" )" );
769     MESSAGE ( "DIST: " << sqrt( sqDist ) );
770 #endif
771
772     if ( sqDist < sqDistance ) { // get better
773       sqDistance = sqDist;
774       solution   = params;
775       nbGetWorst = 0;
776       if ( sqDistance < sqTolerance ) // a solution found
777         break;
778     }
779
780         // look for a next better solution
781     for ( int iP = 1; iP <= 3; iP++ ) {
782       if ( iP == myFaceIndex )
783         continue;
784       // see where we move with a small (=parDelta) step in this direction
785       gp_XYZ nearParams = params;
786       bool onEdge = ( params.Coord( iP ) + parDelta > 1. );
787       if ( onEdge )
788         nearParams.SetCoord( iP, params.Coord( iP ) - parDelta );
789       else
790         nearParams.SetCoord( iP, params.Coord( iP ) + parDelta );
791       ShellPoint( nearParams, Pi );
792       gp_Vec dPi ( P, Pi );
793       if ( onEdge ) dPi *= -1.;
794       // modify a parameter
795       double mag = dPi.Magnitude();
796       if ( mag < DBL_MIN )
797         continue;
798       gp_Vec dir = dPi / mag; // dir we move modifying the parameter
799       double dist = dir * dP; // where we should get to
800       double dPar = dist / mag * parDelta; // predict parameter change
801       double curPar = params.Coord( iP );
802       double par = curPar - dPar; // new parameter value
803       while ( par > 1 || par < 0 ) {
804         dPar /= 2.;
805         par = curPar - dPar;
806       }
807       params.SetCoord( iP, par );
808     }
809
810     nbLoops++;
811   }
812 #ifdef DEBUG_PARAM_COMPUTE
813   myNbIterations += nbLoops*4; // how many times ShellPoint called
814   mySumDist += sqrt( sqDistance );
815   MESSAGE ( " ------ SOLUTION: ( "<<solution.X()<<" "<<solution.Y()<<" "<<solution.Z()<<" )"<< std::endl
816          << " ------ DIST : " << sqrt( sqDistance ) << "\t Tol=" << myTolerance << "\t Nb LOOPS=" << nbLoops << std::endl
817          << " ------ NB IT: " << myNbIterations << ",  SUM DIST: " << mySumDist );
818 #endif
819
820   theParams = solution;
821
822   if ( myFaceIndex > 0 )
823     theParams.SetCoord( myFaceIndex, myFaceParam );
824
825   return true;
826 }
827
828 //=======================================================================
829 //function : VertexParameters
830 //purpose  : return parameters of a vertex given by TShapeID
831 //=======================================================================
832
833 bool SMESH_Block::VertexParameters(const int theVertexID, gp_XYZ& theParams)
834 {
835   switch ( theVertexID ) {
836   case ID_V000: theParams.SetCoord(0., 0., 0.); return true;
837   case ID_V100: theParams.SetCoord(1., 0., 0.); return true;
838   case ID_V110: theParams.SetCoord(1., 1., 0.); return true;
839   case ID_V010: theParams.SetCoord(0., 1., 0.); return true;
840   default:;
841   }
842   return false;
843 }
844
845 //=======================================================================
846 //function : EdgeParameters
847 //purpose  : return parameters of a point given by theU on edge
848 //=======================================================================
849
850 bool SMESH_Block::EdgeParameters(const int theEdgeID, const double theU, gp_XYZ& theParams)
851 {
852   if ( IsEdgeID( theEdgeID )) {
853     vector< int > vertexVec;
854     GetEdgeVertexIDs( theEdgeID, vertexVec );
855     VertexParameters( vertexVec[0], theParams );
856     TEdge& e = myEdge[ theEdgeID - ID_Ex00 ];
857     double param = ( theU - e.EndParam(0) ) / ( e.EndParam(1) - e.EndParam(0) );
858     theParams.SetCoord( e.CoordInd(), param );
859     return true;
860   }
861   return false;
862 }
863
864 //=======================================================================
865 //function : DumpShapeID
866 //purpose  : debug an id of a block sub-shape
867 //=======================================================================
868
869 #define CASEDUMP(id,strm) case id: strm << #id; break;
870
871 ostream& SMESH_Block::DumpShapeID (const int id, ostream& stream)
872 {
873   switch ( id ) {
874   CASEDUMP( ID_V000, stream );
875   CASEDUMP( ID_V100, stream );
876   CASEDUMP( ID_V010, stream );
877   CASEDUMP( ID_V110, stream );
878   CASEDUMP( ID_V001, stream );
879   CASEDUMP( ID_V101, stream );
880   CASEDUMP( ID_V011, stream );
881   CASEDUMP( ID_V111, stream );
882   CASEDUMP( ID_Ex00, stream );
883   CASEDUMP( ID_Ex10, stream );
884   CASEDUMP( ID_Ex01, stream );
885   CASEDUMP( ID_Ex11, stream );
886   CASEDUMP( ID_E0y0, stream );
887   CASEDUMP( ID_E1y0, stream );
888   CASEDUMP( ID_E0y1, stream );
889   CASEDUMP( ID_E1y1, stream );
890   CASEDUMP( ID_E00z, stream );
891   CASEDUMP( ID_E10z, stream );
892   CASEDUMP( ID_E01z, stream );
893   CASEDUMP( ID_E11z, stream );
894   CASEDUMP( ID_Fxy0, stream );
895   CASEDUMP( ID_Fxy1, stream );
896   CASEDUMP( ID_Fx0z, stream );
897   CASEDUMP( ID_Fx1z, stream );
898   CASEDUMP( ID_F0yz, stream );
899   CASEDUMP( ID_F1yz, stream );
900   CASEDUMP( ID_Shell, stream );
901   default: stream << "ID_INVALID";
902   }
903   return stream;
904 }
905
906 //=======================================================================
907 //function : GetShapeIDByParams
908 //purpose  : define an id of the block sub-shape by normlized point coord
909 //=======================================================================
910
911 int SMESH_Block::GetShapeIDByParams ( const gp_XYZ& theCoord )
912 {
913   //   id ( 0 - 26 ) computation:
914
915   //   vertex     ( 0 - 7 )  : id = 1*x + 2*y + 4*z
916
917   //   edge || X  ( 8 - 11 ) : id = 8   + 1*y + 2*z
918   //   edge || Y  ( 12 - 15 ): id = 1*x + 12  + 2*z
919   //   edge || Z  ( 16 - 19 ): id = 1*x + 2*y + 16 
920
921   //   face || XY ( 20 - 21 ): id = 8   + 12  + 1*z - 0
922   //   face || XZ ( 22 - 23 ): id = 8   + 1*y + 16  - 2
923   //   face || YZ ( 24 - 25 ): id = 1*x + 12  + 16  - 4
924
925   static int iAddBnd[]    = { 1, 2, 4 };
926   static int iAddNotBnd[] = { 8, 12, 16 };
927   static int iFaceSubst[] = { 0, 2, 4 };
928
929   int id = 0;
930   int iOnBoundary = 0;
931   for ( int iCoord = 0; iCoord < 3; iCoord++ )
932   {
933     double val = theCoord.Coord( iCoord + 1 );
934     if ( val == 0.0 )
935       iOnBoundary++;
936     else if ( val == 1.0 )
937       id += iAddBnd[ iOnBoundary++ ];
938     else
939       id += iAddNotBnd[ iCoord ];
940   }
941   if ( iOnBoundary == 1 ) // face
942     id -= iFaceSubst[ (id - 20) / 4 ];
943   else if ( iOnBoundary == 0 ) // shell
944     id = 26;
945
946   if ( id > 26 || id < 0 ) {
947     MESSAGE( "GetShapeIDByParams() = " << id
948             <<" "<< theCoord.X() <<" "<< theCoord.Y() <<" "<< theCoord.Z() );
949   }
950
951   return id + 1; // shape ids start at 1
952 }
953
954 //=======================================================================
955 //function : GetOrderedEdges
956 //purpose  : return nb wires and a list of oredered edges
957 //=======================================================================
958
959 int SMESH_Block::GetOrderedEdges (const TopoDS_Face&   theFace,
960                                   TopoDS_Vertex        theFirstVertex,
961                                   list< TopoDS_Edge >& theEdges,
962                                   list< int >  &       theNbVertexInWires)
963 {
964   // put wires in a list, so that an outer wire comes first
965   list<TopoDS_Wire> aWireList;
966   TopoDS_Wire anOuterWire = BRepTools::OuterWire( theFace );
967   aWireList.push_back( anOuterWire );
968   for ( TopoDS_Iterator wIt (theFace); wIt.More(); wIt.Next() )
969     if ( !anOuterWire.IsSame( wIt.Value() ))
970       aWireList.push_back( TopoDS::Wire( wIt.Value() ));
971
972   // loop on edges of wires
973   theNbVertexInWires.clear();
974   list<TopoDS_Wire>::iterator wlIt = aWireList.begin();
975   for ( ; wlIt != aWireList.end(); wlIt++ )
976   {
977     int iE;
978     BRepTools_WireExplorer wExp( *wlIt, theFace );
979     for ( iE = 0; wExp.More(); wExp.Next(), iE++ )
980     {
981       TopoDS_Edge edge = wExp.Current();
982       edge = TopoDS::Edge( edge.Oriented( wExp.Orientation() ));
983       theEdges.push_back( edge );
984     }
985     theNbVertexInWires.push_back( iE );
986     iE = 0;
987     if ( wlIt == aWireList.begin() && theEdges.size() > 1 ) { // the outer wire
988       // orient closed edges
989       list< TopoDS_Edge >::iterator eIt, eIt2;
990       for ( eIt = theEdges.begin(); eIt != theEdges.end(); eIt++ )
991       {
992         TopoDS_Edge& edge = *eIt;
993         if ( TopExp::FirstVertex( edge ).IsSame( TopExp::LastVertex( edge ) ))
994         {
995           eIt2 = eIt;
996           bool isNext = ( eIt2 == theEdges.begin() );
997           TopoDS_Edge edge2 = isNext ? *(++eIt2) : *(--eIt2);
998           double f1,l1,f2,l2;
999           Handle(Geom2d_Curve) c1 = BRep_Tool::CurveOnSurface( edge, theFace, f1,l1 );
1000           Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( edge2, theFace, f2,l2 );
1001           gp_Pnt2d pf = c1->Value( edge.Orientation() == TopAbs_FORWARD ? f1 : l1 );
1002           gp_Pnt2d pl = c1->Value( edge.Orientation() == TopAbs_FORWARD ? l1 : f1 );
1003           bool isFirst = ( edge2.Orientation() == TopAbs_FORWARD ? isNext : !isNext );
1004           gp_Pnt2d p2 = c2->Value( isFirst ? f2 : l2 );
1005           isFirst = ( p2.SquareDistance( pf ) < p2.SquareDistance( pl ));
1006           if ( isNext ? isFirst : !isFirst )
1007             edge.Reverse();
1008           // to make a seam go first
1009           if ( theFirstVertex.IsNull() )
1010             theFirstVertex = TopExp::FirstVertex( edge, true );
1011         }
1012       }
1013       // rotate theEdges until it begins from theFirstVertex
1014       if ( ! theFirstVertex.IsNull() ) {
1015         TopoDS_Vertex vv[2];
1016         TopExp::Vertices( theEdges.front(), vv[0], vv[1], true );
1017         // on closed face, make seam edge the first in the list
1018         while ( !vv[0].IsSame( theFirstVertex ) || vv[0].IsSame( vv[1] ))
1019         {
1020           theEdges.splice(theEdges.end(), theEdges,
1021                           theEdges.begin(), ++theEdges.begin());
1022           TopExp::Vertices( theEdges.front(), vv[0], vv[1], true );
1023           if ( iE++ > theNbVertexInWires.back() ) {
1024 #ifdef _DEBUG_
1025             gp_Pnt p = BRep_Tool::Pnt( theFirstVertex );
1026             MESSAGE ( " : Warning : vertex "<< theFirstVertex.TShape().operator->()
1027                    << " ( " << p.X() << " " << p.Y() << " " << p.Z() << " )" 
1028                    << " not found in outer wire of face "<< theFace.TShape().operator->()
1029                    << " with vertices: " );
1030             wExp.Init( *wlIt, theFace );
1031             for ( int i = 0; wExp.More(); wExp.Next(), i++ )
1032             {
1033               TopoDS_Edge edge = wExp.Current();
1034               edge = TopoDS::Edge( edge.Oriented( wExp.Orientation() ));
1035               TopoDS_Vertex v = TopExp::FirstVertex( edge, true );
1036               gp_Pnt p = BRep_Tool::Pnt( v );
1037               MESSAGE_ADD ( i << " " << v.TShape().operator->() << " "
1038                             << p.X() << " " << p.Y() << " " << p.Z() << " " << std::endl );
1039             }
1040 #endif
1041             break; // break infinite loop
1042           }
1043         }
1044       }
1045     } // end outer wire
1046   }
1047
1048   return aWireList.size();
1049 }
1050 //================================================================================
1051 /*!
1052  * \brief Call it after geometry initialisation
1053  */
1054 //================================================================================
1055
1056 void SMESH_Block::init()
1057 {
1058   myNbIterations = 0;
1059   mySumDist = 0;
1060   myGridComputed = false;
1061 }
1062
1063 //=======================================================================
1064 //function : LoadMeshBlock
1065 //purpose  : prepare to work with theVolume
1066 //=======================================================================
1067
1068 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
1069
1070 bool SMESH_Block::LoadMeshBlock(const SMDS_MeshVolume*        theVolume,
1071                                 const int                     theNode000Index,
1072                                 const int                     theNode001Index,
1073                                 vector<const SMDS_MeshNode*>& theOrderedNodes)
1074 {
1075   MESSAGE(" ::LoadMeshBlock()");
1076   init();
1077
1078   SMDS_VolumeTool vTool;
1079   if (!vTool.Set( theVolume ) || vTool.NbNodes() != 8 ||
1080       !vTool.IsLinked( theNode000Index, theNode001Index )) {
1081     MESSAGE(" Bad arguments ");
1082     return false;
1083   }
1084   vTool.SetExternalNormal();
1085   // In terms of indices used for access to nodes and faces in SMDS_VolumeTool:
1086   int V000, V100, V010, V110, V001, V101, V011, V111; // 8 vertices
1087   int Fxy0, Fxy1; // bottom and top faces
1088   // vertices of faces
1089   vector<int> vFxy0, vFxy1;
1090
1091   V000 = theNode000Index;
1092   V001 = theNode001Index;
1093
1094   // get faces sharing V000 and V001
1095   list<int> fV000, fV001;
1096   int i, iF, iE, iN;
1097   for ( iF = 0; iF < vTool.NbFaces(); ++iF ) {
1098     const int* nid = vTool.GetFaceNodesIndices( iF );
1099     for ( iN = 0; iN < 4; ++iN )
1100       if ( nid[ iN ] == V000 ) {
1101         fV000.push_back( iF );
1102       } else if ( nid[ iN ] == V001 ) {
1103         fV001.push_back( iF );
1104       }
1105   }
1106
1107   // find the bottom (Fxy0), the top (Fxy1) faces
1108   list<int>::iterator fIt1, fIt2, Fxy0Pos;
1109   for ( fIt1 = fV000.begin(); fIt1 != fV000.end(); fIt1++) {
1110     fIt2 = std::find( fV001.begin(), fV001.end(), *fIt1 );
1111     if ( fIt2 != fV001.end() ) { // *fIt1 is in the both lists
1112       fV001.erase( fIt2 ); // erase Fx0z or F0yz from fV001
1113     } else { // *fIt1 is in fV000 only
1114       Fxy0Pos = fIt1; // points to Fxy0
1115     }
1116   }
1117   Fxy0 = *Fxy0Pos;
1118   Fxy1 = fV001.front();
1119   const SMDS_MeshNode** nn = vTool.GetNodes();
1120
1121   // find bottom veritices, their order is that a face normal is external
1122   vFxy0.resize(4);
1123   const int* nid = vTool.GetFaceNodesIndices( Fxy0 );
1124   for ( i = 0; i < 4; ++i )
1125     if ( nid[ i ] == V000 )
1126       break;
1127   for ( iN = 0; iN < 4; ++iN, ++i ) {
1128     if ( i == 4 ) i = 0;
1129     vFxy0[ iN ] = nid[ i ];
1130   }
1131   // find top veritices, their order is that a face normal is external
1132   vFxy1.resize(4);
1133   nid = vTool.GetFaceNodesIndices( Fxy1 );
1134   for ( i = 0; i < 4; ++i )
1135     if ( nid[ i ] == V001 )
1136       break;
1137   for ( iN = 0; iN < 4; ++iN, ++i ) {
1138     if ( i == 4 ) i = 0;
1139     vFxy1[ iN ] = nid[ i ];
1140   }
1141   // find indices of the rest veritices 
1142   V100 = vFxy0[3];
1143   V010 = vFxy0[1];
1144   V110 = vFxy0[2];
1145   V101 = vFxy1[1];
1146   V011 = vFxy1[3];
1147   V111 = vFxy1[2];
1148
1149   // set points coordinates
1150   myPnt[ ID_V000 - 1 ] = gpXYZ( nn[ V000 ] );
1151   myPnt[ ID_V100 - 1 ] = gpXYZ( nn[ V100 ] );
1152   myPnt[ ID_V010 - 1 ] = gpXYZ( nn[ V010 ] );
1153   myPnt[ ID_V110 - 1 ] = gpXYZ( nn[ V110 ] );
1154   myPnt[ ID_V001 - 1 ] = gpXYZ( nn[ V001 ] );
1155   myPnt[ ID_V101 - 1 ] = gpXYZ( nn[ V101 ] );
1156   myPnt[ ID_V011 - 1 ] = gpXYZ( nn[ V011 ] );
1157   myPnt[ ID_V111 - 1 ] = gpXYZ( nn[ V111 ] );
1158
1159   // fill theOrderedNodes
1160   theOrderedNodes.resize( 8 );
1161   theOrderedNodes[ 0 ] = nn[ V000 ];
1162   theOrderedNodes[ 1 ] = nn[ V100 ];
1163   theOrderedNodes[ 2 ] = nn[ V010 ];
1164   theOrderedNodes[ 3 ] = nn[ V110 ];
1165   theOrderedNodes[ 4 ] = nn[ V001 ];
1166   theOrderedNodes[ 5 ] = nn[ V101 ];
1167   theOrderedNodes[ 6 ] = nn[ V011 ];
1168   theOrderedNodes[ 7 ] = nn[ V111 ];
1169   
1170   // fill edges
1171   vector< int > vertexVec;
1172   for ( iE = 0; iE < NbEdges(); ++iE ) {
1173     GetEdgeVertexIDs(( iE + ID_FirstE ), vertexVec );
1174     myEdge[ iE ].Set(( iE + ID_FirstE ),
1175                      myPnt[ vertexVec[0] - 1 ],
1176                      myPnt[ vertexVec[1] - 1 ]);
1177   }
1178
1179   // fill faces' corners
1180   for ( iF = ID_Fxy0; iF < ID_Shell; ++iF )
1181   {
1182     TFace& tFace = myFace[ iF - ID_FirstF ];
1183     vector< int > edgeIdVec(4, -1);
1184     GetFaceEdgesIDs( iF, edgeIdVec );
1185     tFace.Set( iF, myEdge[ edgeIdVec [ 0 ] - ID_Ex00], myEdge[ edgeIdVec [ 1 ] - ID_Ex00]);
1186   }
1187
1188   return true;
1189 }
1190
1191 //=======================================================================
1192 //function : LoadBlockShapes
1193 //purpose  : Initialize block geometry with theShell,
1194 //           add sub-shapes of theBlock to theShapeIDMap so that they get
1195 //           IDs acoording to enum TShapeID
1196 //=======================================================================
1197
1198 bool SMESH_Block::LoadBlockShapes(const TopoDS_Shell&         theShell,
1199                                   const TopoDS_Vertex&        theVertex000,
1200                                   const TopoDS_Vertex&        theVertex001,
1201                                   TopTools_IndexedMapOfOrientedShape& theShapeIDMap )
1202 {
1203   MESSAGE(" ::LoadBlockShapes()");
1204   return ( FindBlockShapes( theShell, theVertex000, theVertex001, theShapeIDMap ) &&
1205            LoadBlockShapes( theShapeIDMap ));
1206 }
1207
1208 //=======================================================================
1209 //function : LoadBlockShapes
1210 //purpose  : add sub-shapes of theBlock to theShapeIDMap so that they get
1211 //           IDs acoording to enum TShapeID
1212 //=======================================================================
1213
1214 bool SMESH_Block::FindBlockShapes(const TopoDS_Shell&         theShell,
1215                                   const TopoDS_Vertex&        theVertex000,
1216                                   const TopoDS_Vertex&        theVertex001,
1217                                   TopTools_IndexedMapOfOrientedShape& theShapeIDMap )
1218 {
1219   MESSAGE(" ::FindBlockShapes()");
1220
1221   // 8 vertices
1222   TopoDS_Shape V000, V100, V010, V110, V001, V101, V011, V111;
1223   // 12 edges
1224   TopoDS_Shape Ex00, Ex10, Ex01, Ex11;
1225   TopoDS_Shape E0y0, E1y0, E0y1, E1y1;
1226   TopoDS_Shape E00z, E10z, E01z, E11z;
1227   // 6 faces
1228   TopoDS_Shape Fxy0, Fx0z, F0yz, Fxy1, Fx1z, F1yz;
1229
1230   // nb of faces bound to a vertex in TopTools_IndexedDataMapOfShapeListOfShape
1231   // filled by TopExp::MapShapesAndAncestors()
1232   const int NB_FACES_BY_VERTEX = 6;
1233
1234   TopTools_IndexedDataMapOfShapeListOfShape vfMap;
1235   TopExp::MapShapesAndAncestors( theShell, TopAbs_VERTEX, TopAbs_FACE, vfMap );
1236   if ( vfMap.Extent() != 8 ) {
1237     MESSAGE(" Wrong nb of vertices in the block: " << vfMap.Extent() );
1238     return false;
1239   }
1240
1241   V000 = theVertex000;
1242   V001 = theVertex001;
1243
1244   if ( V000.IsNull() ) {
1245     // find vertex 000 - the one with smallest coordinates
1246     double minVal = DBL_MAX, minX, val;
1247     for ( int i = 1; i <= 8; i++ ) {
1248       const TopoDS_Vertex& v = TopoDS::Vertex( vfMap.FindKey( i ));
1249       gp_Pnt P = BRep_Tool::Pnt( v );
1250       val = P.X() + P.Y() + P.Z();
1251       if ( val < minVal || ( val == minVal && P.X() < minX )) {
1252         V000 = v;
1253         minVal = val;
1254         minX = P.X();
1255       }
1256     }
1257     // find vertex 001 - the one on the most vertical edge passing through V000
1258     TopTools_IndexedDataMapOfShapeListOfShape veMap;
1259     TopExp::MapShapesAndAncestors( theShell, TopAbs_VERTEX, TopAbs_EDGE, veMap );
1260     gp_Vec dir001 = gp::DZ();
1261     gp_Pnt p000 = BRep_Tool::Pnt( TopoDS::Vertex( V000 ));
1262     double maxVal = -DBL_MAX;
1263     TopTools_ListIteratorOfListOfShape eIt ( veMap.FindFromKey( V000 ));
1264     for (  ; eIt.More(); eIt.Next() ) {
1265       const TopoDS_Edge& e = TopoDS::Edge( eIt.Value() );
1266       TopoDS_Vertex v = TopExp::FirstVertex( e );
1267       if ( v.IsSame( V000 ))
1268         v = TopExp::LastVertex( e );
1269       val = dir001 * gp_Vec( p000, BRep_Tool::Pnt( v )).Normalized();
1270       if ( val > maxVal ) {
1271         V001 = v;
1272         maxVal = val;
1273       }
1274     }
1275   }
1276
1277   // find the bottom (Fxy0), Fx0z and F0yz faces
1278
1279   const TopTools_ListOfShape& f000List = vfMap.FindFromKey( V000 );
1280   const TopTools_ListOfShape& f001List = vfMap.FindFromKey( V001 );
1281   if (f000List.Extent() != NB_FACES_BY_VERTEX ||
1282       f001List.Extent() != NB_FACES_BY_VERTEX ) {
1283     MESSAGE(" LoadBlockShapes() " << f000List.Extent() << " " << f001List.Extent());
1284     return false;
1285   }
1286   TopTools_ListIteratorOfListOfShape f001It, f000It ( f000List );
1287   int i, j, iFound1, iFound2;
1288   for ( j = 0; f000It.More(); f000It.Next(), j++ )
1289   {
1290     if ( NB_FACES_BY_VERTEX == 6 && j % 2 ) continue; // each face encounters twice
1291     const TopoDS_Shape& F = f000It.Value();
1292     for ( i = 0, f001It.Initialize( f001List ); f001It.More(); f001It.Next(), i++ ) {
1293       if ( NB_FACES_BY_VERTEX == 6 && i % 2 ) continue; // each face encounters twice
1294       if ( F.IsSame( f001It.Value() ))
1295         break;
1296     }
1297     if ( f001It.More() ) // Fx0z or F0yz found
1298       if ( Fx0z.IsNull() ) {
1299         Fx0z = F;
1300         iFound1 = i;
1301       } else {
1302         F0yz = F;
1303         iFound2 = i;
1304       }
1305     else // F is the bottom face
1306       Fxy0 = F;
1307   }
1308   if ( Fxy0.IsNull() || Fx0z.IsNull() || F0yz.IsNull() ) {
1309     MESSAGE( Fxy0.IsNull() <<" "<< Fx0z.IsNull() <<" "<< F0yz.IsNull() );
1310     return false;
1311   }
1312
1313   // choose the top face (Fxy1)
1314   for ( i = 0, f001It.Initialize( f001List ); f001It.More(); f001It.Next(), i++ ) {
1315     if ( NB_FACES_BY_VERTEX == 6 && i % 2 ) continue; // each face encounters twice
1316     if ( i != iFound1 && i != iFound2 )
1317       break;
1318   }
1319   Fxy1 = f001It.Value();
1320   if ( Fxy1.IsNull() ) {
1321     MESSAGE(" LoadBlockShapes() error ");
1322     return false;
1323   }
1324
1325   // find bottom edges and veritices
1326   list< TopoDS_Edge > eList;
1327   list< int >         nbVertexInWires;
1328   GetOrderedEdges( TopoDS::Face( Fxy0 ), TopoDS::Vertex( V000 ), eList, nbVertexInWires );
1329   if ( nbVertexInWires.size() != 1 || nbVertexInWires.front() != 4 ) {
1330     MESSAGE(" LoadBlockShapes() error ");
1331     return false;
1332   }
1333   list< TopoDS_Edge >::iterator elIt = eList.begin();
1334   for ( i = 0; elIt != eList.end(); elIt++, i++ )
1335     switch ( i ) {
1336     case 0: E0y0 = *elIt; V010 = TopExp::LastVertex( *elIt, true ); break;
1337     case 1: Ex10 = *elIt; V110 = TopExp::LastVertex( *elIt, true ); break;
1338     case 2: E1y0 = *elIt; V100 = TopExp::LastVertex( *elIt, true ); break;
1339     case 3: Ex00 = *elIt; break;
1340     default:;
1341     }
1342   if ( i != 4 || E0y0.IsNull() || Ex10.IsNull() || E1y0.IsNull() || Ex00.IsNull() ) {
1343     MESSAGE(" LoadBlockShapes() error, eList.size()=" << eList.size());
1344     return false;
1345   }
1346
1347
1348   // find top edges and veritices
1349   eList.clear();
1350   GetOrderedEdges( TopoDS::Face( Fxy1 ), TopoDS::Vertex( V001 ), eList, nbVertexInWires );
1351   if ( nbVertexInWires.size() != 1 || nbVertexInWires.front() != 4 ) {
1352     MESSAGE(" LoadBlockShapes() error ");
1353     return false;
1354   }
1355   for ( i = 0, elIt = eList.begin(); elIt != eList.end(); elIt++, i++ )
1356     switch ( i ) {
1357     case 0: Ex01 = *elIt; V101 = TopExp::LastVertex( *elIt, true ); break;
1358     case 1: E1y1 = *elIt; V111 = TopExp::LastVertex( *elIt, true ); break;
1359     case 2: Ex11 = *elIt; V011 = TopExp::LastVertex( *elIt, true ); break;
1360     case 3: E0y1 = *elIt; break;
1361     default:;
1362     }
1363   if ( i != 4 || Ex01.IsNull() || E1y1.IsNull() || Ex11.IsNull() || E0y1.IsNull() ) {
1364     MESSAGE(" LoadBlockShapes() error, eList.size()=" << eList.size());
1365     return false;
1366   }
1367
1368   // swap Fx0z and F0yz if necessary
1369   TopExp_Explorer exp( Fx0z, TopAbs_VERTEX );
1370   for ( ; exp.More(); exp.Next() ) // Fx0z shares V101 and V100
1371     if ( V101.IsSame( exp.Current() ) || V100.IsSame( exp.Current() ))
1372       break; // V101 or V100 found
1373   if ( !exp.More() ) { // not found
1374     std::swap( Fx0z, F0yz);
1375   }
1376
1377   // find Fx1z and F1yz faces
1378   const TopTools_ListOfShape& f111List = vfMap.FindFromKey( V111 );
1379   const TopTools_ListOfShape& f110List = vfMap.FindFromKey( V110 );
1380   if (f111List.Extent() != NB_FACES_BY_VERTEX ||
1381       f110List.Extent() != NB_FACES_BY_VERTEX ) {
1382     MESSAGE(" LoadBlockShapes() " << f111List.Extent() << " " << f110List.Extent());
1383     return false;
1384   }
1385   TopTools_ListIteratorOfListOfShape f111It, f110It ( f110List);
1386   for ( j = 0 ; f110It.More(); f110It.Next(), j++ ) {
1387     if ( NB_FACES_BY_VERTEX == 6 && j % 2 ) continue; // each face encounters twice
1388     const TopoDS_Shape& F = f110It.Value();
1389     for ( i = 0, f111It.Initialize( f111List ); f111It.More(); f111It.Next(), i++ ) {
1390       if ( NB_FACES_BY_VERTEX == 6 && i % 2 ) continue; // each face encounters twice
1391       if ( F.IsSame( f111It.Value() )) { // Fx1z or F1yz found
1392         if ( Fx1z.IsNull() )
1393           Fx1z = F;
1394         else
1395           F1yz = F;
1396       }
1397     }
1398   }
1399   if ( Fx1z.IsNull() || F1yz.IsNull() ) {
1400     MESSAGE(" LoadBlockShapes() error ");
1401     return false;
1402   }
1403
1404   // swap Fx1z and F1yz if necessary
1405   for ( exp.Init( Fx1z, TopAbs_VERTEX ); exp.More(); exp.Next() )
1406     if ( V010.IsSame( exp.Current() ) || V011.IsSame( exp.Current() ))
1407       break;
1408   if ( !exp.More() ) {
1409     std::swap( Fx1z, F1yz);
1410   }
1411
1412   // find vertical edges
1413   for ( exp.Init( Fx0z, TopAbs_EDGE ); exp.More(); exp.Next() ) {
1414     const TopoDS_Edge& edge = TopoDS::Edge( exp.Current() );
1415     const TopoDS_Shape& vFirst = TopExp::FirstVertex( edge, true );
1416     if ( vFirst.IsSame( V001 ))
1417       E00z = edge;
1418     else if ( vFirst.IsSame( V100 ))
1419       E10z = edge;
1420   }
1421   if ( E00z.IsNull() || E10z.IsNull() ) {
1422     MESSAGE(" LoadBlockShapes() error ");
1423     return false;
1424   }
1425   for ( exp.Init( Fx1z, TopAbs_EDGE ); exp.More(); exp.Next() ) {
1426     const TopoDS_Edge& edge = TopoDS::Edge( exp.Current() );
1427     const TopoDS_Shape& vFirst = TopExp::FirstVertex( edge, true );
1428     if ( vFirst.IsSame( V111 ))
1429       E11z = edge;
1430     else if ( vFirst.IsSame( V010 ))
1431       E01z = edge;
1432   }
1433   if ( E01z.IsNull() || E11z.IsNull() ) {
1434     MESSAGE(" LoadBlockShapes() error ");
1435     return false;
1436   }
1437
1438   // load shapes in theShapeIDMap
1439
1440   theShapeIDMap.Clear();
1441   
1442   theShapeIDMap.Add(V000.Oriented( TopAbs_FORWARD ));
1443   theShapeIDMap.Add(V100.Oriented( TopAbs_FORWARD ));
1444   theShapeIDMap.Add(V010.Oriented( TopAbs_FORWARD ));
1445   theShapeIDMap.Add(V110.Oriented( TopAbs_FORWARD ));
1446   theShapeIDMap.Add(V001.Oriented( TopAbs_FORWARD ));
1447   theShapeIDMap.Add(V101.Oriented( TopAbs_FORWARD ));
1448   theShapeIDMap.Add(V011.Oriented( TopAbs_FORWARD ));
1449   theShapeIDMap.Add(V111.Oriented( TopAbs_FORWARD ));
1450
1451   theShapeIDMap.Add(Ex00);
1452   theShapeIDMap.Add(Ex10);
1453   theShapeIDMap.Add(Ex01);
1454   theShapeIDMap.Add(Ex11);
1455
1456   theShapeIDMap.Add(E0y0);
1457   theShapeIDMap.Add(E1y0);
1458   theShapeIDMap.Add(E0y1);
1459   theShapeIDMap.Add(E1y1);
1460
1461   theShapeIDMap.Add(E00z);
1462   theShapeIDMap.Add(E10z);
1463   theShapeIDMap.Add(E01z);
1464   theShapeIDMap.Add(E11z);
1465
1466   theShapeIDMap.Add(Fxy0);
1467   theShapeIDMap.Add(Fxy1);
1468   theShapeIDMap.Add(Fx0z);
1469   theShapeIDMap.Add(Fx1z);
1470   theShapeIDMap.Add(F0yz);
1471   theShapeIDMap.Add(F1yz);
1472   
1473   theShapeIDMap.Add(theShell);
1474
1475   return true;
1476 }
1477
1478 //================================================================================
1479 /*!
1480  * \brief Initialize block geometry with shapes from theShapeIDMap
1481   * \param theShapeIDMap - map of block subshapes
1482   * \retval bool - is a success
1483  */
1484 //================================================================================
1485
1486 bool SMESH_Block::LoadBlockShapes(const TopTools_IndexedMapOfOrientedShape& theShapeIDMap)
1487 {
1488   init();
1489
1490   // store shapes geometry
1491   for ( int shapeID = 1; shapeID < theShapeIDMap.Extent(); shapeID++ )
1492   {
1493     const TopoDS_Shape& S = theShapeIDMap( shapeID );
1494     switch ( S.ShapeType() )
1495     {
1496     case TopAbs_VERTEX: {
1497
1498       if ( !IsVertexID( ID_V111 )) return false;
1499       myPnt[ shapeID - ID_V000 ] = BRep_Tool::Pnt( TopoDS::Vertex( S )).XYZ();
1500       break;
1501     }
1502     case TopAbs_EDGE: {
1503
1504       if ( !IsEdgeID( shapeID )) return false;
1505       const TopoDS_Edge& edge = TopoDS::Edge( S );
1506       TEdge& tEdge = myEdge[ shapeID - ID_FirstE ];
1507       tEdge.Set( shapeID,
1508                  new BRepAdaptor_Curve( edge ),
1509                  IsForwardEdge( edge, theShapeIDMap ));
1510       break;
1511     }
1512     case TopAbs_FACE: {
1513
1514       if ( !LoadFace( TopoDS::Face( S ), shapeID, theShapeIDMap ))
1515         return false;
1516       break;
1517     }
1518     default: break;
1519     }
1520   } // loop on shapes in theShapeIDMap
1521
1522   return true;
1523 }
1524
1525 //================================================================================
1526 /*!
1527  * \brief Load face geometry
1528   * \param theFace - face
1529   * \param theFaceID - face in-block ID
1530   * \param theShapeIDMap - map of block subshapes
1531   * \retval bool - is a success
1532  * 
1533  * It is enough to compute params or coordinates on the face.
1534  * Face subshapes must be loaded into theShapeIDMap before
1535  */
1536 //================================================================================
1537
1538 bool SMESH_Block::LoadFace(const TopoDS_Face& theFace,
1539                            const int          theFaceID,
1540                            const TopTools_IndexedMapOfOrientedShape& theShapeIDMap)
1541 {
1542   if ( !IsFaceID( theFaceID ) ) return false;
1543   // pcurves
1544   Adaptor2d_Curve2d* c2d[4];
1545   bool isForward[4];
1546   vector< int > edgeIdVec;
1547   GetFaceEdgesIDs( theFaceID, edgeIdVec );
1548   for ( int iE = 0; iE < edgeIdVec.size(); iE++ ) // loop on 4 edges
1549   {
1550     if ( edgeIdVec[ iE ] > theShapeIDMap.Extent() )
1551       return false;
1552     const TopoDS_Edge& edge = TopoDS::Edge( theShapeIDMap( edgeIdVec[ iE ]));
1553     c2d[ iE ] = new BRepAdaptor_Curve2d( edge, theFace );
1554     isForward[ iE ] = IsForwardEdge( edge, theShapeIDMap );
1555   }
1556   TFace& tFace = myFace[ theFaceID - ID_FirstF ];
1557   tFace.Set( theFaceID, new BRepAdaptor_Surface( theFace ), c2d, isForward );
1558   return true;
1559 }
1560
1561 //================================================================================
1562 /*!
1563  * \brief/ Insert theShape into theShapeIDMap with theShapeID
1564   * \param theShape - shape to insert
1565   * \param theShapeID - shape in-block ID
1566   * \param theShapeIDMap - map of block subshapes
1567  */
1568 //================================================================================
1569
1570 bool SMESH_Block::Insert(const TopoDS_Shape& theShape,
1571                          const int           theShapeID,
1572                          TopTools_IndexedMapOfOrientedShape& theShapeIDMap)
1573 {
1574   if ( !theShape.IsNull() && theShapeID > 0 )
1575   {
1576     if ( theShapeIDMap.Contains( theShape ))
1577       return ( theShapeIDMap.FindIndex( theShape ) == theShapeID );
1578
1579     if ( theShapeID <= theShapeIDMap.Extent() ) {
1580         theShapeIDMap.Substitute( theShapeID, theShape );
1581     }
1582     else {
1583       while ( theShapeIDMap.Extent() < theShapeID - 1 ) {
1584         TopoDS_Compound comp;
1585         BRep_Builder().MakeCompound( comp );
1586         theShapeIDMap.Add( comp );
1587       }
1588       theShapeIDMap.Add( theShape );
1589     }
1590     return true;
1591   }
1592   return false;
1593 }
1594
1595 //=======================================================================
1596 //function : GetFaceEdgesIDs
1597 //purpose  : return edges IDs in the order u0, u1, 0v, 1v
1598 //           u0 means "|| u, v == 0"
1599 //=======================================================================
1600
1601 void SMESH_Block::GetFaceEdgesIDs (const int faceID, vector< int >& edgeVec )
1602 {
1603   edgeVec.resize( 4 );
1604   switch ( faceID ) {
1605   case ID_Fxy0:
1606     edgeVec[ 0 ] = ID_Ex00;
1607     edgeVec[ 1 ] = ID_Ex10;
1608     edgeVec[ 2 ] = ID_E0y0;
1609     edgeVec[ 3 ] = ID_E1y0;
1610     break;
1611   case ID_Fxy1:
1612     edgeVec[ 0 ] = ID_Ex01;
1613     edgeVec[ 1 ] = ID_Ex11;
1614     edgeVec[ 2 ] = ID_E0y1;
1615     edgeVec[ 3 ] = ID_E1y1;
1616     break;
1617   case ID_Fx0z:
1618     edgeVec[ 0 ] = ID_Ex00;
1619     edgeVec[ 1 ] = ID_Ex01;
1620     edgeVec[ 2 ] = ID_E00z;
1621     edgeVec[ 3 ] = ID_E10z;
1622     break;
1623   case ID_Fx1z:
1624     edgeVec[ 0 ] = ID_Ex10;
1625     edgeVec[ 1 ] = ID_Ex11;
1626     edgeVec[ 2 ] = ID_E01z;
1627     edgeVec[ 3 ] = ID_E11z;
1628     break;
1629   case ID_F0yz:
1630     edgeVec[ 0 ] = ID_E0y0;
1631     edgeVec[ 1 ] = ID_E0y1;
1632     edgeVec[ 2 ] = ID_E00z;
1633     edgeVec[ 3 ] = ID_E01z;
1634     break;
1635   case ID_F1yz:
1636     edgeVec[ 0 ] = ID_E1y0;
1637     edgeVec[ 1 ] = ID_E1y1;
1638     edgeVec[ 2 ] = ID_E10z;
1639     edgeVec[ 3 ] = ID_E11z;
1640     break;
1641   default:
1642     MESSAGE(" GetFaceEdgesIDs(), wrong face ID: " << faceID );
1643   }
1644 }
1645
1646 //=======================================================================
1647 //function : GetEdgeVertexIDs
1648 //purpose  : return vertex IDs of an edge
1649 //=======================================================================
1650
1651 void SMESH_Block::GetEdgeVertexIDs (const int edgeID, vector< int >& vertexVec )
1652 {
1653   vertexVec.resize( 2 );
1654   switch ( edgeID ) {
1655
1656   case ID_Ex00:
1657     vertexVec[ 0 ] = ID_V000;
1658     vertexVec[ 1 ] = ID_V100;
1659     break;
1660   case ID_Ex10:
1661     vertexVec[ 0 ] = ID_V010;
1662     vertexVec[ 1 ] = ID_V110;
1663     break;
1664   case ID_Ex01:
1665     vertexVec[ 0 ] = ID_V001;
1666     vertexVec[ 1 ] = ID_V101;
1667     break;
1668   case ID_Ex11:
1669     vertexVec[ 0 ] = ID_V011;
1670     vertexVec[ 1 ] = ID_V111;
1671     break;
1672
1673   case ID_E0y0:
1674     vertexVec[ 0 ] = ID_V000;
1675     vertexVec[ 1 ] = ID_V010;
1676     break;
1677   case ID_E1y0:
1678     vertexVec[ 0 ] = ID_V100;
1679     vertexVec[ 1 ] = ID_V110;
1680     break;
1681   case ID_E0y1:
1682     vertexVec[ 0 ] = ID_V001;
1683     vertexVec[ 1 ] = ID_V011;
1684     break;
1685   case ID_E1y1:
1686     vertexVec[ 0 ] = ID_V101;
1687     vertexVec[ 1 ] = ID_V111;
1688     break;
1689
1690   case ID_E00z:
1691     vertexVec[ 0 ] = ID_V000;
1692     vertexVec[ 1 ] = ID_V001;
1693     break;
1694   case ID_E10z:
1695     vertexVec[ 0 ] = ID_V100;
1696     vertexVec[ 1 ] = ID_V101;
1697     break;
1698   case ID_E01z:
1699     vertexVec[ 0 ] = ID_V010;
1700     vertexVec[ 1 ] = ID_V011;
1701     break;
1702   case ID_E11z:
1703     vertexVec[ 0 ] = ID_V110;
1704     vertexVec[ 1 ] = ID_V111;
1705     break;
1706   default:
1707     vertexVec.resize(0);
1708     MESSAGE(" GetEdgeVertexIDs(), wrong edge ID: " << edgeID );
1709   }
1710 }