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