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