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