Salome HOME
0020128: EDF SMESH 926 : Quadratic conversion of BLSURF mesh
[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(), 10*MaxTolerance( 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(), 10*MaxTolerance( F ),/*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 /*!
1526  * \brief Return maximal tolerance of shape
1527  */
1528 //================================================================================
1529
1530 double SMESH_MesherHelper::MaxTolerance( const TopoDS_Shape& shape )
1531 {
1532   double tol = Precision::Confusion();
1533   TopExp_Explorer exp;
1534   for ( exp.Init( shape, TopAbs_FACE ); exp.More(); exp.Next() )
1535     tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Face( exp.Current())));
1536   for ( exp.Init( shape, TopAbs_EDGE ); exp.More(); exp.Next() )
1537     tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Edge( exp.Current())));
1538   for ( exp.Init( shape, TopAbs_VERTEX ); exp.More(); exp.Next() )
1539     tol = Max( tol, BRep_Tool::Tolerance( TopoDS::Vertex( exp.Current())));
1540
1541   return tol;
1542 }
1543
1544 //=======================================================================
1545 //function : IsQuadraticMesh
1546 //purpose  : Check mesh without geometry for: if all elements on this shape are quadratic,
1547 //           quadratic elements will be created.
1548 //           Used then generated 3D mesh without geometry.
1549 //=======================================================================
1550
1551 SMESH_MesherHelper:: MType SMESH_MesherHelper::IsQuadraticMesh()
1552 {
1553   int NbAllEdgsAndFaces=0;
1554   int NbQuadFacesAndEdgs=0;
1555   int NbFacesAndEdges=0;
1556   //All faces and edges
1557   NbAllEdgsAndFaces = myMesh->NbEdges() + myMesh->NbFaces();
1558   
1559   //Quadratic faces and edges
1560   NbQuadFacesAndEdgs = myMesh->NbEdges(ORDER_QUADRATIC) + myMesh->NbFaces(ORDER_QUADRATIC);
1561
1562   //Linear faces and edges
1563   NbFacesAndEdges = myMesh->NbEdges(ORDER_LINEAR) + myMesh->NbFaces(ORDER_LINEAR);
1564   
1565   if (NbAllEdgsAndFaces == NbQuadFacesAndEdgs) {
1566     //Quadratic mesh
1567     return SMESH_MesherHelper::QUADRATIC;
1568   }
1569   else if (NbAllEdgsAndFaces == NbFacesAndEdges) {
1570     //Linear mesh
1571     return SMESH_MesherHelper::LINEAR;
1572   }
1573   else
1574     //Mesh with both type of elements
1575     return SMESH_MesherHelper::COMP;
1576 }
1577
1578 //=======================================================================
1579 //function : GetOtherParam
1580 //purpose  : Return an alternative parameter for a node on seam
1581 //=======================================================================
1582
1583 double SMESH_MesherHelper::GetOtherParam(const double param) const
1584 {
1585   int i = myParIndex & U_periodic ? 0 : 1;
1586   return fabs(param-myPar1[i]) < fabs(param-myPar2[i]) ? myPar2[i] : myPar1[i];
1587 }
1588
1589 //#include <Perf_Meter.hxx>
1590
1591 //=======================================================================
1592 namespace { // Structures used by FixQuadraticElements()
1593 //=======================================================================
1594
1595 #define __DMP__(txt) \
1596   //cout << txt
1597 #define MSG(txt) __DMP__(txt<<endl)
1598 #define MSGBEG(txt) __DMP__(txt)
1599
1600   //const double straightTol2 = 1e-33; // to detect straing links
1601   bool isStraightLink(double linkLen2, double middleNodeMove2)
1602   {
1603     // straight if <node move> < 1/15 * <link length>
1604     return middleNodeMove2 < 1/15./15. * linkLen2;
1605   }
1606
1607   struct QFace;
1608   // ---------------------------------------
1609   /*!
1610    * \brief Quadratic link knowing its faces
1611    */
1612   struct QLink: public SMESH_TLink
1613   {
1614     const SMDS_MeshNode*          _mediumNode;
1615     mutable vector<const QFace* > _faces;
1616     mutable gp_Vec                _nodeMove;
1617     mutable int                   _nbMoves;
1618
1619     QLink(const SMDS_MeshNode* n1, const SMDS_MeshNode* n2, const SMDS_MeshNode* nm):
1620       SMESH_TLink( n1,n2 ), _mediumNode(nm), _nodeMove(0,0,0), _nbMoves(0) {
1621       _faces.reserve(4);
1622       //if ( MediumPos() != SMDS_TOP_3DSPACE )
1623         _nodeMove = MediumPnt() - MiddlePnt();
1624     }
1625     void SetContinuesFaces() const;
1626     const QFace* GetContinuesFace( const QFace* face ) const;
1627     bool OnBoundary() const;
1628     gp_XYZ MiddlePnt() const { return ( XYZ( node1() ) + XYZ( node2() )) / 2.; }
1629     gp_XYZ MediumPnt() const { return XYZ( _mediumNode ); }
1630
1631     SMDS_TypeOfPosition MediumPos() const
1632     { return _mediumNode->GetPosition()->GetTypeOfPosition(); }
1633     SMDS_TypeOfPosition EndPos(bool isSecond) const
1634     { return (isSecond ? node2() : node1())->GetPosition()->GetTypeOfPosition(); }
1635     const SMDS_MeshNode* EndPosNode(SMDS_TypeOfPosition pos) const
1636     { return EndPos(0) == pos ? node1() : EndPos(1) == pos ? node2() : 0; }
1637
1638     void Move(const gp_Vec& move, bool sum=false) const
1639     { _nodeMove += move; _nbMoves += sum ? (_nbMoves==0) : 1; }
1640     gp_XYZ Move() const { return _nodeMove.XYZ() / _nbMoves; }
1641     bool IsMoved() const { return (_nbMoves > 0 && !IsStraight()); }
1642     bool IsStraight() const
1643     { return isStraightLink( (XYZ(node1())-XYZ(node2())).SquareModulus(),
1644                              _nodeMove.SquareMagnitude());
1645     }
1646     bool operator<(const QLink& other) const {
1647       return (node1()->GetID() == other.node1()->GetID() ?
1648               node2()->GetID() < other.node2()->GetID() :
1649               node1()->GetID() < other.node1()->GetID());
1650     }
1651     struct PtrComparator {
1652       bool operator() (const QLink* l1, const QLink* l2 ) const { return *l1 < *l2; }
1653     };
1654   };
1655   // ---------------------------------------------------------
1656   /*!
1657    * \brief Link in the chain of links; it connects two faces
1658    */
1659   struct TChainLink
1660   {
1661     const QLink*         _qlink;
1662     mutable const QFace* _qfaces[2];
1663
1664     TChainLink(const QLink* qlink=0):_qlink(qlink) {
1665       _qfaces[0] = _qfaces[1] = 0;
1666     }
1667     void SetFace(const QFace* face) const { int iF = _qfaces[0] ? 1 : 0; _qfaces[iF]=face; }
1668
1669     bool IsBoundary() const { return !_qfaces[1]; }
1670
1671     void RemoveFace( const QFace* face ) const
1672     { _qfaces[(face == _qfaces[1])] = 0; if (!_qfaces[0]) std::swap(_qfaces[0],_qfaces[1]); }
1673
1674     const QFace* NextFace( const QFace* f ) const
1675     { return _qfaces[0]==f ? _qfaces[1] : _qfaces[0]; }
1676
1677     const SMDS_MeshNode* NextNode( const SMDS_MeshNode* n ) const
1678     { return n == _qlink->node1() ? _qlink->node2() : _qlink->node1(); }
1679
1680     bool operator<(const TChainLink& other) const { return *_qlink < *other._qlink; }
1681
1682     operator bool() const { return (_qlink); }
1683
1684     const QLink* operator->() const { return _qlink; }
1685
1686     gp_Vec Normal() const;
1687   };
1688   // --------------------------------------------------------------------
1689   typedef list< TChainLink > TChain;
1690   typedef set < TChainLink > TLinkSet;
1691   typedef TLinkSet::const_iterator TLinkInSet;
1692
1693   const int theFirstStep = 5;
1694
1695   enum { ERR_OK, ERR_TRI, ERR_PRISM, ERR_UNKNOWN }; // errors of QFace::GetLinkChain()
1696   // --------------------------------------------------------------------
1697   /*!
1698    * \brief Face shared by two volumes and bound by QLinks
1699    */
1700   struct QFace: public TIDSortedElemSet
1701   {
1702     mutable const SMDS_MeshElement* _volumes[2];
1703     mutable vector< const QLink* >  _sides;
1704     mutable bool                    _sideIsAdded[4]; // added in chain of links
1705     gp_Vec                          _normal;
1706 #ifdef _DEBUG_
1707     mutable const SMDS_MeshElement* _face;
1708 #endif
1709
1710     QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face=0 );
1711
1712     void SetVolume(const SMDS_MeshElement* v) const { _volumes[ _volumes[0] ? 1 : 0 ] = v; }
1713
1714     int NbVolumes() const { return !_volumes[0] ? 0 : !_volumes[1] ? 1 : 2; }
1715
1716     void AddSelfToLinks() const {
1717       for ( int i = 0; i < _sides.size(); ++i )
1718         _sides[i]->_faces.push_back( this );
1719     }
1720     int LinkIndex( const QLink* side ) const {
1721       for (int i=0; i<_sides.size(); ++i ) if ( _sides[i] == side ) return i;
1722       return -1;
1723     }
1724     bool GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& err) const;
1725
1726     bool GetLinkChain( TChainLink& link, TChain& chain, SMDS_TypeOfPosition pos, int& err) const
1727     {
1728       int i = LinkIndex( link._qlink );
1729       if ( i < 0 ) return true;
1730       _sideIsAdded[i] = true;
1731       link.SetFace( this );
1732       // continue from opposite link
1733       return GetLinkChain( (i+2)%_sides.size(), chain, pos, err );
1734     }
1735     bool IsBoundary() const { return !_volumes[1]; }
1736
1737     bool Contains( const SMDS_MeshNode* node ) const { return count(node); }
1738
1739     TLinkInSet GetBoundaryLink( const TLinkSet&      links,
1740                                 const TChainLink&    avoidLink,
1741                                 TLinkInSet *         notBoundaryLink = 0,
1742                                 const SMDS_MeshNode* nodeToContain = 0,
1743                                 bool *               isAdjacentUsed = 0,
1744                                 int                  nbRecursionsLeft = -1) const;
1745
1746     TLinkInSet GetLinkByNode( const TLinkSet&      links,
1747                               const TChainLink&    avoidLink,
1748                               const SMDS_MeshNode* nodeToContain) const;
1749
1750     const SMDS_MeshNode* GetNodeInFace() const {
1751       for ( int iL = 0; iL < _sides.size(); ++iL )
1752         if ( _sides[iL]->MediumPos() == SMDS_TOP_FACE ) return _sides[iL]->_mediumNode;
1753       return 0;
1754     }
1755
1756     gp_Vec LinkNorm(const int i, SMESH_MesherHelper* theFaceHelper=0) const;
1757
1758     double MoveByBoundary( const TChainLink&   theLink,
1759                            const gp_Vec&       theRefVec,
1760                            const TLinkSet&     theLinks,
1761                            SMESH_MesherHelper* theFaceHelper=0,
1762                            const double        thePrevLen=0,
1763                            const int           theStep=theFirstStep,
1764                            gp_Vec*             theLinkNorm=0,
1765                            double              theSign=1.0) const;
1766   };
1767
1768   //================================================================================
1769   /*!
1770    * \brief Dump QLink and QFace
1771    */
1772   ostream& operator << (ostream& out, const QLink& l)
1773   {
1774     out <<"QLink nodes: "
1775         << l.node1()->GetID() << " - "
1776         << l._mediumNode->GetID() << " - "
1777         << l.node2()->GetID() << endl;
1778     return out;
1779   }
1780   ostream& operator << (ostream& out, const QFace& f)
1781   {
1782     out <<"QFace nodes: "/*<< &f << "  "*/;
1783     for ( TIDSortedElemSet::const_iterator n = f.begin(); n != f.end(); ++n )
1784       out << (*n)->GetID() << " ";
1785     out << " \tvolumes: "
1786         << (f._volumes[0] ? f._volumes[0]->GetID() : 0) << " "
1787         << (f._volumes[1] ? f._volumes[1]->GetID() : 0);
1788     out << "  \tNormal: "<< f._normal.X() <<", "<<f._normal.Y() <<", "<<f._normal.Z() << endl;
1789     return out;
1790   }
1791
1792   //================================================================================
1793   /*!
1794    * \brief Construct QFace from QLinks 
1795    */
1796   //================================================================================
1797
1798   QFace::QFace( const vector< const QLink*>& links, const SMDS_MeshElement* face )
1799   {
1800     _volumes[0] = _volumes[1] = 0;
1801     _sides = links;
1802     _sideIsAdded[0]=_sideIsAdded[1]=_sideIsAdded[2]=_sideIsAdded[3]=false;
1803     _normal.SetCoord(0,0,0);
1804     for ( int i = 1; i < _sides.size(); ++i ) {
1805       const QLink *l1 = _sides[i-1], *l2 = _sides[i];
1806       insert( l1->node1() ); insert( l1->node2() );
1807       // compute normal
1808       gp_Vec v1( XYZ( l1->node2()), XYZ( l1->node1()));
1809       gp_Vec v2( XYZ( l2->node1()), XYZ( l2->node2()));
1810       if ( l1->node1() != l2->node1() && l1->node2() != l2->node2() )
1811         v1.Reverse(); 
1812       _normal += v1 ^ v2;
1813     }
1814     double normSqSize = _normal.SquareMagnitude();
1815     if ( normSqSize > numeric_limits<double>::min() )
1816       _normal /= sqrt( normSqSize );
1817     else
1818       _normal.SetCoord(1e-33,0,0);
1819
1820 #ifdef _DEBUG_
1821     _face = face;
1822 #endif
1823   }
1824   //================================================================================
1825   /*!
1826    * \brief Make up a chain of links
1827    *  \param iSide - link to add first
1828    *  \param chain - chain to fill in
1829    *  \param pos   - postion of medium nodes the links should have
1830    *  \param error - out, specifies what is wrong
1831    *  \retval bool - false if valid chain can't be built; "valid" means that links
1832    *                 of the chain belongs to rectangles bounding hexahedrons
1833    */
1834   //================================================================================
1835
1836   bool QFace::GetLinkChain( int iSide, TChain& chain, SMDS_TypeOfPosition pos, int& error) const
1837   {
1838     if ( iSide >= _sides.size() ) // wrong argument iSide
1839       return false;
1840     if ( _sideIsAdded[ iSide ]) // already in chain
1841       return true;
1842
1843     if ( _sides.size() != 4 ) { // triangle - visit all my continous faces
1844       MSGBEG( *this );
1845       TLinkSet links;
1846       list< const QFace* > faces( 1, this );
1847       while ( !faces.empty() ) {
1848         const QFace* face = faces.front();
1849         for ( int i = 0; i < face->_sides.size(); ++i ) {
1850           if ( !face->_sideIsAdded[i] && face->_sides[i] ) {
1851             face->_sideIsAdded[i] = true;
1852             // find a face side in the chain
1853             TLinkInSet chLink = links.insert( TChainLink(face->_sides[i])).first;
1854 //             TChain::iterator chLink = chain.begin();
1855 //             for ( ; chLink != chain.end(); ++chLink )
1856 //               if ( chLink->_qlink == face->_sides[i] )
1857 //                 break;
1858 //             if ( chLink == chain.end() )
1859 //               chLink = chain.insert( chain.begin(), TChainLink(face->_sides[i]));
1860             // add a face to a chained link and put a continues face in the queue
1861             chLink->SetFace( face );
1862             if ( face->_sides[i]->MediumPos() >= pos )
1863               if ( const QFace* contFace = face->_sides[i]->GetContinuesFace( face ))
1864                 faces.push_back( contFace );
1865           }
1866         }
1867         faces.pop_front();
1868       }
1869       if ( error < ERR_TRI )
1870         error = ERR_TRI;
1871       chain.insert( chain.end(), links.begin(),links.end() );
1872       return false;
1873     }
1874     _sideIsAdded[iSide] = true; // not to add this link to chain again
1875     const QLink* link = _sides[iSide];
1876     if ( !link)
1877       return true;
1878
1879     // add link into chain
1880     TChain::iterator chLink = chain.insert( chain.begin(), TChainLink(link));
1881     chLink->SetFace( this );
1882     MSGBEG( *this );
1883
1884     // propagate from quadrangle to neighbour faces
1885     if ( link->MediumPos() >= pos ) {
1886       int nbLinkFaces = link->_faces.size();
1887       if ( nbLinkFaces == 4 || nbLinkFaces < 4 && link->OnBoundary()) {
1888         // hexahedral mesh or boundary quadrangles - goto a continous face
1889         if ( const QFace* f = link->GetContinuesFace( this ))
1890           return f->GetLinkChain( *chLink, chain, pos, error );
1891       }
1892       else {
1893         TChainLink chLink(link); // side face of prismatic mesh - visit all faces of iSide
1894         for ( int i = 0; i < nbLinkFaces; ++i )
1895           if ( link->_faces[i] )
1896             link->_faces[i]->GetLinkChain( chLink, chain, pos, error );
1897         if ( error < ERR_PRISM )
1898           error = ERR_PRISM;
1899         return false;
1900       }
1901     }
1902     return true;
1903   }
1904
1905   //================================================================================
1906   /*!
1907    * \brief Return a boundary link of the triangle face
1908    *  \param links - set of all links
1909    *  \param avoidLink - link not to return
1910    *  \param notBoundaryLink - out, neither the returned link nor avoidLink
1911    *  \param nodeToContain - node the returned link must contain; if provided, search
1912    *                         also performed on adjacent faces
1913    *  \param isAdjacentUsed - returns true if link is found in adjacent faces
1914    *  \param nbRecursionsLeft - to limit recursion
1915    */
1916   //================================================================================
1917
1918   TLinkInSet QFace::GetBoundaryLink( const TLinkSet&      links,
1919                                      const TChainLink&    avoidLink,
1920                                      TLinkInSet *         notBoundaryLink,
1921                                      const SMDS_MeshNode* nodeToContain,
1922                                      bool *               isAdjacentUsed,
1923                                      int                  nbRecursionsLeft) const
1924   {
1925     TLinkInSet linksEnd = links.end(), boundaryLink = linksEnd;
1926
1927     typedef list< pair< const QFace*, TLinkInSet > > TFaceLinkList;
1928     TFaceLinkList adjacentFaces;
1929
1930     for ( int iL = 0; iL < _sides.size(); ++iL )
1931     {
1932       if ( avoidLink._qlink == _sides[iL] )
1933         continue;
1934       TLinkInSet link = links.find( _sides[iL] );
1935       if ( link == linksEnd ) continue;
1936       if ( (*link)->MediumPos() > SMDS_TOP_FACE )
1937         continue; // We work on faces here, don't go inside a solid
1938
1939       // check link
1940       if ( link->IsBoundary() ) {
1941         if ( !nodeToContain ||
1942              (*link)->node1() == nodeToContain ||
1943              (*link)->node2() == nodeToContain )
1944         {
1945           boundaryLink = link;
1946           if ( !notBoundaryLink ) break;
1947         }
1948       }
1949       else if ( notBoundaryLink ) {
1950         *notBoundaryLink = link;
1951         if ( boundaryLink != linksEnd ) break;
1952       }
1953
1954       if ( boundaryLink == linksEnd && nodeToContain ) // collect adjacent faces
1955         if ( const QFace* adj = link->NextFace( this ))
1956           if ( adj->Contains( nodeToContain ))
1957             adjacentFaces.push_back( make_pair( adj, link ));
1958     }
1959
1960     if ( isAdjacentUsed ) *isAdjacentUsed = false;
1961     if ( boundaryLink == linksEnd && nodeToContain && nbRecursionsLeft) // check adjacent faces
1962     {
1963       if ( nbRecursionsLeft < 0 )
1964         nbRecursionsLeft = nodeToContain->NbInverseElements();
1965       TFaceLinkList::iterator adj = adjacentFaces.begin();
1966       for ( ; boundaryLink == linksEnd && adj != adjacentFaces.end(); ++adj )
1967         boundaryLink = adj->first->GetBoundaryLink( links, *(adj->second), 0, nodeToContain,
1968                                                     isAdjacentUsed, nbRecursionsLeft-1);
1969       if ( isAdjacentUsed ) *isAdjacentUsed = true;
1970     }
1971     return boundaryLink;
1972   }
1973   //================================================================================
1974   /*!
1975    * \brief Return a link ending at the given node but not avoidLink
1976    */
1977   //================================================================================
1978
1979   TLinkInSet QFace::GetLinkByNode( const TLinkSet&      links,
1980                                    const TChainLink&    avoidLink,
1981                                    const SMDS_MeshNode* nodeToContain) const
1982   {
1983     for ( int i = 0; i < _sides.size(); ++i )
1984       if ( avoidLink._qlink != _sides[i] &&
1985            (_sides[i]->node1() == nodeToContain || _sides[i]->node2() == nodeToContain ))
1986         return links.find( _sides[ i ]);
1987     return links.end();
1988   }
1989
1990   //================================================================================
1991   /*!
1992    * \brief Return normal to the i-th side pointing outside the face
1993    */
1994   //================================================================================
1995
1996   gp_Vec QFace::LinkNorm(const int i, SMESH_MesherHelper* /*uvHelper*/) const
1997   {
1998     gp_Vec norm, vecOut;
1999 //     if ( uvHelper ) {
2000 //       TopoDS_Face face = TopoDS::Face( uvHelper->GetSubShape());
2001 //       const SMDS_MeshNode* inFaceNode = uvHelper->GetNodeUVneedInFaceNode() ? GetNodeInFace() : 0;
2002 //       gp_XY uv1 = uvHelper->GetNodeUV( face, _sides[i]->node1(), inFaceNode );
2003 //       gp_XY uv2 = uvHelper->GetNodeUV( face, _sides[i]->node2(), inFaceNode );
2004 //       norm.SetCoord( uv1.Y() - uv2.Y(), uv2.X() - uv1.X(), 0 );
2005
2006 //       const QLink* otherLink = _sides[(i + 1) % _sides.size()];
2007 //       const SMDS_MeshNode* otherNode =
2008 //         otherLink->node1() == _sides[i]->node1() ? otherLink->node2() : otherLink->node1();
2009 //       gp_XY pIn = uvHelper->GetNodeUV( face, otherNode, inFaceNode );
2010 //       vecOut.SetCoord( uv1.X() - pIn.X(), uv1.Y() - pIn.Y(), 0 );
2011 //     }
2012 //     else {
2013       norm = _normal ^ gp_Vec( XYZ(_sides[i]->node1()), XYZ(_sides[i]->node2()));
2014       gp_XYZ pIn = ( XYZ( _sides[0]->node1() ) +
2015                      XYZ( _sides[0]->node2() ) +
2016                      XYZ( _sides[1]->node1() )) / 3.;
2017       vecOut.SetXYZ( _sides[i]->MiddlePnt() - pIn );
2018       //}
2019     if ( norm * vecOut < 0 )
2020       norm.Reverse();
2021     double mag2 = norm.SquareMagnitude();
2022     if ( mag2 > numeric_limits<double>::min() )
2023       norm /= sqrt( mag2 );
2024     return norm;
2025   }
2026   //================================================================================
2027   /*!
2028    * \brief Move medium node of theLink according to its distance from boundary
2029    *  \param theLink - link to fix
2030    *  \param theRefVec - movement of boundary
2031    *  \param theLinks - all adjacent links of continous triangles
2032    *  \param theFaceHelper - helper is not used so far
2033    *  \param thePrevLen - distance from the boundary
2034    *  \param theStep - number of steps till movement propagation limit
2035    *  \param theLinkNorm - out normal to theLink
2036    *  \param theSign - 1 or -1 depending on movement of boundary
2037    *  \retval double - distance from boundary to propagation limit or other boundary
2038    */
2039   //================================================================================
2040
2041   double QFace::MoveByBoundary( const TChainLink&   theLink,
2042                                 const gp_Vec&       theRefVec,
2043                                 const TLinkSet&     theLinks,
2044                                 SMESH_MesherHelper* theFaceHelper,
2045                                 const double        thePrevLen,
2046                                 const int           theStep,
2047                                 gp_Vec*             theLinkNorm,
2048                                 double              theSign) const
2049   {
2050     if ( !theStep )
2051       return thePrevLen; // propagation limit reached
2052
2053     int iL; // index of theLink
2054     for ( iL = 0; iL < _sides.size(); ++iL )
2055       if ( theLink._qlink == _sides[ iL ])
2056         break;
2057
2058     MSG(string(theStep,'.')<<" Ref( "<<theRefVec.X()<<","<<theRefVec.Y()<<","<<theRefVec.Z()<<" )"
2059         <<" thePrevLen " << thePrevLen);
2060     MSG(string(theStep,'.')<<" "<<*theLink._qlink);
2061
2062     gp_Vec linkNorm = -LinkNorm( iL/*, theFaceHelper*/ ); // normal to theLink
2063     double refProj = theRefVec * linkNorm; // project movement vector to normal of theLink
2064     if ( theStep == theFirstStep )
2065       theSign = refProj < 0. ? -1. : 1.;
2066     else if ( theSign * refProj < 0.4 * theRefVec.Magnitude())
2067       return thePrevLen; // to propagate movement forward only, not in side dir or backward
2068
2069     int iL1 = (iL + 1) % 3, iL2 = (iL + 2) % 3; // indices of the two other links of triangle
2070     TLinkInSet link1 = theLinks.find( _sides[iL1] );
2071     TLinkInSet link2 = theLinks.find( _sides[iL2] );
2072     if ( link1 == theLinks.end() || link2 == theLinks.end() )
2073       return thePrevLen;
2074     const QFace* f1 = link1->NextFace( this ); // adjacent faces
2075     const QFace* f2 = link2->NextFace( this );
2076
2077     // propagate to adjacent faces till limit step or boundary
2078     double len1 = thePrevLen + (theLink->MiddlePnt() - _sides[iL1]->MiddlePnt()).Modulus();
2079     double len2 = thePrevLen + (theLink->MiddlePnt() - _sides[iL2]->MiddlePnt()).Modulus();
2080     gp_Vec linkDir1(0,0,0); // initialize to avoid valgrind error ("Conditional jump...")
2081     gp_Vec linkDir2(0,0,0);
2082     try {
2083       OCC_CATCH_SIGNALS;
2084       if ( f1 )
2085         len1 = f1->MoveByBoundary
2086           ( *link1, theRefVec, theLinks, theFaceHelper, len1, theStep-1, &linkDir1, theSign);
2087       else
2088         linkDir1 = LinkNorm( iL1/*, theFaceHelper*/ );
2089     } catch (...) {
2090       MSG( " --------------- EXCEPTION");
2091       return thePrevLen;
2092     }
2093     try {
2094       OCC_CATCH_SIGNALS;
2095       if ( f2 )
2096         len2 = f2->MoveByBoundary
2097           ( *link2, theRefVec, theLinks, theFaceHelper, len2, theStep-1, &linkDir2, theSign);
2098       else
2099         linkDir2 = LinkNorm( iL2/*, theFaceHelper*/ );
2100     } catch (...) {
2101       MSG( " --------------- EXCEPTION");
2102       return thePrevLen;
2103     }
2104
2105     double fullLen = 0;
2106     if ( theStep != theFirstStep )
2107     {
2108       // choose chain length by direction of propagation most codirected with theRefVec
2109       bool choose1 = ( theRefVec * linkDir1 * theSign > theRefVec * linkDir2 * theSign );
2110       fullLen = choose1 ? len1 : len2;
2111       double r = thePrevLen / fullLen;
2112
2113       gp_Vec move = linkNorm * refProj * ( 1 - r );
2114       theLink->Move( move, true );
2115
2116       MSG(string(theStep,'.')<<" Move "<< theLink->_mediumNode->GetID()<<
2117           " by " << refProj * ( 1 - r ) << " following " <<
2118           (choose1 ? *link1->_qlink : *link2->_qlink));
2119
2120       if ( theLinkNorm ) *theLinkNorm = linkNorm;
2121     }
2122     return fullLen;
2123   }
2124
2125   //================================================================================
2126   /*!
2127    * \brief Find pairs of continues faces 
2128    */
2129   //================================================================================
2130
2131   void QLink::SetContinuesFaces() const
2132   {
2133     //       x0         x - QLink, [-|] - QFace, v - volume
2134     //   v0  |   v1   
2135     //       |          Between _faces of link x2 two vertical faces are continues
2136     // x1----x2-----x3  and two horizontal faces are continues. We set vertical faces
2137     //       |          to _faces[0] and _faces[1] and horizontal faces to
2138     //   v2  |   v3     _faces[2] and _faces[3] (or vise versa).
2139     //       x4
2140
2141     if ( _faces.empty() )
2142       return;
2143     int iFaceCont = -1;
2144     for ( int iF = 1; iFaceCont < 0 && iF < _faces.size(); ++iF )
2145     {
2146       // look for a face bounding none of volumes bound by _faces[0]
2147       bool sameVol = false;
2148       int nbVol = _faces[iF]->NbVolumes();
2149       for ( int iV = 0; !sameVol && iV < nbVol; ++iV )
2150         sameVol = ( _faces[iF]->_volumes[iV] == _faces[0]->_volumes[0] ||
2151                     _faces[iF]->_volumes[iV] == _faces[0]->_volumes[1]);
2152       if ( !sameVol )
2153         iFaceCont = iF;
2154     }
2155     if ( iFaceCont > 0 ) // continues faces found, set one by the other
2156     {
2157       if ( iFaceCont != 1 )
2158         std::swap( _faces[1], _faces[iFaceCont] );
2159     }
2160     else if ( _faces.size() > 1 ) // not found, set NULL by the first face
2161     {
2162       _faces.insert( ++_faces.begin(), 0 );
2163     }
2164   }
2165   //================================================================================
2166   /*!
2167    * \brief Return a face continues to the given one
2168    */
2169   //================================================================================
2170
2171   const QFace* QLink::GetContinuesFace( const QFace* face ) const
2172   {
2173     for ( int i = 0; i < _faces.size(); ++i ) {
2174       if ( _faces[i] == face ) {
2175         int iF = i < 2 ? 1-i : 5-i;
2176         return iF < _faces.size() ? _faces[iF] : 0;
2177       }
2178     }
2179     return 0;
2180   }
2181   //================================================================================
2182   /*!
2183    * \brief True if link is on mesh boundary
2184    */
2185   //================================================================================
2186
2187   bool QLink::OnBoundary() const
2188   {
2189     for ( int i = 0; i < _faces.size(); ++i )
2190       if (_faces[i] && _faces[i]->IsBoundary()) return true;
2191     return false;
2192   }
2193   //================================================================================
2194   /*!
2195    * \brief Return normal of link of the chain
2196    */
2197   //================================================================================
2198
2199   gp_Vec TChainLink::Normal() const {
2200     gp_Vec norm;
2201     if (_qfaces[0]) norm  = _qfaces[0]->_normal;
2202     if (_qfaces[1]) norm += _qfaces[1]->_normal;
2203     return norm;
2204   }
2205   //================================================================================
2206   /*!
2207    * \brief Move medium nodes of vertical links of pentahedrons adjacent by side faces
2208    */
2209   //================================================================================
2210
2211   void fixPrism( TChain& allLinks )
2212   {
2213     // separate boundary links from internal ones
2214     typedef set<const QLink*/*, QLink::PtrComparator*/> QLinkSet;
2215     QLinkSet interLinks, bndLinks1, bndLink2;
2216
2217     bool isCurved = false;
2218     for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2219       if ( (*lnk)->OnBoundary() )
2220         bndLinks1.insert( lnk->_qlink );
2221       else
2222         interLinks.insert( lnk->_qlink );
2223       isCurved = isCurved || !(*lnk)->IsStraight();
2224     }
2225     if ( !isCurved )
2226       return; // no need to move
2227
2228     QLinkSet *curBndLinks = &bndLinks1, *newBndLinks = &bndLink2;
2229
2230     while ( !interLinks.empty() && !curBndLinks->empty() )
2231     {
2232       // propagate movement from boundary links to connected internal links
2233       QLinkSet::iterator bnd = curBndLinks->begin(), bndEnd = curBndLinks->end();
2234       for ( ; bnd != bndEnd; ++bnd )
2235       {
2236         const QLink* bndLink = *bnd;
2237         for ( int i = 0; i < bndLink->_faces.size(); ++i ) // loop on faces of bndLink
2238         {
2239           const QFace* face = bndLink->_faces[i]; // quadrange lateral face of a prism
2240           if ( !face ) continue;
2241           // find and move internal link opposite to bndLink within the face
2242           int interInd = ( face->LinkIndex( bndLink ) + 2 ) % face->_sides.size();
2243           const QLink* interLink = face->_sides[ interInd ];
2244           QLinkSet::iterator pInterLink = interLinks.find( interLink );
2245           if ( pInterLink == interLinks.end() ) continue; // not internal link
2246           interLink->Move( bndLink->_nodeMove );
2247           // treated internal links become new boundary ones
2248           interLinks. erase( pInterLink );
2249           newBndLinks->insert( interLink );
2250         }
2251       }
2252       curBndLinks->clear();
2253       std::swap( curBndLinks, newBndLinks );
2254     }
2255   }
2256
2257   //================================================================================
2258   /*!
2259    * \brief Fix links of continues triangles near curved boundary
2260    */
2261   //================================================================================
2262
2263   void fixTriaNearBoundary( TChain & allLinks, SMESH_MesherHelper& /*helper*/)
2264   {
2265     if ( allLinks.empty() ) return;
2266
2267     TLinkSet linkSet( allLinks.begin(), allLinks.end());
2268     TLinkInSet linkIt = linkSet.begin(), linksEnd = linkSet.end();
2269
2270     // move in 2d if we are on geom face
2271 //     TopoDS_Face face;
2272 //     TopLoc_Location loc;
2273 //     SMESH_MesherHelper faceHelper( *helper.GetMesh());
2274 //     while ( linkIt->IsBoundary()) ++linkIt;
2275 //     if ( linkIt == linksEnd ) return;
2276 //     if ( (*linkIt)->MediumPos() == SMDS_TOP_FACE ) {
2277 //       bool checkPos = true;
2278 //       TopoDS_Shape f = helper.GetSubShapeByNode( (*linkIt)->_mediumNode, helper.GetMeshDS() );
2279 //       if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE ) {
2280 //         face = TopoDS::Face( f );
2281 //         helper.GetNodeUV( face, (*linkIt)->_mediumNode, 0, &checkPos);
2282 //         if (checkPos)
2283 //           face.Nullify();
2284 //         else
2285 //           faceHelper.SetSubShape( face );
2286 //       }
2287 //     }
2288     for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt)
2289     {
2290       if ( linkIt->IsBoundary() && !(*linkIt)->IsStraight() && linkIt->_qfaces[0])
2291       {
2292 //         if ( !face.IsNull() ) {
2293 //           const SMDS_MeshNode* inFaceNode =
2294 //             faceHelper.GetNodeUVneedInFaceNode() ? linkIt->_qfaces[0]->GetNodeInFace() : 0;
2295 //           gp_XY uvm = helper.GetNodeUV( face, (*linkIt)->_mediumNode, inFaceNode );
2296 //           gp_XY uv1 = helper.GetNodeUV( face, (*linkIt)->node1(), inFaceNode);
2297 //           gp_XY uv2 = helper.GetNodeUV( face, (*linkIt)->node2(), inFaceNode);
2298 //           gp_XY uvMove = uvm - helper.GetMiddleUV( BRep_Tool::Surface(face,loc), uv1, uv2);
2299 //           gp_Vec move( uvMove.X(), uvMove.Y(), 0 );
2300 //           linkIt->_qfaces[0]->MoveByBoundary( *linkIt, move, linkSet, &faceHelper );
2301 //         }
2302 //         else {
2303           linkIt->_qfaces[0]->MoveByBoundary( *linkIt, (*linkIt)->_nodeMove, linkSet );
2304           //}
2305       }
2306     }
2307   }
2308
2309   //================================================================================
2310   /*!
2311    * \brief Detect rectangular structure of links and build chains from them
2312    */
2313   //================================================================================
2314
2315   enum TSplitTriaResult {
2316     _OK, _NO_CORNERS, _FEW_ROWS, _MANY_ROWS, _NO_SIDELINK, _BAD_MIDQUAD, _NOT_RECT,
2317     _NO_MIDQUAD, _NO_UPTRIA, _BAD_SET_SIZE, _BAD_CORNER, _BAD_START, _NO_BOTLINK, _TWISTED_CHAIN };
2318
2319   TSplitTriaResult splitTrianglesIntoChains( TChain &            allLinks,
2320                                              vector< TChain> &   resultChains,
2321                                              SMDS_TypeOfPosition pos )
2322   {
2323     // put links in the set and evalute number of result chains by number of boundary links
2324     TLinkSet linkSet;
2325     int nbBndLinks = 0;
2326     for ( TChain::iterator lnk = allLinks.begin(); lnk != allLinks.end(); ++lnk ) {
2327       linkSet.insert( *lnk );
2328       nbBndLinks += lnk->IsBoundary();
2329     }
2330     resultChains.clear();
2331     resultChains.reserve( nbBndLinks / 2 );
2332
2333     TLinkInSet linkIt, linksEnd = linkSet.end();
2334
2335     // find a boundary link with corner node; corner node has position pos-2
2336     // i.e. SMDS_TOP_VERTEX for links on faces and SMDS_TOP_EDGE for
2337     // links in volume
2338     SMDS_TypeOfPosition cornerPos = SMDS_TypeOfPosition(pos-2);
2339     const SMDS_MeshNode* corner = 0;
2340     for ( linkIt = linkSet.begin(); linkIt != linksEnd; ++linkIt )
2341       if ( linkIt->IsBoundary() && (corner = (*linkIt)->EndPosNode(cornerPos)))
2342         break;
2343     if ( !corner)
2344       return _NO_CORNERS;
2345
2346     TLinkInSet           startLink = linkIt;
2347     const SMDS_MeshNode* startCorner = corner;
2348     vector< TChain* >    rowChains;
2349     int iCol = 0;
2350
2351     while ( startLink != linksEnd) // loop on columns
2352     {
2353       // We suppose we have a rectangular structure like shown here. We have found a
2354       //               corner of the rectangle (startCorner) and a boundary link sharing  
2355       //    |/  |/  |  the startCorner (startLink). We are going to loop on rows of the   
2356       //  --o---o---o  structure making several chains at once. One chain (columnChain)   
2357       //    |\  |  /|  starts at startLink and continues upward (we look at the structure 
2358       //  \ | \ | / |  from such point that startLink is on the bottom of the structure). 
2359       //   \|  \|/  |  While going upward we also fill horizontal chains (rowChains) we   
2360       //  --o---o---o  encounter.                                                         
2361       //   /|\  |\  |
2362       //  / | \ | \ |  startCorner
2363       //    |  \|  \|,'
2364       //  --o---o---o
2365       //          `.startLink
2366
2367       if ( resultChains.size() == nbBndLinks / 2 )
2368         return _NOT_RECT;
2369       resultChains.push_back( TChain() );
2370       TChain& columnChain = resultChains.back();
2371
2372       TLinkInSet botLink = startLink; // current horizontal link to go up from
2373       corner = startCorner; // current corner the botLink ends at
2374       int iRow = 0;
2375       while ( botLink != linksEnd ) // loop on rows
2376       {
2377         // add botLink to the columnChain
2378         columnChain.push_back( *botLink );
2379
2380         const QFace* botTria = botLink->_qfaces[0]; // bottom triangle bound by botLink
2381         if ( !botTria )
2382         { // the column ends
2383           if ( botLink == startLink )
2384             return _TWISTED_CHAIN; // issue 0020951
2385           linkSet.erase( botLink );
2386           if ( iRow != rowChains.size() )
2387             return _FEW_ROWS; // different nb of rows in columns
2388           break;
2389         }
2390         // find the link dividing the quadrangle (midQuadLink) and vertical boundary
2391         // link ending at <corner> (sideLink); there are two cases:
2392         // 1) midQuadLink does not end at <corner>, then we easily find it by botTria,
2393         //   since midQuadLink is not at boundary while sideLink is.
2394         // 2) midQuadLink ends at <corner>
2395         bool isCase2;
2396         TLinkInSet midQuadLink = linksEnd;
2397         TLinkInSet sideLink = botTria->GetBoundaryLink( linkSet, *botLink, &midQuadLink,
2398                                                         corner, &isCase2 );
2399         if ( isCase2 ) { // find midQuadLink among links of botTria
2400           midQuadLink = botTria->GetLinkByNode( linkSet, *botLink, corner );
2401           if ( midQuadLink->IsBoundary() )
2402             return _BAD_MIDQUAD;
2403         }
2404         if ( sideLink == linksEnd || midQuadLink == linksEnd || sideLink == midQuadLink )
2405           return sideLink == linksEnd ? _NO_SIDELINK : _NO_MIDQUAD;
2406
2407         // fill chains
2408         columnChain.push_back( *midQuadLink );
2409         if ( iRow >= rowChains.size() ) {
2410           if ( iCol > 0 )
2411             return _MANY_ROWS; // different nb of rows in columns
2412           if ( resultChains.size() == nbBndLinks / 2 )
2413             return _NOT_RECT;
2414           resultChains.push_back( TChain() );
2415           rowChains.push_back( & resultChains.back() );
2416         }
2417         rowChains[iRow]->push_back( *sideLink );
2418         rowChains[iRow]->push_back( *midQuadLink );
2419
2420         const QFace* upTria = midQuadLink->NextFace( botTria ); // upper tria of the rectangle
2421         if ( !upTria)
2422           return _NO_UPTRIA;
2423         if ( iRow == 0 ) {
2424           // prepare startCorner and startLink for the next column
2425           startCorner = startLink->NextNode( startCorner );
2426           if (isCase2)
2427             startLink = botTria->GetBoundaryLink( linkSet, *botLink, 0, startCorner );
2428           else
2429             startLink = upTria->GetBoundaryLink( linkSet, *midQuadLink, 0, startCorner );
2430           // check if no more columns remains
2431           if ( startLink != linksEnd ) {
2432             const SMDS_MeshNode* botNode = startLink->NextNode( startCorner );
2433             if ( (isCase2 ? botTria : upTria)->Contains( botNode ))
2434               startLink = linksEnd; // startLink bounds upTria or botTria
2435             else if ( startLink == botLink || startLink == midQuadLink || startLink == sideLink )
2436               return _BAD_START;
2437           }
2438         }
2439         // find bottom link and corner for the next row
2440         corner = sideLink->NextNode( corner );
2441         // next bottom link ends at the new corner
2442         linkSet.erase( botLink );
2443         botLink = upTria->GetLinkByNode( linkSet, (isCase2 ? *sideLink : *midQuadLink), corner );
2444         if ( botLink == linksEnd || botLink == midQuadLink || botLink == sideLink)
2445           return _NO_BOTLINK;
2446         if ( midQuadLink == startLink || sideLink == startLink )
2447           return _TWISTED_CHAIN; // issue 0020951
2448         linkSet.erase( midQuadLink );
2449         linkSet.erase( sideLink );
2450
2451         // make faces neighboring the found ones be boundary
2452         if ( startLink != linksEnd ) {
2453           const QFace* tria = isCase2 ? botTria : upTria;
2454           for ( int iL = 0; iL < 3; ++iL ) {
2455             linkIt = linkSet.find( tria->_sides[iL] );
2456             if ( linkIt != linksEnd )
2457               linkIt->RemoveFace( tria );
2458           }
2459         }
2460         if ( botLink->_qfaces[0] == upTria || botLink->_qfaces[1] == upTria )
2461           botLink->RemoveFace( upTria ); // make next botTria first in vector
2462
2463         iRow++;
2464       } // loop on rows
2465
2466       iCol++;
2467     }
2468     // In the linkSet, there must remain the last links of rowChains; add them
2469     if ( linkSet.size() != rowChains.size() )
2470       return _BAD_SET_SIZE;
2471     for ( int iRow = 0; iRow < rowChains.size(); ++iRow ) {
2472       // find the link (startLink) ending at startCorner
2473       corner = 0;
2474       for ( startLink = linkSet.begin(); startLink != linksEnd; ++startLink ) {
2475         if ( (*startLink)->node1() == startCorner ) {
2476           corner = (*startLink)->node2(); break;
2477         }
2478         else if ( (*startLink)->node2() == startCorner) {
2479           corner = (*startLink)->node1(); break;
2480         }
2481       }
2482       if ( startLink == linksEnd )
2483         return _BAD_CORNER;
2484       rowChains[ iRow ]->push_back( *startLink );
2485       linkSet.erase( startLink );
2486       startCorner = corner;
2487     }
2488
2489     return _OK;
2490   }
2491 } //namespace
2492
2493 //=======================================================================
2494 /*!
2495  * \brief Move medium nodes of faces and volumes to fix distorted elements
2496  * \param volumeOnly - to fix nodes on faces or not, if the shape is solid
2497  * 
2498  * Issue 0020307: EDF 992 SMESH : Linea/Quadratic with Medium Node on Geometry
2499  */
2500 //=======================================================================
2501
2502 void SMESH_MesherHelper::FixQuadraticElements(bool volumeOnly)
2503 {
2504   // 0. Apply algorithm to solids or geom faces
2505   // ----------------------------------------------
2506   if ( myShape.IsNull() ) {
2507     if ( !myMesh->HasShapeToMesh() ) return;
2508     SetSubShape( myMesh->GetShapeToMesh() );
2509
2510 #ifdef _DEBUG_
2511     int nbSolids = 0;
2512     TopTools_IndexedMapOfShape solids;
2513     TopExp::MapShapes(myShape,TopAbs_SOLID,solids);
2514     nbSolids = solids.Extent();
2515 #endif
2516     TopTools_MapOfShape faces; // faces not in solid or in not meshed solid
2517     for ( TopExp_Explorer f(myShape,TopAbs_FACE,TopAbs_SOLID); f.More(); f.Next() ) {
2518       faces.Add( f.Current() ); // not in solid
2519     }
2520     for ( TopExp_Explorer s(myShape,TopAbs_SOLID); s.More(); s.Next() ) {
2521       if ( myMesh->GetSubMesh( s.Current() )->IsEmpty() ) { // get faces of solid
2522         for ( TopExp_Explorer f( s.Current(), TopAbs_FACE); f.More(); f.Next() )
2523           faces.Add( f.Current() ); // in not meshed solid
2524       }
2525       else { // fix nodes in the solid and its faces
2526         MSG("FIX SOLID " << nbSolids-- << " #" << GetMeshDS()->ShapeToIndex(s.Current()));
2527         SMESH_MesherHelper h(*myMesh);
2528         h.SetSubShape( s.Current() );
2529         h.FixQuadraticElements(false);
2530       }
2531     }
2532     // fix nodes on geom faces
2533 #ifdef _DEBUG_
2534     int nbfaces = faces.Extent();
2535 #endif
2536     for ( TopTools_MapIteratorOfMapOfShape fIt( faces ); fIt.More(); fIt.Next() ) {
2537       MSG("FIX FACE " << nbfaces-- << " #" << GetMeshDS()->ShapeToIndex(fIt.Key()));
2538       SMESH_MesherHelper h(*myMesh);
2539       h.SetSubShape( fIt.Key() );
2540       h.FixQuadraticElements(true);
2541     }
2542     //perf_print_all_meters(1);
2543     return;
2544   }
2545
2546   // 1. Find out type of elements and get iterator on them
2547   // ---------------------------------------------------
2548
2549   SMDS_ElemIteratorPtr elemIt;
2550   SMDSAbs_ElementType elemType = SMDSAbs_All;
2551
2552   SMESH_subMesh* submesh = myMesh->GetSubMeshContaining( myShapeID );
2553   if ( !submesh )
2554     return;
2555   if ( SMESHDS_SubMesh* smDS = submesh->GetSubMeshDS() ) {
2556     elemIt = smDS->GetElements();
2557     if ( elemIt->more() ) {
2558       elemType = elemIt->next()->GetType();
2559       elemIt = smDS->GetElements();
2560     }
2561   }
2562   if ( !elemIt || !elemIt->more() || elemType < SMDSAbs_Face )
2563     return;
2564
2565   // 2. Fill in auxiliary data structures
2566   // ----------------------------------
2567
2568   set< QLink > links;
2569   set< QFace > faces;
2570   set< QLink >::iterator pLink;
2571   set< QFace >::iterator pFace;
2572
2573   bool isCurved = false;
2574   //bool hasRectFaces = false;
2575   set<int> nbElemNodeSet;
2576
2577   if ( elemType == SMDSAbs_Volume )
2578   {
2579     SMDS_VolumeTool volTool;
2580     while ( elemIt->more() ) // loop on volumes
2581     {
2582       const SMDS_MeshElement* vol = elemIt->next();
2583       if ( !vol->IsQuadratic() || !volTool.Set( vol ))
2584         return; //continue;
2585       for ( int iF = 0; iF < volTool.NbFaces(); ++iF ) // loop on faces of volume
2586       {
2587         int nbN = volTool.NbFaceNodes( iF );
2588         nbElemNodeSet.insert( nbN );
2589         const SMDS_MeshNode** faceNodes = volTool.GetFaceNodes( iF );
2590         vector< const QLink* > faceLinks( nbN/2 );
2591         for ( int iN = 0; iN < nbN; iN += 2 ) // loop on links of a face
2592         {
2593           // store QLink
2594           QLink link( faceNodes[iN], faceNodes[iN+2], faceNodes[iN+1] );
2595           pLink = links.insert( link ).first;
2596           faceLinks[ iN/2 ] = & *pLink;
2597           if ( !isCurved )
2598             isCurved = !link.IsStraight();
2599           if ( link.MediumPos() == SMDS_TOP_3DSPACE && !link.IsStraight() )
2600             return; // already fixed
2601         }
2602         // store QFace
2603         pFace = faces.insert( QFace( faceLinks )).first;
2604         if ( pFace->NbVolumes() == 0 )
2605           pFace->AddSelfToLinks();
2606         pFace->SetVolume( vol );
2607 //         hasRectFaces = hasRectFaces ||
2608 //           ( volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_HEXA ||
2609 //             volTool.GetVolumeType() == SMDS_VolumeTool::QUAD_PENTA );
2610 #ifdef _DEBUG_
2611         if ( nbN == 6 )
2612           pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],faceNodes[4]);
2613         else
2614           pFace->_face = GetMeshDS()->FindFace(faceNodes[0],faceNodes[2],
2615                                                faceNodes[4],faceNodes[6] );
2616 #endif
2617       }
2618     }
2619     set< QLink >::iterator pLink = links.begin();
2620     for ( ; pLink != links.end(); ++pLink )
2621       pLink->SetContinuesFaces();
2622   }
2623   else
2624   {
2625     while ( elemIt->more() ) // loop on faces
2626     {
2627       const SMDS_MeshElement* face = elemIt->next();
2628       if ( !face->IsQuadratic() )
2629         continue;
2630       nbElemNodeSet.insert( face->NbNodes() );
2631       int nbN = face->NbNodes()/2;
2632       vector< const QLink* > faceLinks( nbN );
2633       for ( int iN = 0; iN < nbN; ++iN ) // loop on links of a face
2634       {
2635         // store QLink
2636         QLink link( face->GetNode(iN), face->GetNode((iN+1)%nbN), face->GetNode(iN+nbN) );
2637         pLink = links.insert( link ).first;
2638         faceLinks[ iN ] = & *pLink;
2639         if ( !isCurved )
2640           isCurved = !link.IsStraight();
2641       }
2642       // store QFace
2643       pFace = faces.insert( QFace( faceLinks )).first;
2644       pFace->AddSelfToLinks();
2645       //hasRectFaces = ( hasRectFaces || nbN == 4 );
2646     }
2647   }
2648   if ( !isCurved )
2649     return; // no curved edges of faces
2650
2651   // 3. Compute displacement of medium nodes
2652   // -------------------------------------
2653
2654   // two loops on faces: the first is to treat boundary links, the second is for internal ones
2655   TopLoc_Location loc;
2656   // not treat boundary of volumic submesh
2657   int isInside = ( elemType == SMDSAbs_Volume && volumeOnly ) ? 1 : 0;
2658   for ( ; isInside < 2; ++isInside ) {
2659     MSG( "--------------- LOOP (inside=" << isInside << ") ------------------");
2660     SMDS_TypeOfPosition pos = isInside ? SMDS_TOP_3DSPACE : SMDS_TOP_FACE;
2661     SMDS_TypeOfPosition bndPos = isInside ? SMDS_TOP_FACE : SMDS_TOP_EDGE;
2662
2663     for ( pFace = faces.begin(); pFace != faces.end(); ++pFace ) {
2664       if ( bool(isInside) == pFace->IsBoundary() )
2665         continue;
2666       for ( int dir = 0; dir < 2; ++dir ) // 2 directions of propagation from quadrangle
2667       {
2668         MSG( "CHAIN");
2669         // make chain of links connected via continues faces
2670         int error = ERR_OK;
2671         TChain rawChain;
2672         if ( !pFace->GetLinkChain( dir, rawChain, pos, error) && error ==ERR_UNKNOWN ) continue;
2673         rawChain.reverse();
2674         if ( !pFace->GetLinkChain( dir+2, rawChain, pos, error ) && error ==ERR_UNKNOWN ) continue;
2675
2676         vector< TChain > chains;
2677         if ( error == ERR_OK ) { // chain contains continues rectangles
2678           chains.resize(1);
2679           chains[0].splice( chains[0].begin(), rawChain );
2680         }
2681         else if ( error == ERR_TRI ) {  // chain contains continues triangles
2682           TSplitTriaResult res = splitTrianglesIntoChains( rawChain, chains, pos );
2683           if ( res != _OK ) { // not quadrangles split into triangles
2684             fixTriaNearBoundary( rawChain, *this );
2685             break;
2686           }
2687         }
2688         else if ( error == ERR_PRISM ) { // quadrangle side faces of prisms
2689           fixPrism( rawChain );
2690           break;
2691         }
2692         else {
2693           continue;
2694         }
2695         for ( int iC = 0; iC < chains.size(); ++iC )
2696         {
2697           TChain& chain = chains[iC];
2698           if ( chain.empty() ) continue;
2699           if ( chain.front()->IsStraight() && chain.back()->IsStraight() ) {
2700             MSG("3D straight - ignore");
2701             continue;
2702           }
2703           if ( chain.front()->MediumPos() > bndPos ||
2704                chain.back()->MediumPos() > bndPos ) {
2705             MSG("Internal chain - ignore");
2706             continue;
2707           }
2708           // mesure chain length and compute link position along the chain
2709           double chainLen = 0;
2710           vector< double > linkPos;
2711           MSGBEG( "Link medium nodes: ");
2712           TChain::iterator link0 = chain.begin(), link1 = chain.begin(), link2;
2713           for ( ++link1; link1 != chain.end(); ++link1, ++link0 ) {
2714             MSGBEG( (*link0)->_mediumNode->GetID() << "-" <<(*link1)->_mediumNode->GetID()<<" ");
2715             double len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2716             while ( len < numeric_limits<double>::min() ) { // remove degenerated link
2717               link1 = chain.erase( link1 );
2718               if ( link1 == chain.end() )
2719                 break;
2720               len = ((*link0)->MiddlePnt() - (*link1)->MiddlePnt()).Modulus();
2721             }
2722             chainLen += len;
2723             linkPos.push_back( chainLen );
2724           }
2725           MSG("");
2726           if ( linkPos.size() < 2 )
2727             continue;
2728
2729           gp_Vec move0 = chain.front()->_nodeMove;
2730           gp_Vec move1 = chain.back ()->_nodeMove;
2731
2732           TopoDS_Face face;
2733           bool checkUV = true;
2734           if ( !isInside ) {
2735             // compute node displacement of end links in parametric space of face
2736             const SMDS_MeshNode* nodeOnFace = (*(++chain.begin()))->_mediumNode;
2737             TopoDS_Shape f = GetSubShapeByNode( nodeOnFace, GetMeshDS() );
2738             if ( !f.IsNull() && f.ShapeType() == TopAbs_FACE )
2739             {
2740               face = TopoDS::Face( f );
2741               Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
2742               bool isStraight[2];
2743               for ( int is1 = 0; is1 < 2; ++is1 ) // move0 or move1
2744               {
2745                 TChainLink& link = is1 ? chain.back() : chain.front();
2746                 gp_XY uvm = GetNodeUV( face, link->_mediumNode, nodeOnFace, &checkUV);
2747                 gp_XY uv1 = GetNodeUV( face, link->node1(), nodeOnFace, &checkUV);
2748                 gp_XY uv2 = GetNodeUV( face, link->node2(), nodeOnFace, &checkUV);
2749                 gp_XY uv12 = GetMiddleUV( surf, uv1, uv2);
2750                 // uvMove = uvm - uv12
2751                 gp_XY uvMove = applyIn2D(surf, uvm, uv12, gp_XY_Subtracted, /*inPeriod=*/false);
2752                 ( is1 ? move1 : move0 ).SetCoord( uvMove.X(), uvMove.Y(), 0 );
2753                 if ( !is1 ) // correct nodeOnFace for move1 (issue 0020919)
2754                   nodeOnFace = (*(++chain.rbegin()))->_mediumNode;
2755                 isStraight[is1] = isStraightLink( (uv2-uv1).SquareModulus(),uvMove.SquareModulus());
2756               }
2757 //               if ( move0.SquareMagnitude() < straightTol2 &&
2758 //                    move1.SquareMagnitude() < straightTol2 ) {
2759               if ( isStraight[0] && isStraight[1] ) {
2760                 MSG("2D straight - ignore");
2761                 continue; // straight - no need to move nodes of internal links
2762               }
2763             }
2764           }
2765           gp_Trsf trsf;
2766           if ( isInside || face.IsNull() )
2767           {
2768             // compute node displacement of end links in their local coord systems
2769             {
2770               TChainLink& ln0 = chain.front(), ln1 = *(++chain.begin());
2771               trsf.SetTransformation( gp_Ax3( gp::Origin(), ln0.Normal(),
2772                                               gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2773               move0.Transform(trsf);
2774             }
2775             {
2776               TChainLink& ln0 = *(++chain.rbegin()), ln1 = chain.back();
2777               trsf.SetTransformation( gp_Ax3( gp::Origin(), ln1.Normal(),
2778                                               gp_Vec( ln0->MiddlePnt(), ln1->MiddlePnt() )));
2779               move1.Transform(trsf);
2780             }
2781           }
2782           // compute displacement of medium nodes
2783           link2 = chain.begin();
2784           link0 = link2++;
2785           link1 = link2++;
2786           for ( int i = 0; link2 != chain.end(); ++link0, ++link1, ++link2, ++i )
2787           {
2788             double r = linkPos[i] / chainLen;
2789             // displacement in local coord system
2790             gp_Vec move = (1. - r) * move0 + r * move1;
2791             if ( isInside || face.IsNull()) {
2792               // transform to global
2793               gp_Vec x01( (*link0)->MiddlePnt(), (*link1)->MiddlePnt() );
2794               gp_Vec x12( (*link1)->MiddlePnt(), (*link2)->MiddlePnt() );
2795               gp_Vec x = x01.Normalized() + x12.Normalized();
2796               trsf.SetTransformation( gp_Ax3( gp::Origin(), link1->Normal(), x), gp_Ax3() );
2797               move.Transform(trsf);
2798             }
2799             else {
2800               // compute 3D displacement by 2D one
2801               Handle(Geom_Surface) s = BRep_Tool::Surface(face,loc);
2802               gp_XY oldUV   = GetNodeUV( face, (*link1)->_mediumNode, 0, &checkUV);
2803               gp_XY newUV   = applyIn2D( s, oldUV, gp_XY( move.X(),move.Y()), gp_XY_Added);
2804               gp_Pnt newPnt = s->Value( newUV.X(), newUV.Y());
2805               move = gp_Vec( XYZ((*link1)->_mediumNode), newPnt.Transformed(loc) );
2806 #ifdef _DEBUG_
2807               if ( (XYZ((*link1)->node1()) - XYZ((*link1)->node2())).SquareModulus() <
2808                    move.SquareMagnitude())
2809               {
2810                 gp_XY uv0 = GetNodeUV( face, (*link0)->_mediumNode, 0, &checkUV);
2811                 gp_XY uv2 = GetNodeUV( face, (*link2)->_mediumNode, 0, &checkUV);
2812                 MSG( "TOO LONG MOVE \t" <<
2813                      "uv0: "<<uv0.X()<<", "<<uv0.Y()<<" \t" <<
2814                      "uv2: "<<uv2.X()<<", "<<uv2.Y()<<" \t" <<
2815                      "uvOld: "<<oldUV.X()<<", "<<oldUV.Y()<<" \t" <<
2816                      "newUV: "<<newUV.X()<<", "<<newUV.Y()<<" \t");
2817               }
2818 #endif
2819             }
2820             (*link1)->Move( move );
2821             MSG( "Move " << (*link1)->_mediumNode->GetID() << " following "
2822                  << chain.front()->_mediumNode->GetID() <<"-"
2823                  << chain.back ()->_mediumNode->GetID() <<
2824                  " by " << move.Magnitude());
2825           }
2826         } // loop on chains of links
2827       } // loop on 2 directions of propagation from quadrangle
2828     } // loop on faces
2829   }
2830
2831   // 4. Move nodes
2832   // -----------
2833
2834   for ( pLink = links.begin(); pLink != links.end(); ++pLink ) {
2835     if ( pLink->IsMoved() ) {
2836       //gp_Pnt p = pLink->MediumPnt() + pLink->Move();
2837       gp_Pnt p = pLink->MiddlePnt() + pLink->Move();
2838       GetMeshDS()->MoveNode( pLink->_mediumNode, p.X(), p.Y(), p.Z());
2839     }
2840   }
2841 }
2842
2843 //=======================================================================
2844 /*!
2845  * \brief Iterator on ancestors of the given type
2846  */
2847 //=======================================================================
2848
2849 struct TAncestorsIterator : public SMDS_Iterator<const TopoDS_Shape*>
2850 {
2851   TopTools_ListIteratorOfListOfShape _ancIter;
2852   TopAbs_ShapeEnum                   _type;
2853   TAncestorsIterator( const TopTools_ListOfShape& ancestors, TopAbs_ShapeEnum type)
2854     : _ancIter( ancestors ), _type( type )
2855   {
2856     if ( _ancIter.More() && _ancIter.Value().ShapeType() != _type ) next();
2857   }
2858   virtual bool more()
2859   {
2860     return _ancIter.More();
2861   }
2862   virtual const TopoDS_Shape* next()
2863   {
2864     const TopoDS_Shape* s = _ancIter.More() ? & _ancIter.Value() : 0;
2865     if ( _ancIter.More() )
2866       for ( _ancIter.Next();  _ancIter.More(); _ancIter.Next())
2867         if ( _ancIter.Value().ShapeType() == _type )
2868           break;
2869     return s;
2870   }
2871 };
2872
2873 //=======================================================================
2874 /*!
2875  * \brief Return iterator on ancestors of the given type
2876  */
2877 //=======================================================================
2878
2879 PShapeIteratorPtr SMESH_MesherHelper::GetAncestors(const TopoDS_Shape& shape,
2880                                                    const SMESH_Mesh&   mesh,
2881                                                    TopAbs_ShapeEnum    ancestorType)
2882 {
2883   return PShapeIteratorPtr( new TAncestorsIterator( mesh.GetAncestors(shape), ancestorType));
2884 }