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