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