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