Salome HOME
PR: merged from V5_1_4rc1
[modules/smesh.git] / src / SMESH / SMESH_MesherHelper.cxx
1 //  Copyright (C) 2007-2010  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
23 // File:      SMESH_MesherHelper.cxx
24 // Created:   15.02.06 15:22:41
25 // Author:    Sergey KUUL
26 //
27 #include "SMESH_MesherHelper.hxx"
28
29 #include "SMDS_FacePosition.hxx" 
30 #include "SMDS_EdgePosition.hxx"
31 #include "SMDS_VolumeTool.hxx"
32 #include "SMESH_subMesh.hxx"
33
34 #include <BRepAdaptor_Surface.hxx>
35 #include <BRepTools.hxx>
36 #include <BRepTools_WireExplorer.hxx>
37 #include <BRep_Tool.hxx>
38 #include <Geom2d_Curve.hxx>
39 #include <GeomAPI_ProjectPointOnCurve.hxx>
40 #include <GeomAPI_ProjectPointOnSurf.hxx>
41 #include <Geom_Curve.hxx>
42 #include <Geom_Surface.hxx>
43 #include <ShapeAnalysis.hxx>
44 #include <TopExp.hxx>
45 #include <TopExp_Explorer.hxx>
46 #include <TopTools_ListIteratorOfListOfShape.hxx>
47 #include <TopTools_MapIteratorOfMapOfShape.hxx>
48 #include <TopTools_MapOfShape.hxx>
49 #include <TopoDS.hxx>
50 #include <gp_Ax3.hxx>
51 #include <gp_Pnt2d.hxx>
52 #include <gp_Trsf.hxx>
53
54 #include <Standard_Failure.hxx>
55 #include <Standard_ErrorHandler.hxx>
56
57 #include <utilities.h>
58
59 #include <limits>
60
61 using namespace std;
62
63 #define RETURN_BAD_RESULT(msg) { MESSAGE(msg); return false; }
64
65 namespace {
66
67   gp_XYZ XYZ(const SMDS_MeshNode* n) { return gp_XYZ(n->X(), n->Y(), n->Z()); }
68
69   enum { U_periodic = 1, V_periodic = 2 };
70 }
71
72 //================================================================================
73 /*!
74  * \brief Constructor
75  */
76 //================================================================================
77
78 SMESH_MesherHelper::SMESH_MesherHelper(SMESH_Mesh& theMesh)
79   : myParIndex(0), myMesh(&theMesh), myShapeID(0), myCreateQuadratic(false)
80 {
81   myPar1[0] = myPar2[0] = myPar1[1] = myPar2[1] = 0;
82   mySetElemOnShape = ( ! myMesh->HasShapeToMesh() );
83 }
84
85 //=======================================================================
86 //function : ~SMESH_MesherHelper
87 //purpose  : 
88 //=======================================================================
89
90 SMESH_MesherHelper::~SMESH_MesherHelper()
91 {
92   TID2Projector::iterator i_proj = myFace2Projector.begin();
93   for ( ; i_proj != myFace2Projector.end(); ++i_proj )
94     delete i_proj->second;
95 }
96
97 //=======================================================================
98 //function : IsQuadraticSubMesh
99 //purpose  : Check submesh for given shape: if all elements on this shape 
100 //           are quadratic, quadratic elements will be created.
101 //           Also fill myTLinkNodeMap
102 //=======================================================================
103
104 bool SMESH_MesherHelper::IsQuadraticSubMesh(const TopoDS_Shape& aSh)
105 {
106   SMESHDS_Mesh* meshDS = GetMeshDS();
107   // we can create quadratic elements only if all elements
108   // created on subshapes of given shape are quadratic
109   // also we have to fill myTLinkNodeMap
110   myCreateQuadratic = true;
111   mySeamShapeIds.clear();
112   myDegenShapeIds.clear();
113   TopAbs_ShapeEnum subType( aSh.ShapeType()==TopAbs_FACE ? TopAbs_EDGE : TopAbs_FACE );
114   SMDSAbs_ElementType elemType( subType==TopAbs_FACE ? SMDSAbs_Face : SMDSAbs_Edge );
115
116   int nbOldLinks = myTLinkNodeMap.size();
117
118   TopExp_Explorer exp( aSh, subType );
119   for (; exp.More() && myCreateQuadratic; exp.Next()) {
120     if ( SMESHDS_SubMesh * subMesh = meshDS->MeshElements( exp.Current() )) {
121       if ( SMDS_ElemIteratorPtr it = subMesh->GetElements() ) {
122         while(it->more()) {
123           const SMDS_MeshElement* e = it->next();
124           if ( e->GetType() != elemType || !e->IsQuadratic() ) {
125             myCreateQuadratic = false;
126             break;
127           }
128           else {
129             // fill TLinkNodeMap
130             switch ( e->NbNodes() ) {
131             case 3:
132               AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(2)); break;
133             case 6:
134               AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(3));
135               AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(4));
136               AddTLinkNode(e->GetNode(2),e->GetNode(0),e->GetNode(5)); break;
137             case 8:
138               AddTLinkNode(e->GetNode(0),e->GetNode(1),e->GetNode(4));
139               AddTLinkNode(e->GetNode(1),e->GetNode(2),e->GetNode(5));
140               AddTLinkNode(e->GetNode(2),e->GetNode(3),e->GetNode(6));
141               AddTLinkNode(e->GetNode(3),e->GetNode(0),e->GetNode(7));
142               break;
143             default:
144               myCreateQuadratic = false;
145               break;
146             }
147           }
148         }
149       }
150     }
151   }
152
153   if ( nbOldLinks == myTLinkNodeMap.size() )
154     myCreateQuadratic = false;
155
156   if(!myCreateQuadratic) {
157     myTLinkNodeMap.clear();
158   }
159   SetSubShape( aSh );
160
161   return myCreateQuadratic;
162 }
163
164 //=======================================================================
165 //function : SetSubShape
166 //purpose  : Set geomerty to make elements on
167 //=======================================================================
168
169 void SMESH_MesherHelper::SetSubShape(const int aShID)
170 {
171   if ( aShID == myShapeID )
172     return;
173   if ( aShID > 1 )
174     SetSubShape( GetMeshDS()->IndexToShape( aShID ));
175   else
176     SetSubShape( TopoDS_Shape() );
177 }
178
179 //=======================================================================
180 //function : SetSubShape
181 //purpose  : Set geomerty to create elements on
182 //=======================================================================
183
184 void SMESH_MesherHelper::SetSubShape(const TopoDS_Shape& aSh)
185 {
186   if ( myShape.IsSame( aSh ))
187     return;
188
189   myShape = aSh;
190   mySeamShapeIds.clear();
191   myDegenShapeIds.clear();
192
193   if ( myShape.IsNull() ) {
194     myShapeID  = 0;
195     return;
196   }
197   SMESHDS_Mesh* meshDS = GetMeshDS();
198   myShapeID = meshDS->ShapeToIndex(aSh);
199   myParIndex = 0;
200
201   // treatment of periodic faces
202   for ( TopExp_Explorer eF( aSh, TopAbs_FACE ); eF.More(); eF.Next() )
203   {
204     const TopoDS_Face& face = TopoDS::Face( eF.Current() );
205     BRepAdaptor_Surface surface( face );
206     if ( surface.IsUPeriodic() || surface.IsVPeriodic() )
207     {
208       for (TopExp_Explorer exp( face, TopAbs_EDGE ); exp.More(); exp.Next())
209       {
210         // look for a seam edge
211         const TopoDS_Edge& edge = TopoDS::Edge( exp.Current() );
212         if ( BRep_Tool::IsClosed( edge, face )) {
213           // initialize myPar1, myPar2 and myParIndex
214           gp_Pnt2d uv1, uv2;
215           BRep_Tool::UVPoints( edge, face, uv1, uv2 );
216           if ( Abs( uv1.Coord(1) - uv2.Coord(1) ) < Abs( uv1.Coord(2) - uv2.Coord(2) ))
217           {
218             myParIndex |= U_periodic;
219             myPar1[0] = surface.FirstUParameter();
220             myPar2[0] = surface.LastUParameter();
221           }
222           else {
223             myParIndex |= V_periodic;
224             myPar1[1] = surface.FirstVParameter();
225             myPar2[1] = surface.LastVParameter();
226           }
227           // store seam shape indices, negative if shape encounters twice
228           int edgeID = meshDS->ShapeToIndex( edge );
229           mySeamShapeIds.insert( IsSeamShape( edgeID ) ? -edgeID : edgeID );
230           for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() ) {
231             int vertexID = meshDS->ShapeToIndex( v.Current() );
232             mySeamShapeIds.insert( IsSeamShape( vertexID ) ? -vertexID : vertexID );
233           }
234         }
235
236         // look for a degenerated edge
237         if ( BRep_Tool::Degenerated( edge )) {
238           myDegenShapeIds.insert( meshDS->ShapeToIndex( edge ));
239           for ( TopExp_Explorer v( edge, TopAbs_VERTEX ); v.More(); v.Next() )
240             myDegenShapeIds.insert( meshDS->ShapeToIndex( v.Current() ));
241         }
242       }
243     }
244   }
245 }
246
247 //=======================================================================
248 //function : GetNodeUVneedInFaceNode
249 //purpose  : Check if inFaceNode argument is necessary for call GetNodeUV(F,..)
250 //           Return true if the face is periodic.
251 //           If F is Null, answer about subshape set through IsQuadraticSubMesh() or
252 //           * SetSubShape()
253 //=======================================================================
254
255 bool SMESH_MesherHelper::GetNodeUVneedInFaceNode(const TopoDS_Face& F) const
256 {
257   if ( F.IsNull() ) return !mySeamShapeIds.empty();
258
259   if ( !F.IsNull() && !myShape.IsNull() && myShape.IsSame( F ))
260     return !mySeamShapeIds.empty();
261
262   TopLoc_Location loc;
263   Handle(Geom_Surface) aSurface = BRep_Tool::Surface( F,loc );
264   if ( !aSurface.IsNull() )
265     return ( aSurface->IsUPeriodic() || aSurface->IsVPeriodic() );
266
267   return false;
268 }
269
270 //=======================================================================
271 //function : IsMedium
272 //purpose  : 
273 //=======================================================================
274
275 bool SMESH_MesherHelper::IsMedium(const SMDS_MeshNode*      node,
276                                   const SMDSAbs_ElementType typeToCheck)
277 {
278   return SMESH_MeshEditor::IsMedium( node, typeToCheck );
279 }
280
281 //=======================================================================
282 //function : GetSubShapeByNode
283 //purpose  : Return support shape of a node
284 //=======================================================================
285
286 TopoDS_Shape SMESH_MesherHelper::GetSubShapeByNode(const SMDS_MeshNode* node,
287                                                    SMESHDS_Mesh*        meshDS)
288 {
289   int shapeID = node->GetPosition()->GetShapeId();
290   if ( 0 < shapeID && shapeID <= meshDS->MaxShapeIndex() )
291     return meshDS->IndexToShape( shapeID );
292   else
293     return TopoDS_Shape();
294 }
295
296
297 //=======================================================================
298 //function : AddTLinkNode
299 //purpose  : add a link in my data structure
300 //=======================================================================
301
302 void SMESH_MesherHelper::AddTLinkNode(const SMDS_MeshNode* n1,
303                                       const SMDS_MeshNode* n2,
304                                       const SMDS_MeshNode* n12)
305 {
306   // add new record to map
307   SMESH_TLink link( n1, n2 );
308   myTLinkNodeMap.insert( make_pair(link,n12));
309 }
310
311 //=======================================================================
312 //function : GetUVOnSeam
313 //purpose  : Select UV on either of 2 pcurves of a seam edge, closest to the given UV
314 //=======================================================================
315
316 gp_Pnt2d SMESH_MesherHelper::GetUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const
317 {
318   gp_Pnt2d result = uv1;
319   for ( int i = U_periodic; i <= V_periodic ; ++i )
320   {
321     if ( myParIndex & i )
322     {
323       double p1 = uv1.Coord( i );
324       double dp1 = Abs( p1-myPar1[i-1]), dp2 = Abs( p1-myPar2[i-1]);
325       if ( myParIndex == i ||
326            dp1 < ( myPar2[i-1] - myPar2[i-1] ) / 100. ||
327            dp2 < ( myPar2[i-1] - myPar2[i-1] ) / 100. )
328       {
329         double p2 = uv2.Coord( i );
330         double p1Alt = ( dp1 < dp2 ) ? myPar2[i-1] : myPar1[i-1];
331         if ( Abs( p2 - p1 ) > Abs( p2 - p1Alt ))
332           result.SetCoord( i, p1Alt );
333       }
334     }
335   }
336   return result;
337 }
338
339 //=======================================================================
340 //function : GetNodeUV
341 //purpose  : Return node UV on face
342 //=======================================================================
343
344 gp_XY SMESH_MesherHelper::GetNodeUV(const TopoDS_Face&   F,
345                                     const SMDS_MeshNode* n,
346                                     const SMDS_MeshNode* n2,
347                                     bool*                check) const
348 {
349   gp_Pnt2d uv( Precision::Infinite(), Precision::Infinite() );
350   const SMDS_PositionPtr Pos = n->GetPosition();
351   bool uvOK = false;
352   if(Pos->GetTypeOfPosition()==SMDS_TOP_FACE)
353   {
354     // node has position on face
355     const SMDS_FacePosition* fpos =
356       static_cast<const SMDS_FacePosition*>(n->GetPosition());
357     uv.SetCoord(fpos->GetUParameter(),fpos->GetVParameter());
358     if ( check )
359       uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( F ));
360   }
361   else if(Pos->GetTypeOfPosition()==SMDS_TOP_EDGE)
362   {
363     // node has position on edge => it is needed to find
364     // corresponding edge from face, get pcurve for this
365     // edge and retrieve value from this pcurve
366     const SMDS_EdgePosition* epos =
367       static_cast<const SMDS_EdgePosition*>(n->GetPosition());
368     int edgeID = Pos->GetShapeId();
369     TopoDS_Edge E = TopoDS::Edge(GetMeshDS()->IndexToShape(edgeID));
370     double f, l, u = epos->GetUParameter();
371     Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l);
372     bool validU = ( f < u && u < l );
373     if ( validU )
374       uv = C2d->Value( u );
375     else
376       uv.SetCoord(0.,0.);
377     if ( check || !validU )
378       uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), BRep_Tool::Tolerance( E ),/*force=*/ !validU );
379
380     // for a node on a seam edge select one of UVs on 2 pcurves
381     if ( n2 && IsSeamShape( edgeID ) )
382     {
383       uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
384     }
385     else
386     { // adjust uv to period
387       TopLoc_Location loc;
388       Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
389       Standard_Boolean isUPeriodic = S->IsUPeriodic();
390       Standard_Boolean isVPeriodic = S->IsVPeriodic();
391       if ( isUPeriodic || isVPeriodic ) {
392         Standard_Real UF,UL,VF,VL;
393         S->Bounds(UF,UL,VF,VL);
394         if(isUPeriodic)
395           uv.SetX( uv.X() + ShapeAnalysis::AdjustToPeriod(uv.X(),UF,UL));
396         if(isVPeriodic)
397           uv.SetY( uv.Y() + ShapeAnalysis::AdjustToPeriod(uv.Y(),VF,VL));
398       }
399     }
400   }
401   else if(Pos->GetTypeOfPosition()==SMDS_TOP_VERTEX)
402   {
403     if ( int vertexID = n->GetPosition()->GetShapeId() ) {
404       const TopoDS_Vertex& V = TopoDS::Vertex(GetMeshDS()->IndexToShape(vertexID));
405       try {
406         uv = BRep_Tool::Parameters( V, F );
407         uvOK = true;
408       }
409       catch (Standard_Failure& exc) {
410       }
411       if ( !uvOK ) {
412         for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() )
413           uvOK = ( V == vert.Current() );
414         if ( !uvOK ) {
415 #ifdef _DEBUG_
416           MESSAGE ( "SMESH_MesherHelper::GetNodeUV(); Vertex " << vertexID
417                << " not in face " << GetMeshDS()->ShapeToIndex( F ) );
418 #endif
419           // get UV of a vertex closest to the node
420           double dist = 1e100;
421           gp_Pnt pn = XYZ( n );
422           for ( TopExp_Explorer vert(F,TopAbs_VERTEX); !uvOK && vert.More(); vert.Next() ) {
423             TopoDS_Vertex curV = TopoDS::Vertex( vert.Current() );
424             gp_Pnt p = BRep_Tool::Pnt( curV );
425             double curDist = p.SquareDistance( pn );
426             if ( curDist < dist ) {
427               dist = curDist;
428               uv = BRep_Tool::Parameters( curV, F );
429               uvOK = ( dist < DBL_MIN );
430             }
431           }
432         }
433         else {
434           uvOK = false;
435           TopTools_ListIteratorOfListOfShape it( myMesh->GetAncestors( V ));
436           for ( ; it.More(); it.Next() ) {
437             if ( it.Value().ShapeType() == TopAbs_EDGE ) {
438               const TopoDS_Edge & edge = TopoDS::Edge( it.Value() );
439               double f,l;
440               Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(edge, F, f, l);
441               if ( !C2d.IsNull() ) {
442                 double u = ( V == TopExp::FirstVertex( edge ) ) ?  f : l;
443                 uv = C2d->Value( u );
444                 uvOK = true;
445                 break;
446               }
447             }
448           }
449         }
450       }
451       if ( n2 && IsSeamShape( vertexID ) )
452         uv = GetUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
453     }
454   }
455
456   if ( check )
457     *check = uvOK;
458
459   return uv.XY();
460 }
461
462 //=======================================================================
463 //function : CheckNodeUV
464 //purpose  : Check and fix node UV on a face
465 //=======================================================================
466
467 bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face&   F,
468                                      const SMDS_MeshNode* n,
469                                      gp_XY&               uv,
470                                      const double         tol,
471                                      const bool           force) const
472 {
473   if ( force || !myOkNodePosShapes.count( n->GetPosition()->GetShapeId() ))
474   {
475     // check that uv is correct
476     TopLoc_Location loc;
477     Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
478     gp_Pnt nodePnt = XYZ( n );
479     if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
480     if ( Precision::IsInfinite( uv.X() ) ||
481          Precision::IsInfinite( uv.Y() ) ||
482          nodePnt.Distance( surface->Value( uv.X(), uv.Y() )) > tol )
483     {
484       // uv incorrect, project the node to surface
485       GeomAPI_ProjectPointOnSurf& projector = GetProjector( F, loc, tol );
486       projector.Perform( nodePnt );
487       if ( !projector.IsDone() || projector.NbPoints() < 1 )
488       {
489         MESSAGE( "SMESH_MesherHelper::CheckNodeUV() failed to project" );
490         return false;
491       }
492       Quantity_Parameter U,V;
493       projector.LowerDistanceParameters(U,V);
494       if ( nodePnt.Distance( surface->Value( U, V )) > tol )
495       {
496         MESSAGE( "SMESH_MesherHelper::CheckNodeUV(), invalid projection" );
497         return false;
498       }
499       uv.SetCoord( U,V );
500     }
501     else if ( uv.Modulus() > numeric_limits<double>::min() )
502     {
503       ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->GetPosition()->GetShapeId() );
504     }
505   }
506   return true;
507 }
508
509 //=======================================================================
510 //function : GetProjector
511 //purpose  : Return projector intitialized by given face without location, which is returned
512 //=======================================================================
513
514 GeomAPI_ProjectPointOnSurf& SMESH_MesherHelper::GetProjector(const TopoDS_Face& F,
515                                                              TopLoc_Location&   loc,
516                                                              double             tol ) const
517 {
518   Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
519   int faceID = GetMeshDS()->ShapeToIndex( F );
520   TID2Projector& i2proj = const_cast< TID2Projector&>( myFace2Projector );
521   TID2Projector::iterator i_proj = i2proj.find( faceID );
522   if ( i_proj == i2proj.end() )
523   {
524     if ( tol == 0 ) tol = BRep_Tool::Tolerance( F );
525     double U1, U2, V1, V2;
526     surface->Bounds(U1, U2, V1, V2);
527     GeomAPI_ProjectPointOnSurf* proj = new GeomAPI_ProjectPointOnSurf();
528     proj->Init( surface, U1, U2, V1, V2, tol );
529     i_proj = i2proj.insert( make_pair( faceID, proj )).first;
530   }
531   return *( i_proj->second );
532 }
533
534 namespace
535 {
536   gp_XY AverageUV(const gp_XY& uv1, const gp_XY& uv2) { return ( uv1 + uv2 ) / 2.; }
537   gp_XY_FunPtr(Added); // define gp_XY_Added pointer to function calling gp_XY::Added(gp_XY)
538   gp_XY_FunPtr(Subtracted); 
539 }
540
541 //=======================================================================
542 //function : applyIn2D
543 //purpose  : Perform given operation on two 2d points in parameric space of given surface.
544 //           It takes into account period of the surface. Use gp_XY_FunPtr macro
545 //           to easily define pointer to function of gp_XY class.
546 //=======================================================================
547
548 gp_XY SMESH_MesherHelper::applyIn2D(const Handle(Geom_Surface)& surface,
549                                     const gp_XY&                uv1,
550                                     const gp_XY&                uv2,
551                                     xyFunPtr                    fun,
552                                     const bool                  resultInPeriod)
553 {
554   Standard_Boolean isUPeriodic = surface.IsNull() ? false : surface->IsUPeriodic();
555   Standard_Boolean isVPeriodic = surface.IsNull() ? false : surface->IsVPeriodic();
556   if ( !isUPeriodic && !isVPeriodic )
557     return fun(uv1,uv2);
558
559   // move uv2 not far than half-period from uv1
560   double u2 = 
561     uv2.X()+(isUPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.X(),uv1.X(),surface->UPeriod()) :0);
562   double v2 = 
563     uv2.Y()+(isVPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.Y(),uv1.Y(),surface->VPeriod()) :0);
564
565   // execute operation
566   gp_XY res = fun( uv1, gp_XY(u2,v2) );
567
568   // move result within period
569   if ( resultInPeriod )
570   {
571     Standard_Real UF,UL,VF,VL;
572     surface->Bounds(UF,UL,VF,VL);
573     if ( isUPeriodic )
574       res.SetX( res.X() + ShapeAnalysis::AdjustToPeriod(res.X(),UF,UL));
575     if ( isVPeriodic )
576       res.SetY( res.Y() + ShapeAnalysis::AdjustToPeriod(res.Y(),VF,VL));
577   }
578
579   return res;
580 }
581 //=======================================================================
582 //function : GetMiddleUV
583 //purpose  : Return middle UV taking in account surface period
584 //=======================================================================
585
586 gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
587                                       const gp_XY&                p1,
588                                       const gp_XY&                p2)
589 {
590   return applyIn2D( surface, p1, p2, & AverageUV );
591 }
592
593 //=======================================================================
594 //function : GetNodeU
595 //purpose  : Return node U on edge
596 //=======================================================================
597
598 double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge&   E,
599                                     const SMDS_MeshNode* n,
600                                     const SMDS_MeshNode* inEdgeNode,
601                                     bool*                check)
602 {
603   double param = 0;
604   const SMDS_PositionPtr Pos = n->GetPosition();
605   if ( Pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
606   {
607     const SMDS_EdgePosition* epos =
608       static_cast<const SMDS_EdgePosition*>(n->GetPosition());
609     param =  epos->GetUParameter();
610   }
611   else if( Pos->GetTypeOfPosition() == SMDS_TOP_VERTEX )
612   {
613     if ( inEdgeNode && TopExp::FirstVertex( E ).IsSame( TopExp::LastVertex( E ))) // issue 0020128
614     {
615       Standard_Real f,l;
616       BRep_Tool::Range( E, f,l );
617       double uInEdge = GetNodeU( E, inEdgeNode );
618       param = ( fabs( uInEdge - f ) < fabs( l - uInEdge )) ? f : l;
619     }
620     else
621     {
622       SMESHDS_Mesh * meshDS = GetMeshDS();
623       int vertexID = n->GetPosition()->GetShapeId();
624       const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
625       param =  BRep_Tool::Parameter( V, E );
626     }
627   }
628   if ( check )
629     *check = CheckNodeU( E, n, param, BRep_Tool::Tolerance( E ));
630   return param;
631 }
632
633 //=======================================================================
634 //function : CheckNodeU
635 //purpose  : Check and fix node U on an edge
636 //           Return false if U is bad and could not be fixed
637 //=======================================================================
638
639 bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge&   E,
640                                     const SMDS_MeshNode* n,
641                                     double&              u,
642                                     const double         tol,
643                                     const bool           force) const
644 {
645   if ( force || !myOkNodePosShapes.count( n->GetPosition()->GetShapeId() ))
646   {
647     // check that u is correct
648     TopLoc_Location loc; double f,l;
649     Handle(Geom_Curve) curve = BRep_Tool::Curve( E,loc,f,l );
650     if ( curve.IsNull() ) // degenerated edge
651     {
652       if ( u+tol < f || u-tol > l )
653       {
654         double r = Max( 0.5, 1 - tol*n->GetID()); // to get a unique u on edge
655         u =  f*r + l*(1-r);
656       }
657     }
658     else
659     {
660       gp_Pnt nodePnt = SMESH_MeshEditor::TNodeXYZ( n );
661       if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
662       if ( nodePnt.Distance( curve->Value( u )) > tol )
663       {
664         // u incorrect, project the node to the curve
665         GeomAPI_ProjectPointOnCurve projector( nodePnt, curve, f, l );
666         if ( projector.NbPoints() < 1 )
667         {
668           MESSAGE( "SMESH_MesherHelper::CheckNodeU() failed to project" );
669           return false;
670         }
671         Quantity_Parameter U = projector.LowerDistanceParameter();
672         if ( nodePnt.Distance( curve->Value( U )) > tol )
673         {
674           MESSAGE( "SMESH_MesherHelper::CheckNodeU(), invalid projection" );
675           return false;
676         }
677         u = double( U );
678       }
679       else if ( fabs( u ) > numeric_limits<double>::min() )
680       {
681         ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->GetPosition()->GetShapeId() );
682       }
683     }
684   }
685   return true;
686 }
687
688 //=======================================================================
689 //function : GetMediumNode
690 //purpose  : Return existing or create new medium nodes between given ones
691 //=======================================================================
692
693 const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
694                                                        const SMDS_MeshNode* n2,
695                                                        bool                 force3d)
696 {
697   // Find existing node
698
699   SMESH_TLink link(n1,n2);
700   ItTLinkNode itLN = myTLinkNodeMap.find( link );
701   if ( itLN != myTLinkNodeMap.end() ) {
702     return (*itLN).second;
703   }
704
705   // Create medium node
706
707   SMDS_MeshNode* n12;
708   SMESHDS_Mesh* meshDS = GetMeshDS();
709
710   // get type of shape for the new medium node
711   int faceID = -1, edgeID = -1;
712   const SMDS_PositionPtr Pos1 = n1->GetPosition();
713   const SMDS_PositionPtr Pos2 = n2->GetPosition();
714
715   if( myShape.IsNull() )
716   {
717     if( Pos1->GetTypeOfPosition()==SMDS_TOP_FACE ) {
718       faceID = Pos1->GetShapeId();
719     }
720     else if( Pos2->GetTypeOfPosition()==SMDS_TOP_FACE ) {
721       faceID = Pos2->GetShapeId();
722     }
723
724     if( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
725       edgeID = Pos1->GetShapeId();
726     }
727     if( Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
728       edgeID = Pos2->GetShapeId();
729     }
730   }
731   // get positions of the given nodes on shapes
732   TopoDS_Edge E; double u [2];
733   TopoDS_Face F; gp_XY  uv[2];
734   bool uvOK[2] = { false, false };
735   TopAbs_ShapeEnum shapeType = myShape.IsNull() ? TopAbs_SHAPE : myShape.ShapeType();
736   if ( faceID>0 || shapeType == TopAbs_FACE)
737   {
738     if( myShape.IsNull() )
739       F = TopoDS::Face(meshDS->IndexToShape(faceID));
740     else {
741       F = TopoDS::Face(myShape);
742       faceID = myShapeID;
743     }
744     uv[0] = GetNodeUV(F,n1,n2, force3d ? 0 : &uvOK[0]);
745     uv[1] = GetNodeUV(F,n2,n1, force3d ? 0 : &uvOK[1]);
746   }
747   else if (edgeID>0 || shapeType == TopAbs_EDGE)
748   {
749     if( myShape.IsNull() )
750       E = TopoDS::Edge(meshDS->IndexToShape(edgeID));
751     else {
752       E = TopoDS::Edge(myShape);
753       edgeID = myShapeID;
754     }
755     u[0] = GetNodeU(E,n1,n2, force3d ? 0 : &uvOK[0]);
756     u[1] = GetNodeU(E,n2,n1, force3d ? 0 : &uvOK[1]);
757   }
758   if(!force3d)
759   {
760     // we try to create medium node using UV parameters of
761     // nodes, else - medium between corresponding 3d points
762     if( ! F.IsNull() )
763     {
764       if ( uvOK[0] && uvOK[1] )
765       {
766         if ( IsDegenShape( Pos1->GetShapeId() ))
767           if ( myParIndex & U_periodic ) uv[0].SetCoord( 1, uv[1].Coord( 1 ));
768           else                           uv[0].SetCoord( 2, uv[1].Coord( 2 ));
769         else if ( IsDegenShape( Pos2->GetShapeId() ))
770           if ( myParIndex & U_periodic ) uv[1].SetCoord( 1, uv[0].Coord( 1 ));
771           else                           uv[1].SetCoord( 2, uv[0].Coord( 2 ));
772
773         TopLoc_Location loc;
774         Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
775         gp_XY UV = GetMiddleUV( S, uv[0], uv[1] );
776         gp_Pnt P = S->Value( UV.X(), UV.Y() ).Transformed(loc);
777         n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
778         meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y());
779         myTLinkNodeMap.insert(make_pair(link,n12));
780         return n12;
781       }
782     }
783     else if ( !E.IsNull() )
784     {
785       double f,l;
786       Handle(Geom_Curve) C = BRep_Tool::Curve(E, f, l);
787       if(!C.IsNull())
788       {
789         Standard_Boolean isPeriodic = C->IsPeriodic();
790         double U;
791         if(isPeriodic) {
792           Standard_Real Period = C->Period();
793           Standard_Real p = u[1]+ShapeAnalysis::AdjustByPeriod(u[1],u[0],Period);
794           Standard_Real pmid = (u[0]+p)/2.;
795           U = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
796         }
797         else
798           U = (u[0]+u[1])/2.;
799
800         gp_Pnt P = C->Value( U );
801         n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
802         meshDS->SetNodeOnEdge(n12, edgeID, U);
803         myTLinkNodeMap.insert(make_pair(link,n12));
804         return n12;
805       }
806     }
807   }
808   // 3d variant
809   double x = ( n1->X() + n2->X() )/2.;
810   double y = ( n1->Y() + n2->Y() )/2.;
811   double z = ( n1->Z() + n2->Z() )/2.;
812   n12 = meshDS->AddNode(x,y,z);
813   if ( !F.IsNull() )
814   {
815     gp_XY UV = ( uv[0] + uv[1] ) / 2.;
816     CheckNodeUV( F, n12, UV, BRep_Tool::Tolerance( F ), /*force=*/true);
817     meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y() );
818   }
819   else if ( !E.IsNull() )
820   {
821     double U = ( u[0] + u[1] ) / 2.;
822     CheckNodeU( E, n12, U, BRep_Tool::Tolerance( E ), /*force=*/true);
823     meshDS->SetNodeOnEdge(n12, edgeID, U);
824   }
825   else if ( myShapeID > 1 )
826   {
827     meshDS->SetNodeInVolume(n12, myShapeID);
828   }
829   myTLinkNodeMap.insert( make_pair( link, n12 ));
830   return n12;
831 }
832
833 //=======================================================================
834 //function : AddNode
835 //purpose  : Creates a node
836 //=======================================================================
837
838 SMDS_MeshNode* SMESH_MesherHelper::AddNode(double x, double y, double z, int ID)
839 {
840   SMESHDS_Mesh * meshDS = GetMeshDS();
841   SMDS_MeshNode* node = 0;
842   if ( ID )
843     node = meshDS->AddNodeWithID( x, y, z, ID );
844   else
845     node = meshDS->AddNode( x, y, z );
846   if ( mySetElemOnShape && myShapeID > 0 ) {
847     switch ( myShape.ShapeType() ) {
848     case TopAbs_SOLID:  meshDS->SetNodeInVolume( node, myShapeID); break;
849     case TopAbs_SHELL:  meshDS->SetNodeInVolume( node, myShapeID); break;
850     case TopAbs_FACE:   meshDS->SetNodeOnFace(   node, myShapeID); break;
851     case TopAbs_EDGE:   meshDS->SetNodeOnEdge(   node, myShapeID); break;
852     case TopAbs_VERTEX: meshDS->SetNodeOnVertex( node, myShapeID); break;
853     default: ;
854     }
855   }
856   return node;
857 }
858
859 //=======================================================================
860 //function : AddEdge
861 //purpose  : Creates quadratic or linear edge
862 //=======================================================================
863
864 SMDS_MeshEdge* SMESH_MesherHelper::AddEdge(const SMDS_MeshNode* n1,
865                                            const SMDS_MeshNode* n2,
866                                            const int            id,
867                                            const bool           force3d)
868 {
869   SMESHDS_Mesh * meshDS = GetMeshDS();
870   
871   SMDS_MeshEdge* edge = 0;
872   if (myCreateQuadratic) {
873     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
874     if(id)
875       edge = meshDS->AddEdgeWithID(n1, n2, n12, id);
876     else
877       edge = meshDS->AddEdge(n1, n2, n12);
878   }
879   else {
880     if(id)
881       edge = meshDS->AddEdgeWithID(n1, n2, id);
882     else
883       edge = meshDS->AddEdge(n1, n2);
884   }
885
886   if ( mySetElemOnShape && myShapeID > 0 )
887     meshDS->SetMeshElementOnShape( edge, myShapeID );
888
889   return edge;
890 }
891
892 //=======================================================================
893 //function : AddFace
894 //purpose  : Creates quadratic or linear triangle
895 //=======================================================================
896
897 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
898                                            const SMDS_MeshNode* n2,
899                                            const SMDS_MeshNode* n3,
900                                            const int id,
901                                            const bool force3d)
902 {
903   SMESHDS_Mesh * meshDS = GetMeshDS();
904   SMDS_MeshFace* elem = 0;
905
906   if( n1==n2 || n2==n3 || n3==n1 )
907     return elem;
908
909   if(!myCreateQuadratic) {
910     if(id)
911       elem = meshDS->AddFaceWithID(n1, n2, n3, id);
912     else
913       elem = meshDS->AddFace(n1, n2, n3);
914   }
915   else {
916     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
917     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
918     const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
919
920     if(id)
921       elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, id);
922     else
923       elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31);
924   }
925   if ( mySetElemOnShape && myShapeID > 0 )
926     meshDS->SetMeshElementOnShape( elem, myShapeID );
927
928   return elem;
929 }
930
931 //=======================================================================
932 //function : AddFace
933 //purpose  : Creates quadratic or linear quadrangle
934 //=======================================================================
935
936 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
937                                            const SMDS_MeshNode* n2,
938                                            const SMDS_MeshNode* n3,
939                                            const SMDS_MeshNode* n4,
940                                            const int            id,
941                                            const bool           force3d)
942 {
943   SMESHDS_Mesh * meshDS = GetMeshDS();
944   SMDS_MeshFace* elem = 0;
945
946   if( n1==n2 ) {
947     return AddFace(n1,n3,n4,id,force3d);
948   }
949   if( n1==n3 ) {
950     return AddFace(n1,n2,n4,id,force3d);
951   }
952   if( n1==n4 ) {
953     return AddFace(n1,n2,n3,id,force3d);
954   }
955   if( n2==n3 ) {
956     return AddFace(n1,n2,n4,id,force3d);
957   }
958   if( n2==n4 ) {
959     return AddFace(n1,n2,n3,id,force3d);
960   }
961   if( n3==n4 ) {
962     return AddFace(n1,n2,n3,id,force3d);
963   }
964
965   if(!myCreateQuadratic) {
966     if(id)
967       elem = meshDS->AddFaceWithID(n1, n2, n3, n4, id);
968     else
969       elem = meshDS->AddFace(n1, n2, n3, n4);
970   }
971   else {
972     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
973     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
974     const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
975     const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
976
977     if(id)
978       elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, id);
979     else
980       elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41);
981   }
982   if ( mySetElemOnShape && myShapeID > 0 )
983     meshDS->SetMeshElementOnShape( elem, myShapeID );
984
985   return elem;
986 }
987
988 //=======================================================================
989 //function : AddVolume
990 //purpose  : Creates quadratic or linear prism
991 //=======================================================================
992
993 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
994                                                const SMDS_MeshNode* n2,
995                                                const SMDS_MeshNode* n3,
996                                                const SMDS_MeshNode* n4,
997                                                const SMDS_MeshNode* n5,
998                                                const SMDS_MeshNode* n6,
999                                                const int id,
1000                                                const bool force3d)
1001 {
1002   SMESHDS_Mesh * meshDS = GetMeshDS();
1003   SMDS_MeshVolume* elem = 0;
1004   if(!myCreateQuadratic) {
1005     if(id)
1006       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, id);
1007     else
1008       elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6);
1009   }
1010   else {
1011     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1012     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1013     const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1014
1015     const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1016     const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1017     const SMDS_MeshNode* n64 = GetMediumNode(n6,n4,force3d);
1018
1019     const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1020     const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1021     const SMDS_MeshNode* n36 = GetMediumNode(n3,n6,force3d);
1022
1023     if(id)
1024       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, 
1025                                      n12, n23, n31, n45, n56, n64, n14, n25, n36, id);
1026     else
1027       elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6,
1028                                n12, n23, n31, n45, n56, n64, n14, n25, n36);
1029   }
1030   if ( mySetElemOnShape && myShapeID > 0 )
1031     meshDS->SetMeshElementOnShape( elem, myShapeID );
1032
1033   return elem;
1034 }
1035
1036 //=======================================================================
1037 //function : AddVolume
1038 //purpose  : Creates quadratic or linear tetrahedron
1039 //=======================================================================
1040
1041 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1042                                                const SMDS_MeshNode* n2,
1043                                                const SMDS_MeshNode* n3,
1044                                                const SMDS_MeshNode* n4,
1045                                                const int id, 
1046                                                const bool force3d)
1047 {
1048   SMESHDS_Mesh * meshDS = GetMeshDS();
1049   SMDS_MeshVolume* elem = 0;
1050   if(!myCreateQuadratic) {
1051     if(id)
1052       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, id);
1053     else
1054       elem = meshDS->AddVolume(n1, n2, n3, n4);
1055   }
1056   else {
1057     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1058     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1059     const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1060
1061     const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1062     const SMDS_MeshNode* n24 = GetMediumNode(n2,n4,force3d);
1063     const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1064
1065     if(id)
1066       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34, id);
1067     else
1068       elem = meshDS->AddVolume(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34);
1069   }
1070   if ( mySetElemOnShape && myShapeID > 0 )
1071     meshDS->SetMeshElementOnShape( elem, myShapeID );
1072
1073   return elem;
1074 }
1075
1076 //=======================================================================
1077 //function : AddVolume
1078 //purpose  : Creates quadratic or linear pyramid
1079 //=======================================================================
1080
1081 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1082                                                const SMDS_MeshNode* n2,
1083                                                const SMDS_MeshNode* n3,
1084                                                const SMDS_MeshNode* n4,
1085                                                const SMDS_MeshNode* n5,
1086                                                const int id, 
1087                                                const bool force3d)
1088 {
1089   SMDS_MeshVolume* elem = 0;
1090   if(!myCreateQuadratic) {
1091     if(id)
1092       elem = GetMeshDS()->AddVolumeWithID(n1, n2, n3, n4, n5, id);
1093     else
1094       elem = GetMeshDS()->AddVolume(n1, n2, n3, n4, n5);
1095   }
1096   else {
1097     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1098     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1099     const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1100     const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1101
1102     const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1103     const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1104     const SMDS_MeshNode* n35 = GetMediumNode(n3,n5,force3d);
1105     const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1106
1107     if(id)
1108       elem = GetMeshDS()->AddVolumeWithID ( n1,  n2,  n3,  n4,  n5,
1109                                             n12, n23, n34, n41,
1110                                             n15, n25, n35, n45,
1111                                             id);
1112     else
1113       elem = GetMeshDS()->AddVolume( n1,  n2,  n3,  n4,  n5,
1114                                      n12, n23, n34, n41,
1115                                      n15, n25, n35, n45);
1116   }
1117   if ( mySetElemOnShape && myShapeID > 0 )
1118     GetMeshDS()->SetMeshElementOnShape( elem, myShapeID );
1119
1120   return elem;
1121 }
1122
1123 //=======================================================================
1124 //function : AddVolume
1125 //purpose  : Creates quadratic or linear hexahedron
1126 //=======================================================================
1127
1128 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1129                                                const SMDS_MeshNode* n2,
1130                                                const SMDS_MeshNode* n3,
1131                                                const SMDS_MeshNode* n4,
1132                                                const SMDS_MeshNode* n5,
1133                                                const SMDS_MeshNode* n6,
1134                                                const SMDS_MeshNode* n7,
1135                                                const SMDS_MeshNode* n8,
1136                                                const int id,
1137                                                const bool force3d)
1138 {
1139   SMESHDS_Mesh * meshDS = GetMeshDS();
1140   SMDS_MeshVolume* elem = 0;
1141   if(!myCreateQuadratic) {
1142     if(id)
1143       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8, id);
1144     else
1145       elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8);
1146   }
1147   else {
1148     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1149     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1150     const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1151     const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1152
1153     const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1154     const SMDS_MeshNode* n67 = GetMediumNode(n6,n7,force3d);
1155     const SMDS_MeshNode* n78 = GetMediumNode(n7,n8,force3d);
1156     const SMDS_MeshNode* n85 = GetMediumNode(n8,n5,force3d);
1157
1158     const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1159     const SMDS_MeshNode* n26 = GetMediumNode(n2,n6,force3d);
1160     const SMDS_MeshNode* n37 = GetMediumNode(n3,n7,force3d);
1161     const SMDS_MeshNode* n48 = GetMediumNode(n4,n8,force3d);
1162
1163     if(id)
1164       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
1165                                      n12, n23, n34, n41, n56, n67,
1166                                      n78, n85, n15, n26, n37, n48, id);
1167     else
1168       elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
1169                                n12, n23, n34, n41, n56, n67,
1170                                n78, n85, n15, n26, n37, n48);
1171   }
1172   if ( mySetElemOnShape && myShapeID > 0 )
1173     meshDS->SetMeshElementOnShape( elem, myShapeID );
1174
1175   return elem;
1176 }
1177
1178 //=======================================================================
1179 //function : LoadNodeColumns
1180 //purpose  : Load nodes bound to face into a map of node columns
1181 //=======================================================================
1182
1183 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
1184                                          const TopoDS_Face& theFace,
1185                                          const TopoDS_Edge& theBaseEdge,
1186                                          SMESHDS_Mesh*      theMesh)
1187 {
1188   SMESHDS_SubMesh* faceSubMesh = theMesh->MeshElements( theFace );
1189   if ( !faceSubMesh || faceSubMesh->NbElements() == 0 )
1190     return false;
1191
1192   // get nodes on theBaseEdge sorted by param on edge and initialize theParam2ColumnMap with them
1193
1194   map< double, const SMDS_MeshNode*> sortedBaseNodes;
1195   if ( !SMESH_Algo::GetSortedNodesOnEdge( theMesh, theBaseEdge,/*noMedium=*/true, sortedBaseNodes)
1196        || sortedBaseNodes.size() < 2 )
1197     return false;
1198
1199   int nbRows = faceSubMesh->NbElements() / ( sortedBaseNodes.size()-1 ) + 1;
1200   map< double, const SMDS_MeshNode*>::iterator u_n = sortedBaseNodes.begin();
1201   double f = u_n->first, range = sortedBaseNodes.rbegin()->first - f;
1202   for ( ; u_n != sortedBaseNodes.end(); u_n++ )
1203   {
1204     double par = ( u_n->first - f ) / range;
1205     vector<const SMDS_MeshNode*>& nCol = theParam2ColumnMap[ par ];
1206     nCol.resize( nbRows );
1207     nCol[0] = u_n->second;
1208   }
1209
1210   // fill theParam2ColumnMap column by column by passing from nodes on
1211   // theBaseEdge up via mesh faces on theFace
1212
1213   TParam2ColumnMap::iterator par_nVec_2 = theParam2ColumnMap.begin();
1214   TParam2ColumnMap::iterator par_nVec_1 = par_nVec_2++;
1215   TIDSortedElemSet emptySet, avoidSet;
1216   for ( ; par_nVec_2 != theParam2ColumnMap.end(); ++par_nVec_1, ++par_nVec_2 )
1217   {
1218     vector<const SMDS_MeshNode*>& nCol1 = par_nVec_1->second;
1219     vector<const SMDS_MeshNode*>& nCol2 = par_nVec_2->second;
1220
1221     int i1, i2, iRow = 0;
1222     const SMDS_MeshNode *n1 = nCol1[0], *n2 = nCol2[0];
1223     // find face sharing node n1 and n2 and belonging to faceSubMesh
1224     while ( const SMDS_MeshElement* face =
1225             SMESH_MeshEditor::FindFaceInSet( n1, n2, emptySet, avoidSet, &i1, &i2))
1226     {
1227       if ( faceSubMesh->Contains( face ))
1228       {
1229         int nbNodes = face->IsQuadratic() ? face->NbNodes()/2 : face->NbNodes();
1230         if ( nbNodes != 4 )
1231           return false;
1232         n1 = face->GetNode( (i2+2) % 4 ); // opposite corner of quadrangle face
1233         n2 = face->GetNode( (i1+2) % 4 );
1234         if ( ++iRow >= nbRows )
1235           return false;
1236         nCol1[ iRow ] = n1;
1237         nCol2[ iRow ] = n2;
1238         avoidSet.clear();
1239       }
1240       avoidSet.insert( face );
1241     }
1242     if ( iRow + 1 < nbRows ) // compact if necessary
1243       nCol1.resize( iRow + 1 ), nCol2.resize( iRow + 1 );
1244   }
1245   return true;
1246 }
1247
1248 //=======================================================================
1249 //function : NbAncestors
1250 //purpose  : Return number of unique ancestors of the shape
1251 //=======================================================================
1252
1253 int SMESH_MesherHelper::NbAncestors(const TopoDS_Shape& shape,
1254                                     const SMESH_Mesh&   mesh,
1255                                     TopAbs_ShapeEnum    ancestorType/*=TopAbs_SHAPE*/)
1256 {
1257   TopTools_MapOfShape ancestors;
1258   TopTools_ListIteratorOfListOfShape ansIt( mesh.GetAncestors(shape) );
1259   for ( ; ansIt.More(); ansIt.Next() ) {
1260     if ( ancestorType == TopAbs_SHAPE || ansIt.Value().ShapeType() == ancestorType )
1261       ancestors.Add( ansIt.Value() );
1262   }
1263   return ancestors.Extent();
1264 }
1265
1266 //=======================================================================
1267 //function : GetSubShapeOri
1268 //purpose  : Return orientation of sub-shape in the main shape
1269 //=======================================================================
1270
1271 TopAbs_Orientation SMESH_MesherHelper::GetSubShapeOri(const TopoDS_Shape& shape,
1272                                                       const TopoDS_Shape& subShape)
1273 {
1274   TopAbs_Orientation ori = TopAbs_Orientation(-1);
1275   if ( !shape.IsNull() && !subShape.IsNull() )
1276   {
1277     TopExp_Explorer e( shape, subShape.ShapeType() );
1278     if ( shape.Orientation() >= TopAbs_INTERNAL ) // TopAbs_INTERNAL or TopAbs_EXTERNAL
1279       e.Init( shape.Oriented(TopAbs_FORWARD), subShape.ShapeType() );
1280     for ( ; e.More(); e.Next())
1281       if ( subShape.IsSame( e.Current() ))
1282         break;
1283     if ( e.More() )
1284       ori = e.Current().Orientation();
1285   }
1286   return ori;
1287 }
1288
1289 //=======================================================================
1290 //function : IsSubShape
1291 //purpose  : 
1292 //=======================================================================
1293
1294 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape,
1295                                      const TopoDS_Shape& mainShape )
1296 {
1297   if ( !shape.IsNull() && !mainShape.IsNull() )
1298   {
1299     for ( TopExp_Explorer exp( mainShape, shape.ShapeType());
1300           exp.More();
1301           exp.Next() )
1302       if ( shape.IsSame( exp.Current() ))
1303         return true;
1304   }
1305   SCRUTE((shape.IsNull()));
1306   SCRUTE((mainShape.IsNull()));
1307   return false;
1308 }
1309
1310 //=======================================================================
1311 //function : IsSubShape
1312 //purpose  : 
1313 //=======================================================================
1314
1315 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape, SMESH_Mesh* aMesh )
1316 {
1317   if ( shape.IsNull() || !aMesh )
1318     return false;
1319   return
1320     aMesh->GetMeshDS()->ShapeToIndex( shape ) ||
1321     // PAL16202
1322     shape.ShapeType() == TopAbs_COMPOUND && aMesh->GetMeshDS()->IsGroupOfSubShapes( shape );
1323 }
1324
1325 //=======================================================================
1326 //function : IsQuadraticMesh
1327 //purpose  : Check mesh without geometry for: if all elements on this shape are quadratic,
1328 //           quadratic elements will be created.
1329 //           Used then generated 3D mesh without geometry.
1330 //=======================================================================
1331
1332 SMESH_MesherHelper:: MType SMESH_MesherHelper::IsQuadraticMesh()
1333 {
1334   int NbAllEdgsAndFaces=0;
1335   int NbQuadFacesAndEdgs=0;
1336   int NbFacesAndEdges=0;
1337   //All faces and edges
1338   NbAllEdgsAndFaces = myMesh->NbEdges() + myMesh->NbFaces();
1339   
1340   //Quadratic faces and edges
1341   NbQuadFacesAndEdgs = myMesh->NbEdges(ORDER_QUADRATIC) + myMesh->NbFaces(ORDER_QUADRATIC);
1342
1343   //Linear faces and edges
1344   NbFacesAndEdges = myMesh->NbEdges(ORDER_LINEAR) + myMesh->NbFaces(ORDER_LINEAR);
1345   
1346   if (NbAllEdgsAndFaces == NbQuadFacesAndEdgs) {
1347     //Quadratic mesh
1348     return SMESH_MesherHelper::QUADRATIC;
1349   }
1350   else if (NbAllEdgsAndFaces == NbFacesAndEdges) {
1351     //Linear mesh
1352     return SMESH_MesherHelper::LINEAR;
1353   }
1354   else
1355     //Mesh with both type of elements
1356     return SMESH_MesherHelper::COMP;
1357 }
1358
1359 //=======================================================================
1360 //function : GetOtherParam
1361 //purpose  : Return an alternative parameter for a node on seam
1362 //=======================================================================
1363
1364 double SMESH_MesherHelper::GetOtherParam(const double param) const
1365 {
1366   int i = myParIndex & U_periodic ? 0 : 1;
1367   return fabs(param-myPar1[i]) < fabs(param-myPar2[i]) ? myPar2[i] : myPar1[i];
1368 }
1369
1370 //=======================================================================
1371 namespace { // Structures used by FixQuadraticElements()
1372 //=======================================================================
1373
1374 #define __DMP__(txt) \
1375 //cout << txt
1376 #define MSG(txt) __DMP__(txt<<endl)
1377 #define MSGBEG(txt) __DMP__(txt)
1378
1379   const double straightTol2 = 1e-33; // to detect straing links
1380
1381   struct QFace;
1382   // ---------------------------------------
1383   /*!
1384    * \brief Quadratic link knowing its faces
1385    */
1386   struct QLink: public SMESH_TLink
1387   {
1388     const SMDS_MeshNode*          _mediumNode;
1389     mutable vector<const QFace* > _faces;
1390     mutable gp_Vec                _nodeMove;
1391     mutable int                   _nbMoves;
1392
1393     QLink(const SMDS_MeshNode* n1, const SMDS_MeshNode* n2, const SMDS_MeshNode* nm):
1394       SMESH_TLink( n1,n2 ), _mediumNode(nm), _nodeMove(0,0,0), _nbMoves(0) {
1395       _faces.reserve(4);
1396       //if ( MediumPos() != SMDS_TOP_3DSPACE )
1397         _nodeMove = MediumPnt() - MiddlePnt();
1398     }
1399     void SetContinuesFaces() const;
1400     const QFace* GetContinuesFace( const QFace* face ) const;
1401     bool OnBoundary() const;
1402     gp_XYZ MiddlePnt() const { return ( XYZ( node1() ) + XYZ( node2() )) / 2.; }
1403     gp_XYZ MediumPnt() const { return XYZ( _mediumNode ); }
1404
1405     SMDS_TypeOfPosition MediumPos() const
1406     { return _mediumNode->GetPosition()->GetTypeOfPosition(); }
1407     SMDS_TypeOfPosition EndPos(bool isSecond) const
1408     { return (isSecond ? node2() : node1())->GetPosition()->GetTypeOfPosition(); }
1409     const SMDS_MeshNode* EndPosNode(SMDS_TypeOfPosition pos) const
1410     { return EndPos(0) == pos ? node1() : EndPos(1) == pos ? node2() : 0; }
1411
1412     void Move(const gp_Vec& move, bool sum=false) const
1413     { _nodeMove += move; _nbMoves += sum ? (_nbMoves==0) : 1; }
1414     gp_XYZ Move() const { return _nodeMove.XYZ() / _nbMoves; }
1415     bool IsMoved() const { return (_nbMoves > 0 && !IsStraight()); }
1416     bool IsStraight() const { return _nodeMove.SquareMagnitude() <= straightTol2; }
1417
1418     bool operator<(const QLink& other) const {
1419       return (node1()->GetID() == other.node1()->GetID() ?
1420               node2()->GetID() < other.node2()->GetID() :
1421               node1()->GetID() < other.node1()->GetID());
1422     }
1423     struct PtrComparator {
1424       bool operator() (const QLink* l1, const QLink* l2 ) const { return *l1 < *l2; }
1425     };
1426   };
1427   // ---------------------------------------------------------
1428   /*!
1429    * \brief Link in the chain of links; it connects two faces
1430    */
1431   struct TChainLink
1432   {
1433     const QLink*         _qlink;
1434     mutable const QFace* _qfaces[2];
1435
1436     TChainLink(const QLink* qlink=0):_qlink(qlink) {
1437       _qfaces[0] = _qfaces[1] = 0;
1438     }
1439     void SetFace(const QFace* face) { int iF = _qfaces[0] ? 1 : 0; _qfaces[iF]=face; }
1440
1441     bool IsBoundary() const { return !_qfaces[1]; }
1442
1443     void RemoveFace( const QFace* face ) const
1444     { _qfaces[(face == _qfaces[1])] = 0; if (!_qfaces[0]) std::swap(_qfaces[0],_qfaces[1]); }
1445
1446     const QFace* NextFace( const QFace* f ) const
1447     { return _qfaces[0]==f ? _qfaces[1] : _qfaces[0]; }
1448
1449     const SMDS_MeshNode* NextNode( const SMDS_MeshNode* n ) const
1450     { return n == _qlink->node1() ? _qlink->node2() : _qlink->node1(); }
1451
1452     bool operator<(const TChainLink& other) const { return *_qlink < *other._qlink; }
1453
1454     operator bool() const { return (_qlink); }
1455
1456     const QLink* operator->() const { return _qlink; }
1457
1458     gp_Vec Normal() const;
1459   };
1460   // --------------------------------------------------------------------
1461   typedef list< TChainLink > TChain;
1462   typedef set < TChainLink > TLinkSet;
1463   typedef TLinkSet::const_iterator TLinkInSet;
1464
1465   const int theFirstStep = 5;
1466
1467   enum { ERR_OK, ERR_TRI, ERR_PRISM, ERR_UNKNOWN }; // errors of QFace::GetLinkChain()
1468   // --------------------------------------------------------------------
1469   /*!
1470    * \brief Face shared by two volumes and bound by QLinks
1471    */
1472   struct QFace: public TIDSortedElemSet
1473   {
1474     mutable const SMDS_MeshElement* _volumes[2];
1475     mutable vector< const QLink* >  _sides;
1476     mutable bool                    _sideIsAdded[4]; // added in chain of links
1477     gp_Vec                          _normal;
1478 #ifdef _DEBUG_
1479     mutable const SMDS_MeshElement* _face;
1480 #endif
1481
1482     QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face=0 );
1483
1484     void SetVolume(const SMDS_MeshElement* v) const { _volumes[ _volumes[0] ? 1 : 0 ] = v; }
1485
1486     int NbVolumes() const { return !_volumes[0] ? 0 : !_volumes[1] ? 1 : 2; }
1487
1488     void AddSelfToLinks() const {
1489       for ( int i = 0; i < _sides.size(); ++i )
1490         _sides[i]->_faces.push_back( this );
1491     }
1492     int LinkIndex( const QLink* side ) const {
1493       for (int i=0; i<_sides.size(); ++i ) if ( _sides[i] == side ) return i;
1494       return -1;
1495     }
1496     bool GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& err) const;
1497
1498     bool GetLinkChain( TChainLink& link, TChain& chain, SMDS_TypeOfPosition pos, int& err) const
1499     {
1500       int i = LinkIndex( link._qlink );
1501       if ( i < 0 ) return true;
1502       _sideIsAdded[i] = true;
1503       link.SetFace( this );
1504       // continue from opposite link
1505       return GetLinkChain( (i+2)%_sides.size(), chain, pos, err );
1506     }
1507     bool IsBoundary() const { return !_volumes[1]; }
1508
1509     bool Contains( const SMDS_MeshNode* node ) const { return count(node); }
1510
1511     TLinkInSet GetBoundaryLink( const TLinkSet&      links,
1512                                 const TChainLink&    avoidLink,
1513                                 TLinkInSet *         notBoundaryLink = 0,
1514                                 const SMDS_MeshNode* nodeToContain = 0,
1515                                 bool *               isAdjacentUsed = 0,
1516                                 int                  nbRecursionsLeft = -1) const;
1517
1518     TLinkInSet GetLinkByNode( const TLinkSet&      links,
1519                               const TChainLink&    avoidLink,
1520                               const SMDS_MeshNode* nodeToContain) const;
1521
1522     const SMDS_MeshNode* GetNodeInFace() const {
1523       for ( int iL = 0; iL < _sides.size(); ++iL )
1524         if ( _sides[iL]->MediumPos() == SMDS_TOP_FACE ) return _sides[iL]->_mediumNode;
1525       return 0;
1526     }
1527
1528     gp_Vec LinkNorm(const int i, SMESH_MesherHelper* theFaceHelper=0) const;
1529
1530     double MoveByBoundary( const TChainLink&   theLink,
1531                            const gp_Vec&       theRefVec,
1532                            const TLinkSet&     theLinks,
1533                            SMESH_MesherHelper* theFaceHelper=0,
1534                            const double        thePrevLen=0,
1535                            const int           theStep=theFirstStep,
1536                            gp_Vec*             theLinkNorm=0,
1537                            double              theSign=1.0) const;
1538   };
1539
1540   //================================================================================
1541   /*!
1542    * \brief Dump QLink and QFace
1543    */
1544   ostream& operator << (ostream& out, const QLink& l)
1545   {
1546     out <<"QLink nodes: "
1547         << l.node1()->GetID() << " - "
1548         << l._mediumNode->GetID() << " - "
1549         << l.node2()->GetID() << endl;
1550     return out;
1551   }
1552   ostream& operator << (ostream& out, const QFace& f)
1553   {
1554     out <<"QFace nodes: "/*<< &f << "  "*/;
1555     for ( TIDSortedElemSet::const_iterator n = f.begin(); n != f.end(); ++n )
1556       out << (*n)->GetID() << " ";
1557     out << " \tvolumes: "
1558         << (f._volumes[0] ? f._volumes[0]->GetID() : 0) << " "
1559         << (f._volumes[1] ? f._volumes[1]->GetID() : 0);
1560     out << "  \tNormal: "<< f._normal.X() <<", "<<f._normal.Y() <<", "<<f._normal.Z() << endl;
1561     return out;
1562   }
1563
1564   //================================================================================
1565   /*!
1566    * \brief Construct QFace from QLinks 
1567    */
1568   //================================================================================
1569
1570   QFace::QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face )
1571   {
1572     _volumes[0] = _volumes[1] = 0;
1573     _sides = links;
1574     _sideIsAdded[0]=_sideIsAdded[1]=_sideIsAdded[2]=_sideIsAdded[3]=false;
1575     _normal.SetCoord(0,0,0);
1576     for ( int i = 1; i < _sides.size(); ++i ) {
1577       const QLink *l1 = _sides[i-1], *l2 = _sides[i];
1578       insert( l1->node1() ); insert( l1->node2() );
1579       // compute normal
1580       gp_Vec v1( XYZ( l1->node2()), XYZ( l1->node1()));
1581       gp_Vec v2( XYZ( l2->node1()), XYZ( l2->node2()));
1582       if ( l1->node1() != l2->node1() && l1->node2() != l2->node2() )
1583         v1.Reverse(); 
1584       _normal += v1 ^ v2;
1585     }
1586     double normSqSize = _normal.SquareMagnitude();
1587     if ( normSqSize > numeric_limits<double>::min() )
1588       _normal /= sqrt( normSqSize );
1589     else
1590       _normal.SetCoord(1e-33,0,0);
1591
1592 #ifdef _DEBUG_
1593     _face = face;
1594 #endif
1595   }
1596   //================================================================================
1597   /*!
1598    * \brief Make up chain of links
1599    *  \param iSide - link to add first
1600    *  \param chain - chain to fill in
1601    *  \param pos   - postion of medium nodes the links should have
1602    *  \param error - out, specifies what is wrong
1603    *  \retval bool - false if valid chain can't be built; "valid" means that links
1604    *                 of the chain belongs to rectangles bounding hexahedrons
1605    */
1606   //================================================================================
1607
1608   bool QFace::GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
1609   {
1610     if ( iSide >= _sides.size() ) // wrong argument iSide
1611       return false;
1612     if ( _sideIsAdded[ iSide ]) // already in chain
1613       return true;
1614
1615     if ( _sides.size() != 4 ) { // triangle - visit all my continous faces
1616       MSGBEG( *this );
1617       list< const QFace* > faces( 1, this );
1618       for (list< const QFace* >::iterator fIt = faces.begin(); fIt != faces.end(); ++fIt ) {
1619         const QFace* face = *fIt;
1620         for ( int i = 0; i < face->_sides.size(); ++i ) {
1621           if ( !face->_sideIsAdded[i] && face->_sides[i] ) {
1622             face->_sideIsAdded[i] = true;
1623             TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(face->_sides[i]));
1624             chLink->SetFace( face );
1625             if ( face->_sides[i]->MediumPos() >= pos )
1626               if ( const QFace* contFace = face->_sides[i]->GetContinuesFace( face ))
1627                 faces.push_back( contFace );
1628           }
1629         }
1630       }
1631       if ( error < ERR_TRI )
1632         error = ERR_TRI;
1633       return false;
1634     }
1635     _sideIsAdded[iSide] = true; // not to add this link to chain again
1636     const QLink* link = _sides[iSide];
1637     if ( !link)
1638       return true;
1639
1640     // add link into chain
1641     TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(link));
1642     chLink->SetFace( this );
1643     MSGBEG( *this );
1644
1645     // propagate from quadrangle to neighbour faces
1646     if ( link->MediumPos() >= pos ) {
1647       int nbLinkFaces = link->_faces.size();
1648       if ( nbLinkFaces == 4 || nbLinkFaces < 4 && link->OnBoundary()) {
1649         // hexahedral mesh or boundary quadrangles - goto a continous face
1650         if ( const QFace* f = link->GetContinuesFace( this ))
1651           return f->GetLinkChain( *chLink, chain, pos, error );
1652       }
1653       else {
1654         TChainLink chLink(link); // side face of prismatic mesh - visit all faces of iSide
1655         for ( int i = 0; i < nbLinkFaces; ++i )
1656           if ( link->_faces[i] )
1657             link->_faces[i]->GetLinkChain( chLink, chain, pos, error );
1658         if ( error < ERR_PRISM )
1659           error = ERR_PRISM;
1660         return false;
1661       }
1662     }
1663     return true;
1664   }
1665
1666   //================================================================================
1667   /*!
1668    * \brief Return a boundary link of the triangle face
1669    *  \param links - set of all links
1670    *  \param avoidLink - link not to return
1671    *  \param notBoundaryLink - out, neither the returned link nor avoidLink
1672    *  \param nodeToContain - node the returned link must contain; if provided, search
1673    *                         also performed on adjacent faces
1674    *  \param isAdjacentUsed - returns true if link is found in adjacent faces
1675    *  \param nbRecursionsLeft - to limit recursion
1676    */
1677   //================================================================================
1678
1679   TLinkInSet QFace::GetBoundaryLink( const TLinkSet&      links,
1680                                      const TChainLink&    avoidLink,
1681                                      TLinkInSet *         notBoundaryLink,
1682                                      const SMDS_MeshNode* nodeToContain,
1683                                      bool *               isAdjacentUsed,
1684                                      int                  nbRecursionsLeft) const
1685   {
1686     TLinkInSet linksEnd = links.end(), boundaryLink = linksEnd;
1687
1688     typedef list< pair< const QFace*, TLinkInSet > > TFaceLinkList;
1689     TFaceLinkList adjacentFaces;
1690
1691     for ( int iL = 0; iL < _sides.size(); ++iL )
1692     {
1693       if ( avoidLink._qlink == _sides[iL] )
1694         continue;
1695       TLinkInSet link = links.find( _sides[iL] );
1696       if ( link == linksEnd ) continue;
1697       if ( (*link)->MediumPos() > SMDS_TOP_FACE )
1698         continue; // We work on faces here, don't go into a volume
1699
1700       // check link
1701       if ( link->IsBoundary() ) {
1702         if ( !nodeToContain ||
1703              (*link)->node1() == nodeToContain ||
1704              (*link)->node2() == nodeToContain )
1705         {
1706           boundaryLink = link;
1707           if ( !notBoundaryLink ) break;
1708         }
1709       }
1710       else if ( notBoundaryLink ) {
1711         *notBoundaryLink = link;
1712         if ( boundaryLink != linksEnd ) break;
1713       }
1714
1715       if ( boundaryLink == linksEnd && nodeToContain ) // collect adjacent faces
1716         if ( const QFace* adj = link->NextFace( this ))
1717           if ( adj->Contains( nodeToContain ))
1718             adjacentFaces.push_back( make_pair( adj, link ));
1719     }
1720
1721     if ( isAdjacentUsed ) *isAdjacentUsed = false;
1722     if ( boundaryLink == linksEnd && nodeToContain && nbRecursionsLeft) // check adjacent faces
1723     {
1724       if ( nbRecursionsLeft < 0 )
1725         nbRecursionsLeft = nodeToContain->NbInverseElements();
1726       TFaceLinkList::iterator adj = adjacentFaces.begin();
1727       for ( ; boundaryLink == linksEnd && adj != adjacentFaces.end(); ++adj )
1728         boundaryLink = adj->first->GetBoundaryLink( links, *(adj->second), 0, nodeToContain,
1729                                                     isAdjacentUsed, nbRecursionsLeft-1);
1730       if ( isAdjacentUsed ) *isAdjacentUsed = true;
1731     }
1732     return boundaryLink;
1733   }
1734   //================================================================================
1735   /*!
1736    * \brief Return a link ending at the given node but not avoidLink
1737    */
1738   //================================================================================
1739
1740   TLinkInSet QFace::GetLinkByNode( const TLinkSet&      links,
1741                                    const TChainLink&    avoidLink,
1742                                    const SMDS_MeshNode* nodeToContain) const
1743   {
1744     for ( int i = 0; i < _sides.size(); ++i )
1745       if ( avoidLink._qlink != _sides[i] &&
1746            (_sides[i]->node1() == nodeToContain || _sides[i]->node2() == nodeToContain ))
1747         return links.find( _sides[ i ]);
1748     return links.end();
1749   }
1750
1751   //================================================================================
1752   /*!
1753    * \brief Return normal to the i-th side pointing outside the face
1754    */
1755   //================================================================================
1756
1757   gp_Vec QFace::LinkNorm(const int i, SMESH_MesherHelper* /*uvHelper*/) const
1758   {
1759     gp_Vec norm, vecOut;
1760 //     if ( uvHelper ) {
1761 //       TopoDS_Face face = TopoDS::Face( uvHelper->GetSubShape());
1762 //       const SMDS_MeshNode* inFaceNode = uvHelper->GetNodeUVneedInFaceNode() ? GetNodeInFace() : 0;
1763 //       gp_XY uv1 = uvHelper->GetNodeUV( face, _sides[i]->node1(), inFaceNode );
1764 //       gp_XY uv2 = uvHelper->GetNodeUV( face, _sides[i]->node2(), inFaceNode );
1765 //       norm.SetCoord( uv1.Y() - uv2.Y(), uv2.X() - uv1.X(), 0 );
1766
1767 //       const QLink* otherLink = _sides[(i + 1) % _sides.size()];
1768 //       const SMDS_MeshNode* otherNode =
1769 //         otherLink->node1() == _sides[i]->node1() ? otherLink->node2() : otherLink->node1();
1770 //       gp_XY pIn = uvHelper->GetNodeUV( face, otherNode, inFaceNode );
1771 //       vecOut.SetCoord( uv1.X() - pIn.X(), uv1.Y() - pIn.Y(), 0 );
1772 //     }
1773 //     else {
1774       norm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
1775       gp_XYZ pIn = ( XYZ( _sides[0]->node1() ) +
1776                      XYZ( _sides[0]->node2() ) +
1777                      XYZ( _sides[1]->node1() )) / 3.;
1778       vecOut.SetXYZ( _sides[i]->MiddlePnt() - pIn );
1779       //}
1780     if ( norm * vecOut < 0 )
1781       norm.Reverse();
1782     double mag2 = norm.SquareMagnitude();
1783     if ( mag2 > numeric_limits<double>::min() )
1784       norm /= sqrt( mag2 );
1785     return norm;
1786   }
1787   //================================================================================
1788   /*!
1789    * \brief Move medium node of theLink according to its distance from boundary
1790    *  \param theLink - link to fix
1791    *  \param theRefVec - movement of boundary
1792    *  \param theLinks - all adjacent links of continous triangles
1793    *  \param theFaceHelper - helper is not used so far
1794    *  \param thePrevLen - distance from the boundary
1795    *  \param theStep - number of steps till movement propagation limit
1796    *  \param theLinkNorm - out normal to theLink
1797    *  \param theSign - 1 or -1 depending on movement of boundary
1798    *  \retval double - distance from boundary to propagation limit or other boundary
1799    */
1800   //================================================================================
1801
1802   double QFace::MoveByBoundary( const TChainLink&   theLink,
1803                                 const gp_Vec&       theRefVec,
1804                                 const TLinkSet&     theLinks,
1805                                 SMESH_MesherHelper* theFaceHelper,
1806                                 const double        thePrevLen,
1807                                 const int           theStep,
1808                                 gp_Vec*             theLinkNorm,
1809                                 double              theSign) const
1810   {
1811     if ( !theStep )
1812       return thePrevLen; // propagation limit reached
1813
1814     int iL; // index of theLink
1815     for ( iL = 0; iL < _sides.size(); ++iL )
1816       if ( theLink._qlink == _sides[ iL ])
1817         break;
1818
1819     MSG(string(theStep,'.')<<" Ref( "<<theRefVec.X()<<","<<theRefVec.Y()<<","<<theRefVec.Z()<<" )"
1820         <<" thePrevLen " << thePrevLen);
1821     MSG(string(theStep,'.')<<" "<<*theLink._qlink);
1822
1823     gp_Vec linkNorm = -LinkNorm( iL/*, theFaceHelper*/ ); // normal to theLink
1824     double refProj = theRefVec * linkNorm; // project movement vector to normal of theLink
1825     if ( theStep == theFirstStep )
1826       theSign = refProj < 0. ? -1. : 1.;
1827     else if ( theSign * refProj < 0.4 * theRefVec.Magnitude())
1828       return thePrevLen; // to propagate movement forward only, not in side dir or backward
1829
1830     int iL1 = (iL + 1) % 3, iL2 = (iL + 2) % 3; // indices of the two other links of triangle
1831     TLinkInSet link1 = theLinks.find( _sides[iL1] );
1832     TLinkInSet link2 = theLinks.find( _sides[iL2] );
1833     if ( link1 == theLinks.end() || link2 == theLinks.end() )
1834       return thePrevLen;
1835     const QFace* f1 = link1->NextFace( this ); // adjacent faces
1836     const QFace* f2 = link2->NextFace( this );
1837
1838     // propagate to adjacent faces till limit step or boundary
1839     double len1 = thePrevLen + (theLink->MiddlePnt() - _sides[iL1]->MiddlePnt()).Modulus();
1840     double len2 = thePrevLen + (theLink->MiddlePnt() - _sides[iL2]->MiddlePnt()).Modulus();
1841     gp_Vec linkDir1, linkDir2;
1842     try {
1843       OCC_CATCH_SIGNALS;
1844       if ( f1 )
1845         len1 = f1->MoveByBoundary
1846           ( *link1, theRefVec, theLinks, theFaceHelper, len1, theStep-1, &linkDir1, theSign);
1847       else
1848         linkDir1 = LinkNorm( iL1/*, theFaceHelper*/ );
1849     } catch (...) {
1850       MSG( " --------------- EXCEPTION");
1851       return thePrevLen;
1852     }
1853     try {
1854       OCC_CATCH_SIGNALS;
1855       if ( f2 )
1856         len2 = f2->MoveByBoundary
1857           ( *link2, theRefVec, theLinks, theFaceHelper, len2, theStep-1, &linkDir2, theSign);
1858       else
1859         linkDir2 = LinkNorm( iL2/*, theFaceHelper*/ );
1860     } catch (...) {
1861       MSG( " --------------- EXCEPTION");
1862       return thePrevLen;
1863     }
1864
1865     double fullLen = 0;
1866     if ( theStep != theFirstStep )
1867     {
1868       // choose chain length by direction of propagation most codirected with theRefVec
1869       bool choose1 = ( theRefVec * linkDir1 * theSign > theRefVec * linkDir2 * theSign );
1870       fullLen = choose1 ? len1 : len2;
1871       double r = thePrevLen / fullLen;
1872
1873       gp_Vec move = linkNorm * refProj * ( 1 - r );
1874       theLink->Move( move, true );
1875
1876       MSG(string(theStep,'.')<<" Move "<< theLink->_mediumNode->GetID()<<
1877           " by " << refProj * ( 1 - r ) << " following " <<
1878           (choose1 ? *link1->_qlink : *link2->_qlink));
1879
1880       if ( theLinkNorm ) *theLinkNorm = linkNorm;
1881     }
1882     return fullLen;
1883   }
1884
1885   //================================================================================
1886   /*!
1887    * \brief Find pairs of continues faces 
1888    */
1889   //================================================================================
1890
1891   void QLink::SetContinuesFaces() const
1892   {
1893     //       x0         x - QLink, [-|] - QFace, v - volume
1894     //   v0  |   v1   
1895     //       |          Between _faces of link x2 two vertical faces are continues
1896     // x1----x2-----x3  and two horizontal faces are continues. We set vertical faces
1897     //       |          to _faces[0] and _faces[1] and horizontal faces to
1898     //   v2  |   v3     _faces[2] and _faces[3] (or vise versa).
1899     //       x4
1900
1901     if ( _faces.empty() )
1902       return;
1903     int iFaceCont = -1;
1904     for ( int iF = 1; iFaceCont < 0 && iF < _faces.size(); ++iF )
1905     {
1906       // look for a face bounding none of volumes bound by _faces[0]
1907       bool sameVol = false;
1908       int nbVol = _faces[iF]->NbVolumes();
1909       for ( int iV = 0; !sameVol && iV < nbVol; ++iV )
1910         sameVol = ( _faces[iF]->_volumes[iV] == _faces[0]->_volumes[0] ||
1911                     _faces[iF]->_volumes[iV] == _faces[0]->_volumes[1]);
1912       if ( !sameVol )
1913         iFaceCont = iF;
1914     }
1915     if ( iFaceCont > 0 ) // continues faces found, set one by the other
1916     {
1917       if ( iFaceCont != 1 )
1918         std::swap( _faces[1], _faces[iFaceCont] );
1919     }
1920     else if ( _faces.size() > 1 ) // not found, set NULL by the first face
1921     {
1922       _faces.insert( ++_faces.begin(), 0 );
1923     }
1924   }
1925   //================================================================================
1926   /*!
1927    * \brief Return a face continues to the given one
1928    */
1929   //================================================================================
1930
1931   const QFace* QLink::GetContinuesFace( const QFace* face ) const
1932   {
1933     for ( int i = 0; i < _faces.size(); ++i ) {
1934       if ( _faces[i] == face ) {
1935         int iF = i < 2 ? 1-i : 5-i;
1936         return iF < _faces.size() ? _faces[iF] : 0;
1937       }
1938     }
1939     return 0;
1940   }
1941   //================================================================================
1942   /*!
1943    * \brief True if link is on mesh boundary
1944    */
1945   //================================================================================
1946
1947   bool QLink::OnBoundary() const
1948   {
1949     for ( int i = 0; i < _faces.size(); ++i )
1950       if (_faces[i] && _faces[i]->IsBoundary()) return true;
1951     return false;
1952   }
1953   //================================================================================
1954   /*!
1955    * \brief Return normal of link of the chain
1956    */
1957   //================================================================================
1958
1959   gp_Vec TChainLink::Normal() const {
1960     gp_Vec norm;
1961     if (_qfaces[0]) norm  = _qfaces[0]->_normal;
1962     if (_qfaces[1]) norm += _qfaces[1]->_normal;
1963     return norm;
1964   }
1965   //================================================================================
1966   /*!
1967    * \brief Move medium nodes of vertical links of pentahedrons adjacent by side faces
1968    */
1969   //================================================================================
1970
1971   void fixPrism( TChain& allLinks )
1972   {
1973     // separate boundary links from internal ones
1974     typedef set<const QLink*/*, QLink::PtrComparator*/> QLinkSet;
1975     QLinkSet interLinks, bndLinks1, bndLink2;
1976
1977     bool isCurved = false;
1978     for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
1979       if ( (*lnk)->OnBoundary() )
1980         bndLinks1.insert( lnk->_qlink );
1981       else
1982         interLinks.insert( lnk->_qlink );
1983       isCurved = isCurved || !(*lnk)->IsStraight();
1984     }
1985     if ( !isCurved )
1986       return; // no need to move
1987
1988     QLinkSet *curBndLinks = &bndLinks1, *newBndLinks = &bndLink2;
1989
1990     while ( !interLinks.empty() && !curBndLinks->empty() )
1991     {
1992       // propagate movement from boundary links to connected internal links
1993       QLinkSet::iterator bnd = curBndLinks->begin(), bndEnd = curBndLinks->end();
1994       for ( ; bnd != bndEnd; ++bnd )
1995       {
1996         const QLink* bndLink = *bnd;
1997         for ( int i = 0; i < bndLink->_faces.size(); ++i ) // loop on faces of bndLink
1998         {
1999           const QFace* face = bndLink->_faces[i]; // quadrange lateral face of a prism
2000           if ( !face ) continue;
2001           // find and move internal link opposite to bndLink within the face
2002           int interInd = ( face->LinkIndex( bndLink ) + 2 ) % face->_sides.size();
2003           const QLink* interLink = face->_sides[ interInd ];
2004           QLinkSet::iterator pInterLink = interLinks.find( interLink );
2005           if ( pInterLink == interLinks.end() ) continue; // not internal link
2006           interLink->Move( bndLink->_nodeMove );
2007           // treated internal links become new boundary ones
2008           interLinks. erase( pInterLink );
2009           newBndLinks->insert( interLink );
2010         }
2011       }
2012       curBndLinks->clear();
2013       std::swap( curBndLinks, newBndLinks );
2014     }
2015   }
2016
2017   //================================================================================
2018   /*!
2019    * \brief Fix links of continues triangles near curved boundary
2020    */
2021   //================================================================================
2022
2023   void fixTriaNearBoundary( TChain & allLinks, SMESH_MesherHelper& /*helper*/)
2024   {
2025     if ( allLinks.empty() ) return;
2026
2027     TLinkSet linkSet( allLinks.begin(), allLinks.end());
2028     TLinkInSet linkIt = linkSet.begin(), linksEnd = linkSet.end();
2029
2030     // move in 2d if we are on geom face
2031 //     TopoDS_Face face;
2032 //     TopLoc_Location loc;
2033 //     SMESH_MesherHelper faceHelper( *helper.GetMesh());
2034 //     while ( linkIt->IsBoundary()) ++linkIt;
2035 //     if ( linkIt == linksEnd ) return;
2036 //     if ( (*linkIt)->MediumPos() == SMDS_TOP_FACE ) {
2037 //       bool checkPos = true;
2038 //       TopoDS_Shape f = helper.GetSubShapeByNode( (*linkIt)->_mediumNode, helper.GetMeshDS() );
2039 //       if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) {
2040 //         face = TopoDS::Face( f );
2041 //         helper.GetNodeUV( face, (*linkIt)->_mediumNode, 0, &checkPos);
2042 //         if (checkPos)
2043 //           face.Nullify();
2044 //         else
2045 //           faceHelper.SetSubShape( face );
2046 //       }
2047 //     }
2048     for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt)
2049     {
2050       if ( linkIt->IsBoundary() && !(*linkIt)->IsStraight() && linkIt->_qfaces[0])
2051       {
2052 //         if ( !face.IsNull() ) {
2053 //           const SMDS_MeshNode* inFaceNode =
2054 //             faceHelper.GetNodeUVneedInFaceNode() ? linkIt->_qfaces[0]->GetNodeInFace() : 0;
2055 //           gp_XY uvm = helper.GetNodeUV( face, (*linkIt)->_mediumNode, inFaceNode );
2056 //           gp_XY uv1 = helper.GetNodeUV( face, (*linkIt)->node1(), inFaceNode);
2057 //           gp_XY uv2 = helper.GetNodeUV( face, (*linkIt)->node2(), inFaceNode);
2058 //           gp_XY uvMove = uvm - helper.GetMiddleUV( BRep_Tool::Surface(face,loc), uv1, uv2);
2059 //           gp_Vec move( uvMove.X(), uvMove.Y(), 0 );
2060 //           linkIt->_qfaces[0]->MoveByBoundary( *linkIt, move, linkSet, &faceHelper );
2061 //         }
2062 //         else {
2063           linkIt->_qfaces[0]->MoveByBoundary( *linkIt, (*linkIt)->_nodeMove, linkSet );
2064           //}
2065       }
2066     }
2067   }
2068
2069   //================================================================================
2070   /*!
2071    * \brief Detect rectangular structure of links and build chains from them
2072    */
2073   //================================================================================
2074
2075   enum TSplitTriaResult {
2076     _OK, _NO_CORNERS, _FEW_ROWS, _MANY_ROWS, _NO_SIDELINK, _BAD_MIDQUAD, _NOT_RECT,
2077     _NO_MIDQUAD, _NO_UPTRIA, _BAD_SET_SIZE, _BAD_CORNER, _BAD_START, _NO_BOTLINK };
2078
2079   TSplitTriaResult splitTrianglesIntoChains( TChain &            allLinks,
2080                                              vector< TChain> &   resultChains,
2081                                              SMDS_TypeOfPosition pos )
2082   {
2083     // put links in the set and evalute number of result chains by number of boundary links
2084     TLinkSet linkSet;
2085     int nbBndLinks = 0;
2086     for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2087       linkSet.insert( *lnk );
2088       nbBndLinks += lnk->IsBoundary();
2089     }
2090     resultChains.clear();
2091     resultChains.reserve( nbBndLinks / 2 );
2092
2093     TLinkInSet linkIt, linksEnd = linkSet.end();
2094
2095     // find a boundary link with corner node; corner node has position pos-2
2096     // i.e. SMDS_TOP_VERTEX for links on faces and SMDS_TOP_EDGE for
2097     // links in volume
2098     SMDS_TypeOfPosition cornerPos = SMDS_TypeOfPosition(pos-2);
2099     const SMDS_MeshNode* corner = 0;
2100     for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt )
2101       if ( linkIt->IsBoundary() && (corner = (*linkIt)->EndPosNode(cornerPos)))
2102         break;
2103     if ( !corner)
2104       return _NO_CORNERS;
2105
2106     TLinkInSet           startLink = linkIt;
2107     const SMDS_MeshNode* startCorner = corner;
2108     vector< TChain* >    rowChains;
2109     int iCol = 0;
2110
2111     while ( startLink != linksEnd) // loop on columns
2112     {
2113       // We suppose we have a rectangular structure like shown here. We have found a
2114       //               corner of the rectangle (startCorner) and a boundary link sharing  
2115       //    |/  |/  |  the startCorner (startLink). We are going to loop on rows of the   
2116       //  --o---o---o  structure making several chains at once. One chain (columnChain)   
2117       //    |\  |  /|  starts at startLink and continues upward (we look at the structure 
2118       //  \ | \ | / |  from such point that startLink is on the bottom of the structure). 
2119       //   \|  \|/  |  While going upward we also fill horizontal chains (rowChains) we   
2120       //  --o---o---o  encounter.                                                         
2121       //   /|\  |\  |
2122       //  / | \ | \ |  startCorner
2123       //    |  \|  \|,'
2124       //  --o---o---o
2125       //          `.startLink
2126
2127       if ( resultChains.size() == nbBndLinks / 2 )
2128         return _NOT_RECT;
2129       resultChains.push_back( TChain() );
2130       TChain& columnChain = resultChains.back();
2131
2132       TLinkInSet botLink = startLink; // current horizontal link to go up from
2133       corner = startCorner; // current corner the botLink ends at
2134       int iRow = 0;
2135       while ( botLink != linksEnd ) // loop on rows
2136       {
2137         // add botLink to the columnChain
2138         columnChain.push_back( *botLink );
2139
2140         const QFace* botTria = botLink->_qfaces[0]; // bottom triangle bound by botLink
2141         if ( !botTria )
2142         { // the column ends
2143           linkSet.erase( botLink );
2144           if ( iRow != rowChains.size() )
2145             return _FEW_ROWS; // different nb of rows in columns
2146           break;
2147         }
2148         // find the link dividing the quadrangle (midQuadLink) and vertical boundary
2149         // link ending at <corner> (sideLink); there are two cases:
2150         // 1) midQuadLink does not end at <corner>, then we easily find it by botTria,
2151         //   since midQuadLink is not at boundary while sideLink is.
2152         // 2) midQuadLink ends at <corner>
2153         bool isCase2;
2154         TLinkInSet midQuadLink = linksEnd;
2155         TLinkInSet sideLink = botTria->GetBoundaryLink( linkSet, *botLink, &midQuadLink,
2156                                                         corner, &isCase2 );
2157         if ( isCase2 ) { // find midQuadLink among links of botTria
2158           midQuadLink = botTria->GetLinkByNode( linkSet, *botLink, corner );
2159           if ( midQuadLink->IsBoundary() )
2160             return _BAD_MIDQUAD;
2161         }
2162         if ( sideLink == linksEnd || midQuadLink == linksEnd || sideLink == midQuadLink )
2163           return sideLink == linksEnd ? _NO_SIDELINK : _NO_MIDQUAD;
2164
2165         // fill chains
2166         columnChain.push_back( *midQuadLink );
2167         if ( iRow >= rowChains.size() ) {
2168           if ( iCol > 0 )
2169             return _MANY_ROWS; // different nb of rows in columns
2170           if ( resultChains.size() == nbBndLinks / 2 )
2171             return _NOT_RECT;
2172           resultChains.push_back( TChain() );
2173           rowChains.push_back( & resultChains.back() );
2174         }
2175         rowChains[iRow]->push_back( *sideLink );
2176         rowChains[iRow]->push_back( *midQuadLink );
2177
2178         const QFace* upTria = midQuadLink->NextFace( botTria ); // upper tria of the rectangle
2179         if ( !upTria)
2180           return _NO_UPTRIA;
2181         if ( iRow == 0 ) {
2182           // prepare startCorner and startLink for the next column
2183           startCorner = startLink->NextNode( startCorner );
2184           if (isCase2)
2185             startLink = botTria->GetBoundaryLink( linkSet, *botLink, 0, startCorner );
2186           else
2187             startLink = upTria->GetBoundaryLink( linkSet, *midQuadLink, 0, startCorner );
2188           // check if no more columns remains
2189           if ( startLink != linksEnd ) {
2190             const SMDS_MeshNode* botNode = startLink->NextNode( startCorner );
2191             if ( (isCase2 ? botTria : upTria)->Contains( botNode ))
2192               startLink = linksEnd; // startLink bounds upTria or botTria
2193             else if ( startLink == botLink || startLink == midQuadLink || startLink == sideLink )
2194               return _BAD_START;
2195           }
2196         }
2197         // find bottom link and corner for the next row
2198         corner = sideLink->NextNode( corner );
2199         // next bottom link ends at the new corner
2200         linkSet.erase( botLink );
2201         botLink = upTria->GetLinkByNode( linkSet, (isCase2 ? *sideLink : *midQuadLink), corner );
2202         if ( botLink == linksEnd || botLink == (isCase2 ? midQuadLink : sideLink))
2203           return _NO_BOTLINK;
2204         linkSet.erase( midQuadLink );
2205         linkSet.erase( sideLink );
2206
2207         // make faces neighboring the found ones be boundary
2208         if ( startLink != linksEnd ) {
2209           const QFace* tria = isCase2 ? botTria : upTria;
2210           for ( int iL = 0; iL < 3; ++iL ) {
2211             linkIt = linkSet.find( tria->_sides[iL] );
2212             if ( linkIt != linksEnd )
2213               linkIt->RemoveFace( tria );
2214           }
2215         }
2216         if ( botLink->_qfaces[0] == upTria || botLink->_qfaces[1] == upTria )
2217           botLink->RemoveFace( upTria ); // make next botTria first in vector
2218
2219         iRow++;
2220       } // loop on rows
2221
2222       iCol++;
2223     }
2224     // In the linkSet, there must remain the last links of rowChains; add them
2225     if ( linkSet.size() != rowChains.size() )
2226       return _BAD_SET_SIZE;
2227     for ( int iRow = 0; iRow < rowChains.size(); ++iRow ) {
2228       // find the link (startLink) ending at startCorner
2229       corner = 0;
2230       for ( startLink = linkSet.begin(); startLink != linksEnd; ++startLink ) {
2231         if ( (*startLink)->node1() == startCorner ) {
2232           corner = (*startLink)->node2(); break;
2233         }
2234         else if ( (*startLink)->node2() == startCorner) {
2235           corner = (*startLink)->node1(); break;
2236         }
2237       }
2238       if ( startLink == linksEnd )
2239         return _BAD_CORNER;
2240       rowChains[ iRow ]->push_back( *startLink );
2241       linkSet.erase( startLink );
2242       startCorner = corner;
2243     }
2244
2245     return _OK;
2246   }
2247 }
2248
2249 //=======================================================================
2250 /*!
2251  * \brief Move medium nodes of faces and volumes to fix distorted elements
2252  * \param volumeOnly - to fix nodes on faces or not, if the shape is solid
2253  * 
2254  * Issue 0020307: EDF 992 SMESH : Linea/Quadratic with Medium Node on Geometry
2255  */
2256 //=======================================================================
2257
2258 void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
2259 {
2260   // apply algorithm to solids or geom faces
2261   // ----------------------------------------------
2262   if ( myShape.IsNull() ) {
2263     if ( !myMesh->HasShapeToMesh() ) return;
2264     SetSubShape( myMesh->GetShapeToMesh() );
2265
2266     TopTools_MapOfShape faces; // faces not in solid or in not meshed solid
2267     for ( TopExp_Explorer f(myShape,TopAbs_FACE,TopAbs_SOLID); f.More(); f.Next() ) {
2268       faces.Add( f.Current() );
2269     }
2270     for ( TopExp_Explorer s(myShape,TopAbs_SOLID); s.More(); s.Next() ) {
2271       if ( myMesh->GetSubMesh( s.Current() )->IsEmpty() ) { // get faces of solid
2272         for ( TopExp_Explorer f( s.Current(), TopAbs_FACE); f.More(); f.Next() )
2273           faces.Add( f.Current() );
2274       }
2275       else { // fix nodes in the solid and its faces
2276         SMESH_MesherHelper h(*myMesh);
2277         h.SetSubShape( s.Current() );
2278         h.FixQuadraticElements(false);
2279       }
2280     }
2281     // fix nodes on geom faces
2282     for ( TopTools_MapIteratorOfMapOfShape fIt( faces ); fIt.More(); fIt.Next() ) {
2283       SMESH_MesherHelper h(*myMesh);
2284       h.SetSubShape( fIt.Key() );
2285       h.FixQuadraticElements(true);
2286     }
2287     return;
2288   }
2289
2290   // Find out type of elements and get iterator on them
2291   // ---------------------------------------------------
2292
2293   SMDS_ElemIteratorPtr elemIt;
2294   SMDSAbs_ElementType elemType = SMDSAbs_All;
2295
2296   SMESH_subMesh* submesh = myMesh->GetSubMeshContaining( myShapeID );
2297   if ( !submesh )
2298     return;
2299   if ( SMESHDS_SubMesh* smDS = submesh->GetSubMeshDS() ) {
2300     elemIt = smDS->GetElements();
2301     if ( elemIt->more() ) {
2302       elemType = elemIt->next()->GetType();
2303       elemIt = smDS->GetElements();
2304     }
2305   }
2306   if ( !elemIt || !elemIt->more() || elemType < SMDSAbs_Face )
2307     return;
2308
2309   // Fill in auxiliary data structures
2310   // ----------------------------------
2311
2312   set< QLink > links;
2313   set< QFace > faces;
2314   set< QLink >::iterator pLink;
2315   set< QFace >::iterator pFace;
2316
2317   bool isCurved = false;
2318   bool hasRectFaces = false;
2319   set<int> nbElemNodeSet;
2320
2321   if ( elemType == SMDSAbs_Volume )
2322   {
2323     SMDS_VolumeTool volTool;
2324     while ( elemIt->more() ) // loop on volumes
2325     {
2326       const SMDS_MeshElement* vol = elemIt->next();
2327       if ( !vol->IsQuadratic() || !volTool.Set( vol ))
2328         return; //continue;
2329       for ( int iF = 0; iF < volTool.NbFaces(); ++iF ) // loop on faces of volume
2330       {
2331         int nbN = volTool.NbFaceNodes( iF );
2332         nbElemNodeSet.insert( nbN );
2333         const SMDS_MeshNode** faceNodes = volTool.GetFaceNodes( iF );
2334         vector< const QLink* > faceLinks( nbN/2 );
2335         for ( int iN = 0; iN < nbN; iN += 2 ) // loop on links of a face
2336         {
2337           // store QLink
2338           QLink link( faceNodes[iN], faceNodes[iN+2], faceNodes[iN+1] );
2339           pLink = links.insert( link ).first;
2340           faceLinks[ iN/2 ] = & *pLink;
2341           if ( !isCurved )
2342             isCurved = !link.IsStraight();
2343           if ( link.MediumPos() == SMDS_TOP_3DSPACE && !link.IsStraight() )
2344             return; // already fixed
2345         }
2346         // store QFace
2347         pFace = faces.insert( QFace( faceLinks )).first;
2348         if ( pFace->NbVolumes() == 0 )
2349           pFace->AddSelfToLinks();
2350         pFace->SetVolume( vol );
2351         hasRectFaces = hasRectFaces ||
2352           ( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
2353             volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
2354 #ifdef _DEBUG_
2355         if ( nbN == 6 )
2356           pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],faceNodes[4]);
2357         else
2358           pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],
2359                                                faceNodes[4],faceNodes[6] );
2360 #endif
2361       }
2362     }
2363     set< QLink >::iterator pLink = links.begin();
2364     for ( ; pLink != links.end(); ++pLink )
2365       pLink->SetContinuesFaces();
2366   }
2367   else
2368   {
2369     while ( elemIt->more() ) // loop on faces
2370     {
2371       const SMDS_MeshElement* face = elemIt->next();
2372       if ( !face->IsQuadratic() )
2373         continue;
2374       nbElemNodeSet.insert( face->NbNodes() );
2375       int nbN = face->NbNodes()/2;
2376       vector< const QLink* > faceLinks( nbN );
2377       for ( int iN = 0; iN < nbN; ++iN ) // loop on links of a face
2378       {
2379         // store QLink
2380         QLink link( face->GetNode(iN), face->GetNode((iN+1)%nbN), face->GetNode(iN+nbN) );
2381         pLink = links.insert( link ).first;
2382         faceLinks[ iN ] = & *pLink;
2383         if ( !isCurved )
2384           isCurved = !link.IsStraight();
2385       }
2386       // store QFace
2387       pFace = faces.insert( QFace( faceLinks )).first;
2388       pFace->AddSelfToLinks();
2389       hasRectFaces = ( hasRectFaces || nbN == 4 );
2390     }
2391   }
2392   if ( !isCurved )
2393     return; // no curved edges of faces
2394
2395   // Compute displacement of medium nodes
2396   // -------------------------------------
2397
2398   // two loops on faces: the first is to treat boundary links, the second is for internal ones
2399   TopLoc_Location loc;
2400   // not treat boundary of volumic submesh
2401   int isInside = ( elemType == SMDSAbs_Volume && volumeOnly ) ? 1 : 0;
2402   for ( ; isInside < 2; ++isInside ) {
2403     MSG( "--------------- LOOP (inside=" << isInside << ") ------------------");
2404     SMDS_TypeOfPosition pos = isInside ? SMDS_TOP_3DSPACE : SMDS_TOP_FACE;
2405
2406     for ( pFace = faces.begin(); pFace != faces.end(); ++pFace ) {
2407       if ( bool(isInside) == pFace->IsBoundary() )
2408         continue;
2409       for ( int dir = 0; dir < 2; ++dir ) // 2 directions of propagation from quadrangle
2410       {
2411         MSG( "CHAIN");
2412         // make chain of links connected via continues faces
2413         int error = ERR_OK;
2414         TChain rawChain;
2415         if ( !pFace->GetLinkChain( dir, rawChain, pos, error) && error ==ERR_UNKNOWN ) continue;
2416         rawChain.reverse();
2417         if ( !pFace->GetLinkChain( dir+2, rawChain, pos, error ) && error ==ERR_UNKNOWN ) continue;
2418
2419         vector< TChain > chains;
2420         if ( error == ERR_OK ) { // chain contains continues rectangles
2421           chains.resize(1);
2422           chains[0].splice( chains[0].begin(), rawChain );
2423         }
2424         else if ( error == ERR_TRI ) {  // chain contains continues triangles
2425           TSplitTriaResult res = splitTrianglesIntoChains( rawChain, chains, pos );
2426           if ( res != _OK ) { // not quadrangles split into triangles
2427             fixTriaNearBoundary( rawChain, *this );
2428             break;
2429           }
2430         }
2431         else if ( error == ERR_PRISM ) { // quadrangle side faces of prisms
2432           fixPrism( rawChain );
2433           break;
2434         }
2435         else {
2436           continue;
2437         }
2438         for ( int iC = 0; iC < chains.size(); ++iC )
2439         {
2440           TChain& chain = chains[iC];
2441           if ( chain.empty() ) continue;
2442           if ( chain.front()->IsStraight() && chain.back()->IsStraight() ) {
2443             MSG("3D straight");
2444             continue;
2445           }
2446           // mesure chain length and compute link position along the chain
2447           double chainLen = 0;
2448           vector< double > linkPos;
2449           MSGBEG( "Link medium nodes: ");
2450           TChain::iterator link0 = chain.begin(), link1 = chain.begin(), link2;
2451           for ( ++link1; link1 != chain.end(); ++link1, ++link0 ) {
2452             MSGBEG( (*link0)->_mediumNode->GetID() << "-" <<(*link1)->_mediumNode->GetID()<<" ");
2453             double len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2454             while ( len < numeric_limits<double>::min() ) { // remove degenerated link
2455               link1 = chain.erase( link1 );
2456               if ( link1 == chain.end() )
2457                 break;
2458               len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2459             }
2460             chainLen += len;
2461             linkPos.push_back( chainLen );
2462           }
2463           MSG("");
2464           if ( linkPos.size() < 2 )
2465             continue;
2466
2467           gp_Vec move0 = chain.front()->_nodeMove;
2468           gp_Vec move1 = chain.back ()->_nodeMove;
2469
2470           TopoDS_Face face;
2471           bool checkUV = true;
2472           if ( !isInside ) {
2473             // compute node displacement of end links in parametric space of face
2474             const SMDS_MeshNode* nodeOnFace = (*(++chain.begin()))->_mediumNode;
2475             TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() );
2476             if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE )
2477             {
2478               face = TopoDS::Face( f );
2479               Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
2480               for ( int is1 = 0; is1 < 2; ++is1 ) // move0 or move1
2481               {
2482                 TChainLink& link = is1 ? chain.back() : chain.front();
2483                 gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV);
2484                 gp_XY uv1 = GetNodeUV( face, link->node1(), nodeOnFace, &checkUV);
2485                 gp_XY uv2 = GetNodeUV( face, link->node2(), nodeOnFace, &checkUV);
2486                 gp_XY uv12 = GetMiddleUV( surf, uv1, uv2);
2487                 // uvMove = uvm - uv12
2488                 gp_XY uvMove = applyIn2D(surf, uvm, uv12, gp_XY_Subtracted, /*inPeriod=*/false);
2489                 ( is1 ? move1 : move0 ).SetCoord( uvMove.X(), uvMove.Y(), 0 );
2490               }
2491               if ( move0.SquareMagnitude() < straightTol2 &&
2492                    move1.SquareMagnitude() < straightTol2 ) {
2493                 MSG("2D straight");
2494                 continue; // straight - no need to move nodes of internal links
2495               }
2496             }
2497           }
2498           gp_Trsf trsf;
2499           if ( isInside || face.IsNull() )
2500           {
2501             // compute node displacement of end links in their local coord systems
2502             {
2503               TChainLink& ln0 = chain.front(), ln1 = *(++chain.begin());
2504               trsf.SetTransformation( gp_Ax3( gp::Origin(), ln0.Normal(),
2505                                               gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2506               move0.Transform(trsf);
2507             }
2508             {
2509               TChainLink& ln0 = *(++chain.rbegin()), ln1 = chain.back();
2510               trsf.SetTransformation( gp_Ax3( gp::Origin(), ln1.Normal(),
2511                                               gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2512               move1.Transform(trsf);
2513             }
2514           }
2515           // compute displacement of medium nodes
2516           link2 = chain.begin();
2517           link0 = link2++;
2518           link1 = link2++;
2519           for ( int i = 0; link2 != chain.end(); ++link0, ++link1, ++link2, ++i )
2520           {
2521             double r = linkPos[i] / chainLen;
2522             // displacement in local coord system
2523             gp_Vec move = (1. - r) * move0 + r * move1;
2524             if ( isInside || face.IsNull()) {
2525               // transform to global
2526               gp_Vec x01( (*link0)->MiddlePnt(), (*link1)->MiddlePnt() );
2527               gp_Vec x12( (*link1)->MiddlePnt(), (*link2)->MiddlePnt() );
2528               gp_Vec x = x01.Normalized() + x12.Normalized();
2529               trsf.SetTransformation( gp_Ax3( gp::Origin(), link1->Normal(), x), gp_Ax3() );
2530               move.Transform(trsf);
2531             }
2532             else {
2533               // compute 3D displacement by 2D one
2534               Handle(Geom_Surface) s = BRep_Tool::Surface(face,loc);
2535               gp_XY oldUV   = GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV);
2536               gp_XY newUV   = applyIn2D( s, oldUV, gp_XY( move.X(),move.Y()), gp_XY_Added);
2537               gp_Pnt newPnt = s->Value( newUV.X(), newUV.Y());
2538               move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) );
2539 #ifdef _DEBUG_
2540               if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() <
2541                    move.SquareMagnitude())
2542               {
2543                 gp_XY uv0 = GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV);
2544                 gp_XY uv2 = GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV);
2545                 MSG( "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
2546                      "uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" <<
2547                      "uvOld: "<<oldUV.X()<<", "<<oldUV.Y()<<" \t" <<
2548                      "newUV: "<<newUV.X()<<", "<<newUV.Y()<<" \t");
2549               }
2550 #endif
2551             }
2552             (*link1)->Move( move );
2553             MSG( "Move " << (*link1)->_mediumNode->GetID() << " following "
2554                  << chain.front()->_mediumNode->GetID() <<"-"
2555                  << chain.back ()->_mediumNode->GetID() <<
2556                  " by " << move.Magnitude());
2557           }
2558         } // loop on chains of links
2559       } // loop on 2 directions of propagation from quadrangle
2560     } // loop on faces
2561   }
2562
2563   // Move nodes
2564   // -----------
2565
2566   for ( pLink = links.begin(); pLink != links.end(); ++pLink ) {
2567     if ( pLink->IsMoved() ) {
2568       //gp_Pnt p = pLink->MediumPnt() + pLink->Move();
2569       gp_Pnt p = pLink->MiddlePnt() + pLink->Move();
2570       GetMeshDS()->MoveNode( pLink->_mediumNode, p.X(), p.Y(), p.Z());
2571     }
2572   }
2573 }
2574
2575 //=======================================================================
2576 /*!
2577  * \brief Iterator on ancestors of the given type
2578  */
2579 //=======================================================================
2580
2581 struct TAncestorsIterator : public SMDS_Iterator<const TopoDS_Shape*>
2582 {
2583   TopTools_ListIteratorOfListOfShape _ancIter;
2584   TopAbs_ShapeEnum                   _type;
2585   TAncestorsIterator( const TopTools_ListOfShape& ancestors, TopAbs_ShapeEnum type)
2586     : _ancIter( ancestors ), _type( type )
2587   {
2588     if ( _ancIter.More() && _ancIter.Value().ShapeType() != _type ) next();
2589   }
2590   virtual bool more()
2591   {
2592     return _ancIter.More();
2593   }
2594   virtual const TopoDS_Shape* next()
2595   {
2596     const TopoDS_Shape* s = _ancIter.More() ? & _ancIter.Value() : 0;
2597     if ( _ancIter.More() )
2598       for ( _ancIter.Next();  _ancIter.More(); _ancIter.Next())
2599         if ( _ancIter.Value().ShapeType() == _type )
2600           break;
2601     return s;
2602   }
2603 };
2604
2605 //=======================================================================
2606 /*!
2607  * \brief Return iterator on ancestors of the given type
2608  */
2609 //=======================================================================
2610
2611 PShapeIteratorPtr SMESH_MesherHelper::GetAncestors(const TopoDS_Shape& shape,
2612                                                    const SMESH_Mesh&   mesh,
2613                                                    TopAbs_ShapeEnum    ancestorType)
2614 {
2615   return PShapeIteratorPtr( new TAncestorsIterator( mesh.GetAncestors(shape), ancestorType));
2616 }