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