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