]> SALOME platform Git repositories - modules/smesh.git/blob - src/SMESH/SMESH_MesherHelper.cxx
Salome HOME
Rever incorrect integration!!! from 13.09.2009 into the V5_1_main
[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().get());
357     uv.SetCoord(fpos->GetUParameter(),fpos->GetVParameter());
358     if ( check )
359       uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 2*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().get());
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(), 2*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, check ));
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       uv.SetCoord( U,V );
495       if ( nodePnt.Distance( surface->Value( U, V )) > tol )
496       {
497         MESSAGE( "SMESH_MesherHelper::CheckNodeUV(), invalid projection" );
498         return false;
499       }
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 = static_cast<const SMDS_EdgePosition*>( pos.get() );
608     param =  epos->GetUParameter();
609   }
610   else if( pos->GetTypeOfPosition() == SMDS_TOP_VERTEX )
611   {
612     if ( inEdgeNode && TopExp::FirstVertex( E ).IsSame( TopExp::LastVertex( E ))) // issue 0020128
613     {
614       Standard_Real f,l;
615       BRep_Tool::Range( E, f,l );
616       double uInEdge = GetNodeU( E, inEdgeNode );
617       param = ( fabs( uInEdge - f ) < fabs( l - uInEdge )) ? f : l;
618     }
619     else
620     {
621       SMESHDS_Mesh * meshDS = GetMeshDS();
622       int vertexID = pos->GetShapeId();
623       const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
624       param =  BRep_Tool::Parameter( V, E );
625     }
626   }
627   if ( check )
628   {
629     double tol = BRep_Tool::Tolerance( E );
630     double f,l;  BRep_Tool::Range( E, f,l );
631     bool force = ( param < f-tol || param > l+tol );
632     if ( !force && pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
633       force = ( GetMeshDS()->ShapeToIndex( E ) != pos->GetShapeId() );
634
635     *check = CheckNodeU( E, n, param, 2*tol, force );
636   }
637   return param;
638 }
639
640 //=======================================================================
641 //function : CheckNodeU
642 //purpose  : Check and fix node U on an edge
643 //           Return false if U is bad and could not be fixed
644 //=======================================================================
645
646 bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge&   E,
647                                     const SMDS_MeshNode* n,
648                                     double&              u,
649                                     const double         tol,
650                                     const bool           force,
651                                     double*              distance) const
652 {
653   if ( force || !myOkNodePosShapes.count( n->GetPosition()->GetShapeId() ))
654   {
655     // check that u is correct
656     TopLoc_Location loc; double f,l;
657     Handle(Geom_Curve) curve = BRep_Tool::Curve( E,loc,f,l );
658     if ( curve.IsNull() ) // degenerated edge
659     {
660       if ( u+tol < f || u-tol > l )
661       {
662         double r = Max( 0.5, 1 - tol*n->GetID()); // to get a unique u on edge
663         u =  f*r + l*(1-r);
664       }
665     }
666     else
667     {
668       gp_Pnt nodePnt = SMESH_MeshEditor::TNodeXYZ( n );
669       if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
670       double dist = nodePnt.Distance( curve->Value( u ));
671       if ( distance ) *distance = dist;
672       if ( dist > tol )
673       {
674         // u incorrect, project the node to the curve
675         GeomAPI_ProjectPointOnCurve projector( nodePnt, curve, f, l );
676         if ( projector.NbPoints() < 1 )
677         {
678           MESSAGE( "SMESH_MesherHelper::CheckNodeU() failed to project" );
679           return false;
680         }
681         Quantity_Parameter U = projector.LowerDistanceParameter();
682         u = double( U );
683         dist = nodePnt.Distance( curve->Value( U ));
684         if ( distance ) *distance = dist;
685         if ( dist > tol )
686         {
687           MESSAGE( "SMESH_MesherHelper::CheckNodeU(), invalid projection" );
688           return false;
689         }
690         //u = double( U );
691       }
692       else if ( fabs( u ) > numeric_limits<double>::min() )
693       {
694         ((SMESH_MesherHelper*) this)->myOkNodePosShapes.insert( n->GetPosition()->GetShapeId() );
695       }
696       if (( u < f-tol || u > l+tol ) && force )
697       {
698         // node is on vertex but is set on periodic but trimmed edge (issue 0020890)
699         try
700         {
701           // do not use IsPeriodic() as Geom_TrimmedCurve::IsPeriodic () returns false
702           double period = curve->Period();
703           u = ( u < f ) ? u + period : u - period;
704         }
705         catch (Standard_Failure& exc)
706         {
707           return false;
708         }
709       }
710     }
711   }
712   return true;
713 }
714
715 //=======================================================================
716 //function : GetMediumNode
717 //purpose  : Return existing or create new medium nodes between given ones
718 //=======================================================================
719
720 const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
721                                                        const SMDS_MeshNode* n2,
722                                                        bool                 force3d)
723 {
724   // Find existing node
725
726   SMESH_TLink link(n1,n2);
727   ItTLinkNode itLN = myTLinkNodeMap.find( link );
728   if ( itLN != myTLinkNodeMap.end() ) {
729     return (*itLN).second;
730   }
731
732   // Create medium node
733
734   SMDS_MeshNode* n12;
735   SMESHDS_Mesh* meshDS = GetMeshDS();
736
737   // get type of shape for the new medium node
738   int faceID = -1, edgeID = -1;
739   const SMDS_PositionPtr Pos1 = n1->GetPosition();
740   const SMDS_PositionPtr Pos2 = n2->GetPosition();
741
742   if( myShape.IsNull() )
743   {
744     if( Pos1->GetTypeOfPosition()==SMDS_TOP_FACE ) {
745       faceID = Pos1->GetShapeId();
746     }
747     else if( Pos2->GetTypeOfPosition()==SMDS_TOP_FACE ) {
748       faceID = Pos2->GetShapeId();
749     }
750
751     if( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
752       edgeID = Pos1->GetShapeId();
753     }
754     if( Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE ) {
755       edgeID = Pos2->GetShapeId();
756     }
757   }
758   // get positions of the given nodes on shapes
759   TopoDS_Edge E; double u [2];
760   TopoDS_Face F; gp_XY  uv[2];
761   bool uvOK[2] = { false, false };
762   TopAbs_ShapeEnum shapeType = myShape.IsNull() ? TopAbs_SHAPE : myShape.ShapeType();
763   if ( faceID>0 || shapeType == TopAbs_FACE)
764   {
765     if( myShape.IsNull() )
766       F = TopoDS::Face(meshDS->IndexToShape(faceID));
767     else {
768       F = TopoDS::Face(myShape);
769       faceID = myShapeID;
770     }
771     uv[0] = GetNodeUV(F,n1,n2, force3d ? 0 : &uvOK[0]);
772     uv[1] = GetNodeUV(F,n2,n1, force3d ? 0 : &uvOK[1]);
773   }
774   else if (edgeID>0 || shapeType == TopAbs_EDGE)
775   {
776     if ( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE &&
777          Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE &&
778          Pos1->GetShapeId() != Pos2->GetShapeId() ) // issue 0021006
779     return getMediumNodeOnComposedWire(n1,n2,force3d);
780
781     if( myShape.IsNull() )
782       E = TopoDS::Edge(meshDS->IndexToShape(edgeID));
783     else {
784       E = TopoDS::Edge(myShape);
785       edgeID = myShapeID;
786     }
787     u[0] = GetNodeU(E,n1,n2, force3d ? 0 : &uvOK[0]);
788     u[1] = GetNodeU(E,n2,n1, force3d ? 0 : &uvOK[1]);
789   }
790   if(!force3d)
791   {
792     // we try to create medium node using UV parameters of
793     // nodes, else - medium between corresponding 3d points
794     if( ! F.IsNull() )
795     {
796       if ( uvOK[0] && uvOK[1] )
797       {
798         if ( IsDegenShape( Pos1->GetShapeId() ))
799           if ( myParIndex & U_periodic ) uv[0].SetCoord( 1, uv[1].Coord( 1 ));
800           else                           uv[0].SetCoord( 2, uv[1].Coord( 2 ));
801         else if ( IsDegenShape( Pos2->GetShapeId() ))
802           if ( myParIndex & U_periodic ) uv[1].SetCoord( 1, uv[0].Coord( 1 ));
803           else                           uv[1].SetCoord( 2, uv[0].Coord( 2 ));
804
805         TopLoc_Location loc;
806         Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
807         gp_XY UV = GetMiddleUV( S, uv[0], uv[1] );
808         gp_Pnt P = S->Value( UV.X(), UV.Y() ).Transformed(loc);
809         n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
810         meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y());
811         myTLinkNodeMap.insert(make_pair(link,n12));
812         return n12;
813       }
814     }
815     else if ( !E.IsNull() )
816     {
817       double f,l;
818       Handle(Geom_Curve) C = BRep_Tool::Curve(E, f, l);
819       if(!C.IsNull())
820       {
821         Standard_Boolean isPeriodic = C->IsPeriodic();
822         double U;
823         if(isPeriodic) {
824           Standard_Real Period = C->Period();
825           Standard_Real p = u[1]+ShapeAnalysis::AdjustByPeriod(u[1],u[0],Period);
826           Standard_Real pmid = (u[0]+p)/2.;
827           U = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
828         }
829         else
830           U = (u[0]+u[1])/2.;
831
832         gp_Pnt P = C->Value( U );
833         n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
834         meshDS->SetNodeOnEdge(n12, edgeID, U);
835         myTLinkNodeMap.insert(make_pair(link,n12));
836         return n12;
837       }
838     }
839   }
840   // 3d variant
841   double x = ( n1->X() + n2->X() )/2.;
842   double y = ( n1->Y() + n2->Y() )/2.;
843   double z = ( n1->Z() + n2->Z() )/2.;
844   n12 = meshDS->AddNode(x,y,z);
845   if ( !F.IsNull() )
846   {
847     gp_XY UV = ( uv[0] + uv[1] ) / 2.;
848     CheckNodeUV( F, n12, UV, 2*BRep_Tool::Tolerance( F ), /*force=*/true);
849     meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y() );
850   }
851   else if ( !E.IsNull() )
852   {
853     double U = ( u[0] + u[1] ) / 2.;
854     CheckNodeU( E, n12, U, 2*BRep_Tool::Tolerance( E ), /*force=*/true);
855     meshDS->SetNodeOnEdge(n12, edgeID, U);
856   }
857   else if ( myShapeID > 0 )
858   {
859     meshDS->SetNodeInVolume(n12, myShapeID);
860   }
861   myTLinkNodeMap.insert( make_pair( link, n12 ));
862   return n12;
863 }
864
865 //================================================================================
866 /*!
867  * \brief Makes a medium node if nodes reside different edges
868  */
869 //================================================================================
870
871 const SMDS_MeshNode* SMESH_MesherHelper::getMediumNodeOnComposedWire(const SMDS_MeshNode* n1,
872                                                                      const SMDS_MeshNode* n2,
873                                                                      bool                 force3d)
874 {
875   gp_Pnt middle = 0.5 * XYZ(n1) + 0.5 * XYZ(n2);
876   SMDS_MeshNode* n12 = AddNode( middle.X(), middle.Y(), middle.Z() );
877
878   // To find position on edge and 3D position for n12,
879   // project <middle> to 2 edges and select projection most close to <middle>
880
881   double u = 0, distMiddleProj = Precision::Infinite();
882   int iOkEdge = 0;
883   TopoDS_Edge edges[2];
884   for ( int is2nd = 0; is2nd < 2; ++is2nd )
885   {
886     // get an edge
887     const SMDS_MeshNode* n = is2nd ? n2 : n1;
888     TopoDS_Shape shape = GetSubShapeByNode( n, GetMeshDS() );
889     if ( shape.IsNull() || shape.ShapeType() != TopAbs_EDGE )
890       continue;
891
892     // project to get U of projection and distance from middle to projection
893     TopoDS_Edge edge = edges[ is2nd ] = TopoDS::Edge( shape );
894     double node2MiddleDist = middle.Distance( XYZ(n) );
895     double foundU = GetNodeU( edge, n ), foundDist = node2MiddleDist;
896     CheckNodeU( edge, n12, foundU, 2*BRep_Tool::Tolerance(edge), /*force=*/true, &foundDist );
897     if ( foundDist < node2MiddleDist )
898     {
899       distMiddleProj = foundDist;
900       u = foundU;
901       iOkEdge = is2nd;
902     }
903   }
904   if ( Precision::IsInfinite( distMiddleProj ))
905   {
906     // both projections failed; set n12 on the edge of n1 with U of a common vertex
907     TopoDS_Vertex vCommon;
908     if ( TopExp::CommonVertex( edges[0], edges[1], vCommon ))
909       u = BRep_Tool::Parameter( vCommon, edges[0] );
910     else
911     {
912       double f,l, u0 = GetNodeU( edges[0], n1 );
913       BRep_Tool::Range( edges[0],f,l );
914       u = ( fabs(u0-f) < fabs(u0-l) ) ? f : l;
915     }
916     iOkEdge = 0;
917     distMiddleProj = 0;
918   }
919
920   // move n12 to position of a successfull projection
921   double tol = BRep_Tool::Tolerance(edges[ iOkEdge ]);
922   if ( !force3d && distMiddleProj > 2*tol )
923   {
924     TopLoc_Location loc; double f,l;
925     Handle(Geom_Curve) curve = BRep_Tool::Curve( edges[iOkEdge],loc,f,l );
926     gp_Pnt p = curve->Value( u );
927     GetMeshDS()->MoveNode( n12, p.X(), p.Y(), p.Z() );
928   }
929
930   GetMeshDS()->SetNodeOnEdge(n12, edges[iOkEdge], u);
931
932   myTLinkNodeMap.insert( make_pair( SMESH_TLink(n1,n2), n12 ));
933
934   return n12;
935 }
936
937 //=======================================================================
938 //function : AddNode
939 //purpose  : Creates a node
940 //=======================================================================
941
942 SMDS_MeshNode* SMESH_MesherHelper::AddNode(double x, double y, double z, int ID)
943 {
944   SMESHDS_Mesh * meshDS = GetMeshDS();
945   SMDS_MeshNode* node = 0;
946   if ( ID )
947     node = meshDS->AddNodeWithID( x, y, z, ID );
948   else
949     node = meshDS->AddNode( x, y, z );
950   if ( mySetElemOnShape && myShapeID > 0 ) {
951     switch ( myShape.ShapeType() ) {
952     case TopAbs_SOLID:  meshDS->SetNodeInVolume( node, myShapeID); break;
953     case TopAbs_SHELL:  meshDS->SetNodeInVolume( node, myShapeID); break;
954     case TopAbs_FACE:   meshDS->SetNodeOnFace(   node, myShapeID); break;
955     case TopAbs_EDGE:   meshDS->SetNodeOnEdge(   node, myShapeID); break;
956     case TopAbs_VERTEX: meshDS->SetNodeOnVertex( node, myShapeID); break;
957     default: ;
958     }
959   }
960   return node;
961 }
962
963 //=======================================================================
964 //function : AddEdge
965 //purpose  : Creates quadratic or linear edge
966 //=======================================================================
967
968 SMDS_MeshEdge* SMESH_MesherHelper::AddEdge(const SMDS_MeshNode* n1,
969                                            const SMDS_MeshNode* n2,
970                                            const int            id,
971                                            const bool           force3d)
972 {
973   SMESHDS_Mesh * meshDS = GetMeshDS();
974   
975   SMDS_MeshEdge* edge = 0;
976   if (myCreateQuadratic) {
977     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
978     if(id)
979       edge = meshDS->AddEdgeWithID(n1, n2, n12, id);
980     else
981       edge = meshDS->AddEdge(n1, n2, n12);
982   }
983   else {
984     if(id)
985       edge = meshDS->AddEdgeWithID(n1, n2, id);
986     else
987       edge = meshDS->AddEdge(n1, n2);
988   }
989
990   if ( mySetElemOnShape && myShapeID > 0 )
991     meshDS->SetMeshElementOnShape( edge, myShapeID );
992
993   return edge;
994 }
995
996 //=======================================================================
997 //function : AddFace
998 //purpose  : Creates quadratic or linear triangle
999 //=======================================================================
1000
1001 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
1002                                            const SMDS_MeshNode* n2,
1003                                            const SMDS_MeshNode* n3,
1004                                            const int id,
1005                                            const bool force3d)
1006 {
1007   SMESHDS_Mesh * meshDS = GetMeshDS();
1008   SMDS_MeshFace* elem = 0;
1009
1010   if( n1==n2 || n2==n3 || n3==n1 )
1011     return elem;
1012
1013   if(!myCreateQuadratic) {
1014     if(id)
1015       elem = meshDS->AddFaceWithID(n1, n2, n3, id);
1016     else
1017       elem = meshDS->AddFace(n1, n2, n3);
1018   }
1019   else {
1020     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1021     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1022     const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1023
1024     if(id)
1025       elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, id);
1026     else
1027       elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31);
1028   }
1029   if ( mySetElemOnShape && myShapeID > 0 )
1030     meshDS->SetMeshElementOnShape( elem, myShapeID );
1031
1032   return elem;
1033 }
1034
1035 //=======================================================================
1036 //function : AddFace
1037 //purpose  : Creates quadratic or linear quadrangle
1038 //=======================================================================
1039
1040 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
1041                                            const SMDS_MeshNode* n2,
1042                                            const SMDS_MeshNode* n3,
1043                                            const SMDS_MeshNode* n4,
1044                                            const int            id,
1045                                            const bool           force3d)
1046 {
1047   SMESHDS_Mesh * meshDS = GetMeshDS();
1048   SMDS_MeshFace* elem = 0;
1049
1050   if( n1==n2 ) {
1051     return AddFace(n1,n3,n4,id,force3d);
1052   }
1053   if( n1==n3 ) {
1054     return AddFace(n1,n2,n4,id,force3d);
1055   }
1056   if( n1==n4 ) {
1057     return AddFace(n1,n2,n3,id,force3d);
1058   }
1059   if( n2==n3 ) {
1060     return AddFace(n1,n2,n4,id,force3d);
1061   }
1062   if( n2==n4 ) {
1063     return AddFace(n1,n2,n3,id,force3d);
1064   }
1065   if( n3==n4 ) {
1066     return AddFace(n1,n2,n3,id,force3d);
1067   }
1068
1069   if(!myCreateQuadratic) {
1070     if(id)
1071       elem = meshDS->AddFaceWithID(n1, n2, n3, n4, id);
1072     else
1073       elem = meshDS->AddFace(n1, n2, n3, n4);
1074   }
1075   else {
1076     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1077     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1078     const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1079     const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1080
1081     if(id)
1082       elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, id);
1083     else
1084       elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41);
1085   }
1086   if ( mySetElemOnShape && myShapeID > 0 )
1087     meshDS->SetMeshElementOnShape( elem, myShapeID );
1088
1089   return elem;
1090 }
1091
1092 //=======================================================================
1093 //function : AddPolygonalFace
1094 //purpose  : Creates polygon, with additional nodes in quadratic mesh
1095 //=======================================================================
1096
1097 SMDS_MeshFace* SMESH_MesherHelper::AddPolygonalFace (const vector<const SMDS_MeshNode*>& nodes,
1098                                                      const int                           id,
1099                                                      const bool                          force3d)
1100 {
1101   SMESHDS_Mesh * meshDS = GetMeshDS();
1102   SMDS_MeshFace* elem = 0;
1103
1104   if(!myCreateQuadratic) {
1105     if(id)
1106       elem = meshDS->AddPolygonalFaceWithID(nodes, id);
1107     else
1108       elem = meshDS->AddPolygonalFace(nodes);
1109   }
1110   else {
1111     vector<const SMDS_MeshNode*> newNodes;
1112     for ( int i = 0; i < nodes.size(); ++i )
1113     {
1114       const SMDS_MeshNode* n1 = nodes[i];
1115       const SMDS_MeshNode* n2 = nodes[(i+1)/nodes.size()];
1116       const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1117       newNodes.push_back( n1 );
1118       newNodes.push_back( n12 );
1119     }
1120     if(id)
1121       elem = meshDS->AddPolygonalFaceWithID(newNodes, id);
1122     else
1123       elem = meshDS->AddPolygonalFace(newNodes);
1124   }
1125   if ( mySetElemOnShape && myShapeID > 0 )
1126     meshDS->SetMeshElementOnShape( elem, myShapeID );
1127
1128   return elem;
1129 }
1130
1131 //=======================================================================
1132 //function : AddVolume
1133 //purpose  : Creates quadratic or linear prism
1134 //=======================================================================
1135
1136 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1137                                                const SMDS_MeshNode* n2,
1138                                                const SMDS_MeshNode* n3,
1139                                                const SMDS_MeshNode* n4,
1140                                                const SMDS_MeshNode* n5,
1141                                                const SMDS_MeshNode* n6,
1142                                                const int id,
1143                                                const bool force3d)
1144 {
1145   SMESHDS_Mesh * meshDS = GetMeshDS();
1146   SMDS_MeshVolume* elem = 0;
1147   if(!myCreateQuadratic) {
1148     if(id)
1149       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, id);
1150     else
1151       elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6);
1152   }
1153   else {
1154     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1155     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1156     const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1157
1158     const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1159     const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1160     const SMDS_MeshNode* n64 = GetMediumNode(n6,n4,force3d);
1161
1162     const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1163     const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1164     const SMDS_MeshNode* n36 = GetMediumNode(n3,n6,force3d);
1165
1166     if(id)
1167       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, 
1168                                      n12, n23, n31, n45, n56, n64, n14, n25, n36, id);
1169     else
1170       elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6,
1171                                n12, n23, n31, n45, n56, n64, n14, n25, n36);
1172   }
1173   if ( mySetElemOnShape && myShapeID > 0 )
1174     meshDS->SetMeshElementOnShape( elem, myShapeID );
1175
1176   return elem;
1177 }
1178
1179 //=======================================================================
1180 //function : AddVolume
1181 //purpose  : Creates quadratic or linear tetrahedron
1182 //=======================================================================
1183
1184 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1185                                                const SMDS_MeshNode* n2,
1186                                                const SMDS_MeshNode* n3,
1187                                                const SMDS_MeshNode* n4,
1188                                                const int id, 
1189                                                const bool force3d)
1190 {
1191   SMESHDS_Mesh * meshDS = GetMeshDS();
1192   SMDS_MeshVolume* elem = 0;
1193   if(!myCreateQuadratic) {
1194     if(id)
1195       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, id);
1196     else
1197       elem = meshDS->AddVolume(n1, n2, n3, n4);
1198   }
1199   else {
1200     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1201     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1202     const SMDS_MeshNode* n31 = GetMediumNode(n3,n1,force3d);
1203
1204     const SMDS_MeshNode* n14 = GetMediumNode(n1,n4,force3d);
1205     const SMDS_MeshNode* n24 = GetMediumNode(n2,n4,force3d);
1206     const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1207
1208     if(id)
1209       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34, id);
1210     else
1211       elem = meshDS->AddVolume(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34);
1212   }
1213   if ( mySetElemOnShape && myShapeID > 0 )
1214     meshDS->SetMeshElementOnShape( elem, myShapeID );
1215
1216   return elem;
1217 }
1218
1219 //=======================================================================
1220 //function : AddVolume
1221 //purpose  : Creates quadratic or linear pyramid
1222 //=======================================================================
1223
1224 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1225                                                const SMDS_MeshNode* n2,
1226                                                const SMDS_MeshNode* n3,
1227                                                const SMDS_MeshNode* n4,
1228                                                const SMDS_MeshNode* n5,
1229                                                const int id, 
1230                                                const bool force3d)
1231 {
1232   SMDS_MeshVolume* elem = 0;
1233   if(!myCreateQuadratic) {
1234     if(id)
1235       elem = GetMeshDS()->AddVolumeWithID(n1, n2, n3, n4, n5, id);
1236     else
1237       elem = GetMeshDS()->AddVolume(n1, n2, n3, n4, n5);
1238   }
1239   else {
1240     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1241     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1242     const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1243     const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1244
1245     const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1246     const SMDS_MeshNode* n25 = GetMediumNode(n2,n5,force3d);
1247     const SMDS_MeshNode* n35 = GetMediumNode(n3,n5,force3d);
1248     const SMDS_MeshNode* n45 = GetMediumNode(n4,n5,force3d);
1249
1250     if(id)
1251       elem = GetMeshDS()->AddVolumeWithID ( n1,  n2,  n3,  n4,  n5,
1252                                             n12, n23, n34, n41,
1253                                             n15, n25, n35, n45,
1254                                             id);
1255     else
1256       elem = GetMeshDS()->AddVolume( n1,  n2,  n3,  n4,  n5,
1257                                      n12, n23, n34, n41,
1258                                      n15, n25, n35, n45);
1259   }
1260   if ( mySetElemOnShape && myShapeID > 0 )
1261     GetMeshDS()->SetMeshElementOnShape( elem, myShapeID );
1262
1263   return elem;
1264 }
1265
1266 //=======================================================================
1267 //function : AddVolume
1268 //purpose  : Creates quadratic or linear hexahedron
1269 //=======================================================================
1270
1271 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
1272                                                const SMDS_MeshNode* n2,
1273                                                const SMDS_MeshNode* n3,
1274                                                const SMDS_MeshNode* n4,
1275                                                const SMDS_MeshNode* n5,
1276                                                const SMDS_MeshNode* n6,
1277                                                const SMDS_MeshNode* n7,
1278                                                const SMDS_MeshNode* n8,
1279                                                const int id,
1280                                                const bool force3d)
1281 {
1282   SMESHDS_Mesh * meshDS = GetMeshDS();
1283   SMDS_MeshVolume* elem = 0;
1284   if(!myCreateQuadratic) {
1285     if(id)
1286       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8, id);
1287     else
1288       elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8);
1289   }
1290   else {
1291     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1292     const SMDS_MeshNode* n23 = GetMediumNode(n2,n3,force3d);
1293     const SMDS_MeshNode* n34 = GetMediumNode(n3,n4,force3d);
1294     const SMDS_MeshNode* n41 = GetMediumNode(n4,n1,force3d);
1295
1296     const SMDS_MeshNode* n56 = GetMediumNode(n5,n6,force3d);
1297     const SMDS_MeshNode* n67 = GetMediumNode(n6,n7,force3d);
1298     const SMDS_MeshNode* n78 = GetMediumNode(n7,n8,force3d);
1299     const SMDS_MeshNode* n85 = GetMediumNode(n8,n5,force3d);
1300
1301     const SMDS_MeshNode* n15 = GetMediumNode(n1,n5,force3d);
1302     const SMDS_MeshNode* n26 = GetMediumNode(n2,n6,force3d);
1303     const SMDS_MeshNode* n37 = GetMediumNode(n3,n7,force3d);
1304     const SMDS_MeshNode* n48 = GetMediumNode(n4,n8,force3d);
1305
1306     if(id)
1307       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
1308                                      n12, n23, n34, n41, n56, n67,
1309                                      n78, n85, n15, n26, n37, n48, id);
1310     else
1311       elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
1312                                n12, n23, n34, n41, n56, n67,
1313                                n78, n85, n15, n26, n37, n48);
1314   }
1315   if ( mySetElemOnShape && myShapeID > 0 )
1316     meshDS->SetMeshElementOnShape( elem, myShapeID );
1317
1318   return elem;
1319 }
1320
1321 //=======================================================================
1322 //function : AddPolyhedralVolume
1323 //purpose  : Creates polyhedron. In quadratic mesh, adds medium nodes
1324 //=======================================================================
1325
1326 SMDS_MeshVolume*
1327 SMESH_MesherHelper::AddPolyhedralVolume (const std::vector<const SMDS_MeshNode*>& nodes,
1328                                          const std::vector<int>&                  quantities,
1329                                          const int                                id,
1330                                          const bool                               force3d)
1331 {
1332   SMESHDS_Mesh * meshDS = GetMeshDS();
1333   SMDS_MeshVolume* elem = 0;
1334   if(!myCreateQuadratic)
1335   {
1336     if(id)
1337       elem = meshDS->AddPolyhedralVolumeWithID(nodes, quantities, id);
1338     else
1339       elem = meshDS->AddPolyhedralVolume(nodes, quantities);
1340   }
1341   else
1342   {
1343     vector<const SMDS_MeshNode*> newNodes;
1344     vector<int> newQuantities;
1345     for ( int iFace=0, iN=0; iFace < quantities.size(); ++iFace)
1346     {
1347       int nbNodesInFace = quantities[iFace];
1348       newQuantities.push_back(0);
1349       for ( int i = 0; i < nbNodesInFace; ++i )
1350       {
1351         const SMDS_MeshNode* n1 = nodes[ iN + i ];
1352         newNodes.push_back( n1 );
1353         newQuantities.back()++;
1354         
1355         const SMDS_MeshNode* n2 = nodes[ iN + ( i+1==nbNodesInFace ? 0 : i+1 )];
1356 //         if ( n1->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE &&
1357 //              n2->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE )
1358         {
1359           const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1360           newNodes.push_back( n12 );
1361           newQuantities.back()++;
1362         }
1363       }
1364       iN += nbNodesInFace;
1365     }
1366     if(id)
1367       elem = meshDS->AddPolyhedralVolumeWithID( newNodes, newQuantities, id );
1368     else
1369       elem = meshDS->AddPolyhedralVolume( newNodes, newQuantities );
1370   }
1371   if ( mySetElemOnShape && myShapeID > 0 )
1372     meshDS->SetMeshElementOnShape( elem, myShapeID );
1373
1374   return elem;
1375 }
1376
1377 //=======================================================================
1378 //function : LoadNodeColumns
1379 //purpose  : Load nodes bound to face into a map of node columns
1380 //=======================================================================
1381
1382 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
1383                                          const TopoDS_Face& theFace,
1384                                          const TopoDS_Edge& theBaseEdge,
1385                                          SMESHDS_Mesh*      theMesh)
1386 {
1387   SMESHDS_SubMesh* faceSubMesh = theMesh->MeshElements( theFace );
1388   if ( !faceSubMesh || faceSubMesh->NbElements() == 0 )
1389     return false;
1390
1391   // get nodes on theBaseEdge sorted by param on edge and initialize theParam2ColumnMap with them
1392
1393   map< double, const SMDS_MeshNode*> sortedBaseNodes;
1394   if ( !SMESH_Algo::GetSortedNodesOnEdge( theMesh, theBaseEdge,/*noMedium=*/true, sortedBaseNodes)
1395        || sortedBaseNodes.size() < 2 )
1396     return false;
1397
1398   int nbRows = faceSubMesh->NbElements() / ( sortedBaseNodes.size()-1 ) + 1;
1399   map< double, const SMDS_MeshNode*>::iterator u_n = sortedBaseNodes.begin();
1400   double f = u_n->first, range = sortedBaseNodes.rbegin()->first - f;
1401   for ( ; u_n != sortedBaseNodes.end(); u_n++ )
1402   {
1403     double par = ( u_n->first - f ) / range;
1404     vector<const SMDS_MeshNode*>& nCol = theParam2ColumnMap[ par ];
1405     nCol.resize( nbRows );
1406     nCol[0] = u_n->second;
1407   }
1408
1409   // fill theParam2ColumnMap column by column by passing from nodes on
1410   // theBaseEdge up via mesh faces on theFace
1411
1412   TParam2ColumnMap::iterator par_nVec_2 = theParam2ColumnMap.begin();
1413   TParam2ColumnMap::iterator par_nVec_1 = par_nVec_2++;
1414   TIDSortedElemSet emptySet, avoidSet;
1415   for ( ; par_nVec_2 != theParam2ColumnMap.end(); ++par_nVec_1, ++par_nVec_2 )
1416   {
1417     vector<const SMDS_MeshNode*>& nCol1 = par_nVec_1->second;
1418     vector<const SMDS_MeshNode*>& nCol2 = par_nVec_2->second;
1419
1420     int i1, i2, iRow = 0;
1421     const SMDS_MeshNode *n1 = nCol1[0], *n2 = nCol2[0];
1422     // find face sharing node n1 and n2 and belonging to faceSubMesh
1423     while ( const SMDS_MeshElement* face =
1424             SMESH_MeshEditor::FindFaceInSet( n1, n2, emptySet, avoidSet, &i1, &i2))
1425     {
1426       if ( faceSubMesh->Contains( face ))
1427       {
1428         int nbNodes = face->IsQuadratic() ? face->NbNodes()/2 : face->NbNodes();
1429         if ( nbNodes != 4 )
1430           return false;
1431         n1 = face->GetNode( (i2+2) % 4 ); // opposite corner of quadrangle face
1432         n2 = face->GetNode( (i1+2) % 4 );
1433         if ( ++iRow >= nbRows )
1434           return false;
1435         nCol1[ iRow ] = n1;
1436         nCol2[ iRow ] = n2;
1437         avoidSet.clear();
1438       }
1439       avoidSet.insert( face );
1440     }
1441     if ( iRow + 1 < nbRows ) // compact if necessary
1442       nCol1.resize( iRow + 1 ), nCol2.resize( iRow + 1 );
1443   }
1444   return true;
1445 }
1446
1447 //=======================================================================
1448 //function : NbAncestors
1449 //purpose  : Return number of unique ancestors of the shape
1450 //=======================================================================
1451
1452 int SMESH_MesherHelper::NbAncestors(const TopoDS_Shape& shape,
1453                                     const SMESH_Mesh&   mesh,
1454                                     TopAbs_ShapeEnum    ancestorType/*=TopAbs_SHAPE*/)
1455 {
1456   TopTools_MapOfShape ancestors;
1457   TopTools_ListIteratorOfListOfShape ansIt( mesh.GetAncestors(shape) );
1458   for ( ; ansIt.More(); ansIt.Next() ) {
1459     if ( ancestorType == TopAbs_SHAPE || ansIt.Value().ShapeType() == ancestorType )
1460       ancestors.Add( ansIt.Value() );
1461   }
1462   return ancestors.Extent();
1463 }
1464
1465 //=======================================================================
1466 //function : GetSubShapeOri
1467 //purpose  : Return orientation of sub-shape in the main shape
1468 //=======================================================================
1469
1470 TopAbs_Orientation SMESH_MesherHelper::GetSubShapeOri(const TopoDS_Shape& shape,
1471                                                       const TopoDS_Shape& subShape)
1472 {
1473   TopAbs_Orientation ori = TopAbs_Orientation(-1);
1474   if ( !shape.IsNull() && !subShape.IsNull() )
1475   {
1476     TopExp_Explorer e( shape, subShape.ShapeType() );
1477     if ( shape.Orientation() >= TopAbs_INTERNAL ) // TopAbs_INTERNAL or TopAbs_EXTERNAL
1478       e.Init( shape.Oriented(TopAbs_FORWARD), subShape.ShapeType() );
1479     for ( ; e.More(); e.Next())
1480       if ( subShape.IsSame( e.Current() ))
1481         break;
1482     if ( e.More() )
1483       ori = e.Current().Orientation();
1484   }
1485   return ori;
1486 }
1487
1488 //=======================================================================
1489 //function : IsSubShape
1490 //purpose  : 
1491 //=======================================================================
1492
1493 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape,
1494                                      const TopoDS_Shape& mainShape )
1495 {
1496   if ( !shape.IsNull() && !mainShape.IsNull() )
1497   {
1498     for ( TopExp_Explorer exp( mainShape, shape.ShapeType());
1499           exp.More();
1500           exp.Next() )
1501       if ( shape.IsSame( exp.Current() ))
1502         return true;
1503   }
1504   SCRUTE((shape.IsNull()));
1505   SCRUTE((mainShape.IsNull()));
1506   return false;
1507 }
1508
1509 //=======================================================================
1510 //function : IsSubShape
1511 //purpose  : 
1512 //=======================================================================
1513
1514 bool SMESH_MesherHelper::IsSubShape( const TopoDS_Shape& shape, SMESH_Mesh* aMesh )
1515 {
1516   if ( shape.IsNull() || !aMesh )
1517     return false;
1518   return
1519     aMesh->GetMeshDS()->ShapeToIndex( shape ) ||
1520     // PAL16202
1521     shape.ShapeType() == TopAbs_COMPOUND && aMesh->GetMeshDS()->IsGroupOfSubShapes( shape );
1522 }
1523
1524 //=======================================================================
1525 //function : IsQuadraticMesh
1526 //purpose  : Check mesh without geometry for: if all elements on this shape are quadratic,
1527 //           quadratic elements will be created.
1528 //           Used then generated 3D mesh without geometry.
1529 //=======================================================================
1530
1531 SMESH_MesherHelper:: MType SMESH_MesherHelper::IsQuadraticMesh()
1532 {
1533   int NbAllEdgsAndFaces=0;
1534   int NbQuadFacesAndEdgs=0;
1535   int NbFacesAndEdges=0;
1536   //All faces and edges
1537   NbAllEdgsAndFaces = myMesh->NbEdges() + myMesh->NbFaces();
1538   
1539   //Quadratic faces and edges
1540   NbQuadFacesAndEdgs = myMesh->NbEdges(ORDER_QUADRATIC) + myMesh->NbFaces(ORDER_QUADRATIC);
1541
1542   //Linear faces and edges
1543   NbFacesAndEdges = myMesh->NbEdges(ORDER_LINEAR) + myMesh->NbFaces(ORDER_LINEAR);
1544   
1545   if (NbAllEdgsAndFaces == NbQuadFacesAndEdgs) {
1546     //Quadratic mesh
1547     return SMESH_MesherHelper::QUADRATIC;
1548   }
1549   else if (NbAllEdgsAndFaces == NbFacesAndEdges) {
1550     //Linear mesh
1551     return SMESH_MesherHelper::LINEAR;
1552   }
1553   else
1554     //Mesh with both type of elements
1555     return SMESH_MesherHelper::COMP;
1556 }
1557
1558 //=======================================================================
1559 //function : GetOtherParam
1560 //purpose  : Return an alternative parameter for a node on seam
1561 //=======================================================================
1562
1563 double SMESH_MesherHelper::GetOtherParam(const double param) const
1564 {
1565   int i = myParIndex & U_periodic ? 0 : 1;
1566   return fabs(param-myPar1[i]) < fabs(param-myPar2[i]) ? myPar2[i] : myPar1[i];
1567 }
1568
1569 //#include <Perf_Meter.hxx>
1570
1571 //=======================================================================
1572 namespace { // Structures used by FixQuadraticElements()
1573 //=======================================================================
1574
1575 #define __DMP__(txt) \
1576   //cout << txt
1577 #define MSG(txt) __DMP__(txt<<endl)
1578 #define MSGBEG(txt) __DMP__(txt)
1579
1580   //const double straightTol2 = 1e-33; // to detect straing links
1581   bool isStraightLink(double linkLen2, double middleNodeMove2)
1582   {
1583     // straight if <node move> < 1/15 * <link length>
1584     return middleNodeMove2 < 1/15./15. * linkLen2;
1585   }
1586
1587   struct QFace;
1588   // ---------------------------------------
1589   /*!
1590    * \brief Quadratic link knowing its faces
1591    */
1592   struct QLink: public SMESH_TLink
1593   {
1594     const SMDS_MeshNode*          _mediumNode;
1595     mutable vector<const QFace* > _faces;
1596     mutable gp_Vec                _nodeMove;
1597     mutable int                   _nbMoves;
1598
1599     QLink(const SMDS_MeshNode* n1, const SMDS_MeshNode* n2, const SMDS_MeshNode* nm):
1600       SMESH_TLink( n1,n2 ), _mediumNode(nm), _nodeMove(0,0,0), _nbMoves(0) {
1601       _faces.reserve(4);
1602       //if ( MediumPos() != SMDS_TOP_3DSPACE )
1603         _nodeMove = MediumPnt() - MiddlePnt();
1604     }
1605     void SetContinuesFaces() const;
1606     const QFace* GetContinuesFace( const QFace* face ) const;
1607     bool OnBoundary() const;
1608     gp_XYZ MiddlePnt() const { return ( XYZ( node1() ) + XYZ( node2() )) / 2.; }
1609     gp_XYZ MediumPnt() const { return XYZ( _mediumNode ); }
1610
1611     SMDS_TypeOfPosition MediumPos() const
1612     { return _mediumNode->GetPosition()->GetTypeOfPosition(); }
1613     SMDS_TypeOfPosition EndPos(bool isSecond) const
1614     { return (isSecond ? node2() : node1())->GetPosition()->GetTypeOfPosition(); }
1615     const SMDS_MeshNode* EndPosNode(SMDS_TypeOfPosition pos) const
1616     { return EndPos(0) == pos ? node1() : EndPos(1) == pos ? node2() : 0; }
1617
1618     void Move(const gp_Vec& move, bool sum=false) const
1619     { _nodeMove += move; _nbMoves += sum ? (_nbMoves==0) : 1; }
1620     gp_XYZ Move() const { return _nodeMove.XYZ() / _nbMoves; }
1621     bool IsMoved() const { return (_nbMoves > 0 && !IsStraight()); }
1622     bool IsStraight() const
1623     { return isStraightLink( (XYZ(node1())-XYZ(node2())).SquareModulus(),
1624                              _nodeMove.SquareMagnitude());
1625     }
1626     bool operator<(const QLink& other) const {
1627       return (node1()->GetID() == other.node1()->GetID() ?
1628               node2()->GetID() < other.node2()->GetID() :
1629               node1()->GetID() < other.node1()->GetID());
1630     }
1631     struct PtrComparator {
1632       bool operator() (const QLink* l1, const QLink* l2 ) const { return *l1 < *l2; }
1633     };
1634   };
1635   // ---------------------------------------------------------
1636   /*!
1637    * \brief Link in the chain of links; it connects two faces
1638    */
1639   struct TChainLink
1640   {
1641     const QLink*         _qlink;
1642     mutable const QFace* _qfaces[2];
1643
1644     TChainLink(const QLink* qlink=0):_qlink(qlink) {
1645       _qfaces[0] = _qfaces[1] = 0;
1646     }
1647     void SetFace(const QFace* face) const { int iF = _qfaces[0] ? 1 : 0; _qfaces[iF]=face; }
1648
1649     bool IsBoundary() const { return !_qfaces[1]; }
1650
1651     void RemoveFace( const QFace* face ) const
1652     { _qfaces[(face == _qfaces[1])] = 0; if (!_qfaces[0]) std::swap(_qfaces[0],_qfaces[1]); }
1653
1654     const QFace* NextFace( const QFace* f ) const
1655     { return _qfaces[0]==f ? _qfaces[1] : _qfaces[0]; }
1656
1657     const SMDS_MeshNode* NextNode( const SMDS_MeshNode* n ) const
1658     { return n == _qlink->node1() ? _qlink->node2() : _qlink->node1(); }
1659
1660     bool operator<(const TChainLink& other) const { return *_qlink < *other._qlink; }
1661
1662     operator bool() const { return (_qlink); }
1663
1664     const QLink* operator->() const { return _qlink; }
1665
1666     gp_Vec Normal() const;
1667   };
1668   // --------------------------------------------------------------------
1669   typedef list< TChainLink > TChain;
1670   typedef set < TChainLink > TLinkSet;
1671   typedef TLinkSet::const_iterator TLinkInSet;
1672
1673   const int theFirstStep = 5;
1674
1675   enum { ERR_OK, ERR_TRI, ERR_PRISM, ERR_UNKNOWN }; // errors of QFace::GetLinkChain()
1676   // --------------------------------------------------------------------
1677   /*!
1678    * \brief Face shared by two volumes and bound by QLinks
1679    */
1680   struct QFace: public TIDSortedElemSet
1681   {
1682     mutable const SMDS_MeshElement* _volumes[2];
1683     mutable vector< const QLink* >  _sides;
1684     mutable bool                    _sideIsAdded[4]; // added in chain of links
1685     gp_Vec                          _normal;
1686 #ifdef _DEBUG_
1687     mutable const SMDS_MeshElement* _face;
1688 #endif
1689
1690     QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face=0 );
1691
1692     void SetVolume(const SMDS_MeshElement* v) const { _volumes[ _volumes[0] ? 1 : 0 ] = v; }
1693
1694     int NbVolumes() const { return !_volumes[0] ? 0 : !_volumes[1] ? 1 : 2; }
1695
1696     void AddSelfToLinks() const {
1697       for ( int i = 0; i < _sides.size(); ++i )
1698         _sides[i]->_faces.push_back( this );
1699     }
1700     int LinkIndex( const QLink* side ) const {
1701       for (int i=0; i<_sides.size(); ++i ) if ( _sides[i] == side ) return i;
1702       return -1;
1703     }
1704     bool GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& err) const;
1705
1706     bool GetLinkChain( TChainLink& link, TChain& chain, SMDS_TypeOfPosition pos, int& err) const
1707     {
1708       int i = LinkIndex( link._qlink );
1709       if ( i < 0 ) return true;
1710       _sideIsAdded[i] = true;
1711       link.SetFace( this );
1712       // continue from opposite link
1713       return GetLinkChain( (i+2)%_sides.size(), chain, pos, err );
1714     }
1715     bool IsBoundary() const { return !_volumes[1]; }
1716
1717     bool Contains( const SMDS_MeshNode* node ) const { return count(node); }
1718
1719     TLinkInSet GetBoundaryLink( const TLinkSet&      links,
1720                                 const TChainLink&    avoidLink,
1721                                 TLinkInSet *         notBoundaryLink = 0,
1722                                 const SMDS_MeshNode* nodeToContain = 0,
1723                                 bool *               isAdjacentUsed = 0,
1724                                 int                  nbRecursionsLeft = -1) const;
1725
1726     TLinkInSet GetLinkByNode( const TLinkSet&      links,
1727                               const TChainLink&    avoidLink,
1728                               const SMDS_MeshNode* nodeToContain) const;
1729
1730     const SMDS_MeshNode* GetNodeInFace() const {
1731       for ( int iL = 0; iL < _sides.size(); ++iL )
1732         if ( _sides[iL]->MediumPos() == SMDS_TOP_FACE ) return _sides[iL]->_mediumNode;
1733       return 0;
1734     }
1735
1736     gp_Vec LinkNorm(const int i, SMESH_MesherHelper* theFaceHelper=0) const;
1737
1738     double MoveByBoundary( const TChainLink&   theLink,
1739                            const gp_Vec&       theRefVec,
1740                            const TLinkSet&     theLinks,
1741                            SMESH_MesherHelper* theFaceHelper=0,
1742                            const double        thePrevLen=0,
1743                            const int           theStep=theFirstStep,
1744                            gp_Vec*             theLinkNorm=0,
1745                            double              theSign=1.0) const;
1746   };
1747
1748   //================================================================================
1749   /*!
1750    * \brief Dump QLink and QFace
1751    */
1752   ostream& operator << (ostream& out, const QLink& l)
1753   {
1754     out <<"QLink nodes: "
1755         << l.node1()->GetID() << " - "
1756         << l._mediumNode->GetID() << " - "
1757         << l.node2()->GetID() << endl;
1758     return out;
1759   }
1760   ostream& operator << (ostream& out, const QFace& f)
1761   {
1762     out <<"QFace nodes: "/*<< &f << "  "*/;
1763     for ( TIDSortedElemSet::const_iterator n = f.begin(); n != f.end(); ++n )
1764       out << (*n)->GetID() << " ";
1765     out << " \tvolumes: "
1766         << (f._volumes[0] ? f._volumes[0]->GetID() : 0) << " "
1767         << (f._volumes[1] ? f._volumes[1]->GetID() : 0);
1768     out << "  \tNormal: "<< f._normal.X() <<", "<<f._normal.Y() <<", "<<f._normal.Z() << endl;
1769     return out;
1770   }
1771
1772   //================================================================================
1773   /*!
1774    * \brief Construct QFace from QLinks 
1775    */
1776   //================================================================================
1777
1778   QFace::QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face )
1779   {
1780     _volumes[0] = _volumes[1] = 0;
1781     _sides = links;
1782     _sideIsAdded[0]=_sideIsAdded[1]=_sideIsAdded[2]=_sideIsAdded[3]=false;
1783     _normal.SetCoord(0,0,0);
1784     for ( int i = 1; i < _sides.size(); ++i ) {
1785       const QLink *l1 = _sides[i-1], *l2 = _sides[i];
1786       insert( l1->node1() ); insert( l1->node2() );
1787       // compute normal
1788       gp_Vec v1( XYZ( l1->node2()), XYZ( l1->node1()));
1789       gp_Vec v2( XYZ( l2->node1()), XYZ( l2->node2()));
1790       if ( l1->node1() != l2->node1() && l1->node2() != l2->node2() )
1791         v1.Reverse(); 
1792       _normal += v1 ^ v2;
1793     }
1794     double normSqSize = _normal.SquareMagnitude();
1795     if ( normSqSize > numeric_limits<double>::min() )
1796       _normal /= sqrt( normSqSize );
1797     else
1798       _normal.SetCoord(1e-33,0,0);
1799
1800 #ifdef _DEBUG_
1801     _face = face;
1802 #endif
1803   }
1804   //================================================================================
1805   /*!
1806    * \brief Make up a chain of links
1807    *  \param iSide - link to add first
1808    *  \param chain - chain to fill in
1809    *  \param pos   - postion of medium nodes the links should have
1810    *  \param error - out, specifies what is wrong
1811    *  \retval bool - false if valid chain can't be built; "valid" means that links
1812    *                 of the chain belongs to rectangles bounding hexahedrons
1813    */
1814   //================================================================================
1815
1816   bool QFace::GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
1817   {
1818     if ( iSide >= _sides.size() ) // wrong argument iSide
1819       return false;
1820     if ( _sideIsAdded[ iSide ]) // already in chain
1821       return true;
1822
1823     if ( _sides.size() != 4 ) { // triangle - visit all my continous faces
1824       MSGBEG( *this );
1825       TLinkSet links;
1826       list< const QFace* > faces( 1, this );
1827       while ( !faces.empty() ) {
1828         const QFace* face = faces.front();
1829         for ( int i = 0; i < face->_sides.size(); ++i ) {
1830           if ( !face->_sideIsAdded[i] && face->_sides[i] ) {
1831             face->_sideIsAdded[i] = true;
1832             // find a face side in the chain
1833             TLinkInSet chLink = links.insert( TChainLink(face->_sides[i])).first;
1834 //             TChain::iterator chLink = chain.begin();
1835 //             for ( ; chLink != chain.end(); ++chLink )
1836 //               if ( chLink->_qlink == face->_sides[i] )
1837 //                 break;
1838 //             if ( chLink == chain.end() )
1839 //               chLink = chain.insert( chain.begin(), TChainLink(face->_sides[i]));
1840             // add a face to a chained link and put a continues face in the queue
1841             chLink->SetFace( face );
1842             if ( face->_sides[i]->MediumPos() >= pos )
1843               if ( const QFace* contFace = face->_sides[i]->GetContinuesFace( face ))
1844                 faces.push_back( contFace );
1845           }
1846         }
1847         faces.pop_front();
1848       }
1849       if ( error < ERR_TRI )
1850         error = ERR_TRI;
1851       chain.insert( chain.end(), links.begin(),links.end() );
1852       return false;
1853     }
1854     _sideIsAdded[iSide] = true; // not to add this link to chain again
1855     const QLink* link = _sides[iSide];
1856     if ( !link)
1857       return true;
1858
1859     // add link into chain
1860     TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(link));
1861     chLink->SetFace( this );
1862     MSGBEG( *this );
1863
1864     // propagate from quadrangle to neighbour faces
1865     if ( link->MediumPos() >= pos ) {
1866       int nbLinkFaces = link->_faces.size();
1867       if ( nbLinkFaces == 4 || nbLinkFaces < 4 && link->OnBoundary()) {
1868         // hexahedral mesh or boundary quadrangles - goto a continous face
1869         if ( const QFace* f = link->GetContinuesFace( this ))
1870           return f->GetLinkChain( *chLink, chain, pos, error );
1871       }
1872       else {
1873         TChainLink chLink(link); // side face of prismatic mesh - visit all faces of iSide
1874         for ( int i = 0; i < nbLinkFaces; ++i )
1875           if ( link->_faces[i] )
1876             link->_faces[i]->GetLinkChain( chLink, chain, pos, error );
1877         if ( error < ERR_PRISM )
1878           error = ERR_PRISM;
1879         return false;
1880       }
1881     }
1882     return true;
1883   }
1884
1885   //================================================================================
1886   /*!
1887    * \brief Return a boundary link of the triangle face
1888    *  \param links - set of all links
1889    *  \param avoidLink - link not to return
1890    *  \param notBoundaryLink - out, neither the returned link nor avoidLink
1891    *  \param nodeToContain - node the returned link must contain; if provided, search
1892    *                         also performed on adjacent faces
1893    *  \param isAdjacentUsed - returns true if link is found in adjacent faces
1894    *  \param nbRecursionsLeft - to limit recursion
1895    */
1896   //================================================================================
1897
1898   TLinkInSet QFace::GetBoundaryLink( const TLinkSet&      links,
1899                                      const TChainLink&    avoidLink,
1900                                      TLinkInSet *         notBoundaryLink,
1901                                      const SMDS_MeshNode* nodeToContain,
1902                                      bool *               isAdjacentUsed,
1903                                      int                  nbRecursionsLeft) const
1904   {
1905     TLinkInSet linksEnd = links.end(), boundaryLink = linksEnd;
1906
1907     typedef list< pair< const QFace*, TLinkInSet > > TFaceLinkList;
1908     TFaceLinkList adjacentFaces;
1909
1910     for ( int iL = 0; iL < _sides.size(); ++iL )
1911     {
1912       if ( avoidLink._qlink == _sides[iL] )
1913         continue;
1914       TLinkInSet link = links.find( _sides[iL] );
1915       if ( link == linksEnd ) continue;
1916       if ( (*link)->MediumPos() > SMDS_TOP_FACE )
1917         continue; // We work on faces here, don't go inside a solid
1918
1919       // check link
1920       if ( link->IsBoundary() ) {
1921         if ( !nodeToContain ||
1922              (*link)->node1() == nodeToContain ||
1923              (*link)->node2() == nodeToContain )
1924         {
1925           boundaryLink = link;
1926           if ( !notBoundaryLink ) break;
1927         }
1928       }
1929       else if ( notBoundaryLink ) {
1930         *notBoundaryLink = link;
1931         if ( boundaryLink != linksEnd ) break;
1932       }
1933
1934       if ( boundaryLink == linksEnd && nodeToContain ) // collect adjacent faces
1935         if ( const QFace* adj = link->NextFace( this ))
1936           if ( adj->Contains( nodeToContain ))
1937             adjacentFaces.push_back( make_pair( adj, link ));
1938     }
1939
1940     if ( isAdjacentUsed ) *isAdjacentUsed = false;
1941     if ( boundaryLink == linksEnd && nodeToContain && nbRecursionsLeft) // check adjacent faces
1942     {
1943       if ( nbRecursionsLeft < 0 )
1944         nbRecursionsLeft = nodeToContain->NbInverseElements();
1945       TFaceLinkList::iterator adj = adjacentFaces.begin();
1946       for ( ; boundaryLink == linksEnd && adj != adjacentFaces.end(); ++adj )
1947         boundaryLink = adj->first->GetBoundaryLink( links, *(adj->second), 0, nodeToContain,
1948                                                     isAdjacentUsed, nbRecursionsLeft-1);
1949       if ( isAdjacentUsed ) *isAdjacentUsed = true;
1950     }
1951     return boundaryLink;
1952   }
1953   //================================================================================
1954   /*!
1955    * \brief Return a link ending at the given node but not avoidLink
1956    */
1957   //================================================================================
1958
1959   TLinkInSet QFace::GetLinkByNode( const TLinkSet&      links,
1960                                    const TChainLink&    avoidLink,
1961                                    const SMDS_MeshNode* nodeToContain) const
1962   {
1963     for ( int i = 0; i < _sides.size(); ++i )
1964       if ( avoidLink._qlink != _sides[i] &&
1965            (_sides[i]->node1() == nodeToContain || _sides[i]->node2() == nodeToContain ))
1966         return links.find( _sides[ i ]);
1967     return links.end();
1968   }
1969
1970   //================================================================================
1971   /*!
1972    * \brief Return normal to the i-th side pointing outside the face
1973    */
1974   //================================================================================
1975
1976   gp_Vec QFace::LinkNorm(const int i, SMESH_MesherHelper* /*uvHelper*/) const
1977   {
1978     gp_Vec norm, vecOut;
1979 //     if ( uvHelper ) {
1980 //       TopoDS_Face face = TopoDS::Face( uvHelper->GetSubShape());
1981 //       const SMDS_MeshNode* inFaceNode = uvHelper->GetNodeUVneedInFaceNode() ? GetNodeInFace() : 0;
1982 //       gp_XY uv1 = uvHelper->GetNodeUV( face, _sides[i]->node1(), inFaceNode );
1983 //       gp_XY uv2 = uvHelper->GetNodeUV( face, _sides[i]->node2(), inFaceNode );
1984 //       norm.SetCoord( uv1.Y() - uv2.Y(), uv2.X() - uv1.X(), 0 );
1985
1986 //       const QLink* otherLink = _sides[(i + 1) % _sides.size()];
1987 //       const SMDS_MeshNode* otherNode =
1988 //         otherLink->node1() == _sides[i]->node1() ? otherLink->node2() : otherLink->node1();
1989 //       gp_XY pIn = uvHelper->GetNodeUV( face, otherNode, inFaceNode );
1990 //       vecOut.SetCoord( uv1.X() - pIn.X(), uv1.Y() - pIn.Y(), 0 );
1991 //     }
1992 //     else {
1993       norm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
1994       gp_XYZ pIn = ( XYZ( _sides[0]->node1() ) +
1995                      XYZ( _sides[0]->node2() ) +
1996                      XYZ( _sides[1]->node1() )) / 3.;
1997       vecOut.SetXYZ( _sides[i]->MiddlePnt() - pIn );
1998       //}
1999     if ( norm * vecOut < 0 )
2000       norm.Reverse();
2001     double mag2 = norm.SquareMagnitude();
2002     if ( mag2 > numeric_limits<double>::min() )
2003       norm /= sqrt( mag2 );
2004     return norm;
2005   }
2006   //================================================================================
2007   /*!
2008    * \brief Move medium node of theLink according to its distance from boundary
2009    *  \param theLink - link to fix
2010    *  \param theRefVec - movement of boundary
2011    *  \param theLinks - all adjacent links of continous triangles
2012    *  \param theFaceHelper - helper is not used so far
2013    *  \param thePrevLen - distance from the boundary
2014    *  \param theStep - number of steps till movement propagation limit
2015    *  \param theLinkNorm - out normal to theLink
2016    *  \param theSign - 1 or -1 depending on movement of boundary
2017    *  \retval double - distance from boundary to propagation limit or other boundary
2018    */
2019   //================================================================================
2020
2021   double QFace::MoveByBoundary( const TChainLink&   theLink,
2022                                 const gp_Vec&       theRefVec,
2023                                 const TLinkSet&     theLinks,
2024                                 SMESH_MesherHelper* theFaceHelper,
2025                                 const double        thePrevLen,
2026                                 const int           theStep,
2027                                 gp_Vec*             theLinkNorm,
2028                                 double              theSign) const
2029   {
2030     if ( !theStep )
2031       return thePrevLen; // propagation limit reached
2032
2033     int iL; // index of theLink
2034     for ( iL = 0; iL < _sides.size(); ++iL )
2035       if ( theLink._qlink == _sides[ iL ])
2036         break;
2037
2038     MSG(string(theStep,'.')<<" Ref( "<<theRefVec.X()<<","<<theRefVec.Y()<<","<<theRefVec.Z()<<" )"
2039         <<" thePrevLen " << thePrevLen);
2040     MSG(string(theStep,'.')<<" "<<*theLink._qlink);
2041
2042     gp_Vec linkNorm = -LinkNorm( iL/*, theFaceHelper*/ ); // normal to theLink
2043     double refProj = theRefVec * linkNorm; // project movement vector to normal of theLink
2044     if ( theStep == theFirstStep )
2045       theSign = refProj < 0. ? -1. : 1.;
2046     else if ( theSign * refProj < 0.4 * theRefVec.Magnitude())
2047       return thePrevLen; // to propagate movement forward only, not in side dir or backward
2048
2049     int iL1 = (iL + 1) % 3, iL2 = (iL + 2) % 3; // indices of the two other links of triangle
2050     TLinkInSet link1 = theLinks.find( _sides[iL1] );
2051     TLinkInSet link2 = theLinks.find( _sides[iL2] );
2052     if ( link1 == theLinks.end() || link2 == theLinks.end() )
2053       return thePrevLen;
2054     const QFace* f1 = link1->NextFace( this ); // adjacent faces
2055     const QFace* f2 = link2->NextFace( this );
2056
2057     // propagate to adjacent faces till limit step or boundary
2058     double len1 = thePrevLen + (theLink->MiddlePnt() - _sides[iL1]->MiddlePnt()).Modulus();
2059     double len2 = thePrevLen + (theLink->MiddlePnt() - _sides[iL2]->MiddlePnt()).Modulus();
2060     gp_Vec linkDir1(0,0,0); // initialize to avoid valgrind error ("Conditional jump...")
2061     gp_Vec linkDir2(0,0,0);
2062     try {
2063       OCC_CATCH_SIGNALS;
2064       if ( f1 )
2065         len1 = f1->MoveByBoundary
2066           ( *link1, theRefVec, theLinks, theFaceHelper, len1, theStep-1, &linkDir1, theSign);
2067       else
2068         linkDir1 = LinkNorm( iL1/*, theFaceHelper*/ );
2069     } catch (...) {
2070       MSG( " --------------- EXCEPTION");
2071       return thePrevLen;
2072     }
2073     try {
2074       OCC_CATCH_SIGNALS;
2075       if ( f2 )
2076         len2 = f2->MoveByBoundary
2077           ( *link2, theRefVec, theLinks, theFaceHelper, len2, theStep-1, &linkDir2, theSign);
2078       else
2079         linkDir2 = LinkNorm( iL2/*, theFaceHelper*/ );
2080     } catch (...) {
2081       MSG( " --------------- EXCEPTION");
2082       return thePrevLen;
2083     }
2084
2085     double fullLen = 0;
2086     if ( theStep != theFirstStep )
2087     {
2088       // choose chain length by direction of propagation most codirected with theRefVec
2089       bool choose1 = ( theRefVec * linkDir1 * theSign > theRefVec * linkDir2 * theSign );
2090       fullLen = choose1 ? len1 : len2;
2091       double r = thePrevLen / fullLen;
2092
2093       gp_Vec move = linkNorm * refProj * ( 1 - r );
2094       theLink->Move( move, true );
2095
2096       MSG(string(theStep,'.')<<" Move "<< theLink->_mediumNode->GetID()<<
2097           " by " << refProj * ( 1 - r ) << " following " <<
2098           (choose1 ? *link1->_qlink : *link2->_qlink));
2099
2100       if ( theLinkNorm ) *theLinkNorm = linkNorm;
2101     }
2102     return fullLen;
2103   }
2104
2105   //================================================================================
2106   /*!
2107    * \brief Find pairs of continues faces 
2108    */
2109   //================================================================================
2110
2111   void QLink::SetContinuesFaces() const
2112   {
2113     //       x0         x - QLink, [-|] - QFace, v - volume
2114     //   v0  |   v1   
2115     //       |          Between _faces of link x2 two vertical faces are continues
2116     // x1----x2-----x3  and two horizontal faces are continues. We set vertical faces
2117     //       |          to _faces[0] and _faces[1] and horizontal faces to
2118     //   v2  |   v3     _faces[2] and _faces[3] (or vise versa).
2119     //       x4
2120
2121     if ( _faces.empty() )
2122       return;
2123     int iFaceCont = -1;
2124     for ( int iF = 1; iFaceCont < 0 && iF < _faces.size(); ++iF )
2125     {
2126       // look for a face bounding none of volumes bound by _faces[0]
2127       bool sameVol = false;
2128       int nbVol = _faces[iF]->NbVolumes();
2129       for ( int iV = 0; !sameVol && iV < nbVol; ++iV )
2130         sameVol = ( _faces[iF]->_volumes[iV] == _faces[0]->_volumes[0] ||
2131                     _faces[iF]->_volumes[iV] == _faces[0]->_volumes[1]);
2132       if ( !sameVol )
2133         iFaceCont = iF;
2134     }
2135     if ( iFaceCont > 0 ) // continues faces found, set one by the other
2136     {
2137       if ( iFaceCont != 1 )
2138         std::swap( _faces[1], _faces[iFaceCont] );
2139     }
2140     else if ( _faces.size() > 1 ) // not found, set NULL by the first face
2141     {
2142       _faces.insert( ++_faces.begin(), 0 );
2143     }
2144   }
2145   //================================================================================
2146   /*!
2147    * \brief Return a face continues to the given one
2148    */
2149   //================================================================================
2150
2151   const QFace* QLink::GetContinuesFace( const QFace* face ) const
2152   {
2153     for ( int i = 0; i < _faces.size(); ++i ) {
2154       if ( _faces[i] == face ) {
2155         int iF = i < 2 ? 1-i : 5-i;
2156         return iF < _faces.size() ? _faces[iF] : 0;
2157       }
2158     }
2159     return 0;
2160   }
2161   //================================================================================
2162   /*!
2163    * \brief True if link is on mesh boundary
2164    */
2165   //================================================================================
2166
2167   bool QLink::OnBoundary() const
2168   {
2169     for ( int i = 0; i < _faces.size(); ++i )
2170       if (_faces[i] && _faces[i]->IsBoundary()) return true;
2171     return false;
2172   }
2173   //================================================================================
2174   /*!
2175    * \brief Return normal of link of the chain
2176    */
2177   //================================================================================
2178
2179   gp_Vec TChainLink::Normal() const {
2180     gp_Vec norm;
2181     if (_qfaces[0]) norm  = _qfaces[0]->_normal;
2182     if (_qfaces[1]) norm += _qfaces[1]->_normal;
2183     return norm;
2184   }
2185   //================================================================================
2186   /*!
2187    * \brief Move medium nodes of vertical links of pentahedrons adjacent by side faces
2188    */
2189   //================================================================================
2190
2191   void fixPrism( TChain& allLinks )
2192   {
2193     // separate boundary links from internal ones
2194     typedef set<const QLink*/*, QLink::PtrComparator*/> QLinkSet;
2195     QLinkSet interLinks, bndLinks1, bndLink2;
2196
2197     bool isCurved = false;
2198     for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2199       if ( (*lnk)->OnBoundary() )
2200         bndLinks1.insert( lnk->_qlink );
2201       else
2202         interLinks.insert( lnk->_qlink );
2203       isCurved = isCurved || !(*lnk)->IsStraight();
2204     }
2205     if ( !isCurved )
2206       return; // no need to move
2207
2208     QLinkSet *curBndLinks = &bndLinks1, *newBndLinks = &bndLink2;
2209
2210     while ( !interLinks.empty() && !curBndLinks->empty() )
2211     {
2212       // propagate movement from boundary links to connected internal links
2213       QLinkSet::iterator bnd = curBndLinks->begin(), bndEnd = curBndLinks->end();
2214       for ( ; bnd != bndEnd; ++bnd )
2215       {
2216         const QLink* bndLink = *bnd;
2217         for ( int i = 0; i < bndLink->_faces.size(); ++i ) // loop on faces of bndLink
2218         {
2219           const QFace* face = bndLink->_faces[i]; // quadrange lateral face of a prism
2220           if ( !face ) continue;
2221           // find and move internal link opposite to bndLink within the face
2222           int interInd = ( face->LinkIndex( bndLink ) + 2 ) % face->_sides.size();
2223           const QLink* interLink = face->_sides[ interInd ];
2224           QLinkSet::iterator pInterLink = interLinks.find( interLink );
2225           if ( pInterLink == interLinks.end() ) continue; // not internal link
2226           interLink->Move( bndLink->_nodeMove );
2227           // treated internal links become new boundary ones
2228           interLinks. erase( pInterLink );
2229           newBndLinks->insert( interLink );
2230         }
2231       }
2232       curBndLinks->clear();
2233       std::swap( curBndLinks, newBndLinks );
2234     }
2235   }
2236
2237   //================================================================================
2238   /*!
2239    * \brief Fix links of continues triangles near curved boundary
2240    */
2241   //================================================================================
2242
2243   void fixTriaNearBoundary( TChain & allLinks, SMESH_MesherHelper& /*helper*/)
2244   {
2245     if ( allLinks.empty() ) return;
2246
2247     TLinkSet linkSet( allLinks.begin(), allLinks.end());
2248     TLinkInSet linkIt = linkSet.begin(), linksEnd = linkSet.end();
2249
2250     // move in 2d if we are on geom face
2251 //     TopoDS_Face face;
2252 //     TopLoc_Location loc;
2253 //     SMESH_MesherHelper faceHelper( *helper.GetMesh());
2254 //     while ( linkIt->IsBoundary()) ++linkIt;
2255 //     if ( linkIt == linksEnd ) return;
2256 //     if ( (*linkIt)->MediumPos() == SMDS_TOP_FACE ) {
2257 //       bool checkPos = true;
2258 //       TopoDS_Shape f = helper.GetSubShapeByNode( (*linkIt)->_mediumNode, helper.GetMeshDS() );
2259 //       if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) {
2260 //         face = TopoDS::Face( f );
2261 //         helper.GetNodeUV( face, (*linkIt)->_mediumNode, 0, &checkPos);
2262 //         if (checkPos)
2263 //           face.Nullify();
2264 //         else
2265 //           faceHelper.SetSubShape( face );
2266 //       }
2267 //     }
2268     for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt)
2269     {
2270       if ( linkIt->IsBoundary() && !(*linkIt)->IsStraight() && linkIt->_qfaces[0])
2271       {
2272 //         if ( !face.IsNull() ) {
2273 //           const SMDS_MeshNode* inFaceNode =
2274 //             faceHelper.GetNodeUVneedInFaceNode() ? linkIt->_qfaces[0]->GetNodeInFace() : 0;
2275 //           gp_XY uvm = helper.GetNodeUV( face, (*linkIt)->_mediumNode, inFaceNode );
2276 //           gp_XY uv1 = helper.GetNodeUV( face, (*linkIt)->node1(), inFaceNode);
2277 //           gp_XY uv2 = helper.GetNodeUV( face, (*linkIt)->node2(), inFaceNode);
2278 //           gp_XY uvMove = uvm - helper.GetMiddleUV( BRep_Tool::Surface(face,loc), uv1, uv2);
2279 //           gp_Vec move( uvMove.X(), uvMove.Y(), 0 );
2280 //           linkIt->_qfaces[0]->MoveByBoundary( *linkIt, move, linkSet, &faceHelper );
2281 //         }
2282 //         else {
2283           linkIt->_qfaces[0]->MoveByBoundary( *linkIt, (*linkIt)->_nodeMove, linkSet );
2284           //}
2285       }
2286     }
2287   }
2288
2289   //================================================================================
2290   /*!
2291    * \brief Detect rectangular structure of links and build chains from them
2292    */
2293   //================================================================================
2294
2295   enum TSplitTriaResult {
2296     _OK, _NO_CORNERS, _FEW_ROWS, _MANY_ROWS, _NO_SIDELINK, _BAD_MIDQUAD, _NOT_RECT,
2297     _NO_MIDQUAD, _NO_UPTRIA, _BAD_SET_SIZE, _BAD_CORNER, _BAD_START, _NO_BOTLINK, _TWISTED_CHAIN };
2298
2299   TSplitTriaResult splitTrianglesIntoChains( TChain &            allLinks,
2300                                              vector< TChain> &   resultChains,
2301                                              SMDS_TypeOfPosition pos )
2302   {
2303     // put links in the set and evalute number of result chains by number of boundary links
2304     TLinkSet linkSet;
2305     int nbBndLinks = 0;
2306     for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2307       linkSet.insert( *lnk );
2308       nbBndLinks += lnk->IsBoundary();
2309     }
2310     resultChains.clear();
2311     resultChains.reserve( nbBndLinks / 2 );
2312
2313     TLinkInSet linkIt, linksEnd = linkSet.end();
2314
2315     // find a boundary link with corner node; corner node has position pos-2
2316     // i.e. SMDS_TOP_VERTEX for links on faces and SMDS_TOP_EDGE for
2317     // links in volume
2318     SMDS_TypeOfPosition cornerPos = SMDS_TypeOfPosition(pos-2);
2319     const SMDS_MeshNode* corner = 0;
2320     for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt )
2321       if ( linkIt->IsBoundary() && (corner = (*linkIt)->EndPosNode(cornerPos)))
2322         break;
2323     if ( !corner)
2324       return _NO_CORNERS;
2325
2326     TLinkInSet           startLink = linkIt;
2327     const SMDS_MeshNode* startCorner = corner;
2328     vector< TChain* >    rowChains;
2329     int iCol = 0;
2330
2331     while ( startLink != linksEnd) // loop on columns
2332     {
2333       // We suppose we have a rectangular structure like shown here. We have found a
2334       //               corner of the rectangle (startCorner) and a boundary link sharing  
2335       //    |/  |/  |  the startCorner (startLink). We are going to loop on rows of the   
2336       //  --o---o---o  structure making several chains at once. One chain (columnChain)   
2337       //    |\  |  /|  starts at startLink and continues upward (we look at the structure 
2338       //  \ | \ | / |  from such point that startLink is on the bottom of the structure). 
2339       //   \|  \|/  |  While going upward we also fill horizontal chains (rowChains) we   
2340       //  --o---o---o  encounter.                                                         
2341       //   /|\  |\  |
2342       //  / | \ | \ |  startCorner
2343       //    |  \|  \|,'
2344       //  --o---o---o
2345       //          `.startLink
2346
2347       if ( resultChains.size() == nbBndLinks / 2 )
2348         return _NOT_RECT;
2349       resultChains.push_back( TChain() );
2350       TChain& columnChain = resultChains.back();
2351
2352       TLinkInSet botLink = startLink; // current horizontal link to go up from
2353       corner = startCorner; // current corner the botLink ends at
2354       int iRow = 0;
2355       while ( botLink != linksEnd ) // loop on rows
2356       {
2357         // add botLink to the columnChain
2358         columnChain.push_back( *botLink );
2359
2360         const QFace* botTria = botLink->_qfaces[0]; // bottom triangle bound by botLink
2361         if ( !botTria )
2362         { // the column ends
2363           if ( botLink == startLink )
2364             return _TWISTED_CHAIN; // issue 0020951
2365           linkSet.erase( botLink );
2366           if ( iRow != rowChains.size() )
2367             return _FEW_ROWS; // different nb of rows in columns
2368           break;
2369         }
2370         // find the link dividing the quadrangle (midQuadLink) and vertical boundary
2371         // link ending at <corner> (sideLink); there are two cases:
2372         // 1) midQuadLink does not end at <corner>, then we easily find it by botTria,
2373         //   since midQuadLink is not at boundary while sideLink is.
2374         // 2) midQuadLink ends at <corner>
2375         bool isCase2;
2376         TLinkInSet midQuadLink = linksEnd;
2377         TLinkInSet sideLink = botTria->GetBoundaryLink( linkSet, *botLink, &midQuadLink,
2378                                                         corner, &isCase2 );
2379         if ( isCase2 ) { // find midQuadLink among links of botTria
2380           midQuadLink = botTria->GetLinkByNode( linkSet, *botLink, corner );
2381           if ( midQuadLink->IsBoundary() )
2382             return _BAD_MIDQUAD;
2383         }
2384         if ( sideLink == linksEnd || midQuadLink == linksEnd || sideLink == midQuadLink )
2385           return sideLink == linksEnd ? _NO_SIDELINK : _NO_MIDQUAD;
2386
2387         // fill chains
2388         columnChain.push_back( *midQuadLink );
2389         if ( iRow >= rowChains.size() ) {
2390           if ( iCol > 0 )
2391             return _MANY_ROWS; // different nb of rows in columns
2392           if ( resultChains.size() == nbBndLinks / 2 )
2393             return _NOT_RECT;
2394           resultChains.push_back( TChain() );
2395           rowChains.push_back( & resultChains.back() );
2396         }
2397         rowChains[iRow]->push_back( *sideLink );
2398         rowChains[iRow]->push_back( *midQuadLink );
2399
2400         const QFace* upTria = midQuadLink->NextFace( botTria ); // upper tria of the rectangle
2401         if ( !upTria)
2402           return _NO_UPTRIA;
2403         if ( iRow == 0 ) {
2404           // prepare startCorner and startLink for the next column
2405           startCorner = startLink->NextNode( startCorner );
2406           if (isCase2)
2407             startLink = botTria->GetBoundaryLink( linkSet, *botLink, 0, startCorner );
2408           else
2409             startLink = upTria->GetBoundaryLink( linkSet, *midQuadLink, 0, startCorner );
2410           // check if no more columns remains
2411           if ( startLink != linksEnd ) {
2412             const SMDS_MeshNode* botNode = startLink->NextNode( startCorner );
2413             if ( (isCase2 ? botTria : upTria)->Contains( botNode ))
2414               startLink = linksEnd; // startLink bounds upTria or botTria
2415             else if ( startLink == botLink || startLink == midQuadLink || startLink == sideLink )
2416               return _BAD_START;
2417           }
2418         }
2419         // find bottom link and corner for the next row
2420         corner = sideLink->NextNode( corner );
2421         // next bottom link ends at the new corner
2422         linkSet.erase( botLink );
2423         botLink = upTria->GetLinkByNode( linkSet, (isCase2 ? *sideLink : *midQuadLink), corner );
2424         if ( botLink == linksEnd || botLink == midQuadLink || botLink == sideLink)
2425           return _NO_BOTLINK;
2426         if ( midQuadLink == startLink || sideLink == startLink )
2427           return _TWISTED_CHAIN; // issue 0020951
2428         linkSet.erase( midQuadLink );
2429         linkSet.erase( sideLink );
2430
2431         // make faces neighboring the found ones be boundary
2432         if ( startLink != linksEnd ) {
2433           const QFace* tria = isCase2 ? botTria : upTria;
2434           for ( int iL = 0; iL < 3; ++iL ) {
2435             linkIt = linkSet.find( tria->_sides[iL] );
2436             if ( linkIt != linksEnd )
2437               linkIt->RemoveFace( tria );
2438           }
2439         }
2440         if ( botLink->_qfaces[0] == upTria || botLink->_qfaces[1] == upTria )
2441           botLink->RemoveFace( upTria ); // make next botTria first in vector
2442
2443         iRow++;
2444       } // loop on rows
2445
2446       iCol++;
2447     }
2448     // In the linkSet, there must remain the last links of rowChains; add them
2449     if ( linkSet.size() != rowChains.size() )
2450       return _BAD_SET_SIZE;
2451     for ( int iRow = 0; iRow < rowChains.size(); ++iRow ) {
2452       // find the link (startLink) ending at startCorner
2453       corner = 0;
2454       for ( startLink = linkSet.begin(); startLink != linksEnd; ++startLink ) {
2455         if ( (*startLink)->node1() == startCorner ) {
2456           corner = (*startLink)->node2(); break;
2457         }
2458         else if ( (*startLink)->node2() == startCorner) {
2459           corner = (*startLink)->node1(); break;
2460         }
2461       }
2462       if ( startLink == linksEnd )
2463         return _BAD_CORNER;
2464       rowChains[ iRow ]->push_back( *startLink );
2465       linkSet.erase( startLink );
2466       startCorner = corner;
2467     }
2468
2469     return _OK;
2470   }
2471 } //namespace
2472
2473 //=======================================================================
2474 /*!
2475  * \brief Move medium nodes of faces and volumes to fix distorted elements
2476  * \param volumeOnly - to fix nodes on faces or not, if the shape is solid
2477  * 
2478  * Issue 0020307: EDF 992 SMESH : Linea/Quadratic with Medium Node on Geometry
2479  */
2480 //=======================================================================
2481
2482 void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
2483 {
2484   // 0. Apply algorithm to solids or geom faces
2485   // ----------------------------------------------
2486   if ( myShape.IsNull() ) {
2487     if ( !myMesh->HasShapeToMesh() ) return;
2488     SetSubShape( myMesh->GetShapeToMesh() );
2489
2490 #ifdef _DEBUG_
2491     int nbSolids = 0;
2492     TopTools_IndexedMapOfShape solids;
2493     TopExp::MapShapes(myShape,TopAbs_SOLID,solids);
2494     nbSolids = solids.Extent();
2495 #endif
2496     TopTools_MapOfShape faces; // faces not in solid or in not meshed solid
2497     for ( TopExp_Explorer f(myShape,TopAbs_FACE,TopAbs_SOLID); f.More(); f.Next() ) {
2498       faces.Add( f.Current() ); // not in solid
2499     }
2500     for ( TopExp_Explorer s(myShape,TopAbs_SOLID); s.More(); s.Next() ) {
2501       if ( myMesh->GetSubMesh( s.Current() )->IsEmpty() ) { // get faces of solid
2502         for ( TopExp_Explorer f( s.Current(), TopAbs_FACE); f.More(); f.Next() )
2503           faces.Add( f.Current() ); // in not meshed solid
2504       }
2505       else { // fix nodes in the solid and its faces
2506         MSG("FIX SOLID " << nbSolids-- << " #" << GetMeshDS()->ShapeToIndex(s.Current()));
2507         SMESH_MesherHelper h(*myMesh);
2508         h.SetSubShape( s.Current() );
2509         h.FixQuadraticElements(false);
2510       }
2511     }
2512     // fix nodes on geom faces
2513 #ifdef _DEBUG_
2514     int nbfaces = faces.Extent();
2515 #endif
2516     for ( TopTools_MapIteratorOfMapOfShape fIt( faces ); fIt.More(); fIt.Next() ) {
2517       MSG("FIX FACE " << nbfaces-- << " #" << GetMeshDS()->ShapeToIndex(fIt.Key()));
2518       SMESH_MesherHelper h(*myMesh);
2519       h.SetSubShape( fIt.Key() );
2520       h.FixQuadraticElements(true);
2521     }
2522     //perf_print_all_meters(1);
2523     return;
2524   }
2525
2526   // 1. Find out type of elements and get iterator on them
2527   // ---------------------------------------------------
2528
2529   SMDS_ElemIteratorPtr elemIt;
2530   SMDSAbs_ElementType elemType = SMDSAbs_All;
2531
2532   SMESH_subMesh* submesh = myMesh->GetSubMeshContaining( myShapeID );
2533   if ( !submesh )
2534     return;
2535   if ( SMESHDS_SubMesh* smDS = submesh->GetSubMeshDS() ) {
2536     elemIt = smDS->GetElements();
2537     if ( elemIt->more() ) {
2538       elemType = elemIt->next()->GetType();
2539       elemIt = smDS->GetElements();
2540     }
2541   }
2542   if ( !elemIt || !elemIt->more() || elemType < SMDSAbs_Face )
2543     return;
2544
2545   // 2. Fill in auxiliary data structures
2546   // ----------------------------------
2547
2548   set< QLink > links;
2549   set< QFace > faces;
2550   set< QLink >::iterator pLink;
2551   set< QFace >::iterator pFace;
2552
2553   bool isCurved = false;
2554   //bool hasRectFaces = false;
2555   set<int> nbElemNodeSet;
2556
2557   if ( elemType == SMDSAbs_Volume )
2558   {
2559     SMDS_VolumeTool volTool;
2560     while ( elemIt->more() ) // loop on volumes
2561     {
2562       const SMDS_MeshElement* vol = elemIt->next();
2563       if ( !vol->IsQuadratic() || !volTool.Set( vol ))
2564         return; //continue;
2565       for ( int iF = 0; iF < volTool.NbFaces(); ++iF ) // loop on faces of volume
2566       {
2567         int nbN = volTool.NbFaceNodes( iF );
2568         nbElemNodeSet.insert( nbN );
2569         const SMDS_MeshNode** faceNodes = volTool.GetFaceNodes( iF );
2570         vector< const QLink* > faceLinks( nbN/2 );
2571         for ( int iN = 0; iN < nbN; iN += 2 ) // loop on links of a face
2572         {
2573           // store QLink
2574           QLink link( faceNodes[iN], faceNodes[iN+2], faceNodes[iN+1] );
2575           pLink = links.insert( link ).first;
2576           faceLinks[ iN/2 ] = & *pLink;
2577           if ( !isCurved )
2578             isCurved = !link.IsStraight();
2579           if ( link.MediumPos() == SMDS_TOP_3DSPACE && !link.IsStraight() )
2580             return; // already fixed
2581         }
2582         // store QFace
2583         pFace = faces.insert( QFace( faceLinks )).first;
2584         if ( pFace->NbVolumes() == 0 )
2585           pFace->AddSelfToLinks();
2586         pFace->SetVolume( vol );
2587 //         hasRectFaces = hasRectFaces ||
2588 //           ( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
2589 //             volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
2590 #ifdef _DEBUG_
2591         if ( nbN == 6 )
2592           pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],faceNodes[4]);
2593         else
2594           pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],
2595                                                faceNodes[4],faceNodes[6] );
2596 #endif
2597       }
2598     }
2599     set< QLink >::iterator pLink = links.begin();
2600     for ( ; pLink != links.end(); ++pLink )
2601       pLink->SetContinuesFaces();
2602   }
2603   else
2604   {
2605     while ( elemIt->more() ) // loop on faces
2606     {
2607       const SMDS_MeshElement* face = elemIt->next();
2608       if ( !face->IsQuadratic() )
2609         continue;
2610       nbElemNodeSet.insert( face->NbNodes() );
2611       int nbN = face->NbNodes()/2;
2612       vector< const QLink* > faceLinks( nbN );
2613       for ( int iN = 0; iN < nbN; ++iN ) // loop on links of a face
2614       {
2615         // store QLink
2616         QLink link( face->GetNode(iN), face->GetNode((iN+1)%nbN), face->GetNode(iN+nbN) );
2617         pLink = links.insert( link ).first;
2618         faceLinks[ iN ] = & *pLink;
2619         if ( !isCurved )
2620           isCurved = !link.IsStraight();
2621       }
2622       // store QFace
2623       pFace = faces.insert( QFace( faceLinks )).first;
2624       pFace->AddSelfToLinks();
2625       //hasRectFaces = ( hasRectFaces || nbN == 4 );
2626     }
2627   }
2628   if ( !isCurved )
2629     return; // no curved edges of faces
2630
2631   // 3. Compute displacement of medium nodes
2632   // -------------------------------------
2633
2634   // two loops on faces: the first is to treat boundary links, the second is for internal ones
2635   TopLoc_Location loc;
2636   // not treat boundary of volumic submesh
2637   int isInside = ( elemType == SMDSAbs_Volume && volumeOnly ) ? 1 : 0;
2638   for ( ; isInside < 2; ++isInside ) {
2639     MSG( "--------------- LOOP (inside=" << isInside << ") ------------------");
2640     SMDS_TypeOfPosition pos = isInside ? SMDS_TOP_3DSPACE : SMDS_TOP_FACE;
2641     SMDS_TypeOfPosition bndPos = isInside ? SMDS_TOP_FACE : SMDS_TOP_EDGE;
2642
2643     for ( pFace = faces.begin(); pFace != faces.end(); ++pFace ) {
2644       if ( bool(isInside) == pFace->IsBoundary() )
2645         continue;
2646       for ( int dir = 0; dir < 2; ++dir ) // 2 directions of propagation from quadrangle
2647       {
2648         MSG( "CHAIN");
2649         // make chain of links connected via continues faces
2650         int error = ERR_OK;
2651         TChain rawChain;
2652         if ( !pFace->GetLinkChain( dir, rawChain, pos, error) && error ==ERR_UNKNOWN ) continue;
2653         rawChain.reverse();
2654         if ( !pFace->GetLinkChain( dir+2, rawChain, pos, error ) && error ==ERR_UNKNOWN ) continue;
2655
2656         vector< TChain > chains;
2657         if ( error == ERR_OK ) { // chain contains continues rectangles
2658           chains.resize(1);
2659           chains[0].splice( chains[0].begin(), rawChain );
2660         }
2661         else if ( error == ERR_TRI ) {  // chain contains continues triangles
2662           TSplitTriaResult res = splitTrianglesIntoChains( rawChain, chains, pos );
2663           if ( res != _OK ) { // not quadrangles split into triangles
2664             fixTriaNearBoundary( rawChain, *this );
2665             break;
2666           }
2667         }
2668         else if ( error == ERR_PRISM ) { // quadrangle side faces of prisms
2669           fixPrism( rawChain );
2670           break;
2671         }
2672         else {
2673           continue;
2674         }
2675         for ( int iC = 0; iC < chains.size(); ++iC )
2676         {
2677           TChain& chain = chains[iC];
2678           if ( chain.empty() ) continue;
2679           if ( chain.front()->IsStraight() && chain.back()->IsStraight() ) {
2680             MSG("3D straight - ignore");
2681             continue;
2682           }
2683           if ( chain.front()->MediumPos() > bndPos ||
2684                chain.back()->MediumPos() > bndPos ) {
2685             MSG("Internal chain - ignore");
2686             continue;
2687           }
2688           // mesure chain length and compute link position along the chain
2689           double chainLen = 0;
2690           vector< double > linkPos;
2691           MSGBEG( "Link medium nodes: ");
2692           TChain::iterator link0 = chain.begin(), link1 = chain.begin(), link2;
2693           for ( ++link1; link1 != chain.end(); ++link1, ++link0 ) {
2694             MSGBEG( (*link0)->_mediumNode->GetID() << "-" <<(*link1)->_mediumNode->GetID()<<" ");
2695             double len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2696             while ( len < numeric_limits<double>::min() ) { // remove degenerated link
2697               link1 = chain.erase( link1 );
2698               if ( link1 == chain.end() )
2699                 break;
2700               len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2701             }
2702             chainLen += len;
2703             linkPos.push_back( chainLen );
2704           }
2705           MSG("");
2706           if ( linkPos.size() < 2 )
2707             continue;
2708
2709           gp_Vec move0 = chain.front()->_nodeMove;
2710           gp_Vec move1 = chain.back ()->_nodeMove;
2711
2712           TopoDS_Face face;
2713           bool checkUV = true;
2714           if ( !isInside ) {
2715             // compute node displacement of end links in parametric space of face
2716             const SMDS_MeshNode* nodeOnFace = (*(++chain.begin()))->_mediumNode;
2717             TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() );
2718             if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE )
2719             {
2720               face = TopoDS::Face( f );
2721               Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
2722               bool isStraight[2];
2723               for ( int is1 = 0; is1 < 2; ++is1 ) // move0 or move1
2724               {
2725                 TChainLink& link = is1 ? chain.back() : chain.front();
2726                 gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV);
2727                 gp_XY uv1 = GetNodeUV( face, link->node1(), nodeOnFace, &checkUV);
2728                 gp_XY uv2 = GetNodeUV( face, link->node2(), nodeOnFace, &checkUV);
2729                 gp_XY uv12 = GetMiddleUV( surf, uv1, uv2);
2730                 // uvMove = uvm - uv12
2731                 gp_XY uvMove = applyIn2D(surf, uvm, uv12, gp_XY_Subtracted, /*inPeriod=*/false);
2732                 ( is1 ? move1 : move0 ).SetCoord( uvMove.X(), uvMove.Y(), 0 );
2733                 if ( !is1 ) // correct nodeOnFace for move1 (issue 0020919)
2734                   nodeOnFace = (*(++chain.rbegin()))->_mediumNode;
2735                 isStraight[is1] = isStraightLink( (uv2-uv1).SquareModulus(),uvMove.SquareModulus());
2736               }
2737 //               if ( move0.SquareMagnitude() < straightTol2 &&
2738 //                    move1.SquareMagnitude() < straightTol2 ) {
2739               if ( isStraight[0] && isStraight[1] ) {
2740                 MSG("2D straight - ignore");
2741                 continue; // straight - no need to move nodes of internal links
2742               }
2743             }
2744           }
2745           gp_Trsf trsf;
2746           if ( isInside || face.IsNull() )
2747           {
2748             // compute node displacement of end links in their local coord systems
2749             {
2750               TChainLink& ln0 = chain.front(), ln1 = *(++chain.begin());
2751               trsf.SetTransformation( gp_Ax3( gp::Origin(), ln0.Normal(),
2752                                               gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2753               move0.Transform(trsf);
2754             }
2755             {
2756               TChainLink& ln0 = *(++chain.rbegin()), ln1 = chain.back();
2757               trsf.SetTransformation( gp_Ax3( gp::Origin(), ln1.Normal(),
2758                                               gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2759               move1.Transform(trsf);
2760             }
2761           }
2762           // compute displacement of medium nodes
2763           link2 = chain.begin();
2764           link0 = link2++;
2765           link1 = link2++;
2766           for ( int i = 0; link2 != chain.end(); ++link0, ++link1, ++link2, ++i )
2767           {
2768             double r = linkPos[i] / chainLen;
2769             // displacement in local coord system
2770             gp_Vec move = (1. - r) * move0 + r * move1;
2771             if ( isInside || face.IsNull()) {
2772               // transform to global
2773               gp_Vec x01( (*link0)->MiddlePnt(), (*link1)->MiddlePnt() );
2774               gp_Vec x12( (*link1)->MiddlePnt(), (*link2)->MiddlePnt() );
2775               gp_Vec x = x01.Normalized() + x12.Normalized();
2776               trsf.SetTransformation( gp_Ax3( gp::Origin(), link1->Normal(), x), gp_Ax3() );
2777               move.Transform(trsf);
2778             }
2779             else {
2780               // compute 3D displacement by 2D one
2781               Handle(Geom_Surface) s = BRep_Tool::Surface(face,loc);
2782               gp_XY oldUV   = GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV);
2783               gp_XY newUV   = applyIn2D( s, oldUV, gp_XY( move.X(),move.Y()), gp_XY_Added);
2784               gp_Pnt newPnt = s->Value( newUV.X(), newUV.Y());
2785               move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) );
2786 #ifdef _DEBUG_
2787               if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() <
2788                    move.SquareMagnitude())
2789               {
2790                 gp_XY uv0 = GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV);
2791                 gp_XY uv2 = GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV);
2792                 MSG( "TOO LONG MOVE \t" <<
2793                      "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
2794                      "uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" <<
2795                      "uvOld: "<<oldUV.X()<<", "<<oldUV.Y()<<" \t" <<
2796                      "newUV: "<<newUV.X()<<", "<<newUV.Y()<<" \t");
2797               }
2798 #endif
2799             }
2800             (*link1)->Move( move );
2801             MSG( "Move " << (*link1)->_mediumNode->GetID() << " following "
2802                  << chain.front()->_mediumNode->GetID() <<"-"
2803                  << chain.back ()->_mediumNode->GetID() <<
2804                  " by " << move.Magnitude());
2805           }
2806         } // loop on chains of links
2807       } // loop on 2 directions of propagation from quadrangle
2808     } // loop on faces
2809   }
2810
2811   // 4. Move nodes
2812   // -----------
2813
2814   for ( pLink = links.begin(); pLink != links.end(); ++pLink ) {
2815     if ( pLink->IsMoved() ) {
2816       //gp_Pnt p = pLink->MediumPnt() + pLink->Move();
2817       gp_Pnt p = pLink->MiddlePnt() + pLink->Move();
2818       GetMeshDS()->MoveNode( pLink->_mediumNode, p.X(), p.Y(), p.Z());
2819     }
2820   }
2821 }
2822
2823 //=======================================================================
2824 /*!
2825  * \brief Iterator on ancestors of the given type
2826  */
2827 //=======================================================================
2828
2829 struct TAncestorsIterator : public SMDS_Iterator<const TopoDS_Shape*>
2830 {
2831   TopTools_ListIteratorOfListOfShape _ancIter;
2832   TopAbs_ShapeEnum                   _type;
2833   TAncestorsIterator( const TopTools_ListOfShape& ancestors, TopAbs_ShapeEnum type)
2834     : _ancIter( ancestors ), _type( type )
2835   {
2836     if ( _ancIter.More() && _ancIter.Value().ShapeType() != _type ) next();
2837   }
2838   virtual bool more()
2839   {
2840     return _ancIter.More();
2841   }
2842   virtual const TopoDS_Shape* next()
2843   {
2844     const TopoDS_Shape* s = _ancIter.More() ? & _ancIter.Value() : 0;
2845     if ( _ancIter.More() )
2846       for ( _ancIter.Next();  _ancIter.More(); _ancIter.Next())
2847         if ( _ancIter.Value().ShapeType() == _type )
2848           break;
2849     return s;
2850   }
2851 };
2852
2853 //=======================================================================
2854 /*!
2855  * \brief Return iterator on ancestors of the given type
2856  */
2857 //=======================================================================
2858
2859 PShapeIteratorPtr SMESH_MesherHelper::GetAncestors(const TopoDS_Shape& shape,
2860                                                    const SMESH_Mesh&   mesh,
2861                                                    TopAbs_ShapeEnum    ancestorType)
2862 {
2863   return PShapeIteratorPtr( new TAncestorsIterator( mesh.GetAncestors(shape), ancestorType));
2864 }