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