Salome HOME
52863: Convert To Bi-Quadratic incorrectly locate in-face central nodes
[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       {
654         if ( !IsSubShape( V, F ))
655         {
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         {
674           uvOK = false;
675           TopTools_ListIteratorOfListOfShape it( myMesh->GetAncestors( V ));
676           for ( ; it.More(); it.Next() ) {
677             if ( it.Value().ShapeType() == TopAbs_EDGE ) {
678               const TopoDS_Edge & edge = TopoDS::Edge( it.Value() );
679               double f,l;
680               Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(edge, F, f, l);
681               if ( !C2d.IsNull() ) {
682                 double u = ( V == IthVertex( 0, edge )) ?  f : l;
683                 uv = C2d->Value( u );
684                 uvOK = true;
685                 break;
686               }
687             }
688           }
689           if ( !uvOK && V.Orientation() == TopAbs_INTERNAL )
690           {
691             Handle(ShapeAnalysis_Surface) projector = GetSurface( F );
692             if ( n2 ) uv = GetNodeUV( F, n2 );
693             if ( Precision::IsInfinite( uv.X() ))
694               uv = projector->NextValueOfUV( uv, BRep_Tool::Pnt( V ), BRep_Tool::Tolerance( F ));
695             else
696               uv = projector->ValueOfUV( BRep_Tool::Pnt( V ), BRep_Tool::Tolerance( F ));
697             uvOK = ( projector->Gap() < getFaceMaxTol( F ));
698           }
699         }
700       }
701       if ( n2 && IsSeamShape( vertexID ))
702       {
703         bool isSeam = ( myShape.IsSame( F ));
704         if ( !isSeam ) {
705           SMESH_MesherHelper h( *myMesh );
706           h.SetSubShape( F );
707           isSeam = IsSeamShape( vertexID );
708         }
709
710         if ( isSeam )
711           uv = getUVOnSeam( uv, GetNodeUV( F, n2, 0 ));
712       }
713     }
714   }
715   else
716   {
717     uvOK = CheckNodeUV( F, n, uv.ChangeCoord(), 2.*getFaceMaxTol( F ));
718   }
719
720   if ( check && !uvOK )
721     *check = uvOK;
722
723   return uv.XY();
724 }
725
726 //=======================================================================
727 //function : CheckNodeUV
728 //purpose  : Check and fix node UV on a face
729 //=======================================================================
730
731 bool SMESH_MesherHelper::CheckNodeUV(const TopoDS_Face&   F,
732                                      const SMDS_MeshNode* n,
733                                      gp_XY&               uv,
734                                      const double         tol,
735                                      const bool           force,
736                                      double               distXYZ[4]) const
737 {
738   int  shapeID = n->getshapeId();
739   bool infinit;
740   if (( infinit = ( Precision::IsInfinite( uv.X() ) || Precision::IsInfinite( uv.Y() ))) ||
741       ( force ) ||
742       ( uv.X() == 0. && uv.Y() == 0. ) ||
743       ( toCheckPosOnShape( shapeID )))
744   {
745     // check that uv is correct
746     TopLoc_Location loc;
747     Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
748     gp_Pnt nodePnt = XYZ( n ), surfPnt(0,0,0);
749     double dist = 0;
750     if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
751     if ( infinit ||
752          (dist = nodePnt.Distance( surfPnt = surface->Value( uv.X(), uv.Y() ))) > tol )
753     {
754       setPosOnShapeValidity( shapeID, false );
755       if ( !infinit && distXYZ ) {
756         surfPnt.Transform( loc );
757         distXYZ[0] = dist;
758         distXYZ[1] = surfPnt.X(); distXYZ[2] = surfPnt.Y(); distXYZ[3]=surfPnt.Z();
759       }
760       // uv incorrect, project the node to surface
761       GeomAPI_ProjectPointOnSurf& projector = GetProjector( F, loc, tol );
762       projector.Perform( nodePnt );
763       if ( !projector.IsDone() || projector.NbPoints() < 1 )
764       {
765         MESSAGE( "SMESH_MesherHelper::CheckNodeUV() failed to project" );
766         return false;
767       }
768       Quantity_Parameter U,V;
769       projector.LowerDistanceParameters(U,V);
770       uv.SetCoord( U,V );
771       surfPnt = surface->Value( U, V );
772       dist = nodePnt.Distance( surfPnt );
773       if ( distXYZ ) {
774         surfPnt.Transform( loc );
775         distXYZ[0] = dist;
776         distXYZ[1] = surfPnt.X(); distXYZ[2] = surfPnt.Y(); distXYZ[3]=surfPnt.Z();
777       }
778       if ( dist > tol )
779       {
780         MESSAGE( "SMESH_MesherHelper::CheckNodeUV(), invalid projection" );
781         return false;
782       }
783       // store the fixed UV on the face
784       if ( myShape.IsSame(F) && shapeID == myShapeID && myFixNodeParameters )
785         const_cast<SMDS_MeshNode*>(n)->SetPosition
786           ( SMDS_PositionPtr( new SMDS_FacePosition( U, V )));
787     }
788     else if ( myShape.IsSame(F) && uv.Modulus() > numeric_limits<double>::min() )
789     {
790       setPosOnShapeValidity( shapeID, true );
791     }
792   }
793   return true;
794 }
795
796 //=======================================================================
797 //function : GetProjector
798 //purpose  : Return projector intitialized by given face without location, which is returned
799 //=======================================================================
800
801 GeomAPI_ProjectPointOnSurf& SMESH_MesherHelper::GetProjector(const TopoDS_Face& F,
802                                                              TopLoc_Location&   loc,
803                                                              double             tol ) const
804 {
805   Handle(Geom_Surface) surface = BRep_Tool::Surface( F,loc );
806   int faceID = GetMeshDS()->ShapeToIndex( F );
807   TID2ProjectorOnSurf& i2proj = const_cast< TID2ProjectorOnSurf&>( myFace2Projector );
808   TID2ProjectorOnSurf::iterator i_proj = i2proj.find( faceID );
809   if ( i_proj == i2proj.end() )
810   {
811     if ( tol == 0 ) tol = BRep_Tool::Tolerance( F );
812     double U1, U2, V1, V2;
813     surface->Bounds(U1, U2, V1, V2);
814     GeomAPI_ProjectPointOnSurf* proj = new GeomAPI_ProjectPointOnSurf();
815     proj->Init( surface, U1, U2, V1, V2, tol );
816     i_proj = i2proj.insert( make_pair( faceID, proj )).first;
817   }
818   return *( i_proj->second );
819 }
820
821 //=======================================================================
822 //function : GetSurface
823 //purpose  : Return a cached ShapeAnalysis_Surface of a FACE
824 //=======================================================================
825
826 Handle(ShapeAnalysis_Surface) SMESH_MesherHelper::GetSurface(const TopoDS_Face& F ) const
827 {
828   Handle(Geom_Surface) surface = BRep_Tool::Surface( F );
829   int faceID = GetMeshDS()->ShapeToIndex( F );
830   TID2Surface::iterator i_surf = myFace2Surface.find( faceID );
831   if ( i_surf == myFace2Surface.end() && faceID )
832   {
833     Handle(ShapeAnalysis_Surface) surf( new ShapeAnalysis_Surface( surface ));
834     i_surf = myFace2Surface.insert( make_pair( faceID, surf )).first;
835   }
836   return i_surf->second;
837 }
838
839 namespace
840 {
841   gp_XY AverageUV(const gp_XY& uv1, const gp_XY& uv2) { return ( uv1 + uv2 ) / 2.; }
842   gp_XY_FunPtr(Added); // define gp_XY_Added pointer to function calling gp_XY::Added(gp_XY)
843   gp_XY_FunPtr(Subtracted); 
844 }
845
846 //=======================================================================
847 //function : ApplyIn2D
848 //purpose  : Perform given operation on two 2d points in parameric space of given surface.
849 //           It takes into account period of the surface. Use gp_XY_FunPtr macro
850 //           to easily define pointer to function of gp_XY class.
851 //=======================================================================
852
853 gp_XY SMESH_MesherHelper::ApplyIn2D(Handle(Geom_Surface) surface,
854                                     const gp_XY&         uv1,
855                                     const gp_XY&         uv2,
856                                     xyFunPtr             fun,
857                                     const bool           resultInPeriod)
858 {
859   if ( surface->IsKind(STANDARD_TYPE(Geom_RectangularTrimmedSurface )))
860     surface = Handle(Geom_RectangularTrimmedSurface)::DownCast( surface )->BasisSurface();
861   Standard_Boolean isUPeriodic = surface.IsNull() ? false : surface->IsUPeriodic();
862   Standard_Boolean isVPeriodic = surface.IsNull() ? false : surface->IsVPeriodic();
863   if ( !isUPeriodic && !isVPeriodic )
864     return fun(uv1,uv2);
865
866   // move uv2 not far than half-period from uv1
867   double u2 = 
868     uv2.X()+(isUPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.X(),uv1.X(),surface->UPeriod()) :0);
869   double v2 = 
870     uv2.Y()+(isVPeriodic ? ShapeAnalysis::AdjustByPeriod(uv2.Y(),uv1.Y(),surface->VPeriod()) :0);
871
872   // execute operation
873   gp_XY res = fun( uv1, gp_XY(u2,v2) );
874
875   // move result within period
876   if ( resultInPeriod )
877   {
878     Standard_Real UF,UL,VF,VL;
879     surface->Bounds(UF,UL,VF,VL);
880     if ( isUPeriodic )
881       res.SetX( res.X() + ShapeAnalysis::AdjustToPeriod(res.X(),UF,UL));
882     if ( isVPeriodic )
883       res.SetY( res.Y() + ShapeAnalysis::AdjustToPeriod(res.Y(),VF,VL));
884   }
885
886   return res;
887 }
888
889 //=======================================================================
890 //function : AdjustByPeriod
891 //purpose  : Move node positions on a FACE within surface period
892 //=======================================================================
893
894 void SMESH_MesherHelper::AdjustByPeriod( const TopoDS_Face& face, gp_XY uv[], const int nbUV )
895 {
896   SMESH_MesherHelper h( *myMesh ), *ph = face.IsSame( myShape ) ? this : &h;
897   ph->SetSubShape( face );
898
899   for ( int iCoo = U_periodic; iCoo <= V_periodic; ++iCoo )
900     if ( ph->GetPeriodicIndex() & iCoo )
901     {
902       const double period = ( ph->myPar2[iCoo-1] - ph->myPar1[iCoo-1] );
903       const double xRef = uv[0].Coord( iCoo );
904       for ( int i = 1; i < nbUV; ++i )
905       {
906         double x = uv[i].Coord( iCoo );
907         double dx = ShapeAnalysis::AdjustByPeriod( x, xRef, period );
908         uv[i].SetCoord( iCoo, x + dx );
909       }
910     }
911 }
912
913 //=======================================================================
914 //function : GetMiddleUV
915 //purpose  : Return middle UV taking in account surface period
916 //=======================================================================
917
918 gp_XY SMESH_MesherHelper::GetMiddleUV(const Handle(Geom_Surface)& surface,
919                                       const gp_XY&                p1,
920                                       const gp_XY&                p2)
921 {
922   // NOTE:
923   // the proper place of getting basic surface seems to be in ApplyIn2D()
924   // but we put it here to decrease a risk of regressions just before releasing a version
925   // Handle(Geom_Surface) surf = surface;
926   // while ( !surf.IsNull() && surf->IsKind(STANDARD_TYPE(Geom_RectangularTrimmedSurface )))
927   //   surf = Handle(Geom_RectangularTrimmedSurface)::DownCast( surf )->BasisSurface();
928
929   return ApplyIn2D( surface, p1, p2, & AverageUV );
930 }
931
932 //=======================================================================
933 //function : GetCenterUV
934 //purpose  : Return UV for the central node of a biquadratic triangle
935 //=======================================================================
936
937 gp_XY SMESH_MesherHelper::GetCenterUV(const gp_XY& uv1,
938                                       const gp_XY& uv2, 
939                                       const gp_XY& uv3, 
940                                       const gp_XY& uv12,
941                                       const gp_XY& uv23,
942                                       const gp_XY& uv31,
943                                       bool *       isBadTria/*=0*/)
944 {
945   bool badTria;
946   gp_XY uvAvg = ( uv12 + uv23 + uv31 ) / 3.;
947
948   if      (( badTria = (( uvAvg - uv1 ) * ( uvAvg - uv23 ) > 0 )))
949     uvAvg = ( uv1 + uv23 ) / 2.;
950   else if (( badTria = (( uvAvg - uv2 ) * ( uvAvg - uv31 ) > 0 )))
951     uvAvg = ( uv2 + uv31 ) / 2.;
952   else if (( badTria = (( uvAvg - uv3 ) * ( uvAvg - uv12 ) > 0 )))
953     uvAvg = ( uv3 + uv12 ) / 2.;
954
955   if ( isBadTria )
956     *isBadTria = badTria;
957   return uvAvg;
958 }
959
960 //=======================================================================
961 //function : GetNodeU
962 //purpose  : Return node U on edge
963 //=======================================================================
964
965 double SMESH_MesherHelper::GetNodeU(const TopoDS_Edge&   E,
966                                     const SMDS_MeshNode* n,
967                                     const SMDS_MeshNode* inEdgeNode,
968                                     bool*                check) const
969 {
970   double param = Precision::Infinite();
971
972   const SMDS_PositionPtr pos = n->GetPosition();
973   if ( pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
974   {
975     const SMDS_EdgePosition* epos = static_cast<const SMDS_EdgePosition*>( pos );
976     param =  epos->GetUParameter();
977   }
978   else if( pos->GetTypeOfPosition() == SMDS_TOP_VERTEX )
979   {
980     if ( inEdgeNode && TopExp::FirstVertex( E ).IsSame( TopExp::LastVertex( E ))) // issue 0020128
981     {
982       Standard_Real f,l;
983       BRep_Tool::Range( E, f,l );
984       double uInEdge = GetNodeU( E, inEdgeNode );
985       param = ( fabs( uInEdge - f ) < fabs( l - uInEdge )) ? f : l;
986     }
987     else
988     {
989       SMESHDS_Mesh * meshDS = GetMeshDS();
990       int vertexID = n->getshapeId();
991       const TopoDS_Vertex& V = TopoDS::Vertex(meshDS->IndexToShape(vertexID));
992       param =  BRep_Tool::Parameter( V, E );
993     }
994   }
995   if ( check )
996   {
997     double tol = BRep_Tool::Tolerance( E );
998     double f,l;  BRep_Tool::Range( E, f,l );
999     bool force = ( param < f-tol || param > l+tol );
1000     if ( !force && pos->GetTypeOfPosition()==SMDS_TOP_EDGE )
1001       force = ( GetMeshDS()->ShapeToIndex( E ) != n->getshapeId() );
1002
1003     *check = CheckNodeU( E, n, param, 2*tol, force );
1004   }
1005   return param;
1006 }
1007
1008 //=======================================================================
1009 //function : CheckNodeU
1010 //purpose  : Check and fix node U on an edge
1011 //           Return false if U is bad and could not be fixed
1012 //=======================================================================
1013
1014 bool SMESH_MesherHelper::CheckNodeU(const TopoDS_Edge&   E,
1015                                     const SMDS_MeshNode* n,
1016                                     double&              u,
1017                                     const double         tol,
1018                                     const bool           force,
1019                                     double               distXYZ[4]) const
1020 {
1021   int  shapeID = n->getshapeId();
1022   bool infinit;
1023   if (( infinit = Precision::IsInfinite( u )) ||
1024       ( force ) ||
1025       ( u == 0. ) ||
1026       ( toCheckPosOnShape( shapeID )))
1027   {
1028     TopLoc_Location loc; double f,l;
1029     Handle(Geom_Curve) curve = BRep_Tool::Curve( E,loc,f,l );
1030     if ( curve.IsNull() ) // degenerated edge
1031     {
1032       if ( u+tol < f || u-tol > l )
1033       {
1034         double r = Max( 0.5, 1 - tol*n->GetID()); // to get a unique u on edge
1035         u =  f*r + l*(1-r);
1036       }
1037     }
1038     else
1039     {
1040       gp_Pnt nodePnt = SMESH_TNodeXYZ( n );
1041       if ( !loc.IsIdentity() ) nodePnt.Transform( loc.Transformation().Inverted() );
1042       gp_Pnt curvPnt;
1043       double dist = 2*tol;
1044       if ( !infinit )
1045       {
1046         curvPnt = curve->Value( u );
1047         dist    = nodePnt.Distance( curvPnt );
1048         if ( distXYZ ) {
1049           curvPnt.Transform( loc );
1050           distXYZ[0] = dist;
1051           distXYZ[1] = curvPnt.X(); distXYZ[2] = curvPnt.Y(); distXYZ[3]=curvPnt.Z();
1052         }
1053       }
1054       if ( dist > tol )
1055       {
1056         setPosOnShapeValidity( shapeID, false );
1057         // u incorrect, project the node to the curve
1058         int edgeID = GetMeshDS()->ShapeToIndex( E );
1059         TID2ProjectorOnCurve& i2proj = const_cast< TID2ProjectorOnCurve&>( myEdge2Projector );
1060         TID2ProjectorOnCurve::iterator i_proj =
1061           i2proj.insert( make_pair( edgeID, (GeomAPI_ProjectPointOnCurve*) 0 )).first;
1062         if ( !i_proj->second  )
1063         {
1064           i_proj->second = new GeomAPI_ProjectPointOnCurve();
1065           i_proj->second->Init( curve, f, l );
1066         }
1067         GeomAPI_ProjectPointOnCurve* projector = i_proj->second;
1068         projector->Perform( nodePnt );
1069         if ( projector->NbPoints() < 1 )
1070         {
1071           MESSAGE( "SMESH_MesherHelper::CheckNodeU() failed to project" );
1072           return false;
1073         }
1074         Quantity_Parameter U = projector->LowerDistanceParameter();
1075         u = double( U );
1076         MESSAGE(" f " << f << " l " << l << " u " << u);
1077         curvPnt = curve->Value( u );
1078         dist = nodePnt.Distance( curvPnt );
1079         if ( distXYZ ) {
1080           curvPnt.Transform( loc );
1081           distXYZ[0] = dist;
1082           distXYZ[1] = curvPnt.X(); distXYZ[2] = curvPnt.Y(); distXYZ[3]=curvPnt.Z();
1083         }
1084         if ( dist > tol )
1085         {
1086           MESSAGE( "SMESH_MesherHelper::CheckNodeU(), invalid projection" );
1087           MESSAGE("distance " << dist << " " << tol );
1088           return false;
1089         }
1090         // store the fixed U on the edge
1091         if ( myShape.IsSame(E) && shapeID == myShapeID && myFixNodeParameters )
1092           const_cast<SMDS_MeshNode*>(n)->SetPosition
1093             ( SMDS_PositionPtr( new SMDS_EdgePosition( U )));
1094       }
1095       else if ( fabs( u ) > numeric_limits<double>::min() )
1096       {
1097         setPosOnShapeValidity( shapeID, true );
1098       }
1099       if (( u < f-tol || u > l+tol ) && force )
1100       {
1101         MESSAGE("u < f-tol || u > l+tol  ; u " << u << " f " << f << " l " << l);
1102         // node is on vertex but is set on periodic but trimmed edge (issue 0020890)
1103         try
1104         {
1105           // do not use IsPeriodic() as Geom_TrimmedCurve::IsPeriodic () returns false
1106           double period = curve->Period();
1107           u = ( u < f ) ? u + period : u - period;
1108         }
1109         catch (Standard_Failure& exc)
1110         {
1111           return false;
1112         }
1113       }
1114     }
1115   }
1116   return true;
1117 }
1118
1119 //=======================================================================
1120 //function : GetMediumPos
1121 //purpose  : Return index and type of the shape  (EDGE or FACE only) to
1122 //           set a medium node on
1123 //param    : useCurSubShape - if true, returns the shape set via SetSubShape()
1124 //           if any
1125 //param    : expectedSupport - shape type corresponding to element being created,
1126 //                             e.g TopAbs_EDGE if SMDSAbs_Edge is created
1127 //                             basing on \a n1 and \a n2
1128 // Calling GetMediumPos() with useCurSubShape=true is OK only for the
1129 // case where the lower dim mesh is already constructed and converted to quadratic,
1130 // else, nodes on EDGEs are assigned to FACE, for example.
1131 //=======================================================================
1132
1133 std::pair<int, TopAbs_ShapeEnum>
1134 SMESH_MesherHelper::GetMediumPos(const SMDS_MeshNode* n1,
1135                                  const SMDS_MeshNode* n2,
1136                                  const bool           useCurSubShape,
1137                                  TopAbs_ShapeEnum     expectedSupport)
1138 {
1139   if ( useCurSubShape && !myShape.IsNull() )
1140     return std::make_pair( myShapeID, myShape.ShapeType() );
1141
1142   TopAbs_ShapeEnum shapeType = TopAbs_SHAPE;
1143   int              shapeID = -1;
1144   TopoDS_Shape     shape;
1145
1146   if (( myShapeID == n1->getshapeId() || myShapeID == n2->getshapeId() ) && myShapeID > 0 )
1147   {
1148     shapeType = myShape.ShapeType();
1149     shapeID   = myShapeID;
1150   }
1151   else if ( n1->getshapeId() == n2->getshapeId() )
1152   {
1153     shapeID = n2->getshapeId();
1154     shape = GetSubShapeByNode( n1, GetMeshDS() );
1155   }
1156   else // 2 different shapes
1157   {
1158     const SMDS_TypeOfPosition Pos1 = n1->GetPosition()->GetTypeOfPosition();
1159     const SMDS_TypeOfPosition Pos2 = n2->GetPosition()->GetTypeOfPosition();
1160
1161     if ( Pos1 == SMDS_TOP_3DSPACE || Pos2 == SMDS_TOP_3DSPACE )
1162     {
1163       // in SOLID
1164     }
1165     else if ( Pos1 == SMDS_TOP_FACE || Pos2 == SMDS_TOP_FACE )
1166     {
1167       // in FACE or SOLID
1168       if ( Pos1 != SMDS_TOP_FACE || Pos2 != SMDS_TOP_FACE ) // not 2 FACEs
1169       {
1170         if ( Pos1 != SMDS_TOP_FACE ) std::swap( n1,n2 );
1171         TopoDS_Shape F = GetSubShapeByNode( n1, GetMeshDS() );
1172         TopoDS_Shape S = GetSubShapeByNode( n2, GetMeshDS() );
1173         if ( IsSubShape( S, F ))
1174         {
1175           shapeType = TopAbs_FACE;
1176           shapeID   = n1->getshapeId();
1177         }
1178       }
1179     }
1180     else if ( Pos1 == SMDS_TOP_EDGE && Pos2 == SMDS_TOP_EDGE )
1181     {
1182       TopoDS_Shape E1 = GetSubShapeByNode( n1, GetMeshDS() );
1183       TopoDS_Shape E2 = GetSubShapeByNode( n2, GetMeshDS() );
1184       shape = GetCommonAncestor( E1, E2, *myMesh, TopAbs_FACE );
1185     }
1186     else if ( Pos1 == SMDS_TOP_VERTEX && Pos2 == SMDS_TOP_VERTEX )
1187     {
1188       TopoDS_Shape V1 = GetSubShapeByNode( n1, GetMeshDS() );
1189       TopoDS_Shape V2 = GetSubShapeByNode( n2, GetMeshDS() );
1190       shape = GetCommonAncestor( V1, V2, *myMesh, TopAbs_EDGE );
1191       if ( shape.IsNull() ) shape = GetCommonAncestor( V1, V2, *myMesh, TopAbs_FACE );
1192     }
1193     else // on VERTEX and EDGE
1194     {
1195       if ( Pos1 != SMDS_TOP_VERTEX ) std::swap( n1,n2 );
1196       TopoDS_Shape V = GetSubShapeByNode( n1, GetMeshDS() );
1197       TopoDS_Shape E = GetSubShapeByNode( n2, GetMeshDS() );
1198       if ( IsSubShape( V, E ))
1199         shape = E;
1200       else
1201         shape = GetCommonAncestor( V, E, *myMesh, TopAbs_FACE );
1202     }
1203   }
1204
1205   if ( !shape.IsNull() )
1206   {
1207     if ( shapeID < 1 )
1208       shapeID = GetMeshDS()->ShapeToIndex( shape );
1209     shapeType = shape.ShapeType(); // EDGE or FACE
1210
1211     if ( expectedSupport < shapeType &&
1212          expectedSupport != TopAbs_SHAPE &&
1213          !myShape.IsNull() &&
1214          myShape.ShapeType() == expectedSupport )
1215     {
1216       // e.g. a side of triangle connects nodes on the same EDGE but does not
1217       // lie on this EDGE (an arc with a coarse mesh)
1218       // =>  shapeType == TopAbs_EDGE, expectedSupport == TopAbs_FACE;
1219       // hope that myShape is a right shape, return it if the found shape
1220       // has converted elements of corresponding dim (segments in our example)
1221       int nbConvertedElems = 0;
1222       SMDSAbs_ElementType type = ( shapeType == TopAbs_FACE ? SMDSAbs_Face : SMDSAbs_Edge );
1223       for ( int iN = 0; iN < 2; ++iN )
1224       {
1225         const SMDS_MeshNode* n = iN ? n2 : n1;
1226         SMDS_ElemIteratorPtr it = n->GetInverseElementIterator( type );
1227         while ( it->more() )
1228         {
1229           const SMDS_MeshElement* elem = it->next();
1230           if ( elem->getshapeId() == shapeID &&
1231                elem->IsQuadratic() )
1232           {
1233             ++nbConvertedElems;
1234             break;
1235           }
1236         }
1237       }
1238       if ( nbConvertedElems == 2 )
1239       {
1240         shapeType = myShape.ShapeType();
1241         shapeID   = myShapeID;
1242       }
1243     }
1244   }
1245   return make_pair( shapeID, shapeType );
1246 }
1247
1248 //=======================================================================
1249 //function : GetCentralNode
1250 //purpose  : Return existing or create a new central node for a quardilateral
1251 //           quadratic face given its 8 nodes.
1252 //@param   : force3d - true means node creation in between the given nodes,
1253 //           else node position is found on a geometrical face if any.
1254 //=======================================================================
1255
1256 const SMDS_MeshNode* SMESH_MesherHelper::GetCentralNode(const SMDS_MeshNode* n1,
1257                                                         const SMDS_MeshNode* n2,
1258                                                         const SMDS_MeshNode* n3,
1259                                                         const SMDS_MeshNode* n4,
1260                                                         const SMDS_MeshNode* n12,
1261                                                         const SMDS_MeshNode* n23,
1262                                                         const SMDS_MeshNode* n34,
1263                                                         const SMDS_MeshNode* n41,
1264                                                         bool                 force3d)
1265 {
1266   SMDS_MeshNode *centralNode = 0; // central node to return
1267
1268   // Find an existing central node
1269
1270   TBiQuad keyOfMap(n1,n2,n3,n4);
1271   std::map<TBiQuad, const SMDS_MeshNode* >::iterator itMapCentralNode;
1272   itMapCentralNode = myMapWithCentralNode.find( keyOfMap );
1273   if ( itMapCentralNode != myMapWithCentralNode.end() ) 
1274   {
1275     return (*itMapCentralNode).second;
1276   }
1277
1278   // Get type of shape for the new central node
1279
1280   TopAbs_ShapeEnum shapeType = TopAbs_SHAPE;
1281   int              solidID = -1;
1282   int              faceID = -1;
1283   TopoDS_Shape     shape;
1284   TopTools_ListIteratorOfListOfShape it;
1285
1286   std::map< int, int > faceId2nbNodes;
1287   std::map< int, int > ::iterator itMapWithIdFace;
1288   
1289   SMESHDS_Mesh* meshDS = GetMeshDS();
1290   
1291   // check if a face lies on a FACE, i.e. its all corner nodes lie either on the FACE or
1292   // on sub-shapes of the FACE
1293   if ( GetMesh()->HasShapeToMesh() )
1294   {
1295     const SMDS_MeshNode* nodes[] = { n1, n2, n3, n4 };
1296     for(int i = 0; i < 4; i++)
1297     {
1298       shape = GetSubShapeByNode( nodes[i], meshDS );
1299       if ( shape.IsNull() ) break;
1300       if ( shape.ShapeType() == TopAbs_SOLID )
1301       {
1302         solidID   = nodes[i]->getshapeId();
1303         shapeType = TopAbs_SOLID;
1304         break;
1305       }
1306       if ( shape.ShapeType() == TopAbs_FACE )
1307       {
1308         faceID          = nodes[i]->getshapeId();
1309         itMapWithIdFace = faceId2nbNodes.insert( std::make_pair( faceID, 0 ) ).first;
1310         itMapWithIdFace->second++;
1311       }
1312       else
1313       {
1314         PShapeIteratorPtr it = GetAncestors( shape, *GetMesh(), TopAbs_FACE );
1315         while ( const TopoDS_Shape* face = it->next() )
1316         {
1317           faceID = meshDS->ShapeToIndex( *face );
1318           itMapWithIdFace = faceId2nbNodes.insert( std::make_pair( faceID, 0 )).first;
1319           itMapWithIdFace->second++;
1320         }
1321       }
1322     }
1323   }
1324   if ( solidID < 1 && !faceId2nbNodes.empty() ) // SOLID not found
1325   {
1326     // find ID of the FACE the four corner nodes belong to
1327     itMapWithIdFace = faceId2nbNodes.find( myShapeID ); // IPAL52698
1328     if ( itMapWithIdFace != faceId2nbNodes.end() &&
1329          itMapWithIdFace->second == 4 )
1330     {
1331       shapeType = TopAbs_FACE;
1332       faceID = myShapeID;
1333     }
1334     else
1335     {
1336       itMapWithIdFace = faceId2nbNodes.begin();
1337       for ( ; itMapWithIdFace != faceId2nbNodes.end(); ++itMapWithIdFace)
1338       {
1339         if ( itMapWithIdFace->second == 4 )
1340         {
1341           shapeType = TopAbs_FACE;
1342           faceID = (*itMapWithIdFace).first;
1343           break;
1344         }
1345       }
1346     }
1347   }
1348
1349   TopoDS_Face F;
1350   if ( shapeType == TopAbs_FACE )
1351   {
1352     F = TopoDS::Face( meshDS->IndexToShape( faceID ));
1353   }
1354
1355   // Create a node
1356
1357   gp_XY  uvAvg;
1358   gp_Pnt P;
1359   bool toCheck = true;
1360   if ( !F.IsNull() && !force3d )
1361   {
1362     Handle(ShapeAnalysis_Surface) surface = GetSurface( F );
1363     if ( HasDegeneratedEdges() || surface->HasSingularities( 1e-7 ))
1364     {
1365       gp_Pnt center = calcTFI (0.5, 0.5, // IPAL0052863
1366                                SMESH_TNodeXYZ(n1),  SMESH_TNodeXYZ(n2),
1367                                SMESH_TNodeXYZ(n3),  SMESH_TNodeXYZ(n4),
1368                                SMESH_TNodeXYZ(n12), SMESH_TNodeXYZ(n23),
1369                                SMESH_TNodeXYZ(n34), SMESH_TNodeXYZ(n41));
1370       gp_Pnt2d uv12 = GetNodeUV( F, n12, n3, &toCheck );
1371       uvAvg = surface->NextValueOfUV( uv12, center, BRep_Tool::Tolerance( F )).XY();
1372     }
1373     else
1374     {
1375       gp_XY uv[8] = {
1376         GetNodeUV( F,n1,  n3, &toCheck ),
1377         GetNodeUV( F,n2,  n4, &toCheck ),
1378         GetNodeUV( F,n3,  n1, &toCheck ),
1379         GetNodeUV( F,n4,  n2, &toCheck ),
1380         GetNodeUV( F,n12, n3 ),
1381         GetNodeUV( F,n23, n4 ),
1382         GetNodeUV( F,n34, n2 ),
1383         GetNodeUV( F,n41, n2 )
1384       };
1385       AdjustByPeriod( F, uv, 8 ); // put uv[] within a period (IPAL52698)
1386
1387       uvAvg = calcTFI (0.5, 0.5, uv[0],uv[1],uv[2],uv[3], uv[4],uv[5],uv[6],uv[7] );
1388     }
1389     P = surface->Value( uvAvg );
1390     centralNode = meshDS->AddNode( P.X(), P.Y(), P.Z() );
1391     // if ( mySetElemOnShape ) node is not elem!
1392     meshDS->SetNodeOnFace( centralNode, faceID, uvAvg.X(), uvAvg.Y() );
1393   }
1394   else // ( force3d || F.IsNull() )
1395   {
1396     P = calcTFI (0.5, 0.5,
1397                  SMESH_TNodeXYZ(n1),  SMESH_TNodeXYZ(n2),
1398                  SMESH_TNodeXYZ(n3),  SMESH_TNodeXYZ(n4),
1399                  SMESH_TNodeXYZ(n12), SMESH_TNodeXYZ(n23),
1400                  SMESH_TNodeXYZ(n34), SMESH_TNodeXYZ(n41));
1401     centralNode = meshDS->AddNode( P.X(), P.Y(), P.Z() );
1402
1403     if ( !F.IsNull() ) // force3d
1404     {
1405       uvAvg = (GetNodeUV(F,n1,n3,&toCheck) +
1406                GetNodeUV(F,n2,n4,&toCheck) +
1407                GetNodeUV(F,n3,n1,&toCheck) +
1408                GetNodeUV(F,n4,n2,&toCheck)) / 4;
1409       //CheckNodeUV( F, centralNode, uvAvg, 2*BRep_Tool::Tolerance( F ), /*force=*/true);
1410       meshDS->SetNodeOnFace( centralNode, faceID, uvAvg.X(), uvAvg.Y() );
1411     }
1412     else if ( solidID > 0 )
1413     {
1414       meshDS->SetNodeInVolume( centralNode, solidID );
1415     }
1416     else if ( myShapeID > 0 && mySetElemOnShape )
1417     {
1418       meshDS->SetMeshElementOnShape( centralNode, myShapeID );
1419     }
1420   }
1421   myMapWithCentralNode.insert( std::make_pair( keyOfMap, centralNode ) );
1422   return centralNode;
1423 }
1424
1425 //=======================================================================
1426 //function : GetCentralNode
1427 //purpose  : Return existing or create a new central node for a
1428 //           quadratic triangle given its 6 nodes.
1429 //@param   : force3d - true means node creation in between the given nodes,
1430 //           else node position is found on a geometrical face if any.
1431 //=======================================================================
1432
1433 const SMDS_MeshNode* SMESH_MesherHelper::GetCentralNode(const SMDS_MeshNode* n1,
1434                                                         const SMDS_MeshNode* n2,
1435                                                         const SMDS_MeshNode* n3,
1436                                                         const SMDS_MeshNode* n12,
1437                                                         const SMDS_MeshNode* n23,
1438                                                         const SMDS_MeshNode* n31,
1439                                                         bool                 force3d)
1440 {
1441   SMDS_MeshNode *centralNode = 0; // central node to return
1442
1443   // Find an existing central node
1444
1445   TBiQuad keyOfMap(n1,n2,n3);
1446   std::map<TBiQuad, const SMDS_MeshNode* >::iterator itMapCentralNode;
1447   itMapCentralNode = myMapWithCentralNode.find( keyOfMap );
1448   if ( itMapCentralNode != myMapWithCentralNode.end() ) 
1449   {
1450     return (*itMapCentralNode).second;
1451   }
1452
1453   // Get type of shape for the new central node
1454
1455   TopAbs_ShapeEnum shapeType = TopAbs_SHAPE;
1456   int              solidID = -1;
1457   int              faceID = -1;
1458   TopoDS_Shape     shape;
1459   TopTools_ListIteratorOfListOfShape it;
1460
1461   std::map< int, int > faceId2nbNodes;
1462   std::map< int, int > ::iterator itMapWithIdFace;
1463   
1464   SMESHDS_Mesh* meshDS = GetMeshDS();
1465   
1466   // check if a face lies on a FACE, i.e. its all corner nodes lie either on the FACE or
1467   // on sub-shapes of the FACE
1468   if ( GetMesh()->HasShapeToMesh() )
1469   {
1470     const SMDS_MeshNode* nodes[] = { n1, n2, n3 };
1471     for(int i = 0; i < 3; i++)
1472     {
1473       shape = GetSubShapeByNode( nodes[i], meshDS );
1474       if ( shape.IsNull() ) break;
1475       if ( shape.ShapeType() == TopAbs_SOLID )
1476       {
1477         solidID   = nodes[i]->getshapeId();
1478         shapeType = TopAbs_SOLID;
1479         break;
1480       }
1481       if ( shape.ShapeType() == TopAbs_FACE )
1482       {
1483         faceID          = nodes[i]->getshapeId();
1484         itMapWithIdFace = faceId2nbNodes.insert( std::make_pair( faceID, 0 ) ).first;
1485         itMapWithIdFace->second++;
1486       }
1487       else
1488       {
1489         PShapeIteratorPtr it = GetAncestors(shape, *GetMesh(), TopAbs_FACE );
1490         while ( const TopoDS_Shape* face = it->next() )
1491         {
1492           faceID = meshDS->ShapeToIndex( *face );
1493           itMapWithIdFace = faceId2nbNodes.insert( std::make_pair( faceID, 0 ) ).first;
1494           itMapWithIdFace->second++;
1495         }
1496       }
1497     }
1498   }
1499   if ( solidID < 1 && !faceId2nbNodes.empty() ) // SOLID not found
1500   {
1501     // find ID of the FACE the four corner nodes belong to
1502     itMapWithIdFace = faceId2nbNodes.find( myShapeID ); // IPAL52698
1503     if ( itMapWithIdFace != faceId2nbNodes.end() &&
1504          itMapWithIdFace->second == 4 )
1505     {
1506       shapeType = TopAbs_FACE;
1507       faceID = myShapeID;
1508     }
1509     else
1510     {
1511       itMapWithIdFace = faceId2nbNodes.begin();
1512       for ( ; itMapWithIdFace != faceId2nbNodes.end(); ++itMapWithIdFace)
1513       {
1514         if ( itMapWithIdFace->second == 3 )
1515         {
1516           shapeType = TopAbs_FACE;
1517           faceID = (*itMapWithIdFace).first;
1518           break;
1519         }
1520       }
1521     }
1522   }
1523
1524   TopoDS_Face F;
1525   gp_XY       uvAvg;
1526
1527   if ( shapeType == TopAbs_FACE )
1528   {
1529     F = TopoDS::Face( meshDS->IndexToShape( faceID ));
1530     bool checkOK = true, badTria = false;
1531     gp_XY uv[6] = {
1532       GetNodeUV( F, n1, n23, &checkOK ),
1533       GetNodeUV( F, n2, n31, &checkOK ),
1534       GetNodeUV( F, n3, n12, &checkOK ),
1535       GetNodeUV( F, n12, n3, &checkOK ),
1536       GetNodeUV( F, n23, n1, &checkOK ),
1537       GetNodeUV( F, n31, n2, &checkOK )
1538     };
1539     AdjustByPeriod( F, uv, 6 ); // put uv[] within a period (IPAL52698)
1540
1541     uvAvg = GetCenterUV( uv[0],uv[1],uv[2], uv[3],uv[4],uv[5], &badTria );
1542
1543     if ( badTria || !checkOK )
1544       force3d = true;
1545   }
1546
1547   // Create a central node
1548
1549   gp_Pnt P;
1550   if ( !F.IsNull() && !force3d )
1551   {
1552     TopLoc_Location        loc;
1553     Handle( Geom_Surface ) S = BRep_Tool::Surface( F, loc );
1554     P = S->Value( uvAvg.X(), uvAvg.Y() ).Transformed( loc );
1555     centralNode = meshDS->AddNode( P.X(), P.Y(), P.Z() );
1556     // if ( mySetElemOnShape ) node is not elem!
1557     meshDS->SetNodeOnFace( centralNode, faceID, uvAvg.X(), uvAvg.Y() );
1558   }
1559   else // ( force3d || F.IsNull() )
1560   {
1561     P = ( SMESH_TNodeXYZ( n12 ) +
1562           SMESH_TNodeXYZ( n23 ) +
1563           SMESH_TNodeXYZ( n31 ) ) / 3;
1564     centralNode = meshDS->AddNode( P.X(), P.Y(), P.Z() );
1565
1566     if ( !F.IsNull() ) // force3d
1567     {
1568       meshDS->SetNodeOnFace( centralNode, faceID, uvAvg.X(), uvAvg.Y() );
1569     }
1570     else if ( solidID > 0 )
1571     {
1572       meshDS->SetNodeInVolume( centralNode, solidID );
1573     }
1574     else if ( myShapeID > 0 && mySetElemOnShape )
1575     {
1576       meshDS->SetMeshElementOnShape( centralNode, myShapeID );
1577     }
1578   }
1579   myMapWithCentralNode.insert( std::make_pair( keyOfMap, centralNode ) );
1580   return centralNode;
1581 }
1582
1583 //=======================================================================
1584 //function : GetMediumNode
1585 //purpose  : Return existing or create a new medium node between given ones
1586 //=======================================================================
1587
1588 const SMDS_MeshNode* SMESH_MesherHelper::GetMediumNode(const SMDS_MeshNode* n1,
1589                                                        const SMDS_MeshNode* n2,
1590                                                        bool                 force3d,
1591                                                        TopAbs_ShapeEnum     expectedSupport)
1592 {
1593   // Find existing node
1594
1595   SMESH_TLink link(n1,n2);
1596   ItTLinkNode itLN = myTLinkNodeMap.find( link );
1597   if ( itLN != myTLinkNodeMap.end() ) {
1598     return (*itLN).second;
1599   }
1600
1601   // Create medium node
1602
1603   SMDS_MeshNode* n12;
1604   SMESHDS_Mesh* meshDS = GetMeshDS();
1605
1606   if ( IsSeamShape( n1->getshapeId() ))
1607     // to get a correct UV of a node on seam, the second node must have checked UV
1608     std::swap( n1, n2 );
1609
1610   // get type of shape for the new medium node
1611   int faceID = -1, edgeID = -1;
1612   TopoDS_Edge E; double u [2];
1613   TopoDS_Face F; gp_XY  uv[2];
1614   bool uvOK[2] = { true, true };
1615   const bool useCurSubShape = ( !myShape.IsNull() && myShape.ShapeType() == TopAbs_EDGE );
1616
1617   pair<int, TopAbs_ShapeEnum> pos = GetMediumPos( n1, n2, useCurSubShape, expectedSupport );
1618
1619   // get positions of the given nodes on shapes
1620   if ( pos.second == TopAbs_FACE )
1621   {
1622     F = TopoDS::Face(meshDS->IndexToShape( faceID = pos.first ));
1623     uv[0] = GetNodeUV(F,n1,n2, force3d ? 0 : &uvOK[0]);
1624     if (( !force3d ) &&
1625         ( HasDegeneratedEdges() || GetSurface( F )->HasSingularities( 1e-7 )))
1626     {
1627       // IPAL52850 (degen VERTEX not at singularity)
1628       // project middle point to a surface
1629       SMESH_TNodeXYZ p1( n1 ), p2( n2 );
1630       gp_Pnt pMid = 0.5 * ( p1 + p2 );
1631       Handle(ShapeAnalysis_Surface) projector = GetSurface( F );
1632       gp_Pnt2d uvMid;
1633       if ( uvOK[0] )
1634         uvMid = projector->NextValueOfUV( uv[0], pMid, BRep_Tool::Tolerance( F ));
1635       else
1636         uvMid = projector->ValueOfUV( pMid, getFaceMaxTol( F ));
1637       if ( projector->Gap() * projector->Gap() < ( p1 - p2 ).SquareModulus() / 4 )
1638       {
1639         gp_Pnt pProj = projector->Value( uvMid );
1640         n12  = meshDS->AddNode( pProj.X(), pProj.Y(), pProj.Z() );
1641         meshDS->SetNodeOnFace( n12, faceID, uvMid.X(), uvMid.Y() );
1642         myTLinkNodeMap.insert( make_pair ( link, n12 ));
1643         return n12;
1644       }
1645     }
1646     uv[1] = GetNodeUV(F,n2,n1, force3d ? 0 : &uvOK[1]);
1647   }
1648   else if ( pos.second == TopAbs_EDGE )
1649   {
1650     const SMDS_PositionPtr Pos1 = n1->GetPosition();
1651     const SMDS_PositionPtr Pos2 = n2->GetPosition();
1652     if ( Pos1->GetTypeOfPosition()==SMDS_TOP_EDGE &&
1653          Pos2->GetTypeOfPosition()==SMDS_TOP_EDGE &&
1654          n1->getshapeId() != n2->getshapeId() )
1655     {
1656       // issue 0021006
1657       return getMediumNodeOnComposedWire(n1,n2,force3d);
1658     }
1659     E = TopoDS::Edge(meshDS->IndexToShape( edgeID = pos.first ));
1660     try {
1661       u[0] = GetNodeU(E,n1,n2, force3d ? 0 : &uvOK[0]);
1662       u[1] = GetNodeU(E,n2,n1, force3d ? 0 : &uvOK[1]);
1663     }
1664     catch ( Standard_Failure& f )
1665     {
1666       // issue 22502 / a node is on VERTEX not belonging to E
1667       // issue 22568 / both nodes are on non-connected VERTEXes
1668       return getMediumNodeOnComposedWire(n1,n2,force3d);
1669     }
1670   }
1671
1672   if ( !force3d & uvOK[0] && uvOK[1] )
1673   {
1674     // we try to create medium node using UV parameters of
1675     // nodes, else - medium between corresponding 3d points
1676     if( ! F.IsNull() )
1677     {
1678       //if ( uvOK[0] && uvOK[1] )
1679       {
1680         if ( IsDegenShape( n1->getshapeId() )) {
1681           if ( myParIndex & U_periodic ) uv[0].SetCoord( 1, uv[1].Coord( 1 ));
1682           else                           uv[0].SetCoord( 2, uv[1].Coord( 2 ));
1683         }
1684         else if ( IsDegenShape( n2->getshapeId() )) {
1685           if ( myParIndex & U_periodic ) uv[1].SetCoord( 1, uv[0].Coord( 1 ));
1686           else                           uv[1].SetCoord( 2, uv[0].Coord( 2 ));
1687         }
1688         TopLoc_Location loc;
1689         Handle(Geom_Surface) S = BRep_Tool::Surface(F,loc);
1690         gp_XY UV = GetMiddleUV( S, uv[0], uv[1] );
1691         gp_Pnt P = S->Value( UV.X(), UV.Y() ).Transformed(loc);
1692         n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
1693         // if ( mySetElemOnShape ) node is not elem!
1694         meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y());
1695         myTLinkNodeMap.insert(make_pair(link,n12));
1696         return n12;
1697       }
1698     }
1699     else if ( !E.IsNull() )
1700     {
1701       double f,l;
1702       Handle(Geom_Curve) C = BRep_Tool::Curve(E, f, l);
1703       if(!C.IsNull())
1704       {
1705         Standard_Boolean isPeriodic = C->IsPeriodic();
1706         double U;
1707         if(isPeriodic) {
1708           Standard_Real Period = C->Period();
1709           Standard_Real p = u[1]+ShapeAnalysis::AdjustByPeriod(u[1],u[0],Period);
1710           Standard_Real pmid = (u[0]+p)/2.;
1711           U = pmid+ShapeAnalysis::AdjustToPeriod(pmid,C->FirstParameter(),C->LastParameter());
1712         }
1713         else
1714           U = (u[0]+u[1])/2.;
1715
1716         gp_Pnt P = C->Value( U );
1717         n12 = meshDS->AddNode(P.X(), P.Y(), P.Z());
1718         //if ( mySetElemOnShape ) node is not elem!
1719         meshDS->SetNodeOnEdge(n12, edgeID, U);
1720         myTLinkNodeMap.insert(make_pair(link,n12));
1721         return n12;
1722       }
1723     }
1724   }
1725
1726   // 3d variant
1727   double x = ( n1->X() + n2->X() )/2.;
1728   double y = ( n1->Y() + n2->Y() )/2.;
1729   double z = ( n1->Z() + n2->Z() )/2.;
1730   n12 = meshDS->AddNode(x,y,z);
1731
1732   //if ( mySetElemOnShape ) node is not elem!
1733   {
1734     if ( !F.IsNull() )
1735     {
1736       gp_XY UV = ( uv[0] + uv[1] ) / 2.;
1737       CheckNodeUV( F, n12, UV, 2 * BRep_Tool::Tolerance( F ), /*force=*/true);
1738       meshDS->SetNodeOnFace(n12, faceID, UV.X(), UV.Y() );
1739     }
1740     else if ( !E.IsNull() )
1741     {
1742       double U = ( u[0] + u[1] ) / 2.;
1743       CheckNodeU( E, n12, U, 2 * BRep_Tool::Tolerance( E ), /*force=*/true);
1744       meshDS->SetNodeOnEdge(n12, edgeID, U);
1745     }
1746     else if ( myShapeID > 0 && mySetElemOnShape )
1747     {
1748       meshDS->SetMeshElementOnShape(n12, myShapeID);
1749     }
1750   }
1751
1752   myTLinkNodeMap.insert( make_pair( link, n12 ));
1753   return n12;
1754 }
1755
1756 //================================================================================
1757 /*!
1758  * \brief Makes a medium node if nodes reside different edges
1759  */
1760 //================================================================================
1761
1762 const SMDS_MeshNode* SMESH_MesherHelper::getMediumNodeOnComposedWire(const SMDS_MeshNode* n1,
1763                                                                      const SMDS_MeshNode* n2,
1764                                                                      bool                 force3d)
1765 {
1766   SMESH_TNodeXYZ p1( n1 ), p2( n2 );
1767   gp_Pnt      middle = 0.5 * p1 + 0.5 * p2;
1768   SMDS_MeshNode* n12 = AddNode( middle.X(), middle.Y(), middle.Z() );
1769
1770   // To find position on edge and 3D position for n12,
1771   // project <middle> to 2 edges and select projection most close to <middle>
1772
1773   TopoDS_Edge bestEdge;
1774   double u = 0, distMiddleProj = Precision::Infinite(), distXYZ[4], f,l;
1775
1776   // get shapes under the nodes
1777   TopoDS_Shape shape[2];
1778   int nbShapes = 0;
1779   for ( int is2nd = 0; is2nd < 2; ++is2nd )
1780   {
1781     const SMDS_MeshNode* n = is2nd ? n2 : n1;
1782     TopoDS_Shape S = GetSubShapeByNode( n, GetMeshDS() );
1783     if ( !S.IsNull() )
1784       shape[ nbShapes++ ] = S;
1785   }
1786   // get EDGEs
1787   vector< TopoDS_Shape > edges;
1788   for ( int iS = 0; iS < nbShapes; ++iS )
1789   {
1790     switch ( shape[iS].ShapeType() ) {
1791     case TopAbs_EDGE:
1792     {
1793       edges.push_back( shape[iS] );
1794       break;
1795     }
1796     case TopAbs_VERTEX:
1797     {
1798       TopoDS_Shape edge;
1799       if ( nbShapes == 2 && iS==0 && shape[1-iS].ShapeType() == TopAbs_VERTEX )
1800         edge = GetCommonAncestor( shape[iS], shape[1-iS], *myMesh, TopAbs_EDGE );
1801
1802       if ( edge.IsNull() )
1803       {
1804         PShapeIteratorPtr eIt = GetAncestors( shape[iS], *myMesh, TopAbs_EDGE );
1805         while( const TopoDS_Shape* e = eIt->next() )
1806           edges.push_back( *e );
1807       }
1808       break;
1809     }
1810     case TopAbs_FACE:
1811     {
1812       if ( nbShapes == 1 || shape[1-iS].ShapeType() < TopAbs_EDGE )
1813         for ( TopExp_Explorer e( shape[iS], TopAbs_EDGE ); e.More(); e.Next() )
1814           edges.push_back( e.Current() );
1815       break;
1816     }
1817     default:
1818       continue;
1819     }
1820   }
1821   // project to get U of projection and distance from middle to projection
1822   for ( size_t iE = 0; iE < edges.size(); ++iE )
1823   {
1824     const TopoDS_Edge& edge = TopoDS::Edge( edges[ iE ]);
1825     distXYZ[0] = distMiddleProj;
1826     double testU = 0;
1827     CheckNodeU( edge, n12, testU, 2 * BRep_Tool::Tolerance(edge), /*force=*/true, distXYZ );
1828     if ( distXYZ[0] < distMiddleProj )
1829     {
1830       distMiddleProj = distXYZ[0];
1831       u = testU;
1832       bestEdge = edge;
1833     }
1834   }
1835   // {
1836   //   // both projections failed; set n12 on the edge of n1 with U of a common vertex
1837   //   TopoDS_Vertex vCommon;
1838   //   if ( TopExp::CommonVertex( edges[0], edges[1], vCommon ))
1839   //     u = BRep_Tool::Parameter( vCommon, edges[0] );
1840   //   else
1841   //   {
1842   //     double f,l, u0 = GetNodeU( edges[0], n1 );
1843   //     BRep_Tool::Range( edges[0],f,l );
1844   //     u = ( fabs(u0-f) < fabs(u0-l) ) ? f : l;
1845   //   }
1846   //   iOkEdge = 0;
1847   //   distMiddleProj = 0;
1848   // }
1849
1850   if ( !bestEdge.IsNull() )
1851   {
1852     // move n12 to position of a successfull projection
1853     //double tol = BRep_Tool::Tolerance(edges[ iOkEdge ]);
1854     if ( !force3d /*&& distMiddleProj > 2*tol*/ )
1855     {
1856       TopLoc_Location loc;
1857       Handle(Geom_Curve) curve = BRep_Tool::Curve( bestEdge,loc,f,l );
1858       gp_Pnt p = curve->Value( u ).Transformed( loc );
1859       GetMeshDS()->MoveNode( n12, p.X(), p.Y(), p.Z() );
1860     }
1861     //if ( mySetElemOnShape ) node is not elem!
1862     {
1863       int edgeID = GetMeshDS()->ShapeToIndex( bestEdge );
1864       if ( edgeID != n12->getshapeId() )
1865         GetMeshDS()->UnSetNodeOnShape( n12 );
1866       GetMeshDS()->SetNodeOnEdge(n12, edgeID, u);
1867     }
1868   }
1869   myTLinkNodeMap.insert( make_pair( SMESH_TLink(n1,n2), n12 ));
1870
1871   return n12;
1872 }
1873
1874 //=======================================================================
1875 //function : AddNode
1876 //purpose  : Creates a node
1877 //=======================================================================
1878
1879 SMDS_MeshNode* SMESH_MesherHelper::AddNode(double x, double y, double z, int ID,
1880                                            double u, double v)
1881 {
1882   SMESHDS_Mesh * meshDS = GetMeshDS();
1883   SMDS_MeshNode* node = 0;
1884   if ( ID )
1885     node = meshDS->AddNodeWithID( x, y, z, ID );
1886   else
1887     node = meshDS->AddNode( x, y, z );
1888   if ( mySetElemOnShape && myShapeID > 0 ) { // node is not elem ?
1889     switch ( myShape.ShapeType() ) {
1890     case TopAbs_SOLID:  meshDS->SetNodeInVolume( node, myShapeID);       break;
1891     case TopAbs_SHELL:  meshDS->SetNodeInVolume( node, myShapeID);       break;
1892     case TopAbs_FACE:   meshDS->SetNodeOnFace(   node, myShapeID, u, v); break;
1893     case TopAbs_EDGE:   meshDS->SetNodeOnEdge(   node, myShapeID, u);    break;
1894     case TopAbs_VERTEX: meshDS->SetNodeOnVertex( node, myShapeID);       break;
1895     default: ;
1896     }
1897   }
1898   return node;
1899 }
1900
1901 //=======================================================================
1902 //function : AddEdge
1903 //purpose  : Creates quadratic or linear edge
1904 //=======================================================================
1905
1906 SMDS_MeshEdge* SMESH_MesherHelper::AddEdge(const SMDS_MeshNode* n1,
1907                                            const SMDS_MeshNode* n2,
1908                                            const int            id,
1909                                            const bool           force3d)
1910 {
1911   SMESHDS_Mesh * meshDS = GetMeshDS();
1912   
1913   SMDS_MeshEdge* edge = 0;
1914   if (myCreateQuadratic) {
1915     const SMDS_MeshNode* n12 = GetMediumNode(n1,n2,force3d);
1916     if(id)
1917       edge = meshDS->AddEdgeWithID(n1, n2, n12, id);
1918     else
1919       edge = meshDS->AddEdge(n1, n2, n12);
1920   }
1921   else {
1922     if(id)
1923       edge = meshDS->AddEdgeWithID(n1, n2, id);
1924     else
1925       edge = meshDS->AddEdge(n1, n2);
1926   }
1927
1928   if ( mySetElemOnShape && myShapeID > 0 )
1929     meshDS->SetMeshElementOnShape( edge, myShapeID );
1930
1931   return edge;
1932 }
1933
1934 //=======================================================================
1935 //function : AddFace
1936 //purpose  : Creates quadratic or linear triangle
1937 //=======================================================================
1938
1939 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
1940                                            const SMDS_MeshNode* n2,
1941                                            const SMDS_MeshNode* n3,
1942                                            const int id,
1943                                            const bool force3d)
1944 {
1945   SMESHDS_Mesh * meshDS = GetMeshDS();
1946   SMDS_MeshFace* elem = 0;
1947
1948   if( n1==n2 || n2==n3 || n3==n1 )
1949     return elem;
1950
1951   if(!myCreateQuadratic) {
1952     if(id)
1953       elem = meshDS->AddFaceWithID(n1, n2, n3, id);
1954     else
1955       elem = meshDS->AddFace(n1, n2, n3);
1956   }
1957   else {
1958     const SMDS_MeshNode* n12 = GetMediumNode( n1, n2, force3d, TopAbs_FACE );
1959     const SMDS_MeshNode* n23 = GetMediumNode( n2, n3, force3d, TopAbs_FACE );
1960     const SMDS_MeshNode* n31 = GetMediumNode( n3, n1, force3d, TopAbs_FACE );
1961     if(myCreateBiQuadratic)
1962     {
1963      const SMDS_MeshNode* nCenter = GetCentralNode(n1, n2, n3, n12, n23, n31, force3d);
1964      if(id)
1965        elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, nCenter, id);
1966      else
1967        elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31, nCenter);
1968     }
1969     else
1970     {
1971       if(id)
1972         elem = meshDS->AddFaceWithID(n1, n2, n3, n12, n23, n31, id);
1973       else
1974         elem = meshDS->AddFace(n1, n2, n3, n12, n23, n31);
1975     }
1976   }
1977   if ( mySetElemOnShape && myShapeID > 0 )
1978     meshDS->SetMeshElementOnShape( elem, myShapeID );
1979
1980   return elem;
1981 }
1982
1983 //=======================================================================
1984 //function : AddFace
1985 //purpose  : Creates bi-quadratic, quadratic or linear quadrangle
1986 //=======================================================================
1987
1988 SMDS_MeshFace* SMESH_MesherHelper::AddFace(const SMDS_MeshNode* n1,
1989                                            const SMDS_MeshNode* n2,
1990                                            const SMDS_MeshNode* n3,
1991                                            const SMDS_MeshNode* n4,
1992                                            const int            id,
1993                                            const bool           force3d)
1994 {
1995   SMESHDS_Mesh * meshDS = GetMeshDS();
1996   SMDS_MeshFace* elem = 0;
1997
1998   if( n1==n2 ) {
1999     return AddFace(n1,n3,n4,id,force3d);
2000   }
2001   if( n1==n3 ) {
2002     return AddFace(n1,n2,n4,id,force3d);
2003   }
2004   if( n1==n4 ) {
2005     return AddFace(n1,n2,n3,id,force3d);
2006   }
2007   if( n2==n3 ) {
2008     return AddFace(n1,n2,n4,id,force3d);
2009   }
2010   if( n2==n4 ) {
2011     return AddFace(n1,n2,n3,id,force3d);
2012   }
2013   if( n3==n4 ) {
2014     return AddFace(n1,n2,n3,id,force3d);
2015   }
2016
2017   if(!myCreateQuadratic) {
2018     if(id)
2019       elem = meshDS->AddFaceWithID(n1, n2, n3, n4, id);
2020     else
2021       elem = meshDS->AddFace(n1, n2, n3, n4);
2022   }
2023   else {
2024     const SMDS_MeshNode* n12 = GetMediumNode( n1, n2, force3d, TopAbs_FACE );
2025     const SMDS_MeshNode* n23 = GetMediumNode( n2, n3, force3d, TopAbs_FACE );
2026     const SMDS_MeshNode* n34 = GetMediumNode( n3, n4, force3d, TopAbs_FACE );
2027     const SMDS_MeshNode* n41 = GetMediumNode( n4, n1, force3d, TopAbs_FACE );
2028     if(myCreateBiQuadratic)
2029     {
2030      const SMDS_MeshNode* nCenter = GetCentralNode(n1, n2, n3, n4, n12, n23, n34, n41, force3d);
2031      if(id)
2032        elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, nCenter, id);
2033      else
2034        elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41, nCenter);
2035     }
2036     else
2037     {
2038       if(id)
2039         elem = meshDS->AddFaceWithID(n1, n2, n3, n4, n12, n23, n34, n41, id);
2040       else
2041         elem = meshDS->AddFace(n1, n2, n3, n4, n12, n23, n34, n41);
2042     }
2043   }
2044   if ( mySetElemOnShape && myShapeID > 0 )
2045     meshDS->SetMeshElementOnShape( elem, myShapeID );
2046
2047   return elem;
2048 }
2049
2050 //=======================================================================
2051 //function : AddPolygonalFace
2052 //purpose  : Creates polygon, with additional nodes in quadratic mesh
2053 //=======================================================================
2054
2055 SMDS_MeshFace* SMESH_MesherHelper::AddPolygonalFace (const vector<const SMDS_MeshNode*>& nodes,
2056                                                      const int                           id,
2057                                                      const bool                          force3d)
2058 {
2059   SMESHDS_Mesh * meshDS = GetMeshDS();
2060   SMDS_MeshFace* elem = 0;
2061
2062   if(!myCreateQuadratic)
2063   {
2064     if(id)
2065       elem = meshDS->AddPolygonalFaceWithID(nodes, id);
2066     else
2067       elem = meshDS->AddPolygonalFace(nodes);
2068   }
2069   else
2070   {
2071     vector<const SMDS_MeshNode*> newNodes( nodes.size() * 2 );
2072     newNodes = nodes;
2073     for ( int i = 0; i < nodes.size(); ++i )
2074     {
2075       const SMDS_MeshNode* n1 = nodes[i];
2076       const SMDS_MeshNode* n2 = nodes[(i+1)%nodes.size()];
2077       const SMDS_MeshNode* n12 = GetMediumNode( n1, n2, force3d, TopAbs_FACE );
2078       newNodes.push_back( n12 );
2079     }
2080     if(id)
2081       elem = meshDS->AddQuadPolygonalFaceWithID(newNodes, id);
2082     else
2083       elem = meshDS->AddQuadPolygonalFace(newNodes);
2084   }
2085   if ( mySetElemOnShape && myShapeID > 0 )
2086     meshDS->SetMeshElementOnShape( elem, myShapeID );
2087
2088   return elem;
2089 }
2090
2091 //=======================================================================
2092 //function : AddVolume
2093 //purpose  : Creates quadratic or linear prism
2094 //=======================================================================
2095
2096 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
2097                                                const SMDS_MeshNode* n2,
2098                                                const SMDS_MeshNode* n3,
2099                                                const SMDS_MeshNode* n4,
2100                                                const SMDS_MeshNode* n5,
2101                                                const SMDS_MeshNode* n6,
2102                                                const int id,
2103                                                const bool force3d)
2104 {
2105   SMESHDS_Mesh * meshDS = GetMeshDS();
2106   SMDS_MeshVolume* elem = 0;
2107   if(!myCreateQuadratic) {
2108     if(id)
2109       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, id);
2110     else
2111       elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6);
2112   }
2113   else {
2114     const SMDS_MeshNode* n12 = GetMediumNode( n1, n2, force3d, TopAbs_SOLID );
2115     const SMDS_MeshNode* n23 = GetMediumNode( n2, n3, force3d, TopAbs_SOLID );
2116     const SMDS_MeshNode* n31 = GetMediumNode( n3, n1, force3d, TopAbs_SOLID );
2117
2118     const SMDS_MeshNode* n45 = GetMediumNode( n4, n5, force3d, TopAbs_SOLID );
2119     const SMDS_MeshNode* n56 = GetMediumNode( n5, n6, force3d, TopAbs_SOLID );
2120     const SMDS_MeshNode* n64 = GetMediumNode( n6, n4, force3d, TopAbs_SOLID );
2121
2122     const SMDS_MeshNode* n14 = GetMediumNode( n1, n4, force3d, TopAbs_SOLID );
2123     const SMDS_MeshNode* n25 = GetMediumNode( n2, n5, force3d, TopAbs_SOLID );
2124     const SMDS_MeshNode* n36 = GetMediumNode( n3, n6, force3d, TopAbs_SOLID );
2125
2126     if(id)
2127       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6,
2128                                      n12, n23, n31, n45, n56, n64, n14, n25, n36, id);
2129     else
2130       elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6,
2131                                n12, n23, n31, n45, n56, n64, n14, n25, n36);
2132   }
2133   if ( mySetElemOnShape && myShapeID > 0 )
2134     meshDS->SetMeshElementOnShape( elem, myShapeID );
2135
2136   return elem;
2137 }
2138
2139 //=======================================================================
2140 //function : AddVolume
2141 //purpose  : Creates quadratic or linear tetrahedron
2142 //=======================================================================
2143
2144 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
2145                                                const SMDS_MeshNode* n2,
2146                                                const SMDS_MeshNode* n3,
2147                                                const SMDS_MeshNode* n4,
2148                                                const int id,
2149                                                const bool force3d)
2150 {
2151   SMESHDS_Mesh * meshDS = GetMeshDS();
2152   SMDS_MeshVolume* elem = 0;
2153   if(!myCreateQuadratic) {
2154     if(id)
2155       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, id);
2156     else
2157       elem = meshDS->AddVolume(n1, n2, n3, n4);
2158   }
2159   else {
2160     const SMDS_MeshNode* n12 = GetMediumNode( n1, n2, force3d, TopAbs_SOLID );
2161     const SMDS_MeshNode* n23 = GetMediumNode( n2, n3, force3d, TopAbs_SOLID );
2162     const SMDS_MeshNode* n31 = GetMediumNode( n3, n1, force3d, TopAbs_SOLID );
2163
2164     const SMDS_MeshNode* n14 = GetMediumNode( n1, n4, force3d, TopAbs_SOLID );
2165     const SMDS_MeshNode* n24 = GetMediumNode( n2, n4, force3d, TopAbs_SOLID );
2166     const SMDS_MeshNode* n34 = GetMediumNode( n3, n4, force3d, TopAbs_SOLID );
2167
2168     if(id)
2169       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34, id);
2170     else
2171       elem = meshDS->AddVolume(n1, n2, n3, n4, n12, n23, n31, n14, n24, n34);
2172   }
2173   if ( mySetElemOnShape && myShapeID > 0 )
2174     meshDS->SetMeshElementOnShape( elem, myShapeID );
2175
2176   return elem;
2177 }
2178
2179 //=======================================================================
2180 //function : AddVolume
2181 //purpose  : Creates quadratic or linear pyramid
2182 //=======================================================================
2183
2184 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
2185                                                const SMDS_MeshNode* n2,
2186                                                const SMDS_MeshNode* n3,
2187                                                const SMDS_MeshNode* n4,
2188                                                const SMDS_MeshNode* n5,
2189                                                const int id,
2190                                                const bool force3d)
2191 {
2192   SMDS_MeshVolume* elem = 0;
2193   if(!myCreateQuadratic) {
2194     if(id)
2195       elem = GetMeshDS()->AddVolumeWithID(n1, n2, n3, n4, n5, id);
2196     else
2197       elem = GetMeshDS()->AddVolume(n1, n2, n3, n4, n5);
2198   }
2199   else {
2200     const SMDS_MeshNode* n12 = GetMediumNode( n1, n2, force3d, TopAbs_SOLID );
2201     const SMDS_MeshNode* n23 = GetMediumNode( n2, n3, force3d, TopAbs_SOLID );
2202     const SMDS_MeshNode* n34 = GetMediumNode( n3, n4, force3d, TopAbs_SOLID );
2203     const SMDS_MeshNode* n41 = GetMediumNode( n4, n1, force3d, TopAbs_SOLID );
2204
2205     const SMDS_MeshNode* n15 = GetMediumNode( n1, n5, force3d, TopAbs_SOLID );
2206     const SMDS_MeshNode* n25 = GetMediumNode( n2, n5, force3d, TopAbs_SOLID );
2207     const SMDS_MeshNode* n35 = GetMediumNode( n3, n5, force3d, TopAbs_SOLID );
2208     const SMDS_MeshNode* n45 = GetMediumNode( n4, n5, force3d, TopAbs_SOLID );
2209
2210     if(id)
2211       elem = GetMeshDS()->AddVolumeWithID ( n1,  n2,  n3,  n4,  n5,
2212                                             n12, n23, n34, n41,
2213                                             n15, n25, n35, n45,
2214                                             id);
2215     else
2216       elem = GetMeshDS()->AddVolume( n1,  n2,  n3,  n4,  n5,
2217                                      n12, n23, n34, n41,
2218                                      n15, n25, n35, n45);
2219   }
2220   if ( mySetElemOnShape && myShapeID > 0 )
2221     GetMeshDS()->SetMeshElementOnShape( elem, myShapeID );
2222
2223   return elem;
2224 }
2225
2226 //=======================================================================
2227 //function : AddVolume
2228 //purpose  : Creates tri-quadratic, quadratic or linear hexahedron
2229 //=======================================================================
2230
2231 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
2232                                                const SMDS_MeshNode* n2,
2233                                                const SMDS_MeshNode* n3,
2234                                                const SMDS_MeshNode* n4,
2235                                                const SMDS_MeshNode* n5,
2236                                                const SMDS_MeshNode* n6,
2237                                                const SMDS_MeshNode* n7,
2238                                                const SMDS_MeshNode* n8,
2239                                                const int id,
2240                                                const bool force3d)
2241 {
2242   SMESHDS_Mesh * meshDS = GetMeshDS();
2243   SMDS_MeshVolume* elem = 0;
2244   if(!myCreateQuadratic) {
2245     if(id)
2246       elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8, id);
2247     else
2248       elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8);
2249   }
2250   else {
2251     const SMDS_MeshNode* n12 = GetMediumNode( n1, n2, force3d, TopAbs_SOLID );
2252     const SMDS_MeshNode* n23 = GetMediumNode( n2, n3, force3d, TopAbs_SOLID );
2253     const SMDS_MeshNode* n34 = GetMediumNode( n3, n4, force3d, TopAbs_SOLID );
2254     const SMDS_MeshNode* n41 = GetMediumNode( n4, n1, force3d, TopAbs_SOLID );
2255
2256     const SMDS_MeshNode* n56 = GetMediumNode( n5, n6, force3d, TopAbs_SOLID );
2257     const SMDS_MeshNode* n67 = GetMediumNode( n6, n7, force3d, TopAbs_SOLID );
2258     const SMDS_MeshNode* n78 = GetMediumNode( n7, n8, force3d, TopAbs_SOLID );
2259     const SMDS_MeshNode* n85 = GetMediumNode( n8, n5, force3d, TopAbs_SOLID );
2260
2261     const SMDS_MeshNode* n15 = GetMediumNode( n1, n5, force3d, TopAbs_SOLID );
2262     const SMDS_MeshNode* n26 = GetMediumNode( n2, n6, force3d, TopAbs_SOLID );
2263     const SMDS_MeshNode* n37 = GetMediumNode( n3, n7, force3d, TopAbs_SOLID );
2264     const SMDS_MeshNode* n48 = GetMediumNode( n4, n8, force3d, TopAbs_SOLID );
2265     if ( myCreateBiQuadratic )
2266     {
2267       const SMDS_MeshNode* n1234 = GetCentralNode( n1,n2,n3,n4,n12,n23,n34,n41,force3d );
2268       const SMDS_MeshNode* n1256 = GetCentralNode( n1,n2,n5,n6,n12,n26,n56,n15,force3d );
2269       const SMDS_MeshNode* n2367 = GetCentralNode( n2,n3,n6,n7,n23,n37,n67,n26,force3d );
2270       const SMDS_MeshNode* n3478 = GetCentralNode( n3,n4,n7,n8,n34,n48,n78,n37,force3d );
2271       const SMDS_MeshNode* n1458 = GetCentralNode( n1,n4,n5,n8,n41,n48,n15,n85,force3d );
2272       const SMDS_MeshNode* n5678 = GetCentralNode( n5,n6,n7,n8,n56,n67,n78,n85,force3d );
2273
2274       vector<gp_XYZ> pointsOnShapes( SMESH_Block::ID_Shell );
2275
2276       pointsOnShapes[ SMESH_Block::ID_V000 ] = SMESH_TNodeXYZ( n4 );
2277       pointsOnShapes[ SMESH_Block::ID_V100 ] = SMESH_TNodeXYZ( n8 );
2278       pointsOnShapes[ SMESH_Block::ID_V010 ] = SMESH_TNodeXYZ( n3 );
2279       pointsOnShapes[ SMESH_Block::ID_V110 ] = SMESH_TNodeXYZ( n7 );
2280       pointsOnShapes[ SMESH_Block::ID_V001 ] = SMESH_TNodeXYZ( n1 );
2281       pointsOnShapes[ SMESH_Block::ID_V101 ] = SMESH_TNodeXYZ( n5 );
2282       pointsOnShapes[ SMESH_Block::ID_V011 ] = SMESH_TNodeXYZ( n2 );
2283       pointsOnShapes[ SMESH_Block::ID_V111 ] = SMESH_TNodeXYZ( n6 );
2284
2285       pointsOnShapes[ SMESH_Block::ID_Ex00 ] = SMESH_TNodeXYZ( n48 );
2286       pointsOnShapes[ SMESH_Block::ID_Ex10 ] = SMESH_TNodeXYZ( n37 );
2287       pointsOnShapes[ SMESH_Block::ID_E0y0 ] = SMESH_TNodeXYZ( n15 );
2288       pointsOnShapes[ SMESH_Block::ID_E1y0 ] = SMESH_TNodeXYZ( n26 );
2289       pointsOnShapes[ SMESH_Block::ID_Ex01 ] = SMESH_TNodeXYZ( n34 );
2290       pointsOnShapes[ SMESH_Block::ID_Ex11 ] = SMESH_TNodeXYZ( n78 );
2291       pointsOnShapes[ SMESH_Block::ID_E0y1 ] = SMESH_TNodeXYZ( n12 );
2292       pointsOnShapes[ SMESH_Block::ID_E1y1 ] = SMESH_TNodeXYZ( n56 );
2293       pointsOnShapes[ SMESH_Block::ID_E00z ] = SMESH_TNodeXYZ( n41 );
2294       pointsOnShapes[ SMESH_Block::ID_E10z ] = SMESH_TNodeXYZ( n85 );
2295       pointsOnShapes[ SMESH_Block::ID_E01z ] = SMESH_TNodeXYZ( n23 );
2296       pointsOnShapes[ SMESH_Block::ID_E11z ] = SMESH_TNodeXYZ( n67 );
2297
2298       pointsOnShapes[ SMESH_Block::ID_Fxy0 ] = SMESH_TNodeXYZ( n3478 );
2299       pointsOnShapes[ SMESH_Block::ID_Fxy1 ] = SMESH_TNodeXYZ( n1256 );
2300       pointsOnShapes[ SMESH_Block::ID_Fx0z ] = SMESH_TNodeXYZ( n1458 );
2301       pointsOnShapes[ SMESH_Block::ID_Fx1z ] = SMESH_TNodeXYZ( n2367 );
2302       pointsOnShapes[ SMESH_Block::ID_F0yz ] = SMESH_TNodeXYZ( n1234 );
2303       pointsOnShapes[ SMESH_Block::ID_F1yz ] = SMESH_TNodeXYZ( n5678 );
2304
2305       gp_XYZ centerCube(0.5, 0.5, 0.5);
2306       gp_XYZ nCenterElem;
2307       SMESH_Block::ShellPoint( centerCube, pointsOnShapes, nCenterElem );
2308       const SMDS_MeshNode* nCenter =
2309         meshDS->AddNode( nCenterElem.X(), nCenterElem.Y(), nCenterElem.Z() );
2310       meshDS->SetNodeInVolume( nCenter, myShapeID );
2311
2312       if(id)
2313         elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
2314                                        n12, n23, n34, n41, n56, n67,
2315                                        n78, n85, n15, n26, n37, n48,
2316                                        n1234, n1256, n2367, n3478, n1458, n5678, nCenter, id);
2317       else
2318         elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
2319                                  n12, n23, n34, n41, n56, n67,
2320                                  n78, n85, n15, n26, n37, n48,
2321                                  n1234, n1256, n2367, n3478, n1458, n5678, nCenter);
2322     }
2323     else
2324     {
2325       if(id)
2326         elem = meshDS->AddVolumeWithID(n1, n2, n3, n4, n5, n6, n7, n8,
2327                                        n12, n23, n34, n41, n56, n67,
2328                                        n78, n85, n15, n26, n37, n48, id);
2329       else
2330         elem = meshDS->AddVolume(n1, n2, n3, n4, n5, n6, n7, n8,
2331                                  n12, n23, n34, n41, n56, n67,
2332                                  n78, n85, n15, n26, n37, n48);
2333     }
2334   }
2335   if ( mySetElemOnShape && myShapeID > 0 )
2336     meshDS->SetMeshElementOnShape( elem, myShapeID );
2337
2338   return elem;
2339 }
2340
2341 //=======================================================================
2342 //function : AddVolume
2343 //purpose  : Creates LINEAR!!!!!!!!! octahedron
2344 //=======================================================================
2345
2346 SMDS_MeshVolume* SMESH_MesherHelper::AddVolume(const SMDS_MeshNode* n1,
2347                                                const SMDS_MeshNode* n2,
2348                                                const SMDS_MeshNode* n3,
2349                                                const SMDS_MeshNode* n4,
2350                                                const SMDS_MeshNode* n5,
2351                                                const SMDS_MeshNode* n6,
2352                                                const SMDS_MeshNode* n7,
2353                                                const SMDS_MeshNode* n8,
2354                                                const SMDS_MeshNode* n9,
2355                                                const SMDS_MeshNode* n10,
2356                                                const SMDS_MeshNode* n11,
2357                                                const SMDS_MeshNode* n12,
2358                                                const int id, 
2359                                                bool force3d)
2360 {
2361   SMESHDS_Mesh * meshDS = GetMeshDS();
2362   SMDS_MeshVolume* elem = 0;
2363   if(id)
2364     elem = meshDS->AddVolumeWithID(n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,id);
2365   else
2366     elem = meshDS->AddVolume(n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12);
2367   if ( mySetElemOnShape && myShapeID > 0 )
2368     meshDS->SetMeshElementOnShape( elem, myShapeID );
2369   return elem;
2370 }
2371
2372 //=======================================================================
2373 //function : AddPolyhedralVolume
2374 //purpose  : Creates polyhedron. In quadratic mesh, adds medium nodes
2375 //=======================================================================
2376
2377 SMDS_MeshVolume*
2378 SMESH_MesherHelper::AddPolyhedralVolume (const std::vector<const SMDS_MeshNode*>& nodes,
2379                                          const std::vector<int>&                  quantities,
2380                                          const int                                id,
2381                                          const bool                               force3d)
2382 {
2383   SMESHDS_Mesh * meshDS = GetMeshDS();
2384   SMDS_MeshVolume* elem = 0;
2385   if(!myCreateQuadratic)
2386   {
2387     if(id)
2388       elem = meshDS->AddPolyhedralVolumeWithID(nodes, quantities, id);
2389     else
2390       elem = meshDS->AddPolyhedralVolume(nodes, quantities);
2391   }
2392   else
2393   {
2394     vector<const SMDS_MeshNode*> newNodes;
2395     vector<int> newQuantities;
2396     for ( int iFace=0, iN=0; iFace < quantities.size(); ++iFace)
2397     {
2398       int nbNodesInFace = quantities[iFace];
2399       newQuantities.push_back(0);
2400       for ( int i = 0; i < nbNodesInFace; ++i )
2401       {
2402         const SMDS_MeshNode* n1 = nodes[ iN + i ];
2403         newNodes.push_back( n1 );
2404         newQuantities.back()++;
2405         
2406         const SMDS_MeshNode* n2 = nodes[ iN + ( i+1==nbNodesInFace ? 0 : i+1 )];
2407 //         if ( n1->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE &&
2408 //              n2->GetPosition()->GetTypeOfPosition() != SMDS_TOP_3DSPACE )
2409         {
2410           const SMDS_MeshNode* n12 = GetMediumNode( n1, n2, force3d, TopAbs_SOLID );
2411           newNodes.push_back( n12 );
2412           newQuantities.back()++;
2413         }
2414       }
2415       iN += nbNodesInFace;
2416     }
2417     if(id)
2418       elem = meshDS->AddPolyhedralVolumeWithID( newNodes, newQuantities, id );
2419     else
2420       elem = meshDS->AddPolyhedralVolume( newNodes, newQuantities );
2421   }
2422   if ( mySetElemOnShape && myShapeID > 0 )
2423     meshDS->SetMeshElementOnShape( elem, myShapeID );
2424
2425   return elem;
2426 }
2427
2428 namespace
2429 {
2430   //================================================================================
2431   /*!
2432    * \brief Check if a node belongs to any face of sub-mesh
2433    */
2434   //================================================================================
2435
2436   bool isNodeInSubMesh( const SMDS_MeshNode* n, const SMESHDS_SubMesh* sm )
2437   {
2438     SMDS_ElemIteratorPtr fIt = n->GetInverseElementIterator( SMDSAbs_Face );
2439     while ( fIt->more() )
2440       if ( sm->Contains( fIt->next() ))
2441         return true;
2442     return false;
2443   }
2444 }
2445
2446 //=======================================================================
2447 //function : IsSameElemGeometry
2448 //purpose  : Returns true if all elements of a sub-mesh are of same shape
2449 //=======================================================================
2450
2451 bool SMESH_MesherHelper::IsSameElemGeometry(const SMESHDS_SubMesh* smDS,
2452                                             SMDSAbs_GeometryType   shape,
2453                                             const bool             nullSubMeshRes)
2454 {
2455   if ( !smDS ) return nullSubMeshRes;
2456
2457   SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
2458   while ( elemIt->more() ) {
2459     const SMDS_MeshElement* e = elemIt->next();
2460     if ( e->GetGeomType() != shape )
2461       return false;
2462   }
2463   return true;
2464 }
2465
2466 //=======================================================================
2467 //function : LoadNodeColumns
2468 //purpose  : Load nodes bound to face into a map of node columns
2469 //=======================================================================
2470
2471 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
2472                                          const TopoDS_Face& theFace,
2473                                          const TopoDS_Edge& theBaseEdge,
2474                                          SMESHDS_Mesh*      theMesh,
2475                                          SMESH_ProxyMesh*   theProxyMesh)
2476 {
2477   return LoadNodeColumns(theParam2ColumnMap,
2478                          theFace,
2479                          std::list<TopoDS_Edge>(1,theBaseEdge),
2480                          theMesh,
2481                          theProxyMesh);
2482 }
2483
2484 //=======================================================================
2485 //function : LoadNodeColumns
2486 //purpose  : Load nodes bound to face into a map of node columns
2487 //=======================================================================
2488
2489 bool SMESH_MesherHelper::LoadNodeColumns(TParam2ColumnMap &            theParam2ColumnMap,
2490                                          const TopoDS_Face&            theFace,
2491                                          const std::list<TopoDS_Edge>& theBaseSide,
2492                                          SMESHDS_Mesh*                 theMesh,
2493                                          SMESH_ProxyMesh*              theProxyMesh)
2494 {
2495   // get a right sub-mesh of theFace
2496
2497   const SMESHDS_SubMesh* faceSubMesh = 0;
2498   if ( theProxyMesh )
2499   {
2500     faceSubMesh = theProxyMesh->GetSubMesh( theFace );
2501     if ( !faceSubMesh ||
2502          faceSubMesh->NbElements() == 0 ||
2503          theProxyMesh->IsTemporary( faceSubMesh->GetElements()->next() ))
2504     {
2505       // can use a proxy sub-mesh with not temporary elements only
2506       faceSubMesh = 0;
2507       theProxyMesh = 0;
2508     }
2509   }
2510   if ( !faceSubMesh )
2511     faceSubMesh = theMesh->MeshElements( theFace );
2512   if ( !faceSubMesh || faceSubMesh->NbElements() == 0 )
2513     return false;
2514
2515   if ( theParam2ColumnMap.empty() )
2516   {
2517     // get data of edges for normalization of params
2518     vector< double > length;
2519     double fullLen = 0;
2520     list<TopoDS_Edge>::const_iterator edge;
2521     {
2522       for ( edge = theBaseSide.begin(); edge != theBaseSide.end(); ++edge )
2523       {
2524         double len = std::max( 1e-10, SMESH_Algo::EdgeLength( *edge ));
2525         fullLen += len;
2526         length.push_back( len );
2527       }
2528     }
2529
2530     // get nodes on theBaseEdge sorted by param on edge and initialize theParam2ColumnMap with them
2531     edge = theBaseSide.begin();
2532     for ( int iE = 0; edge != theBaseSide.end(); ++edge, ++iE )
2533     {
2534       map< double, const SMDS_MeshNode*> sortedBaseNN;
2535       SMESH_Algo::GetSortedNodesOnEdge( theMesh, *edge,/*noMedium=*/true, sortedBaseNN );
2536
2537       map< double, const SMDS_MeshNode*>::iterator u_n;
2538       // pb with mesh_Projection_2D_00/A1 fixed by adding expectedSupport arg to GetMediumPos()
2539       // so the following solution is commented (hope forever :)
2540       //
2541       // SMESH_Algo::GetSortedNodesOnEdge( theMesh, *edge,/*noMedium=*/true, sortedBaseNN,
2542       // // SMDSAbs_Edge here is needed to be coherent with
2543       // // StdMeshers_FaceSide used by Quadrangle to get nodes
2544       // // on EDGE; else pb in mesh_Projection_2D_00/A1 where a
2545       // // medium node on EDGE is medium in a triangle but not
2546       // // in a segment
2547       // SMDSAbs_Edge );
2548       // if ( faceSubMesh->GetElements()->next()->IsQuadratic() )
2549       //   // filter off nodes medium in faces on theFace (same pb with mesh_Projection_2D_00/A1)
2550       //   for ( u_n = sortedBaseNN.begin(); u_n != sortedBaseNN.end() ;     )
2551       //   {
2552       //     const SMDS_MeshNode* node = u_n->second;
2553       //     SMDS_ElemIteratorPtr faceIt = node->GetInverseElementIterator( SMDSAbs_Face );
2554       //     if ( faceIt->more() && node ) {
2555       //       const SMDS_MeshElement* face = faceIt->next();
2556       //       if ( faceSubMesh->Contains( face ) && face->IsMediumNode( node ))
2557       //         node = 0;
2558       //     }
2559       //     if ( !node )
2560       //       sortedBaseNN.erase( u_n++ );
2561       //     else
2562       //       ++u_n;
2563       //   }
2564       if ( sortedBaseNN.empty() ) continue;
2565
2566       u_n = sortedBaseNN.begin();
2567       if ( theProxyMesh ) // from sortedBaseNN remove nodes not shared by faces of faceSubMesh
2568       {
2569         const SMDS_MeshNode* n1 = (++sortedBaseNN.begin())->second;
2570         const SMDS_MeshNode* n2 = (++sortedBaseNN.rbegin())->second;
2571         bool allNodesAreProxy = ( n1 != theProxyMesh->GetProxyNode( n1 ) &&
2572                                   n2 != theProxyMesh->GetProxyNode( n2 ));
2573         if ( allNodesAreProxy )
2574           for ( u_n = sortedBaseNN.begin(); u_n != sortedBaseNN.end(); u_n++ )
2575             u_n->second = theProxyMesh->GetProxyNode( u_n->second );
2576
2577         if ( u_n = sortedBaseNN.begin(), !isNodeInSubMesh( u_n->second, faceSubMesh ))
2578         {
2579           while ( ++u_n != sortedBaseNN.end() && !isNodeInSubMesh( u_n->second, faceSubMesh ));
2580           sortedBaseNN.erase( sortedBaseNN.begin(), u_n );
2581         }
2582         if ( !sortedBaseNN.empty() )
2583           if ( u_n = --sortedBaseNN.end(), !isNodeInSubMesh( u_n->second, faceSubMesh ))
2584           {
2585             while ( u_n != sortedBaseNN.begin() && !isNodeInSubMesh( (--u_n)->second, faceSubMesh ));
2586             sortedBaseNN.erase( ++u_n, sortedBaseNN.end() );
2587           }
2588         if ( sortedBaseNN.empty() ) continue;
2589       }
2590
2591       double f, l;
2592       BRep_Tool::Range( *edge, f, l );
2593       if ( edge->Orientation() == TopAbs_REVERSED ) std::swap( f, l );
2594       const double coeff = 1. / ( l - f ) * length[iE] / fullLen;
2595       const double prevPar = theParam2ColumnMap.empty() ? 0 : theParam2ColumnMap.rbegin()->first;
2596       for ( u_n = sortedBaseNN.begin(); u_n != sortedBaseNN.end(); u_n++ )
2597       {
2598         double par = prevPar + coeff * ( u_n->first - f );
2599         TParam2ColumnMap::iterator u2nn =
2600           theParam2ColumnMap.insert( theParam2ColumnMap.end(), make_pair( par, TNodeColumn()));
2601         u2nn->second.push_back( u_n->second );
2602       }
2603     }
2604     if ( theParam2ColumnMap.size() < 2 )
2605       return false;
2606   }
2607
2608   // nb rows of nodes
2609   int prevNbRows     = theParam2ColumnMap.begin()->second.size(); // current, at least 1 here
2610   int expectedNbRows = faceSubMesh->NbElements() / ( theParam2ColumnMap.size()-1 ); // to be added
2611
2612   // fill theParam2ColumnMap column by column by passing from nodes on
2613   // theBaseEdge up via mesh faces on theFace
2614
2615   TParam2ColumnMap::iterator par_nVec_1, par_nVec_2;
2616   par_nVec_2 = theParam2ColumnMap.begin();
2617   par_nVec_1 = par_nVec_2++;
2618   TIDSortedElemSet emptySet, avoidSet;
2619   for ( ; par_nVec_2 != theParam2ColumnMap.end(); ++par_nVec_1, ++par_nVec_2 )
2620   {
2621     vector<const SMDS_MeshNode*>& nCol1 = par_nVec_1->second;
2622     vector<const SMDS_MeshNode*>& nCol2 = par_nVec_2->second;
2623     nCol1.resize( prevNbRows + expectedNbRows );
2624     nCol2.resize( prevNbRows + expectedNbRows );
2625
2626     int i1, i2, foundNbRows = 0;
2627     const SMDS_MeshNode *n1 = nCol1[ prevNbRows-1 ];
2628     const SMDS_MeshNode *n2 = nCol2[ prevNbRows-1 ];
2629     // find face sharing node n1 and n2 and belonging to faceSubMesh
2630     while ( const SMDS_MeshElement* face =
2631             SMESH_MeshAlgos::FindFaceInSet( n1, n2, emptySet, avoidSet, &i1, &i2))
2632     {
2633       if ( faceSubMesh->Contains( face ))
2634       {
2635         int nbNodes = face->NbCornerNodes();
2636         if ( nbNodes != 4 )
2637           return false;
2638         if ( foundNbRows + 1 > expectedNbRows )
2639           return false;
2640         n1 = face->GetNode( (i2+2) % 4 ); // opposite corner of quadrangle face
2641         n2 = face->GetNode( (i1+2) % 4 );
2642         nCol1[ prevNbRows + foundNbRows] = n1;
2643         nCol2[ prevNbRows + foundNbRows] = n2;
2644         ++foundNbRows;
2645       }
2646       avoidSet.insert( face );
2647     }
2648     if ( foundNbRows != expectedNbRows )
2649       return false;
2650     avoidSet.clear();
2651   }
2652   return ( theParam2ColumnMap.size() > 1 &&
2653            theParam2ColumnMap.begin()->second.size() == prevNbRows + expectedNbRows );
2654 }
2655
2656 namespace
2657 {
2658   //================================================================================
2659   /*!
2660    * \brief Return true if a node is at a corner of a 2D structured mesh of FACE
2661    */
2662   //================================================================================
2663
2664   bool isCornerOfStructure( const SMDS_MeshNode*   n,
2665                             const SMESHDS_SubMesh* faceSM,
2666                             SMESH_MesherHelper&    faceAnalyser )
2667   {
2668     int nbFacesInSM = 0;
2669     if ( n ) {
2670       SMDS_ElemIteratorPtr fIt = n->GetInverseElementIterator( SMDSAbs_Face );
2671       while ( fIt->more() )
2672         nbFacesInSM += faceSM->Contains( fIt->next() );
2673     }
2674     if ( nbFacesInSM == 1 )
2675       return true;
2676
2677     if ( nbFacesInSM == 2 && n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX )
2678     {
2679       return faceAnalyser.IsRealSeam( n->getshapeId() );
2680     }
2681     return false;
2682   }
2683 }
2684
2685 //=======================================================================
2686 //function : IsStructured
2687 //purpose  : Return true if 2D mesh on FACE is a structured rectangle
2688 //=======================================================================
2689
2690 bool SMESH_MesherHelper::IsStructured( SMESH_subMesh* faceSM )
2691 {
2692   SMESHDS_SubMesh* fSM = faceSM->GetSubMeshDS();
2693   if ( !fSM || fSM->NbElements() == 0 )
2694     return false;
2695
2696   list< TopoDS_Edge > edges;
2697   list< int > nbEdgesInWires;
2698   int nbWires = SMESH_Block::GetOrderedEdges( TopoDS::Face( faceSM->GetSubShape() ),
2699                                               edges, nbEdgesInWires );
2700   if ( nbWires != 1 /*|| nbEdgesInWires.front() != 4*/ ) // allow composite sides
2701     return false;
2702
2703   // algo: find corners of a structure and then analyze nb of faces and
2704   // length of structure sides
2705
2706   SMESHDS_Mesh* meshDS = faceSM->GetFather()->GetMeshDS();
2707   SMESH_MesherHelper faceAnalyser( *faceSM->GetFather() );
2708   faceAnalyser.SetSubShape( faceSM->GetSubShape() );
2709
2710   // rotate edges to get the first node being at corner
2711   // (in principle it's not necessary because so far none SALOME algo can make
2712   //  such a structured mesh that all corner nodes are not on VERTEXes)
2713   bool isCorner     = false;
2714   int nbRemainEdges = nbEdgesInWires.front();
2715   do {
2716     TopoDS_Vertex V = IthVertex( 0, edges.front() );
2717     isCorner = isCornerOfStructure( SMESH_Algo::VertexNode( V, meshDS ),
2718                                     fSM, faceAnalyser);
2719     if ( !isCorner ) {
2720       edges.splice( edges.end(), edges, edges.begin() );
2721       --nbRemainEdges;
2722     }
2723   }
2724   while ( !isCorner && nbRemainEdges > 0 );
2725
2726   if ( !isCorner )
2727     return false;
2728
2729   // get all nodes from EDGEs
2730   list< const SMDS_MeshNode* > nodes;
2731   list< TopoDS_Edge >::iterator edge = edges.begin();
2732   for ( ; edge != edges.end(); ++edge )
2733   {
2734     map< double, const SMDS_MeshNode* > u2Nodes;
2735     if ( !SMESH_Algo::GetSortedNodesOnEdge( meshDS, *edge,
2736                                             /*skipMedium=*/true, u2Nodes ))
2737       return false;
2738
2739     list< const SMDS_MeshNode* > edgeNodes;
2740     map< double, const SMDS_MeshNode* >::iterator u2n = u2Nodes.begin();
2741     for ( ; u2n != u2Nodes.end(); ++u2n )
2742       edgeNodes.push_back( u2n->second );
2743     if ( edge->Orientation() == TopAbs_REVERSED )
2744       edgeNodes.reverse();
2745
2746     if ( !nodes.empty() && nodes.back() == edgeNodes.front() )
2747       edgeNodes.pop_front();
2748     nodes.splice( nodes.end(), edgeNodes, edgeNodes.begin(), edgeNodes.end() );
2749   }
2750
2751   // get length of structured sides
2752   vector<int> nbEdgesInSide;
2753   int nbEdges = 0;
2754   list< const SMDS_MeshNode* >::iterator n = ++nodes.begin();
2755   for ( ; n != nodes.end(); ++n )
2756   {
2757     ++nbEdges;
2758     if ( isCornerOfStructure( *n, fSM, faceAnalyser )) {
2759       nbEdgesInSide.push_back( nbEdges );
2760       nbEdges = 0;
2761     }
2762   }
2763
2764   // checks
2765   if ( nbEdgesInSide.size() != 4 )
2766     return false;
2767   if ( nbEdgesInSide[0] != nbEdgesInSide[2] )
2768     return false;
2769   if ( nbEdgesInSide[1] != nbEdgesInSide[3] )
2770     return false;
2771   if ( nbEdgesInSide[0] * nbEdgesInSide[1] != fSM->NbElements() )
2772     return false;
2773
2774   return true;
2775 }
2776
2777 //=======================================================================
2778 //function : IsDistorted2D
2779 //purpose  : Return true if 2D mesh on FACE is ditorted
2780 //=======================================================================
2781
2782 bool SMESH_MesherHelper::IsDistorted2D( SMESH_subMesh* faceSM,
2783                                         bool           checkUV)
2784 {
2785   if ( !faceSM || faceSM->GetSubShape().ShapeType() != TopAbs_FACE )
2786     return false;
2787
2788   bool haveBadFaces = false;
2789
2790   SMESH_MesherHelper helper( *faceSM->GetFather() );
2791   helper.SetSubShape( faceSM->GetSubShape() );
2792
2793   const TopoDS_Face&  F = TopoDS::Face( faceSM->GetSubShape() );
2794   SMESHDS_SubMesh* smDS = helper.GetMeshDS()->MeshElements( F );
2795   if ( !smDS || smDS->NbElements() == 0 ) return false;
2796
2797   SMDS_ElemIteratorPtr faceIt = smDS->GetElements();
2798   double prevArea = 0;
2799   vector< const SMDS_MeshNode* > nodes;
2800   vector< gp_XY >                uv;
2801   bool* toCheckUV = checkUV ? & checkUV : 0;
2802   while ( faceIt->more() && !haveBadFaces )
2803   {
2804     const SMDS_MeshElement* face = faceIt->next();
2805
2806     // get nodes
2807     nodes.resize( face->NbCornerNodes() );
2808     SMDS_MeshElement::iterator n = face->begin_nodes();
2809     for ( size_t i = 0; i < nodes.size(); ++n, ++i )
2810       nodes[ i ] = *n;
2811
2812     // avoid elems on degenarate shapes as UV on them can be wrong
2813     if ( helper.HasDegeneratedEdges() )
2814     {
2815       bool isOnDegen = false;
2816       for ( size_t i = 0; ( i < nodes.size() && !isOnDegen ); ++i )
2817         isOnDegen = helper.IsDegenShape( nodes[ i ]->getshapeId() );
2818       if ( isOnDegen )
2819         continue;
2820     }
2821     // prepare to getting UVs
2822     const SMDS_MeshNode* inFaceNode = 0;
2823     if ( helper.HasSeam() ) {
2824       for ( size_t i = 0; ( i < nodes.size() && !inFaceNode ); ++i )
2825         if ( !helper.IsSeamShape( nodes[ i ]->getshapeId() ))
2826           inFaceNode = nodes[ i ];
2827       if ( !inFaceNode )
2828         continue;
2829     }
2830     // get UVs
2831     uv.resize( nodes.size() );
2832     for ( size_t i = 0; i < nodes.size(); ++i )
2833       uv[ i ] = helper.GetNodeUV( F, nodes[ i ], inFaceNode, toCheckUV );
2834
2835     // compare orientation of triangles
2836     double faceArea = 0;
2837     for ( int iT = 0, nbT = nodes.size()-2; iT < nbT; ++iT )
2838     {
2839       gp_XY v1 = uv[ iT+1 ] - uv[ 0 ];
2840       gp_XY v2 = uv[ iT+2 ] - uv[ 0 ];
2841       faceArea += v2 ^ v1;
2842     }
2843     haveBadFaces = ( faceArea * prevArea < 0 );
2844     prevArea = faceArea;
2845   }
2846
2847   return haveBadFaces;
2848 }
2849
2850 //================================================================================
2851 /*!
2852  * \brief Find out elements orientation on a geometrical face
2853  * \param theFace - The face correctly oriented in the shape being meshed
2854  * \retval bool - true if the face normal and the normal of first element
2855  *                in the correspoding submesh point in different directions
2856  */
2857 //================================================================================
2858
2859 bool SMESH_MesherHelper::IsReversedSubMesh (const TopoDS_Face& theFace)
2860 {
2861   if ( theFace.IsNull() )
2862     return false;
2863
2864   // find out orientation of a meshed face
2865   int faceID = GetMeshDS()->ShapeToIndex( theFace );
2866   TopoDS_Shape aMeshedFace = GetMeshDS()->IndexToShape( faceID );
2867   bool isReversed = ( theFace.Orientation() != aMeshedFace.Orientation() );
2868
2869   const SMESHDS_SubMesh * aSubMeshDSFace = GetMeshDS()->MeshElements( faceID );
2870   if ( !aSubMeshDSFace )
2871     return isReversed;
2872
2873   // find an element on a bounday of theFace
2874   SMDS_ElemIteratorPtr iteratorElem = aSubMeshDSFace->GetElements();
2875   const SMDS_MeshNode* nn[2];
2876   while ( iteratorElem->more() ) // loop on elements on theFace
2877   {
2878     const SMDS_MeshElement* elem = iteratorElem->next();
2879     if ( ! elem ) continue;
2880
2881     // look for 2 nodes on EDGE
2882     int nbNodes = elem->NbCornerNodes();
2883     nn[0] = elem->GetNode( nbNodes-1 );
2884     for ( int iN = 0; iN < nbNodes; ++iN )
2885     {
2886       nn[1] = elem->GetNode( iN );
2887       if ( nn[0]->GetPosition()->GetDim() < 2 &&
2888            nn[1]->GetPosition()->GetDim() < 2 )
2889       {
2890         TopoDS_Shape s0 = GetSubShapeByNode( nn[0], GetMeshDS() );
2891         TopoDS_Shape s1 = GetSubShapeByNode( nn[1], GetMeshDS() );
2892         TopoDS_Shape  E = GetCommonAncestor( s0, s1, *myMesh, TopAbs_EDGE );
2893         if ( !E.IsNull() && !s0.IsSame( s1 ))
2894         {
2895           // is E seam edge?
2896           int nb = 0;
2897           for ( TopExp_Explorer exp( theFace, TopAbs_EDGE ); exp.More(); exp.Next() )
2898             if ( E.IsSame( exp.Current() )) {
2899               ++nb;
2900               E = exp.Current(); // to know orientation
2901             }
2902           if ( nb == 1 )
2903           {
2904             bool ok = true;
2905             double u0 = GetNodeU( TopoDS::Edge( E ), nn[0], nn[1], &ok );
2906             double u1 = GetNodeU( TopoDS::Edge( E ), nn[1], nn[0], &ok );
2907             if ( ok )
2908             {
2909               isReversed = ( u0 > u1 );
2910               if ( E.Orientation() == TopAbs_REVERSED )
2911                 isReversed = !isReversed;
2912               return isReversed;
2913             }
2914           }
2915         }
2916       }
2917       nn[0] = nn[1];
2918     }
2919   }
2920
2921   // find an element with a good normal
2922   gp_Vec Ne;
2923   bool normalOK = false;
2924   gp_XY uv;
2925   iteratorElem = aSubMeshDSFace->GetElements();
2926   while ( !normalOK && iteratorElem->more() ) // loop on elements on theFace
2927   {
2928     const SMDS_MeshElement* elem = iteratorElem->next();
2929     if ( ! SMESH_MeshAlgos::FaceNormal( elem, const_cast<gp_XYZ&>( Ne.XYZ() ), /*normalized=*/0 ))
2930       continue;
2931     normalOK = true;
2932
2933     // get UV of a node inside theFACE
2934     SMDS_ElemIteratorPtr nodesIt = elem->nodesIterator();
2935     const SMDS_MeshNode* nInFace = 0;
2936     int iPosDim = SMDS_TOP_VERTEX;
2937     while ( nodesIt->more() ) // loop on nodes
2938     {
2939       const SMDS_MeshNode* n = static_cast<const SMDS_MeshNode*>( nodesIt->next() );
2940       if ( n->GetPosition()->GetTypeOfPosition() >= iPosDim )
2941       {
2942         nInFace = n;
2943         iPosDim = n->GetPosition()->GetTypeOfPosition();
2944       }
2945     }
2946     uv = GetNodeUV( theFace, nInFace, 0, &normalOK );
2947   }
2948   if ( !normalOK )
2949     return isReversed;
2950
2951   // face normal at node position
2952   TopLoc_Location loc;
2953   Handle(Geom_Surface) surf = BRep_Tool::Surface( theFace, loc );
2954   // if ( surf.IsNull() || surf->Continuity() < GeomAbs_C1 )
2955   // some surfaces not detected as GeomAbs_C1 are nevertheless correct for meshing
2956   if ( surf.IsNull() || surf->Continuity() < GeomAbs_C0 )
2957     return isReversed;
2958
2959   gp_Vec d1u, d1v; gp_Pnt p;
2960   surf->D1( uv.X(), uv.Y(), p, d1u, d1v );
2961   gp_Vec Nf = (d1u ^ d1v).Transformed( loc );
2962
2963   if ( theFace.Orientation() == TopAbs_REVERSED )
2964     Nf.Reverse();
2965
2966   return Ne * Nf < 0.;
2967 }
2968
2969 //=======================================================================
2970 //function : Count
2971 //purpose  : Count nb of sub-shapes
2972 //=======================================================================
2973
2974 int SMESH_MesherHelper::Count(const TopoDS_Shape&    shape,
2975                               const TopAbs_ShapeEnum type,
2976                               const bool             ignoreSame)
2977 {
2978   if ( ignoreSame ) {
2979     TopTools_IndexedMapOfShape map;
2980     TopExp::MapShapes( shape, type, map );
2981     return map.Extent();
2982   }
2983   else {
2984     int nb = 0;
2985     for ( TopExp_Explorer exp( shape, type ); exp.More(); exp.Next() )
2986       ++nb;
2987     return nb;
2988   }
2989 }
2990
2991 //=======================================================================
2992 //function : NbAncestors
2993 //purpose  : Return number of unique ancestors of the shape
2994 //=======================================================================
2995
2996 int SMESH_MesherHelper::NbAncestors(const TopoDS_Shape& shape,
2997          &