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