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