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