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