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