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