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