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