Salome HOME
Merge changes from 'master' branch.
[modules/smesh.git] / src / StdMeshers / StdMeshers_QuadFromMedialAxis_1D2D.cxx
1 // Copyright (C) 2007-2016  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 // File      : StdMeshers_QuadFromMedialAxis_1D2D.cxx
23 // Created   : Wed Jun  3 17:33:45 2015
24 // Author    : Edward AGAPOV (eap)
25
26 #include "StdMeshers_QuadFromMedialAxis_1D2D.hxx"
27
28 #include "SMESHDS_Mesh.hxx"
29 #include "SMESH_Block.hxx"
30 #include "SMESH_Gen.hxx"
31 #include "SMESH_MAT2d.hxx"
32 #include "SMESH_Mesh.hxx"
33 #include "SMESH_MeshEditor.hxx"
34 #include "SMESH_MesherHelper.hxx"
35 #include "SMESH_ProxyMesh.hxx"
36 #include "SMESH_subMesh.hxx"
37 #include "SMESH_subMeshEventListener.hxx"
38 #include "StdMeshers_FaceSide.hxx"
39 #include "StdMeshers_LayerDistribution.hxx"
40 #include "StdMeshers_NumberOfLayers.hxx"
41 #include "StdMeshers_Regular_1D.hxx"
42 #include "StdMeshers_ViscousLayers2D.hxx"
43
44 #include <BRepAdaptor_Curve.hxx>
45 #include <BRepBuilderAPI_MakeEdge.hxx>
46 #include <BRepTools.hxx>
47 #include <BRep_Tool.hxx>
48 #include <GeomAPI_Interpolate.hxx>
49 #include <Geom_Surface.hxx>
50 #include <Precision.hxx>
51 #include <TColgp_HArray1OfPnt.hxx>
52 #include <TopExp.hxx>
53 #include <TopExp_Explorer.hxx>
54 #include <TopLoc_Location.hxx>
55 #include <TopTools_MapOfShape.hxx>
56 #include <TopoDS.hxx>
57 #include <TopoDS_Edge.hxx>
58 #include <TopoDS_Face.hxx>
59 #include <TopoDS_Vertex.hxx>
60 #include <gp_Pnt.hxx>
61
62 #include <list>
63 #include <vector>
64
65 using namespace std;
66
67 //================================================================================
68 /*!
69  * \brief 1D algo
70  */
71 class StdMeshers_QuadFromMedialAxis_1D2D::Algo1D : public StdMeshers_Regular_1D
72 {
73 public:
74   Algo1D(SMESH_Gen* gen):
75     StdMeshers_Regular_1D( gen->GetANewId(), gen )
76   {
77   }
78   void SetSegmentLength( double len )
79   {
80     SMESH_Algo::_usedHypList.clear();
81     _value[ BEG_LENGTH_IND ] = len;
82     _value[ PRECISION_IND  ] = 1e-7;
83     _hypType = LOCAL_LENGTH;
84   }
85   void SetRadialDistribution( const SMESHDS_Hypothesis* hyp )
86   {
87     SMESH_Algo::_usedHypList.clear();
88     if ( !hyp )
89       return;
90
91     if ( const StdMeshers_NumberOfLayers* nl =
92          dynamic_cast< const StdMeshers_NumberOfLayers* >( hyp ))
93     {
94       _ivalue[ NB_SEGMENTS_IND  ] = nl->GetNumberOfLayers();
95       _ivalue[ DISTR_TYPE_IND ]   = 0;
96       _hypType = NB_SEGMENTS;
97     }
98     if ( const StdMeshers_LayerDistribution* ld =
99          dynamic_cast< const StdMeshers_LayerDistribution* >( hyp ))
100     {
101       if ( SMESH_Hypothesis* h = ld->GetLayerDistribution() )
102       {
103         SMESH_Algo::_usedHypList.clear();
104         SMESH_Algo::_usedHypList.push_back( h );
105       }
106     }
107   }
108   void ComputeDistribution(SMESH_MesherHelper& theHelper,
109                            const gp_Pnt&       thePnt1,
110                            const gp_Pnt&       thePnt2,
111                            list< double >&     theParams)
112   {
113     SMESH_Mesh& mesh = *theHelper.GetMesh();
114     TopoDS_Edge edge = BRepBuilderAPI_MakeEdge( thePnt1, thePnt2 );
115
116     SMESH_Hypothesis::Hypothesis_Status aStatus;
117     CheckHypothesis( mesh, edge, aStatus );
118
119     theParams.clear();
120     BRepAdaptor_Curve C3D(edge);
121     double f = C3D.FirstParameter(), l = C3D.LastParameter(), len = thePnt1.Distance( thePnt2 );
122     if ( !StdMeshers_Regular_1D::computeInternalParameters( mesh, C3D, len, f, l, theParams, false))
123     {
124       for ( size_t i = 1; i < 15; ++i )
125         theParams.push_back( i/15. ); // ????
126     }
127     else
128     {
129       for (list<double>::iterator itU = theParams.begin(); itU != theParams.end(); ++itU )
130         *itU /= len;
131     }
132   }
133   virtual const list <const SMESHDS_Hypothesis *> &
134   GetUsedHypothesis(SMESH_Mesh &, const TopoDS_Shape &, const bool)
135   {
136     return SMESH_Algo::_usedHypList;
137   }
138   virtual bool CheckHypothesis(SMESH_Mesh&                          aMesh,
139                                const TopoDS_Shape&                  aShape,
140                                SMESH_Hypothesis::Hypothesis_Status& aStatus)
141   {
142     if ( !SMESH_Algo::_usedHypList.empty() )
143       return StdMeshers_Regular_1D::CheckHypothesis( aMesh, aShape, aStatus );
144     return true;
145   }
146 };
147  
148 //================================================================================
149 /*!
150  * \brief Constructor sets algo features
151  */
152 //================================================================================
153
154 StdMeshers_QuadFromMedialAxis_1D2D::StdMeshers_QuadFromMedialAxis_1D2D(int        hypId,
155                                                                        SMESH_Gen* gen)
156   : StdMeshers_Quadrangle_2D(hypId, gen),
157     _regular1D( 0 )
158 {
159   _name = "QuadFromMedialAxis_1D2D";
160   _shapeType = (1 << TopAbs_FACE);
161   _onlyUnaryInput          = true;  // FACE by FACE so far
162   _requireDiscreteBoundary = false; // make 1D by myself
163   _supportSubmeshes        = true;  // make 1D by myself
164   _neededLowerHyps[ 1 ]    = true;  // suppress warning on hiding a global 1D algo
165   _neededLowerHyps[ 2 ]    = true;  // suppress warning on hiding a global 2D algo
166   _compatibleHypothesis.clear();
167   _compatibleHypothesis.push_back("ViscousLayers2D");
168   _compatibleHypothesis.push_back("LayerDistribution2D");
169   _compatibleHypothesis.push_back("NumberOfLayers2D");
170 }
171
172 //================================================================================
173 /*!
174  * \brief Destructor
175  */
176 //================================================================================
177
178 StdMeshers_QuadFromMedialAxis_1D2D::~StdMeshers_QuadFromMedialAxis_1D2D()
179 {
180   delete _regular1D;
181   _regular1D = 0;
182 }
183
184 //================================================================================
185 /*!
186  * \brief Check if needed hypotheses are present
187  */
188 //================================================================================
189
190 bool StdMeshers_QuadFromMedialAxis_1D2D::CheckHypothesis(SMESH_Mesh&         aMesh,
191                                                          const TopoDS_Shape& aShape,
192                                                          Hypothesis_Status&  aStatus)
193 {
194   aStatus = HYP_OK;
195
196   // get one main optional hypothesis
197   const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
198   _hyp2D  = hyps.empty() ? 0 : hyps.front();
199
200   return true; // does not require hypothesis
201 }
202
203 namespace
204 {
205   typedef map< const SMDS_MeshNode*, list< const SMDS_MeshNode* > > TMergeMap;
206   
207   //================================================================================
208   /*!
209    * \brief Sinuous face
210    */
211   struct SinuousFace
212   {
213     FaceQuadStruct::Ptr          _quad;
214     vector< TopoDS_Edge >        _edges;
215     vector< TopoDS_Edge >        _sinuSide[2], _shortSide[2];
216     vector< TopoDS_Edge >        _sinuEdges;
217     vector< Handle(Geom_Curve) > _sinuCurves;
218     int                          _nbWires;
219     list< int >                  _nbEdgesInWire;
220     TMergeMap                    _nodesToMerge;
221
222     SinuousFace( const TopoDS_Face& f ): _quad( new FaceQuadStruct )
223     {
224       list< TopoDS_Edge > edges;
225       _nbWires = SMESH_Block::GetOrderedEdges (f, edges, _nbEdgesInWire);
226       _edges.assign( edges.begin(), edges.end() );
227
228       _quad->side.resize( 4 );
229       _quad->face = f;
230     }
231     const TopoDS_Face& Face() const { return _quad->face; }
232     bool IsRing() const { return _shortSide[0].empty() && !_sinuSide[0].empty(); }
233   };
234
235   //================================================================================
236   /*!
237    * \brief Temporary mesh
238    */
239   struct TmpMesh : public SMESH_Mesh
240   {
241     TmpMesh()
242     {
243       _myMeshDS = new SMESHDS_Mesh(/*id=*/0, /*isEmbeddedMode=*/true);
244     }
245   };
246
247   //================================================================================
248   /*!
249    * \brief Event listener which removes mesh from EDGEs when 2D hyps change
250    */
251   struct EdgeCleaner : public SMESH_subMeshEventListener
252   {
253     int _prevAlgoEvent;
254     EdgeCleaner():
255       SMESH_subMeshEventListener( /*isDeletable=*/true,
256                                   "StdMeshers_QuadFromMedialAxis_1D2D::EdgeCleaner")
257     {
258       _prevAlgoEvent = -1;
259     }
260     virtual void ProcessEvent(const int                       event,
261                               const int                       eventType,
262                               SMESH_subMesh*                  faceSubMesh,
263                               SMESH_subMeshEventListenerData* data,
264                               const SMESH_Hypothesis*         hyp)
265     {
266       if ( eventType == SMESH_subMesh::ALGO_EVENT )
267       {
268         _prevAlgoEvent = event;
269         return;
270       }
271       // SMESH_subMesh::COMPUTE_EVENT
272       if ( _prevAlgoEvent == SMESH_subMesh::REMOVE_HYP ||
273            _prevAlgoEvent == SMESH_subMesh::REMOVE_ALGO ||
274            _prevAlgoEvent == SMESH_subMesh::MODIF_HYP )
275       {
276         SMESH_subMeshIteratorPtr smIt = faceSubMesh->getDependsOnIterator(/*includeSelf=*/false);
277         while ( smIt->more() )
278           smIt->next()->ComputeStateEngine( SMESH_subMesh::CLEAN );
279       }
280       _prevAlgoEvent = -1;
281     }
282   };
283
284   //================================================================================
285   /*!
286    * \brief Return a member of a std::pair
287    */
288   //================================================================================
289
290   template< typename T >
291   T& get( std::pair< T, T >& thePair, bool is2nd )
292   {
293     return is2nd ? thePair.second : thePair.first;
294   }
295
296   //================================================================================
297   /*!
298    * \brief Select two EDGEs from a map, either mapped to least values or to max values
299    */
300   //================================================================================
301
302   // template< class TVal2EdgesMap >
303   // void getTwo( bool                 least,
304   //              TVal2EdgesMap&       map,
305   //              vector<TopoDS_Edge>& twoEdges,
306   //              vector<TopoDS_Edge>& otherEdges)
307   // {
308   //   twoEdges.clear();
309   //   otherEdges.clear();
310   //   if ( least )
311   //   {
312   //     TVal2EdgesMap::iterator i = map.begin();
313   //     twoEdges.push_back( i->second );
314   //     twoEdges.push_back( ++i->second );
315   //     for ( ; i != map.end(); ++i )
316   //       otherEdges.push_back( i->second );
317   //   }
318   //   else
319   //   {
320   //     TVal2EdgesMap::reverse_iterator i = map.rbegin();
321   //     twoEdges.push_back( i->second );
322   //     twoEdges.push_back( ++i->second );
323   //     for ( ; i != map.rend(); ++i )
324   //       otherEdges.push_back( i->second );
325   //   }
326   //   TopoDS_Vertex v;
327   //   if ( TopExp::CommonVertex( twoEdges[0], twoEdges[1], v ))
328   //   {
329   //     twoEdges.clear(); // two EDGEs must not be connected
330   //     otherEdges.clear();
331   //   }
332   // }
333
334   //================================================================================
335   /*!
336    * \brief Finds out a minimal segment length given EDGEs will be divided into.
337    *        This length is further used to discretize the Medial Axis
338    */
339   //================================================================================
340
341   double getMinSegLen(SMESH_MesherHelper&        theHelper,
342                       const vector<TopoDS_Edge>& theEdges)
343   {
344     TmpMesh tmpMesh;
345     SMESH_Mesh* mesh = theHelper.GetMesh();
346
347     vector< SMESH_Algo* > algos( theEdges.size() );
348     for ( size_t i = 0; i < theEdges.size(); ++i )
349     {
350       SMESH_subMesh* sm = mesh->GetSubMesh( theEdges[i] );
351       algos[i] = sm->GetAlgo();
352     }
353
354     int nbSegDflt = mesh->GetGen() ? mesh->GetGen()->GetDefaultNbSegments() : 15;
355     double minSegLen = Precision::Infinite();
356
357     for ( size_t i = 0; i < theEdges.size(); ++i )
358     {
359       SMESH_subMesh* sm = mesh->GetSubMesh( theEdges[i] );
360       if ( SMESH_Algo::IsStraight( theEdges[i], /*degenResult=*/true ))
361         continue;
362       // get algo
363       size_t iOpp = ( theEdges.size() == 4 ? (i+2)%4 : i );
364       SMESH_Algo*  algo = sm->GetAlgo();
365       if ( !algo ) algo = algos[ iOpp ];
366       // get hypo
367       SMESH_Hypothesis::Hypothesis_Status status = SMESH_Hypothesis::HYP_MISSING;
368       if ( algo )
369       {
370         if ( !algo->CheckHypothesis( *mesh, theEdges[i], status ))
371           algo->CheckHypothesis( *mesh, theEdges[iOpp], status );
372       }
373       // compute
374       if ( status != SMESH_Hypothesis::HYP_OK )
375       {
376         minSegLen = Min( minSegLen, SMESH_Algo::EdgeLength( theEdges[i] ) / nbSegDflt );
377       }
378       else
379       {
380         tmpMesh.Clear();
381         tmpMesh.ShapeToMesh( TopoDS_Shape());
382         tmpMesh.ShapeToMesh( theEdges[i] );
383         try {
384           if ( !mesh->GetGen() ) continue; // tmp mesh
385           mesh->GetGen()->Compute( tmpMesh, theEdges[i], SMESH_Gen::SHAPE_ONLY_UPWARD ); // make nodes on VERTEXes
386           if ( !algo->Compute( tmpMesh, theEdges[i] ))
387             continue;
388         }
389         catch (...) {
390           continue;
391         }
392         SMDS_EdgeIteratorPtr segIt = tmpMesh.GetMeshDS()->edgesIterator();
393         while ( segIt->more() )
394         {
395           const SMDS_MeshElement* seg = segIt->next();
396           double len = SMESH_TNodeXYZ( seg->GetNode(0) ).Distance( seg->GetNode(1) );
397           minSegLen = Min( minSegLen, len );
398         }
399       }
400     }
401     if ( Precision::IsInfinite( minSegLen ))
402       minSegLen = mesh->GetShapeDiagonalSize() / nbSegDflt;
403
404     return minSegLen;
405   }
406
407   //================================================================================
408   /*!
409    * \brief Returns EDGEs located between two VERTEXes at which given MA branches end
410    *  \param [in] br1 - one MA branch
411    *  \param [in] br2 - one more MA branch
412    *  \param [in] allEdges - all EDGEs of a FACE
413    *  \param [out] shortEdges - the found EDGEs
414    *  \return bool - is OK or not
415    */
416   //================================================================================
417
418   bool getConnectedEdges( const SMESH_MAT2d::Branch* br1,
419                           const SMESH_MAT2d::Branch* br2,
420                           const vector<TopoDS_Edge>& allEdges,
421                           vector<TopoDS_Edge>&       shortEdges)
422   {
423     vector< size_t > edgeIDs[4];
424     br1->getGeomEdges( edgeIDs[0], edgeIDs[1] );
425     br2->getGeomEdges( edgeIDs[2], edgeIDs[3] );
426
427     // EDGEs returned by a Branch form a connected chain with a VERTEX where
428     // the Branch ends at the chain middle. One of end EDGEs of the chain is common
429     // with either end EDGE of the chain of the other Branch, or the chains are connected
430     // at a common VERTEX;
431
432     // Get indices of end EDGEs of the branches
433     bool vAtStart1 = ( br1->getEnd(0)->_type == SMESH_MAT2d::BE_ON_VERTEX );
434     bool vAtStart2 = ( br2->getEnd(0)->_type == SMESH_MAT2d::BE_ON_VERTEX );
435     size_t iEnd[4] = {
436       vAtStart1 ? edgeIDs[0].back() : edgeIDs[0][0],
437       vAtStart1 ? edgeIDs[1].back() : edgeIDs[1][0],
438       vAtStart2 ? edgeIDs[2].back() : edgeIDs[2][0],
439       vAtStart2 ? edgeIDs[3].back() : edgeIDs[3][0]
440     };
441
442     set< size_t > connectedIDs;
443     TopoDS_Vertex vCommon;
444     // look for the same EDGEs
445     for ( int i = 0; i < 2; ++i )
446       for ( int j = 2; j < 4; ++j )
447         if ( iEnd[i] == iEnd[j] )
448         {
449           connectedIDs.insert( edgeIDs[i].begin(), edgeIDs[i].end() );
450           connectedIDs.insert( edgeIDs[j].begin(), edgeIDs[j].end() );
451           i = j = 4;
452         }
453     if ( connectedIDs.empty() )
454       // look for connected EDGEs
455       for ( int i = 0; i < 2; ++i )
456         for ( int j = 2; j < 4; ++j )
457           if ( TopExp::CommonVertex( allEdges[ iEnd[i]], allEdges[ iEnd[j]], vCommon ))
458           {
459             connectedIDs.insert( edgeIDs[i].begin(), edgeIDs[i].end() );
460             connectedIDs.insert( edgeIDs[j].begin(), edgeIDs[j].end() );
461             i = j = 4;
462           }
463     if ( connectedIDs.empty() ||                     // nothing
464          allEdges.size() - connectedIDs.size() < 2 ) // too many
465       return false;
466
467     // set shortEdges in the order as in allEdges
468     if ( connectedIDs.count( 0 ) &&
469          connectedIDs.count( allEdges.size()-1 ))
470     {
471       size_t iE = allEdges.size()-1;
472       while ( connectedIDs.count( iE-1 ))
473         --iE;
474       for ( size_t i = 0; i < connectedIDs.size(); ++i )
475       {
476         shortEdges.push_back( allEdges[ iE ]);
477         iE = ( iE + 1 ) % allEdges.size();
478       }
479     }
480     else
481     {
482       set< size_t >::iterator i = connectedIDs.begin();
483       for ( ; i != connectedIDs.end(); ++i )
484         shortEdges.push_back( allEdges[ *i ]);
485     }
486     return true;
487   }
488
489   //================================================================================
490   /*!
491    * \brief Find EDGEs to discretize using projection from MA
492    *  \param [in,out] theSinuFace - the FACE to be meshed
493    *  \return bool - OK or not
494    *
495    * It separates all EDGEs into four sides of a quadrangle connected in the order:
496    * theSinuEdges[0], theShortEdges[0], theSinuEdges[1], theShortEdges[1]
497    */
498   //================================================================================
499
500   bool getSinuousEdges( SMESH_MesherHelper& theHelper,
501                         SinuousFace&        theSinuFace)
502   {
503     vector<TopoDS_Edge> * theSinuEdges  = & theSinuFace._sinuSide [0];
504     vector<TopoDS_Edge> * theShortEdges = & theSinuFace._shortSide[0];
505     theSinuEdges[0].clear();
506     theSinuEdges[1].clear();
507     theShortEdges[0].clear();
508     theShortEdges[1].clear();
509    
510     vector<TopoDS_Edge> & allEdges = theSinuFace._edges;
511     const size_t nbEdges = allEdges.size();
512     if ( nbEdges < 4 && theSinuFace._nbWires == 1 )
513       return false;
514
515     if ( theSinuFace._nbWires == 2 ) // ring
516     {
517       size_t nbOutEdges = theSinuFace._nbEdgesInWire.front();
518       theSinuEdges[0].assign ( allEdges.begin(), allEdges.begin() + nbOutEdges );
519       theSinuEdges[1].assign ( allEdges.begin() + nbOutEdges, allEdges.end() );
520       theSinuFace._sinuEdges = allEdges;
521       return true;
522     }
523     if ( theSinuFace._nbWires > 2 )
524       return false;
525
526     // create MedialAxis to find short edges by analyzing MA branches
527     double minSegLen = getMinSegLen( theHelper, allEdges );
528     SMESH_MAT2d::MedialAxis ma( theSinuFace.Face(), allEdges, minSegLen * 3 );
529
530     // in an initial request case, theFace represents a part of a river with almost parallel banks
531     // so there should be two branch points
532     using SMESH_MAT2d::BranchEnd;
533     using SMESH_MAT2d::Branch;
534     const vector< const BranchEnd* >& braPoints = ma.getBranchPoints();
535     if ( braPoints.size() < 2 )
536       return false;
537     TopTools_MapOfShape shortMap;
538     size_t nbBranchPoints = 0;
539     for ( size_t i = 0; i < braPoints.size(); ++i )
540     {
541       vector< const Branch* > vertBranches; // branches with an end on VERTEX
542       for ( size_t ib = 0; ib < braPoints[i]->_branches.size(); ++ib )
543       {
544         const Branch* branch = braPoints[i]->_branches[ ib ];
545         if ( branch->hasEndOfType( SMESH_MAT2d::BE_ON_VERTEX ))
546           vertBranches.push_back( branch );
547       }
548       if ( vertBranches.size() != 2 || braPoints[i]->_branches.size() != 3)
549         continue;
550
551       // get common EDGEs of two branches
552       if ( !getConnectedEdges( vertBranches[0], vertBranches[1],
553                                allEdges, theShortEdges[ nbBranchPoints > 0 ] ))
554         return false;
555
556       for ( size_t iS = 0; iS < theShortEdges[ nbBranchPoints > 0 ].size(); ++iS )
557         shortMap.Add( theShortEdges[ nbBranchPoints > 0 ][ iS ]);
558
559       ++nbBranchPoints;
560     }
561
562     if ( nbBranchPoints != 2 )
563       return false;
564
565     // add to theSinuEdges all edges that are not theShortEdges
566     vector< vector<TopoDS_Edge> > sinuEdges(1);
567     TopoDS_Vertex vCommon;
568     for ( size_t i = 0; i < allEdges.size(); ++i )
569     {
570       if ( !shortMap.Contains( allEdges[i] ))
571       {
572         if ( !sinuEdges.back().empty() )
573           if ( !TopExp::CommonVertex( sinuEdges.back().back(), allEdges[ i ], vCommon ))
574             sinuEdges.resize( sinuEdges.size() + 1 );
575
576         sinuEdges.back().push_back( allEdges[i] );
577       }
578     }
579     if ( sinuEdges.size() == 3 )
580     {
581       if ( !TopExp::CommonVertex( sinuEdges.back().back(), sinuEdges[0][0], vCommon ))
582         return false;
583       vector<TopoDS_Edge>& last = sinuEdges.back();
584       last.insert( last.end(), sinuEdges[0].begin(), sinuEdges[0].end() );
585       sinuEdges[0].swap( last );
586       sinuEdges.resize( 2 );
587     }
588     if ( sinuEdges.size() != 2 )
589       return false;
590
591     theSinuEdges[0].swap( sinuEdges[0] );
592     theSinuEdges[1].swap( sinuEdges[1] );
593
594     if ( !TopExp::CommonVertex( theSinuEdges[0].back(), theShortEdges[0][0], vCommon ) ||
595          !vCommon.IsSame( theHelper.IthVertex( 1, theSinuEdges[0].back() )))
596       theShortEdges[0].swap( theShortEdges[1] );
597
598     theSinuFace._sinuEdges = theSinuEdges[0];
599     theSinuFace._sinuEdges.insert( theSinuFace._sinuEdges.end(),
600                                    theSinuEdges[1].begin(), theSinuEdges[1].end() );
601
602     return ( theShortEdges[0].size() > 0 && theShortEdges[1].size() > 0 &&
603              theSinuEdges [0].size() > 0 && theSinuEdges [1].size() > 0 );
604
605     // the sinuous EDGEs can be composite and C0 continuous,
606     // therefore we use a complex criterion to find TWO short non-sinuous EDGEs
607     // and the rest EDGEs will be treated as sinuous.
608     // A short edge should have the following features:
609     // a) straight
610     // b) short
611     // c) with convex corners at ends
612     // d) far from the other short EDGE
613
614     // vector< double > isStraightEdge( nbEdges, 0 ); // criterion value
615
616     // // a0) evaluate continuity
617     // const double contiWgt = 0.5; // weight of continuity in the criterion
618     // multimap< int, TopoDS_Edge > continuity;
619     // for ( size_t i = 0; i < nbEdges; ++I )
620     // {
621     //   BRepAdaptor_Curve curve( allEdges[i] );
622     //   GeomAbs_Shape C = GeomAbs_CN;
623     //   try:
624     //     C = curve.Continuity(); // C0, G1, C1, G2, C2, C3, CN
625     //   catch ( Standard_Failure ) {}
626     //   continuity.insert( make_pair( C, allEdges[i] ));
627     //   isStraight[i] += double( C ) / double( CN ) * contiWgt;
628     // }
629
630     // // try to choose by continuity
631     // int mostStraight = (int) continuity.rbegin()->first;
632     // int lessStraight = (int) continuity.begin()->first;
633     // if ( mostStraight != lessStraight )
634     // {
635     //   int nbStraight = continuity.count( mostStraight );
636     //   if ( nbStraight == 2 )
637     //   {
638     //     getTwo( /*least=*/false, continuity, theShortEdges, theSinuEdges );
639     //   }
640     //   else if ( nbStraight == 3 && nbEdges == 4 )
641     //   {
642     //     theSinuEdges.push_back( continuity.begin()->second );
643     //     vector<TopoDS_Edge>::iterator it =
644     //       std::find( allEdges.begin(), allEdges.end(), theSinuEdges[0] );
645     //     int i = std::distance( allEdges.begin(), it );
646     //     theSinuEdges .push_back( allEdges[( i+2 )%4 ]);
647     //     theShortEdges.push_back( allEdges[( i+1 )%4 ]);
648     //     theShortEdges.push_back( allEdges[( i+3 )%4 ]);
649     //   }
650     //   if ( theShortEdges.size() == 2 )
651     //     return true;
652     // }
653
654     // // a) curvature; evaluate aspect ratio
655     // {
656     //   const double curvWgt = 0.5;
657     //   for ( size_t i = 0; i < nbEdges; ++I )
658     //   {
659     //     BRepAdaptor_Curve curve( allEdges[i] );
660     //     double curvature = 1;
661     //     if ( !curve.IsClosed() )
662     //     {
663     //       const double f = curve.FirstParameter(), l = curve.LastParameter();
664     //       gp_Pnt pf = curve.Value( f ), pl = curve.Value( l );
665     //       gp_Lin line( pf, pl.XYZ() - pf.XYZ() );
666     //       double distMax = 0;
667     //       for ( double u = f; u < l; u += (l-f)/30. )
668     //         distMax = Max( distMax, line.SquareDistance( curve.Value( u )));
669     //       curvature = Sqrt( distMax ) / ( pf.Distance( pl ));
670     //     }
671     //     isStraight[i] += curvWgt / (              curvature + 1e-20 );
672     //   }
673     // }
674     // // b) length
675     // {
676     //   const double lenWgt = 0.5;
677     //   for ( size_t i = 0; i < nbEdges; ++I )
678     //   {
679     //     double length = SMESH_Algo::Length( allEdges[i] );
680     //     if ( length > 0 )
681     //       isStraight[i] += lenWgt / length;
682     //   }
683     // }
684     // // c) with convex corners at ends
685     // {
686     //   const double cornerWgt = 0.25;
687     //   for ( size_t i = 0; i < nbEdges; ++I )
688     //   {
689     //     double convex = 0;
690     //     int iPrev = SMESH_MesherHelper::WrapIndex( int(i)-1, nbEdges );
691     //     int iNext = SMESH_MesherHelper::WrapIndex( int(i)+1, nbEdges );
692     //     TopoDS_Vertex v = helper.IthVertex( 0, allEdges[i] );
693     //     double angle = SMESH_MesherHelper::GetAngle( allEdges[iPrev], allEdges[i], theFace, v );
694     //     if ( angle < M_PI ) // [-PI; PI]
695     //       convex += ( angle + M_PI ) / M_PI / M_PI;
696     //     v = helper.IthVertex( 1, allEdges[i] );
697     //     angle = SMESH_MesherHelper::GetAngle( allEdges[iNext], allEdges[i], theFace, v );
698     //     if ( angle < M_PI ) // [-PI; PI]
699     //       convex += ( angle + M_PI ) / M_PI / M_PI;
700     //     isStraight[i] += cornerWgt * convex;
701     //   }
702     // }
703   }
704
705   //================================================================================
706   /*!
707    * \brief Creates an EDGE from a sole branch of MA
708    */
709   //================================================================================
710
711   TopoDS_Edge makeEdgeFromMA( SMESH_MesherHelper&            theHelper,
712                               const SMESH_MAT2d::MedialAxis& theMA,
713                               const double                   theMinSegLen)
714   {
715     if ( theMA.nbBranches() != 1 )
716       return TopoDS_Edge();
717
718     vector< gp_XY > uv;
719     theMA.getPoints( theMA.getBranch(0), uv );
720     if ( uv.size() < 2 )
721       return TopoDS_Edge();
722
723     TopoDS_Face face = TopoDS::Face( theHelper.GetSubShape() );
724     Handle(Geom_Surface) surface = BRep_Tool::Surface( face );
725
726     vector< gp_Pnt > pnt;
727     pnt.reserve( uv.size() * 2 );
728     pnt.push_back( surface->Value( uv[0].X(), uv[0].Y() ));
729     for ( size_t i = 1; i < uv.size(); ++i )
730     {
731       gp_Pnt p = surface->Value( uv[i].X(), uv[i].Y() );
732       int nbDiv = int( p.Distance( pnt.back() ) / theMinSegLen );
733       for ( int iD = 1; iD < nbDiv; ++iD )
734       {
735         double  R = iD / double( nbDiv );
736         gp_XY uvR = uv[i-1] * (1 - R) + uv[i] * R;
737         pnt.push_back( surface->Value( uvR.X(), uvR.Y() ));
738       }
739       pnt.push_back( p );
740     }
741
742     // cout << "from salome.geom import geomBuilder" << endl;
743     // cout << "geompy = geomBuilder.New()" << endl;
744     Handle(TColgp_HArray1OfPnt) points = new TColgp_HArray1OfPnt(1, pnt.size());
745     for ( size_t i = 0; i < pnt.size(); ++i )
746     {
747       gp_Pnt& p = pnt[i];
748       points->SetValue( i+1, p );
749       // cout << "geompy.MakeVertex( "<< p.X()<<", " << p.Y()<<", " << p.Z()
750       //      <<" theName = 'p_" << i << "')" << endl;
751     }
752
753     GeomAPI_Interpolate interpol( points, /*isClosed=*/false, gp::Resolution());
754     interpol.Perform();
755     if ( !interpol.IsDone())
756       return TopoDS_Edge();
757
758     TopoDS_Edge branchEdge = BRepBuilderAPI_MakeEdge(interpol.Curve());
759     return branchEdge;
760   }
761
762   //================================================================================
763   /*!
764    * \brief Returns a type of shape, to which a hypothesis used to mesh a given edge is assigned
765    */
766   //================================================================================
767
768   TopAbs_ShapeEnum getHypShape( SMESH_Mesh* mesh, const TopoDS_Shape& edge )
769   {
770     TopAbs_ShapeEnum shapeType = TopAbs_SHAPE;
771
772     SMESH_subMesh* sm = mesh->GetSubMesh( edge );
773     SMESH_Algo*  algo = sm->GetAlgo();
774     if ( !algo ) return shapeType;
775
776     const list <const SMESHDS_Hypothesis *> & hyps =
777       algo->GetUsedHypothesis( *mesh, edge, /*ignoreAuxiliary=*/true );
778     if ( hyps.empty() ) return shapeType;
779
780     TopoDS_Shape shapeOfHyp =
781       SMESH_MesherHelper::GetShapeOfHypothesis( hyps.front(), edge, mesh);
782
783     return SMESH_MesherHelper::GetGroupType( shapeOfHyp, /*woCompound=*/true);
784   }
785
786   //================================================================================
787   /*!
788    * \brief Discretize a sole branch of MA an returns parameters of divisions on MA
789    */
790   //================================================================================
791
792   bool divideMA( SMESH_MesherHelper&            theHelper,
793                  const SMESH_MAT2d::MedialAxis& theMA,
794                  const SinuousFace&             theSinuFace,
795                  SMESH_Algo*                    the1dAlgo,
796                  const double                   theMinSegLen,
797                  vector<double>&                theMAParams )
798   {
799     // Check if all EDGEs of one size are meshed, then MA discretization is not needed
800     SMESH_Mesh* mesh = theHelper.GetMesh();
801     size_t nbComputedEdges[2] = { 0, 0 };
802     for ( size_t iS = 0; iS < 2; ++iS )
803       for ( size_t i = 0; i < theSinuFace._sinuSide[iS].size(); ++i )
804       {
805         const TopoDS_Edge& sinuEdge = theSinuFace._sinuSide[iS][i];
806         SMESH_subMesh*           sm = mesh->GetSubMesh( sinuEdge );
807         bool             isComputed = ( !sm->IsEmpty() );
808         if ( isComputed )
809         {
810           TopAbs_ShapeEnum shape = getHypShape( mesh, sinuEdge );
811           if ( shape == TopAbs_SHAPE || shape <= TopAbs_FACE )
812           {
813             // EDGE computed using global hypothesis -> clear it
814             bool hasComputedFace = false;
815             PShapeIteratorPtr faceIt = theHelper.GetAncestors( sinuEdge, *mesh, TopAbs_FACE );
816             while ( const TopoDS_Shape* face = faceIt->next() )
817               if (( !face->IsSame( theSinuFace.Face() )) &&
818                   ( hasComputedFace = !mesh->GetSubMesh( *face )->IsEmpty() ))
819                 break;
820             if ( !hasComputedFace )
821             {
822               sm->ComputeStateEngine( SMESH_subMesh::CLEAN );
823               isComputed = false;
824             }
825           }
826         }
827         nbComputedEdges[ iS ] += isComputed;
828       }
829     if ( nbComputedEdges[0] == theSinuFace._sinuSide[0].size() ||
830          nbComputedEdges[1] == theSinuFace._sinuSide[1].size() )
831       return true; // discretization is not needed
832
833     // Make MA EDGE
834     TopoDS_Edge branchEdge = makeEdgeFromMA( theHelper, theMA, theMinSegLen );
835     if ( branchEdge.IsNull() )
836       return false;
837
838     // const char* file = "/misc/dn25/salome/eap/salome/misc/tmp/MAedge.brep";
839     // BRepTools::Write( branchEdge, file);
840     // cout << "Write " << file << endl;
841
842
843     // Find 1D algo to mesh branchEdge
844   
845     // look for a most local 1D hyp assigned to the FACE
846     int mostSimpleShape = -1, maxShape = TopAbs_EDGE;
847     TopoDS_Edge edge;
848     for ( size_t i = 0; i < theSinuFace._sinuEdges.size(); ++i )
849     {
850       TopAbs_ShapeEnum shapeType = getHypShape( mesh, theSinuFace._sinuEdges[i] );
851       if ( mostSimpleShape < shapeType && shapeType < maxShape )
852       {
853         edge = theSinuFace._sinuEdges[i];
854         mostSimpleShape = shapeType;
855       }
856     }
857
858     SMESH_Algo* algo = the1dAlgo;
859     if ( mostSimpleShape > -1 )
860     {
861       algo = mesh->GetSubMesh( edge )->GetAlgo();
862       SMESH_Hypothesis::Hypothesis_Status status;
863       if ( !algo->CheckHypothesis( *mesh, edge, status ))
864         algo = the1dAlgo;
865     }
866
867     TmpMesh tmpMesh;
868     tmpMesh.ShapeToMesh( branchEdge );
869     try {
870       mesh->GetGen()->Compute( tmpMesh, branchEdge, SMESH_Gen::SHAPE_ONLY_UPWARD ); // make nodes on VERTEXes
871       if ( !algo->Compute( tmpMesh, branchEdge ))
872         return false;
873     }
874     catch (...) {
875       return false;
876     }
877     return SMESH_Algo::GetNodeParamOnEdge( tmpMesh.GetMeshDS(), branchEdge, theMAParams );
878   }
879
880   //================================================================================
881   /*!
882    * \brief Select division parameters on MA and make them coincide at ends with
883    *        projections of VERTEXes to MA for a given pair of opposite EDGEs
884    *  \param [in] theEdgePairInd - index of the EDGE pair
885    *  \param [in] theDivPoints - the BranchPoint's dividing MA into parts each
886    *         corresponding to a unique pair of opposite EDGEs
887    *  \param [in] theMAParams - the MA division parameters
888    *  \param [out] theSelectedMAParams - the selected MA parameters
889    *  \return bool - is OK
890    */
891   //================================================================================
892
893   bool getParamsForEdgePair( const size_t                              theEdgePairInd,
894                              const vector< SMESH_MAT2d::BranchPoint >& theDivPoints,
895                              const vector<double>&                     theMAParams,
896                              vector<double>&                           theSelectedMAParams)
897   {
898     if ( theDivPoints.empty() )
899     {
900       theSelectedMAParams = theMAParams;
901       return true;
902     }
903     if ( theEdgePairInd > theDivPoints.size() || theMAParams.empty() )
904       return false;
905
906     // find a range of params to copy
907
908     double par1 = 0;
909     size_t iPar1 = 0;
910     if ( theEdgePairInd > 0 )
911     {
912       const SMESH_MAT2d::BranchPoint& bp = theDivPoints[ theEdgePairInd-1 ];
913       bp._branch->getParameter( bp, par1 );
914       while ( theMAParams[ iPar1 ] < par1 ) ++iPar1;
915       if ( par1 - theMAParams[ iPar1-1 ] < theMAParams[ iPar1 ] - par1 )
916         --iPar1;
917     }
918
919     double par2 = 1;
920     size_t iPar2 = theMAParams.size() - 1;
921     if ( theEdgePairInd < theDivPoints.size() )
922     {
923       const SMESH_MAT2d::BranchPoint& bp = theDivPoints[ theEdgePairInd ];
924       bp._branch->getParameter( bp, par2 );
925       iPar2 = iPar1;
926       while ( theMAParams[ iPar2 ] < par2 ) ++iPar2;
927       if ( par2 - theMAParams[ iPar2-1 ] < theMAParams[ iPar2 ] - par2 )
928         --iPar2;
929     }
930
931     theSelectedMAParams.assign( theMAParams.begin() + iPar1,
932                                 theMAParams.begin() + iPar2 + 1 );
933
934     // adjust theSelectedMAParams to fit between par1 and par2
935
936     double d = par1 - theSelectedMAParams[0];
937     double f = ( par2 - par1 ) / ( theSelectedMAParams.back() - theSelectedMAParams[0] );
938
939     for ( size_t i = 0; i < theSelectedMAParams.size(); ++i )
940     {
941       theSelectedMAParams[i] += d;
942       theSelectedMAParams[i] = par1 + ( theSelectedMAParams[i] - par1 ) * f;
943     }
944
945     return true;
946   }
947
948   //--------------------------------------------------------------------------------
949   // node or node parameter on EDGE
950   struct NodePoint
951   {
952     const SMDS_MeshNode* _node;
953     double               _u;
954     size_t               _edgeInd; // index in theSinuEdges vector
955
956     NodePoint(): _node(0), _u(0), _edgeInd(-1) {}
957     NodePoint(const SMDS_MeshNode* n, double u, size_t iEdge ): _node(n), _u(u), _edgeInd(iEdge) {}
958     NodePoint(double u, size_t iEdge) : _node(0), _u(u), _edgeInd(iEdge) {}
959     NodePoint(const SMESH_MAT2d::BoundaryPoint& p) : _node(0), _u(p._param), _edgeInd(p._edgeIndex) {}
960     gp_Pnt Point(const vector< Handle(Geom_Curve) >& curves) const
961     {
962       return _node ? SMESH_TNodeXYZ(_node) : curves[ _edgeInd ]->Value( _u );
963     }
964   };
965   typedef multimap< double, pair< NodePoint, NodePoint > > TMAPar2NPoints;
966
967   //================================================================================
968   /*!
969    * \brief Finds a VERTEX corresponding to a point on EDGE, which is also filled
970    *        with a node on the VERTEX, present or created
971    *  \param [in,out] theNodePnt - the node position on the EDGE
972    *  \param [in] theSinuEdges - the sinuous EDGEs
973    *  \param [in] theMeshDS - the mesh
974    *  \return bool - true if the \a theBndPnt is on VERTEX
975    */
976   //================================================================================
977
978   bool findVertexAndNode( NodePoint&                  theNodePnt,
979                           const vector<TopoDS_Edge>&  theSinuEdges,
980                           SMESHDS_Mesh*               theMeshDS = 0,
981                           size_t                      theEdgeIndPrev = 0,
982                           size_t                      theEdgeIndNext = 0)
983   {
984     if ( theNodePnt._edgeInd >= theSinuEdges.size() )
985       return false;
986
987     double f,l;
988     BRep_Tool::Range( theSinuEdges[ theNodePnt._edgeInd ], f,l );
989     const double tol = 1e-3 * ( l - f );
990
991     TopoDS_Vertex V;
992     if      ( Abs( f - theNodePnt._u ) < tol )
993       V = SMESH_MesherHelper::IthVertex( 0, theSinuEdges[ theNodePnt._edgeInd ], /*CumOri=*/false);
994     else if ( Abs( l - theNodePnt._u ) < tol )
995       V = SMESH_MesherHelper::IthVertex( 1, theSinuEdges[ theNodePnt._edgeInd ], /*CumOri=*/false);
996     else if ( theEdgeIndPrev != theEdgeIndNext )
997       TopExp::CommonVertex( theSinuEdges[theEdgeIndPrev], theSinuEdges[theEdgeIndNext], V );
998
999     if ( !V.IsNull() && theMeshDS )
1000     {
1001       theNodePnt._node = SMESH_Algo::VertexNode( V, theMeshDS );
1002       if ( !theNodePnt._node )
1003       {
1004         gp_Pnt p = BRep_Tool::Pnt( V );
1005         theNodePnt._node = theMeshDS->AddNode( p.X(), p.Y(), p.Z() );
1006         theMeshDS->SetNodeOnVertex( theNodePnt._node, V );
1007       }
1008     }
1009     return !V.IsNull();
1010   }
1011
1012   //================================================================================
1013   /*!
1014    * \brief Add to the map of NodePoint's those on VERTEXes
1015    *  \param [in,out] theHelper - the helper
1016    *  \param [in] theMA - Medial Axis
1017    *  \param [in] theMinSegLen - minimal segment length
1018    *  \param [in] theDivPoints - projections of VERTEXes to MA
1019    *  \param [in] theSinuEdges - the sinuous EDGEs
1020    *  \param [in] theSideEdgeIDs - indices of sinuous EDGEs per side
1021    *  \param [in] theIsEdgeComputed - is sinuous EGDE is meshed
1022    *  \param [in,out] thePointsOnE - the map to fill
1023    *  \param [out] theNodes2Merge - the map of nodes to merge
1024    */
1025   //================================================================================
1026
1027   bool projectVertices( SMESH_MesherHelper&                 theHelper,
1028                         const SMESH_MAT2d::MedialAxis&      theMA,
1029                         vector< SMESH_MAT2d::BranchPoint >& theDivPoints,
1030                         const vector< std::size_t > &       theEdgeIDs1,
1031                         const vector< std::size_t > &       theEdgeIDs2,
1032                         const vector< bool >&               theIsEdgeComputed,
1033                         TMAPar2NPoints &                    thePointsOnE,
1034                         SinuousFace&                        theSinuFace)
1035   {
1036     if ( theDivPoints.empty() )
1037       return true;
1038
1039     SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
1040     const vector< TopoDS_Edge >&     theSinuEdges = theSinuFace._sinuEdges;
1041     const vector< Handle(Geom_Curve) >& theCurves = theSinuFace._sinuCurves;
1042
1043     double uMA;
1044     SMESH_MAT2d::BoundaryPoint bp[2]; // 2 sinuous sides
1045     const SMESH_MAT2d::Branch& branch = *theMA.getBranch(0);
1046     {
1047       // add to thePointsOnE NodePoint's of ends of theSinuEdges
1048       if ( !branch.getBoundaryPoints( 0., bp[0], bp[1] ) ||
1049            !theMA.getBoundary().moveToClosestEdgeEnd( bp[0] )) return false;
1050       if ( !theSinuFace.IsRing() &&
1051            !theMA.getBoundary().moveToClosestEdgeEnd( bp[1] )) return false;
1052       NodePoint np0( bp[0] ), np1( bp[1] );
1053       findVertexAndNode( np0, theSinuEdges, meshDS );
1054       findVertexAndNode( np1, theSinuEdges, meshDS );
1055       thePointsOnE.insert( make_pair( -0.1, make_pair( np0, np1 )));
1056     }
1057     if ( !theSinuFace.IsRing() )
1058     {
1059       if ( !branch.getBoundaryPoints( 1., bp[0], bp[1] ) ||
1060            !theMA.getBoundary().moveToClosestEdgeEnd( bp[0] ) ||
1061            !theMA.getBoundary().moveToClosestEdgeEnd( bp[1] )) return false;
1062       NodePoint np0( bp[0] ), np1( bp[1] );
1063       findVertexAndNode( np0, theSinuEdges, meshDS );
1064       findVertexAndNode( np1, theSinuEdges, meshDS );
1065       thePointsOnE.insert( make_pair( 1.1, make_pair( np0, np1)));
1066     }
1067     else
1068     {
1069       // project a VERTEX of outer sinuous side corresponding to branch(0.)
1070       // which is not included into theDivPoints
1071       if ( ! ( theDivPoints[0]._iEdge     == 0 &&
1072                theDivPoints[0]._edgeParam == 0. )) // recursive call
1073       {
1074         SMESH_MAT2d::BranchPoint brp( &branch, 0, 0. );
1075         vector< SMESH_MAT2d::BranchPoint > divPoint( 1, brp );
1076         vector< std::size_t > edgeIDs1(2), edgeIDs2(2);
1077         edgeIDs1[0] = theEdgeIDs1.back();
1078         edgeIDs1[1] = theEdgeIDs1[0];
1079         edgeIDs2[0] = theEdgeIDs2.back();
1080         edgeIDs2[1] = theEdgeIDs2[0];
1081         projectVertices( theHelper, theMA, divPoint, edgeIDs1, edgeIDs2,
1082                          theIsEdgeComputed, thePointsOnE, theSinuFace );
1083       }
1084     }
1085
1086     // project theDivPoints and keep projections to merge
1087
1088     TMAPar2NPoints::iterator u2NP;
1089     vector< TMAPar2NPoints::iterator > projToMerge;
1090     for ( size_t i = 0; i < theDivPoints.size(); ++i )
1091     {
1092       if ( !branch.getParameter( theDivPoints[i], uMA ))
1093         return false;
1094       if ( !branch.getBoundaryPoints( theDivPoints[i], bp[0], bp[1] ))
1095         return false;
1096
1097       NodePoint  np[2] = {
1098         NodePoint( bp[0] ),
1099         NodePoint( bp[1] )
1100       };
1101       bool isVertex[2] = {
1102         findVertexAndNode( np[0], theSinuEdges, meshDS, theEdgeIDs1[i], theEdgeIDs1[i+1] ),
1103         findVertexAndNode( np[1], theSinuEdges, meshDS, theEdgeIDs2[i], theEdgeIDs2[i+1] )
1104       };
1105       const size_t iVert = isVertex[0] ? 0 : 1; // side with a VERTEX
1106       const size_t iNode = 1 - iVert;           // opposite (meshed?) side
1107
1108       if ( isVertex[0] != isVertex[1] ) // try to find an opposite VERTEX
1109       {
1110         theMA.getBoundary().moveToClosestEdgeEnd( bp[iNode] ); // EDGE -> VERTEX
1111         SMESH_MAT2d::BranchPoint brp;
1112         theMA.getBoundary().getBranchPoint( bp[iNode], brp );  // WIRE -> MA
1113         SMESH_MAT2d::BoundaryPoint bp2[2];
1114         branch.getBoundaryPoints( brp, bp2[0], bp2[1] );       // MA -> WIRE
1115         NodePoint np2[2] = { NodePoint( bp2[0]), NodePoint( bp2[1]) };
1116         findVertexAndNode( np2[0], theSinuEdges, meshDS );
1117         findVertexAndNode( np2[1], theSinuEdges, meshDS );
1118         if ( np2[ iVert ]._node == np[ iVert ]._node &&
1119              np2[ iNode ]._node)
1120         {
1121           np[ iNode ] = np2[ iNode ];
1122           isVertex[ iNode ] = true;
1123         }
1124       }
1125
1126       u2NP = thePointsOnE.insert( make_pair( uMA, make_pair( np[0], np[1])));
1127
1128       if ( !isVertex[0] && !isVertex[1] ) return false; // error
1129       if ( isVertex[0] && isVertex[1] )
1130         continue;
1131
1132       // bool isOppComputed = theIsEdgeComputed[ np[ iNode ]._edgeInd ];
1133       // if ( isOppComputed )
1134         projToMerge.push_back( u2NP );
1135     }
1136
1137     // merge projections
1138
1139     for ( size_t i = 0; i < projToMerge.size(); ++i )
1140     {
1141       u2NP = projToMerge[i];
1142       const size_t iVert = get( u2NP->second, 0 )._node ? 0 : 1; // side with a VERTEX
1143       const size_t iNode = 1 - iVert;                            // opposite (meshed?) side
1144
1145       // a VERTEX is projected on a meshed EDGE; there are two options:
1146       // 1) a projected point is joined with a closet node if a strip between this and neighbor
1147       // projection is WIDE enough; joining is done by creating a node coincident with the
1148       //  existing node which will be merged together after all;
1149       // 2) a neighbor projection is merged with this one if it is TOO CLOSE; a node of deleted
1150       // projection is set to the BoundaryPoint of this projection
1151
1152       // evaluate distance to neighbor projections
1153       const double rShort = 0.33;
1154       bool isShortPrev[2], isShortNext[2], isPrevCloser[2];
1155       TMAPar2NPoints::iterator u2NPPrev = u2NP, u2NPNext = u2NP;
1156       --u2NPPrev; ++u2NPNext;
1157       for ( int iS = 0; iS < 2; ++iS ) // side with Vertex and side with Nodes
1158       {
1159         NodePoint np     = get( u2NP->second,     iS );
1160         NodePoint npPrev = get( u2NPPrev->second, iS );
1161         NodePoint npNext = get( u2NPNext->second, iS );
1162         gp_Pnt     p     = np    .Point( theCurves );
1163         gp_Pnt     pPrev = npPrev.Point( theCurves );
1164         gp_Pnt     pNext = npNext.Point( theCurves );
1165         double  distPrev = p.Distance( pPrev );
1166         double  distNext = p.Distance( pNext );
1167         double         r = distPrev / ( distPrev + distNext );
1168         isShortPrev [iS] = ( r < rShort );
1169         isShortNext [iS] = (( 1 - r ) < rShort );
1170         isPrevCloser[iS] = (( r < 0.5 ) && ( theSinuFace.IsRing() || u2NPPrev->first > 0 ));
1171       }
1172
1173       TMAPar2NPoints::iterator u2NPClose;
1174
1175       if (( isShortPrev[0] && isShortPrev[1] ) || // option 2) -> remove a too close projection
1176           ( isShortNext[0] && isShortNext[1] ))
1177       {
1178         u2NPClose = isPrevCloser[0] ? u2NPPrev : u2NPNext;
1179         NodePoint& npProj  = get( u2NP->second,      iNode ); // NP of VERTEX projection
1180         NodePoint& npVert  = get( u2NP->second,      iVert ); // NP of VERTEX
1181         NodePoint npCloseN = get( u2NPClose->second, iNode ); // NP close to npProj
1182         NodePoint npCloseV = get( u2NPClose->second, iVert ); // NP close to npVert
1183         if ( !npCloseV._node || npCloseV._node == npVert._node )
1184         {
1185           npProj = npCloseN;
1186           thePointsOnE.erase( u2NPClose );
1187           continue;
1188         }
1189         else
1190         {
1191           // can't remove the neighbor projection as it is also from VERTEX -> option 1)
1192         }
1193       }
1194       // else: option 1) - wide enough -> "duplicate" existing node
1195       {
1196         u2NPClose = isPrevCloser[ iNode ] ? u2NPPrev : u2NPNext;
1197         NodePoint& npProj   = get( u2NP->second,      iNode ); // NP of VERTEX projection
1198         NodePoint& npCloseN = get( u2NPClose->second, iNode ); // NP close to npProj
1199         npProj = npCloseN;
1200         npProj._node = 0;
1201         //npProj._edgeInd = npCloseN._edgeInd;
1202         // npProj._u       = npCloseN._u + 1e-3 * Abs( get( u2NPPrev->second, iNode )._u -
1203         //                                             get( u2NPNext->second, iNode )._u );
1204         // gp_Pnt        p = npProj.Point( theCurves );
1205         // npProj._node    = meshDS->AddNode( p.X(), p.Y(), p.Z() );
1206         // meshDS->SetNodeOnEdge( npProj._node, theSinuEdges[ npProj._edgeInd ], npProj._u  );
1207
1208         //theNodes2Merge[ npCloseN._node ].push_back( npProj._node );
1209       }
1210     }
1211
1212     // remove auxiliary NodePoint's of ends of theSinuEdges
1213     for ( u2NP = thePointsOnE.begin(); u2NP->first < 0; )
1214       thePointsOnE.erase( u2NP++ );
1215     thePointsOnE.erase( 1.1 );
1216
1217     return true;
1218   }
1219
1220   double getUOnEdgeByPoint( const size_t     iEdge,
1221                             const NodePoint* point,
1222                             SinuousFace&     sinuFace )
1223   {
1224     if ( point->_edgeInd == iEdge )
1225       return point->_u;
1226
1227     TopoDS_Vertex V0 = TopExp::FirstVertex( sinuFace._sinuEdges[ iEdge ]);
1228     TopoDS_Vertex V1 = TopExp::LastVertex ( sinuFace._sinuEdges[ iEdge ]);
1229     gp_Pnt p0 = BRep_Tool::Pnt( V0 );
1230     gp_Pnt p1 = BRep_Tool::Pnt( V1 );
1231     gp_Pnt  p = point->Point( sinuFace._sinuCurves );
1232
1233     double f,l;
1234     BRep_Tool::Range( sinuFace._sinuEdges[ iEdge ], f,l );
1235     return p.SquareDistance( p0 ) < p.SquareDistance( p1 ) ? f : l;
1236   }
1237
1238   //================================================================================
1239   /*!
1240    * \brief Move coincident nodes to make node params on EDGE unique
1241    *  \param [in] theHelper - the helper
1242    *  \param [in] thePointsOnE - nodes on two opposite river sides
1243    *  \param [in] theSinuFace - the sinuous FACE
1244    *  \param [out] theNodes2Merge - the map of nodes to merge
1245    */
1246   //================================================================================
1247
1248   void separateNodes( SMESH_MesherHelper&            theHelper,
1249                       const SMESH_MAT2d::MedialAxis& theMA,
1250                       TMAPar2NPoints &               thePointsOnE,
1251                       SinuousFace&                   theSinuFace,
1252                       const vector< bool >&          theIsComputedEdge)
1253   {
1254     if ( thePointsOnE.size() < 2 )
1255       return;
1256
1257     SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
1258     const vector<TopoDS_Edge>&    theSinuEdges = theSinuFace._sinuEdges;
1259     const vector< Handle(Geom_Curve) >& curves = theSinuFace._sinuCurves;
1260
1261     //SMESH_MAT2d::BoundaryPoint bp[2];
1262     //const SMESH_MAT2d::Branch& branch = *theMA.getBranch(0);
1263
1264     typedef TMAPar2NPoints::iterator TIterator;
1265
1266     for ( int iSide = 0; iSide < 2; ++iSide ) // loop on two sinuous sides
1267     {
1268       // get a tolerance to compare points
1269       double tol = Precision::Confusion();
1270       for ( size_t i = 0; i < theSinuFace._sinuSide[ iSide ].size(); ++i )
1271         tol = Max( tol , BRep_Tool::Tolerance( theSinuFace._sinuSide[ iSide ][ i ]));
1272
1273       // find coincident points
1274       TIterator u2NP = thePointsOnE.begin();
1275       vector< TIterator > sameU2NP( 1, u2NP++ );
1276       while ( u2NP != thePointsOnE.end() )
1277       {
1278         for ( ; u2NP != thePointsOnE.end(); ++u2NP )
1279         {
1280           NodePoint& np1 = get( sameU2NP.back()->second, iSide );
1281           NodePoint& np2 = get( u2NP           ->second, iSide );
1282
1283           if (( !np1._node || !np2._node ) &&
1284               ( np1.Point( curves ).SquareDistance( np2.Point( curves )) < tol*tol ))
1285           {
1286             sameU2NP.push_back( u2NP );
1287           }
1288           else if ( sameU2NP.size() == 1 )
1289           {
1290             sameU2NP[ 0 ] = u2NP;
1291           }
1292           else
1293           {
1294             break;
1295           }
1296         }
1297
1298         if ( sameU2NP.size() > 1 )
1299         {
1300           // find an existing node on VERTEX among sameU2NP and get underlying EDGEs
1301           const SMDS_MeshNode* existingNode = 0;
1302           set< size_t > edgeInds;
1303           NodePoint* np;
1304           for ( size_t i = 0; i < sameU2NP.size(); ++i )
1305           {
1306             np = &get( sameU2NP[i]->second, iSide );
1307             if ( np->_node )
1308               if ( !existingNode || np->_node->GetPosition()->GetDim() == 0 )
1309                 existingNode = np->_node;
1310             edgeInds.insert( np->_edgeInd );
1311           }
1312           list< const SMDS_MeshNode* >& mergeNodes = theSinuFace._nodesToMerge[ existingNode ];
1313
1314           TIterator u2NPprev = sameU2NP.front();
1315           TIterator u2NPnext = sameU2NP.back() ;
1316           if ( u2NPprev->first < 0. ) ++u2NPprev;
1317           if ( u2NPnext->first > 1. ) --u2NPnext;
1318
1319           set< size_t >::iterator edgeID = edgeInds.begin();
1320           for ( ; edgeID != edgeInds.end(); ++edgeID )
1321           {
1322             // get U range on iEdge within which the equal points will be distributed
1323             double u0, u1;
1324             np = &get( u2NPprev->second, iSide );
1325             u0 = getUOnEdgeByPoint( *edgeID, np, theSinuFace );
1326
1327             np = &get( u2NPnext->second, iSide );
1328             u1 = getUOnEdgeByPoint( *edgeID, np, theSinuFace );
1329
1330             if ( u0 == u1 )
1331             {
1332               if ( u2NPprev != thePointsOnE.begin() ) --u2NPprev;
1333               if ( u2NPnext != --thePointsOnE.end() ) ++u2NPnext;
1334               np = &get( u2NPprev->second, iSide );
1335               u0 = getUOnEdgeByPoint( *edgeID, np, theSinuFace );
1336               np = &get( u2NPnext->second, iSide );
1337               u1 = getUOnEdgeByPoint( *edgeID, np, theSinuFace );
1338             }
1339
1340             // distribute points and create nodes
1341             double du = ( u1 - u0 ) / ( sameU2NP.size() + 1 /*!existingNode*/ );
1342             double u  = u0 + du;
1343             for ( size_t i = 0; i < sameU2NP.size(); ++i )
1344             {
1345               np = &get( sameU2NP[i]->second, iSide );
1346               if ( !np->_node && *edgeID == np->_edgeInd )
1347               {
1348                 np->_u = u;
1349                 u += du;
1350                 gp_Pnt p = np->Point( curves );
1351                 np->_node = meshDS->AddNode( p.X(), p.Y(), p.Z() );
1352                 meshDS->SetNodeOnEdge( np->_node, theSinuEdges[ *edgeID ], np->_u  );
1353
1354                 if ( theIsComputedEdge[ *edgeID ])
1355                   mergeNodes.push_back( np->_node );
1356               }
1357             }
1358           }
1359
1360           sameU2NP.resize( 1 );
1361           u2NP = ++sameU2NP.back();
1362           sameU2NP[ 0 ] = u2NP;
1363
1364         } // if ( sameU2NP.size() > 1 )
1365       } // while ( u2NP != thePointsOnE.end() )
1366     } // for ( int iSide = 0; iSide < 2; ++iSide )
1367
1368     return;
1369   } // separateNodes()
1370
1371
1372   //================================================================================
1373   /*!
1374    * \brief Find association of nodes existing on the sinuous sides of a ring
1375    *
1376    * TMAPar2NPoints filled here is used in setQuadSides() only if theSinuFace.IsRing()
1377    * to find most distant nodes of the inner and the outer wires
1378    */
1379   //================================================================================
1380
1381   void assocNodes( SMESH_MesherHelper&            theHelper,
1382                    SinuousFace&                   theSinuFace,
1383                    const SMESH_MAT2d::MedialAxis& theMA,
1384                    TMAPar2NPoints &               thePointsOnE )
1385   {
1386     SMESH_Mesh*     mesh = theHelper.GetMesh();
1387     SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
1388
1389     list< TopoDS_Edge > ee1( theSinuFace._sinuSide [0].begin(), theSinuFace._sinuSide [0].end() );
1390     list< TopoDS_Edge > ee2( theSinuFace._sinuSide [1].begin(), theSinuFace._sinuSide [1].end() );
1391     StdMeshers_FaceSide sideOut( theSinuFace.Face(), ee1, mesh, true, true, &theHelper );
1392     StdMeshers_FaceSide sideIn ( theSinuFace.Face(), ee2, mesh, true, true, &theHelper );
1393     const UVPtStructVec& uvsOut = sideOut.GetUVPtStruct();
1394     const UVPtStructVec& uvsIn  = sideIn.GetUVPtStruct();
1395     // if ( uvs1.size() != uvs2.size() )
1396     //   return;
1397
1398     const SMESH_MAT2d::Branch& branch = *theMA.getBranch(0);
1399     SMESH_MAT2d::BoundaryPoint bp[2];
1400     SMESH_MAT2d::BranchPoint   brp;
1401
1402     map< double, const SMDS_MeshNode* > nodeParams; // params of existing nodes
1403     map< double, const SMDS_MeshNode* >::iterator u2n;
1404
1405     // find a node of sideOut most distant from sideIn
1406
1407     vector< BRepAdaptor_Curve > curvesIn( theSinuFace._sinuSide[1].size() );
1408     for ( size_t iE = 0; iE < theSinuFace._sinuSide[1].size(); ++iE )
1409       curvesIn[ iE ].Initialize( theSinuFace._sinuSide[1][iE] );
1410
1411     double maxDist = 0;
1412     SMESH_MAT2d::BoundaryPoint bpIn; // closest IN point
1413     const SMDS_MeshNode*        nOut = 0;
1414     const size_t              nbEOut = theSinuFace._sinuSide[0].size();
1415     for ( size_t iE = 0; iE < nbEOut; ++iE )
1416     {
1417       const TopoDS_Edge& E = theSinuFace._sinuSide[0][iE];
1418
1419       if ( !SMESH_Algo::GetSortedNodesOnEdge( meshDS, E, /*skipMedium=*/true, nodeParams ))
1420         return;
1421       for ( u2n = nodeParams.begin(); u2n != nodeParams.end(); ++u2n )
1422       {
1423         // point on EDGE (u2n) --> MA point (brp)
1424         if ( !theMA.getBoundary().getBranchPoint( iE, u2n->first, brp ) ||
1425              !branch.getBoundaryPoints( brp, bp[0], bp[1] ))
1426           return;
1427         gp_Pnt pOut = SMESH_TNodeXYZ( u2n->second );
1428         gp_Pnt pIn  = curvesIn[ bp[1]._edgeIndex - nbEOut ].Value( bp[1]._param );
1429         double dist = pOut.SquareDistance( pIn );
1430         if ( dist > maxDist )
1431         {
1432           maxDist = dist;
1433           nOut    = u2n->second;
1434           bpIn    = bp[1];
1435         }
1436       }
1437     }
1438     const SMDS_MeshNode* nIn = 0;
1439     if ( !SMESH_Algo::GetSortedNodesOnEdge( meshDS,
1440                                             theSinuFace._sinuEdges[ bpIn._edgeIndex ],
1441                                             /*skipMedium=*/true,
1442                                             nodeParams ))
1443       return;
1444     u2n = nodeParams.lower_bound( bpIn._param );
1445     if ( u2n == nodeParams.end() )
1446       nIn = nodeParams.rbegin()->second;
1447     else
1448       nIn = u2n->second;
1449     
1450     // find position of distant nodes in uvsOut and uvsIn
1451     size_t iDistOut, iDistIn;
1452     for ( iDistOut = 0; iDistOut < uvsOut.size(); ++iDistOut )
1453     {
1454       if ( uvsOut[iDistOut].node == nOut )
1455         break;
1456     }
1457     for ( iDistIn = 0; iDistIn < uvsIn.size(); ++iDistIn )
1458     {
1459       if ( uvsIn[iDistIn].node == nIn )
1460         break;
1461     }
1462     if ( iDistOut == uvsOut.size() || iDistIn == uvsIn.size() )
1463       return;
1464
1465     // store opposite nodes in thePointsOnE (param and EDGE have no sense)
1466     pair< NodePoint, NodePoint > oppNodes( NodePoint( nOut, 0, 0 ), NodePoint( nIn, 0, 0));
1467     thePointsOnE.insert( make_pair( uvsOut[ iDistOut ].normParam, oppNodes ));
1468     int iOut = iDistOut, iIn = iDistIn;
1469     int i, nbNodes = std::min( uvsOut.size(), uvsIn.size() );
1470     if ( nbNodes > 5 ) nbNodes = 5;
1471     for ( i = 0, ++iOut, --iIn; i < nbNodes; ++iOut, --iIn, ++i )
1472     {
1473       iOut = theHelper.WrapIndex( iOut, uvsOut.size() );
1474       iIn  = theHelper.WrapIndex( iIn,  uvsIn.size()  );
1475       oppNodes.first._node  = uvsOut[ iOut ].node;
1476       oppNodes.second._node = uvsIn[ iIn ].node;
1477       thePointsOnE.insert( make_pair( uvsOut[ iOut ].normParam, oppNodes ));
1478     }
1479
1480     return;
1481   } // assocNodes()
1482
1483   //================================================================================
1484   /*!
1485    * \brief Setup sides of SinuousFace::_quad
1486    *  \param [in] theHelper - helper
1487    *  \param [in] thePointsOnEdges - NodePoint's on sinuous sides
1488    *  \param [in,out] theSinuFace - the FACE
1489    *  \param [in] the1dAlgo - algorithm to use for radial discretization of a ring FACE
1490    *  \return bool - is OK
1491    */
1492   //================================================================================
1493
1494   bool setQuadSides(SMESH_MesherHelper&   theHelper,
1495                     const TMAPar2NPoints& thePointsOnEdges,
1496                     SinuousFace&          theFace,
1497                     SMESH_Algo*           the1dAlgo)
1498   {
1499     SMESH_Mesh*               mesh = theHelper.GetMesh();
1500     const TopoDS_Face&        face = theFace._quad->face;
1501     SMESH_ProxyMesh::Ptr proxyMesh = StdMeshers_ViscousLayers2D::Compute( *mesh, face );
1502     if ( !proxyMesh )
1503       return false;
1504
1505     list< TopoDS_Edge > side[4];
1506     side[0].insert( side[0].end(), theFace._shortSide[0].begin(), theFace._shortSide[0].end() );
1507     side[1].insert( side[1].end(), theFace._sinuSide [1].begin(), theFace._sinuSide [1].end() );
1508     side[2].insert( side[2].end(), theFace._shortSide[1].begin(), theFace._shortSide[1].end() );
1509     side[3].insert( side[3].end(), theFace._sinuSide [0].begin(), theFace._sinuSide [0].end() );
1510
1511     for ( int i = 0; i < 4; ++i )
1512     {
1513       theFace._quad->side[i] = StdMeshers_FaceSide::New( face, side[i], mesh, i < QUAD_TOP_SIDE,
1514                                                          /*skipMediumNodes=*/true,
1515                                                          &theHelper, proxyMesh );
1516     }
1517
1518     if ( theFace.IsRing() )
1519     {
1520       // --------------------------------------
1521       // Discretize a ring in radial direction
1522       // --------------------------------------
1523
1524       if ( thePointsOnEdges.size() < 4 )
1525         return false;
1526
1527       int nbOut = theFace._quad->side[ 1 ].GetUVPtStruct().size();
1528       int nbIn  = theFace._quad->side[ 3 ].GetUVPtStruct().size();
1529       if ( nbOut == 0 || nbIn == 0 )
1530         return false;
1531
1532       // find most distant opposite nodes
1533       double maxDist = 0, dist;
1534       TMAPar2NPoints::const_iterator u2NPdist, u2NP = thePointsOnEdges.begin();
1535       for ( ; u2NP != thePointsOnEdges.end(); ++u2NP )
1536       {
1537         SMESH_TNodeXYZ        xyz( u2NP->second.first._node ); // node out
1538         dist = xyz.SquareDistance( u2NP->second.second._node );// node in
1539         if ( dist > maxDist )
1540         {
1541           u2NPdist = u2NP;
1542           maxDist  = dist;
1543         }
1544       }
1545       // compute distribution of radial nodes
1546       list< double > params; // normalized params
1547       static_cast< StdMeshers_QuadFromMedialAxis_1D2D::Algo1D* >
1548         ( the1dAlgo )->ComputeDistribution( theHelper,
1549                                             SMESH_TNodeXYZ( u2NPdist->second.first._node ),
1550                                             SMESH_TNodeXYZ( u2NPdist->second.second._node ),
1551                                             params );
1552
1553       // add a radial quad side
1554
1555       theHelper.SetElementsOnShape( true );
1556       u2NP = thePointsOnEdges.begin();
1557       const SMDS_MeshNode* nOut = u2NP->second.first._node;
1558       const SMDS_MeshNode*  nIn = u2NP->second.second._node;
1559       nOut = proxyMesh->GetProxyNode( nOut );
1560       nIn  = proxyMesh->GetProxyNode( nIn );
1561       gp_XY uvOut = theHelper.GetNodeUV( face, nOut );
1562       gp_XY uvIn  = theHelper.GetNodeUV( face, nIn );
1563       Handle(Geom_Surface) surface = BRep_Tool::Surface( face );
1564       UVPtStructVec uvsNew; UVPtStruct uvPt;
1565       uvPt.node = nOut;
1566       uvPt.u    = uvOut.X();
1567       uvPt.v    = uvOut.Y();
1568       uvsNew.push_back( uvPt );
1569       for (list<double>::iterator itU = params.begin(); itU != params.end(); ++itU )
1570       {
1571         gp_XY uv  = ( 1 - *itU ) * uvOut + *itU * uvIn; // applied in direction Out -> In
1572         gp_Pnt p  = surface->Value( uv.X(), uv.Y() );
1573         uvPt.node = theHelper.AddNode( p.X(), p.Y(), p.Z(), /*id=*/0, uv.X(), uv.Y() );
1574         uvPt.u    = uv.X();
1575         uvPt.v    = uv.Y();
1576         uvsNew.push_back( uvPt );
1577       }
1578       uvPt.node = nIn;
1579       uvPt.u    = uvIn.X();
1580       uvPt.v    = uvIn.Y();
1581       uvsNew.push_back( uvPt );
1582
1583       theFace._quad->side[ 0 ] = StdMeshers_FaceSide::New( uvsNew );
1584       theFace._quad->side[ 2 ] = theFace._quad->side[ 0 ];
1585       if ( nbIn != nbOut )
1586         theFace._quad->side[ 2 ] = StdMeshers_FaceSide::New( uvsNew );
1587
1588       // assure that the outer sinuous side starts at nOut
1589       {
1590         const UVPtStructVec& uvsOut = theFace._quad->side[ 3 ].GetUVPtStruct(); // _sinuSide[0]
1591         size_t i; // find UVPtStruct holding nOut
1592         for ( i = 0; i < uvsOut.size(); ++i )
1593           if ( nOut == uvsOut[i].node )
1594             break;
1595         if ( i == uvsOut.size() )
1596           return false;
1597
1598         if ( i != 0  &&  i != uvsOut.size()-1 )
1599         {
1600           // create a new OUT quad side
1601           uvsNew.clear();
1602           uvsNew.reserve( uvsOut.size() );
1603           uvsNew.insert( uvsNew.end(), uvsOut.begin() + i, uvsOut.end() );
1604           uvsNew.insert( uvsNew.end(), uvsOut.begin() + 1, uvsOut.begin() + i + 1);
1605           theFace._quad->side[ 3 ] = StdMeshers_FaceSide::New( uvsNew );
1606         }
1607       }
1608
1609       // rotate the IN side if opposite nodes of IN and OUT sides don't match
1610
1611       const SMDS_MeshNode * nIn0 = theFace._quad->side[ 1 ].First().node;
1612       if ( nIn0 != nIn )
1613       {
1614         nIn  = proxyMesh->GetProxyNode( nIn );
1615         const UVPtStructVec& uvsIn = theFace._quad->side[ 1 ].GetUVPtStruct(); // _sinuSide[1]
1616         size_t i; // find UVPtStruct holding nIn
1617         for ( i = 0; i < uvsIn.size(); ++i )
1618           if ( nIn == uvsIn[i].node )
1619             break;
1620         if ( i == uvsIn.size() )
1621           return false;
1622
1623         // create a new IN quad side
1624         uvsNew.clear();
1625         uvsNew.reserve( uvsIn.size() );
1626         uvsNew.insert( uvsNew.end(), uvsIn.begin() + i, uvsIn.end() );
1627         uvsNew.insert( uvsNew.end(), uvsIn.begin() + 1, uvsIn.begin() + i + 1);
1628         theFace._quad->side[ 1 ] = StdMeshers_FaceSide::New( uvsNew );
1629       }
1630
1631       if ( theFace._quad->side[ 1 ].GetUVPtStruct().empty() ||
1632            theFace._quad->side[ 3 ].GetUVPtStruct().empty() )
1633         return false;
1634
1635     } // if ( theFace.IsRing() )
1636
1637     return true;
1638
1639   } // setQuadSides()
1640
1641   //================================================================================
1642   /*!
1643    * \brief Divide the sinuous EDGEs by projecting the division point of Medial
1644    *        Axis to the EGDEs
1645    *  \param [in] theHelper - the helper
1646    *  \param [in] theMinSegLen - minimal segment length
1647    *  \param [in] theMA - the Medial Axis
1648    *  \param [in] theMAParams - parameters of division points of \a theMA
1649    *  \param [in] theSinuEdges - the EDGEs to make nodes on
1650    *  \param [in] theSinuSide0Size - the number of EDGEs in the 1st sinuous side
1651    *  \param [in] the1dAlgo - algorithm to use for radial discretization of a ring FACE
1652    *  \return bool - is OK or not
1653    */
1654   //================================================================================
1655
1656   bool computeSinuEdges( SMESH_MesherHelper&        theHelper,
1657                          double                     /*theMinSegLen*/,
1658                          SMESH_MAT2d::MedialAxis&   theMA,
1659                          vector<double>&            theMAParams,
1660                          SinuousFace&               theSinuFace,
1661                          SMESH_Algo*                the1dAlgo)
1662   {
1663     if ( theMA.nbBranches() != 1 )
1664       return false;
1665
1666     // normalize theMAParams
1667     for ( size_t i = 0; i < theMAParams.size(); ++i )
1668       theMAParams[i] /= theMAParams.back();
1669
1670
1671     SMESH_Mesh*     mesh = theHelper.GetMesh();
1672     SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
1673     double f,l;
1674
1675     // get data of sinuous EDGEs and remove unnecessary nodes
1676     const vector< TopoDS_Edge >& theSinuEdges = theSinuFace._sinuEdges;
1677     vector< Handle(Geom_Curve) >& curves      = theSinuFace._sinuCurves;
1678     vector< int >                edgeIDs   ( theSinuEdges.size() ); // IDs in the main shape
1679     vector< bool >               isComputed( theSinuEdges.size() );
1680     curves.resize( theSinuEdges.size(), 0 );
1681     bool                         allComputed = true;
1682     for ( size_t i = 0; i < theSinuEdges.size(); ++i )
1683     {
1684       curves[i] = BRep_Tool::Curve( theSinuEdges[i], f,l );
1685       if ( !curves[i] )
1686         return false;
1687       SMESH_subMesh* sm = mesh->GetSubMesh( theSinuEdges[i] );
1688       edgeIDs   [i] = sm->GetId();
1689       isComputed[i] = ( !sm->IsEmpty() );
1690       if ( !isComputed[i] )
1691         allComputed = false;
1692     }
1693
1694     const SMESH_MAT2d::Branch& branch = *theMA.getBranch(0);
1695     SMESH_MAT2d::BoundaryPoint bp[2];
1696
1697     TMAPar2NPoints pointsOnE;
1698     // check that computed EDGEs are opposite and equally meshed
1699     if ( allComputed )
1700     {
1701       // int nbNodes[2] = { 0, 0 };
1702       // for ( int iSide = 0; iSide < 2; ++iSide ) // loop on two sinuous sides
1703       //   nbNodes[ iSide ] += meshDS->MeshElements( theSinuFace._sinuSide[ iSide ])->NbNodes() - 1;
1704
1705       // if ( nbNodes[0] != nbNodes[1] )
1706       //   return false;
1707
1708       if ( theSinuFace.IsRing() )
1709         assocNodes( theHelper, theSinuFace, theMA, pointsOnE );
1710     }
1711     else
1712     {
1713       vector< std::size_t > edgeIDs1, edgeIDs2; // indices in theSinuEdges
1714       vector< SMESH_MAT2d::BranchPoint > divPoints;
1715       branch.getOppositeGeomEdges( edgeIDs1, edgeIDs2, divPoints );
1716
1717       for ( size_t i = 0; i < edgeIDs1.size(); ++i )
1718         if ( isComputed[ edgeIDs1[i]] &&
1719              isComputed[ edgeIDs2[i]] )
1720         {
1721           int nbNodes1 = meshDS->MeshElements(edgeIDs[ edgeIDs1[i]] )->NbNodes();
1722           int nbNodes2 = meshDS->MeshElements(edgeIDs[ edgeIDs2[i]] )->NbNodes();
1723           if ( nbNodes1 != nbNodes2 )
1724             return false;
1725           if (( int(i)-1 >= 0 ) &&
1726               ( edgeIDs1[i-1] == edgeIDs1[i] ||
1727                 edgeIDs2[i-1] == edgeIDs2[i] ))
1728             return false;
1729           if (( i+1 < edgeIDs1.size() ) &&
1730               ( edgeIDs1[i+1] == edgeIDs1[i] ||
1731                 edgeIDs2[i+1] == edgeIDs2[i] ))
1732             return false;
1733         }
1734
1735       // map (param on MA) to (parameters of nodes on a pair of theSinuEdges)
1736       vector<double> maParams;
1737       set<int>       projectedEdges; // treated EDGEs which 'isComputed'
1738
1739       // compute params of nodes on EDGEs by projecting division points from MA
1740
1741       for ( size_t iEdgePair = 0; iEdgePair < edgeIDs1.size(); ++iEdgePair )
1742         // loop on pairs of opposite EDGEs
1743       {
1744         if ( projectedEdges.count( edgeIDs1[ iEdgePair ]) ||
1745              projectedEdges.count( edgeIDs2[ iEdgePair ]) )
1746           continue;
1747
1748         // --------------------------------------------------------------------------------
1749         if ( isComputed[ edgeIDs1[ iEdgePair ]] !=                    // one EDGE is meshed
1750              isComputed[ edgeIDs2[ iEdgePair ]])
1751         {
1752           // "projection" from one side to the other
1753
1754           size_t iEdgeComputed = edgeIDs1[iEdgePair], iSideComputed = 0;
1755           if ( !isComputed[ iEdgeComputed ])
1756             ++iSideComputed, iEdgeComputed = edgeIDs2[iEdgePair];
1757
1758           map< double, const SMDS_MeshNode* > nodeParams; // params of existing nodes
1759           if ( !SMESH_Algo::GetSortedNodesOnEdge( meshDS, theSinuEdges[ iEdgeComputed ], /*skipMedium=*/true, nodeParams ))
1760             return false;
1761
1762           projectedEdges.insert( iEdgeComputed );
1763
1764           SMESH_MAT2d::BoundaryPoint& bndPnt = bp[ 1-iSideComputed ];
1765           SMESH_MAT2d::BranchPoint brp;
1766           NodePoint npN, npB; // NodePoint's initialized by node and BoundaryPoint
1767           NodePoint& np0 = iSideComputed ? npB : npN;
1768           NodePoint& np1 = iSideComputed ? npN : npB;
1769
1770           double maParam1st, maParamLast, maParam;
1771           if ( !theMA.getBoundary().getBranchPoint( iEdgeComputed, nodeParams.begin()->first, brp ))
1772             return false;
1773           branch.getParameter( brp, maParam1st );
1774           if ( !theMA.getBoundary().getBranchPoint( iEdgeComputed, nodeParams.rbegin()->first, brp ))
1775             return false;
1776           branch.getParameter( brp, maParamLast );
1777
1778           map< double, const SMDS_MeshNode* >::iterator u2n = nodeParams.begin(), u2nEnd = nodeParams.end();
1779           TMAPar2NPoints::iterator end = pointsOnE.end(), pos = end;
1780           TMAPar2NPoints::iterator & hint = (maParamLast > maParam1st) ? end : pos;
1781           for ( ++u2n, --u2nEnd; u2n != u2nEnd; ++u2n )
1782           {
1783             // point on EDGE (u2n) --> MA point (brp)
1784             if ( !theMA.getBoundary().getBranchPoint( iEdgeComputed, u2n->first, brp ))
1785               return false;
1786             // MA point --> points on 2 EDGEs (bp)
1787             if ( !branch.getBoundaryPoints( brp, bp[0], bp[1] ) ||
1788                  !branch.getParameter( brp, maParam ))
1789               return false;
1790
1791             npN = NodePoint( u2n->second, u2n->first, iEdgeComputed );
1792             npB = NodePoint( bndPnt );
1793             pos = pointsOnE.insert( hint, make_pair( maParam, make_pair( np0, np1 )));
1794           }
1795         }
1796         // --------------------------------------------------------------------------------
1797         else if ( !isComputed[ edgeIDs1[ iEdgePair ]] &&         // none of EDGEs is meshed
1798                   !isComputed[ edgeIDs2[ iEdgePair ]])
1799         {
1800           // "projection" from MA
1801           maParams.clear();
1802           if ( !getParamsForEdgePair( iEdgePair, divPoints, theMAParams, maParams ))
1803             return false;
1804
1805           for ( size_t i = 1; i < maParams.size()-1; ++i )
1806           {
1807             if ( !branch.getBoundaryPoints( maParams[i], bp[0], bp[1] ))
1808               return false;
1809
1810             pointsOnE.insert( pointsOnE.end(), make_pair( maParams[i], make_pair( NodePoint(bp[0]),
1811                                                                                   NodePoint(bp[1]))));
1812           }
1813         }
1814         // --------------------------------------------------------------------------------
1815         else if ( isComputed[ edgeIDs1[ iEdgePair ]] &&             // equally meshed EDGES
1816                   isComputed[ edgeIDs2[ iEdgePair ]])
1817         {
1818           // add existing nodes
1819
1820           size_t iE0 = edgeIDs1[ iEdgePair ];
1821           size_t iE1 = edgeIDs2[ iEdgePair ];
1822           map< double, const SMDS_MeshNode* > nodeParams[2]; // params of existing nodes
1823           if ( !SMESH_Algo::GetSortedNodesOnEdge( meshDS, theSinuEdges[ iE0 ],
1824                                                   /*skipMedium=*/false, nodeParams[0] ) ||
1825                !SMESH_Algo::GetSortedNodesOnEdge( meshDS, theSinuEdges[ iE1 ],
1826                                                   /*skipMedium=*/false, nodeParams[1] ) ||
1827                nodeParams[0].size() != nodeParams[1].size() )
1828             return false;
1829
1830           if ( nodeParams[0].size() <= 2 )
1831             continue; // nodes on VERTEXes only
1832
1833           bool reverse = ( theSinuEdges[0].Orientation() == theSinuEdges[1].Orientation() );
1834           double maParam;
1835           SMESH_MAT2d::BranchPoint brp;
1836           std::pair< NodePoint, NodePoint > npPair;
1837
1838           map< double, const SMDS_MeshNode* >::iterator
1839             u2n0F = ++nodeParams[0].begin(),
1840             u2n1F = ++nodeParams[1].begin();
1841           map< double, const SMDS_MeshNode* >::reverse_iterator
1842             u2n1R = ++nodeParams[1].rbegin();
1843           for ( ; u2n0F != nodeParams[0].end(); ++u2n0F )
1844           {
1845             if ( !theMA.getBoundary().getBranchPoint( iE0, u2n0F->first, brp ) ||
1846                  !branch.getParameter( brp, maParam ))
1847               return false;
1848
1849             npPair.first = NodePoint( u2n0F->second, u2n0F->first, iE0 );
1850             if ( reverse )
1851             {
1852               npPair.second = NodePoint( u2n1R->second, u2n1R->first, iE1 );
1853               ++u2n1R;
1854             }
1855             else
1856             {
1857               npPair.second = NodePoint( u2n1F->second, u2n1F->first, iE1 );
1858               ++u2n1F;
1859             }
1860             pointsOnE.insert( make_pair( maParam, npPair ));
1861           }
1862         }
1863       }  // loop on pairs of opposite EDGEs
1864
1865       if ( !projectVertices( theHelper, theMA, divPoints, edgeIDs1, edgeIDs2,
1866                              isComputed, pointsOnE, theSinuFace ))
1867         return false;
1868
1869       separateNodes( theHelper, theMA, pointsOnE, theSinuFace, isComputed );
1870
1871       // create nodes
1872       TMAPar2NPoints::iterator u2np = pointsOnE.begin();
1873       for ( ; u2np != pointsOnE.end(); ++u2np )
1874       {
1875         NodePoint* np[2] = { & u2np->second.first, & u2np->second.second };
1876         for ( int iSide = 0; iSide < 2; ++iSide )
1877         {
1878           if ( np[ iSide ]->_node ) continue;
1879           size_t       iEdge = np[ iSide ]->_edgeInd;
1880           double           u = np[ iSide ]->_u;
1881           gp_Pnt           p = curves[ iEdge ]->Value( u );
1882           np[ iSide ]->_node = meshDS->AddNode( p.X(), p.Y(), p.Z() );
1883           meshDS->SetNodeOnEdge( np[ iSide ]->_node, edgeIDs[ iEdge ], u );
1884         }
1885       }
1886
1887       // create mesh segments on EDGEs
1888       theHelper.SetElementsOnShape( false );
1889       TopoDS_Face face = TopoDS::Face( theHelper.GetSubShape() );
1890       for ( size_t i = 0; i < theSinuEdges.size(); ++i )
1891       {
1892         SMESH_subMesh* sm = mesh->GetSubMesh( theSinuEdges[i] );
1893         if ( sm->GetSubMeshDS() && sm->GetSubMeshDS()->NbElements() > 0 )
1894           continue;
1895
1896         StdMeshers_FaceSide side( face, theSinuEdges[i], mesh,
1897                                   /*isFwd=*/true, /*skipMediumNodes=*/true, &theHelper );
1898         vector<const SMDS_MeshNode*> nodes = side.GetOrderedNodes();
1899         for ( size_t in = 1; in < nodes.size(); ++in )
1900         {
1901           const SMDS_MeshElement* seg = theHelper.AddEdge( nodes[in-1], nodes[in], 0, false );
1902           meshDS->SetMeshElementOnShape( seg, edgeIDs[ i ] );
1903         }
1904       }
1905
1906       // update sub-meshes on VERTEXes
1907       for ( size_t i = 0; i < theSinuEdges.size(); ++i )
1908       {
1909         mesh->GetSubMesh( theHelper.IthVertex( 0, theSinuEdges[i] ))
1910           ->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1911         mesh->GetSubMesh( theHelper.IthVertex( 1, theSinuEdges[i] ))
1912           ->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1913       }
1914     }
1915
1916     // Setup sides of a quadrangle
1917     if ( !setQuadSides( theHelper, pointsOnE, theSinuFace, the1dAlgo ))
1918       return false;
1919
1920     return true;
1921   }
1922
1923   //================================================================================
1924   /*!
1925    * \brief Mesh short EDGEs
1926    */
1927   //================================================================================
1928
1929   bool computeShortEdges( SMESH_MesherHelper&        theHelper,
1930                           const vector<TopoDS_Edge>& theShortEdges,
1931                           SMESH_Algo*                the1dAlgo,
1932                           const bool                 theHasRadialHyp,
1933                           const bool                 theIs2nd)
1934   {
1935     SMESH_Hypothesis::Hypothesis_Status aStatus;
1936     for ( size_t i = 0; i < theShortEdges.size(); ++i )
1937     {
1938       if ( !theHasRadialHyp )
1939         // use global hyps
1940         theHelper.GetGen()->Compute( *theHelper.GetMesh(), theShortEdges[i],
1941                                      SMESH_Gen::SHAPE_ONLY_UPWARD );
1942
1943       SMESH_subMesh* sm = theHelper.GetMesh()->GetSubMesh(theShortEdges[i] );
1944       if ( sm->IsEmpty() )
1945       {
1946         // use 2D hyp or minSegLen
1947         try {
1948           // compute VERTEXes
1949           SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/false);
1950           while ( smIt->more() )
1951             smIt->next()->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1952
1953           // compute EDGE
1954           the1dAlgo->CheckHypothesis( *theHelper.GetMesh(), theShortEdges[i], aStatus );
1955           if ( !the1dAlgo->Compute( *theHelper.GetMesh(), theShortEdges[i] ))
1956             return false;
1957         }
1958         catch (...) {
1959           return false;
1960         }
1961         sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1962         if ( sm->IsEmpty() )
1963           return false;
1964       }
1965     }
1966     return true;
1967   }
1968
1969   inline double area( const UVPtStruct& p1, const UVPtStruct& p2, const UVPtStruct& p3 )
1970   {
1971     gp_XY v1 = p2.UV() - p1.UV();
1972     gp_XY v2 = p3.UV() - p1.UV();
1973     return v2 ^ v1;
1974   }
1975
1976   bool ellipticSmooth( FaceQuadStruct::Ptr quad, int nbLoops )
1977   {
1978     //nbLoops = 10;
1979     if ( quad->uv_grid.empty() )
1980       return true;
1981
1982     int nbhoriz  = quad->iSize;
1983     int nbvertic = quad->jSize;
1984
1985     const double dksi = 0.5, deta = 0.5;
1986     const double  dksi2 = dksi*dksi, deta2 = deta*deta;
1987     double err = 0., g11, g22, g12;
1988     //int nbErr = 0;
1989
1990     FaceQuadStruct& q = *quad;
1991     UVPtStruct pNew;
1992
1993     //double refArea = area( q.UVPt(0,0), q.UVPt(1,0), q.UVPt(1,1) );
1994
1995     for ( int iLoop = 0; iLoop < nbLoops; ++iLoop )
1996     {
1997       err = 0;
1998       for ( int i = 1; i < nbhoriz - 1; i++ )
1999         for ( int j = 1; j < nbvertic - 1; j++ )
2000         {
2001           g11 = ( (q.U(i,j+1) - q.U(i,j-1))*(q.U(i,j+1) - q.U(i,j-1))/dksi2 +
2002                   (q.V(i,j+1) - q.V(i,j-1))*(q.V(i,j+1) - q.V(i,j-1))/deta2 )/4;
2003
2004           g22 = ( (q.U(i+1,j) - q.U(i-1,j))*(q.U(i+1,j) - q.U(i-1,j))/dksi2 +
2005                   (q.V(i+1,j) - q.V(i-1,j))*(q.V(i+1,j) - q.V(i-1,j))/deta2 )/4;
2006
2007           g12 = ( (q.U(i+1,j) - q.U(i-1,j))*(q.U(i,j+1) - q.U(i,j-1))/dksi2 +
2008                   (q.V(i+1,j) - q.V(i-1,j))*(q.V(i,j+1) - q.V(i,j-1))/deta2 )/(4*dksi*deta);
2009
2010           pNew.u = dksi2/(2*(g11+g22)) * (g11*(q.U(i+1,j) + q.U(i-1,j))/dksi2 +
2011                                           g22*(q.U(i,j+1) + q.U(i,j-1))/dksi2
2012                                           - 0.5*g12*q.U(i+1,j+1) + 0.5*g12*q.U(i-1,j+1) +
2013                                           - 0.5*g12*q.U(i-1,j-1) + 0.5*g12*q.U(i+1,j-1));
2014
2015           pNew.v = deta2/(2*(g11+g22)) * (g11*(q.V(i+1,j) + q.V(i-1,j))/deta2 +
2016                                           g22*(q.V(i,j+1) + q.V(i,j-1))/deta2
2017                                           - 0.5*g12*q.V(i+1,j+1) + 0.5*g12*q.V(i-1,j+1) +
2018                                           - 0.5*g12*q.V(i-1,j-1) + 0.5*g12*q.V(i+1,j-1));
2019
2020           // if (( refArea * area( q.UVPt(i-1,j-1), q.UVPt(i,j-1), pNew ) > 0 ) &&
2021           //     ( refArea * area( q.UVPt(i+1,j-1), q.UVPt(i+1,j), pNew ) > 0 ) &&
2022           //     ( refArea * area( q.UVPt(i+1,j+1), q.UVPt(i,j+1), pNew ) > 0 ) &&
2023           //     ( refArea * area( q.UVPt(i-1,j), q.UVPt(i-1,j-1), pNew ) > 0 ))
2024           {
2025             err += sqrt(( q.U(i,j) - pNew.u ) * ( q.U(i,j) - pNew.u ) +
2026                         ( q.V(i,j) - pNew.v ) * ( q.V(i,j) - pNew.v ));
2027             q.U(i,j) = pNew.u;
2028             q.V(i,j) = pNew.v;
2029           }
2030           // else if ( ++nbErr < 10 )
2031           // {
2032           //   cout << i << ", " << j << endl;
2033           //   cout << "x = ["
2034           //        << "[ " << q.U(i-1,j-1) << ", " <<q.U(i,j-1) << ", " << q.U(i+1,j-1) << " ],"
2035           //        << "[ " << q.U(i-1,j-0) << ", " <<q.U(i,j-0) << ", " << q.U(i+1,j-0) << " ],"
2036           //        << "[ " << q.U(i-1,j+1) << ", " <<q.U(i,j+1) << ", " << q.U(i+1,j+1) << " ]]" << endl;
2037           //   cout << "y = ["
2038           //        << "[ " << q.V(i-1,j-1) << ", " <<q.V(i,j-1) << ", " << q.V(i+1,j-1) << " ],"
2039           //        << "[ " << q.V(i-1,j-0) << ", " <<q.V(i,j-0) << ", " << q.V(i+1,j-0) << " ],"
2040           //        << "[ " << q.V(i-1,j+1) << ", " <<q.V(i,j+1) << ", " << q.V(i+1,j+1) << " ]]" << endl<<endl;
2041           // }
2042         }
2043
2044       if ( err / ( nbhoriz - 2 ) / ( nbvertic - 2 ) < 1e-6 )
2045         break;
2046     }
2047     //cout << " ERR " << err / ( nbhoriz - 2 ) / ( nbvertic - 2 ) << endl;
2048
2049     return true;
2050   }
2051
2052   //================================================================================
2053   /*!
2054    * \brief Remove temporary node
2055    */
2056   //================================================================================
2057
2058   void mergeNodes( SMESH_MesherHelper& theHelper,
2059                    SinuousFace&        theSinuFace )
2060   {
2061     SMESH_MeshEditor editor( theHelper.GetMesh() );
2062     SMESH_MeshEditor::TListOfListOfNodes nodesGroups;
2063
2064     TMergeMap::iterator n2nn = theSinuFace._nodesToMerge.begin();
2065     for ( ; n2nn != theSinuFace._nodesToMerge.end(); ++n2nn )
2066     {
2067       if ( !n2nn->first ) continue;
2068       nodesGroups.push_back( list< const SMDS_MeshNode* >() );
2069       list< const SMDS_MeshNode* > & group = nodesGroups.back();
2070
2071       group.push_back( n2nn->first );
2072       group.splice( group.end(), n2nn->second );
2073     }
2074     editor.MergeNodes( nodesGroups );
2075   }
2076
2077 } // namespace
2078
2079 //================================================================================
2080 /*!
2081  * \brief Sets event listener which removes mesh from EDGEs when 2D hyps change
2082  */
2083 //================================================================================
2084
2085 void StdMeshers_QuadFromMedialAxis_1D2D::SetEventListener(SMESH_subMesh* faceSubMesh)
2086 {
2087   faceSubMesh->SetEventListener( new EdgeCleaner, 0, faceSubMesh );
2088 }
2089
2090 //================================================================================
2091 /*!
2092  * \brief Create quadrangle elements
2093  *  \param [in] theHelper - the helper
2094  *  \param [in] theFace - the face to mesh
2095  *  \param [in] theSinuEdges - the sinuous EDGEs
2096  *  \param [in] theShortEdges - the short EDGEs
2097  *  \return bool - is OK or not
2098  */
2099 //================================================================================
2100
2101 bool StdMeshers_QuadFromMedialAxis_1D2D::computeQuads( SMESH_MesherHelper& theHelper,
2102                                                        FaceQuadStruct::Ptr theQuad)
2103 {
2104   StdMeshers_Quadrangle_2D::myHelper     = &theHelper;
2105   StdMeshers_Quadrangle_2D::myNeedSmooth = false;
2106   StdMeshers_Quadrangle_2D::myCheckOri   = false;
2107   StdMeshers_Quadrangle_2D::myQuadList.clear();
2108
2109   int nbNodesShort0 = theQuad->side[0].NbPoints();
2110   int nbNodesShort1 = theQuad->side[2].NbPoints();
2111   int nbNodesSinu0  = theQuad->side[1].NbPoints();
2112   int nbNodesSinu1  = theQuad->side[3].NbPoints();
2113
2114   // compute UV of internal points
2115   myQuadList.push_back( theQuad );
2116   // if ( !StdMeshers_Quadrangle_2D::setNormalizedGrid( theQuad ))
2117   //   return false;
2118
2119   // elliptic smooth of internal points to get boundary cell normal to the boundary
2120   bool isRing = theQuad->side[0].grid->Edge(0).IsNull();
2121   if ( !isRing ) {
2122     if ( !StdMeshers_Quadrangle_2D::setNormalizedGrid( theQuad ))
2123       return false;
2124     ellipticSmooth( theQuad, 1 );
2125   }
2126   // create quadrangles
2127   bool ok;
2128   theHelper.SetElementsOnShape( true );
2129   if ( nbNodesShort0 == nbNodesShort1 && nbNodesSinu0 == nbNodesSinu1 )
2130     ok = StdMeshers_Quadrangle_2D::computeQuadDominant( *theHelper.GetMesh(),
2131                                                         theQuad->face, theQuad );
2132   else
2133     ok = StdMeshers_Quadrangle_2D::computeTriangles( *theHelper.GetMesh(),
2134                                                      theQuad->face, theQuad );
2135
2136   StdMeshers_Quadrangle_2D::myHelper = 0;
2137
2138   return ok;
2139 }
2140
2141 //================================================================================
2142 /*!
2143  * \brief Generate quadrangle mesh
2144  */
2145 //================================================================================
2146
2147 bool StdMeshers_QuadFromMedialAxis_1D2D::Compute(SMESH_Mesh&         theMesh,
2148                                                  const TopoDS_Shape& theShape)
2149 {
2150   SMESH_MesherHelper helper( theMesh );
2151   helper.SetSubShape( theShape );
2152
2153   TopoDS_Face F = TopoDS::Face( theShape );
2154   if ( F.Orientation() >= TopAbs_INTERNAL ) F.Orientation( TopAbs_FORWARD );
2155
2156   SinuousFace sinuFace( F );
2157
2158   _progress = 0.01;
2159
2160   if ( getSinuousEdges( helper, sinuFace ))
2161   {
2162     _progress = 0.4;
2163
2164     double minSegLen = getMinSegLen( helper, sinuFace._sinuEdges );
2165     SMESH_MAT2d::MedialAxis ma( F, sinuFace._sinuEdges, minSegLen, /*ignoreCorners=*/true );
2166
2167     if ( !_regular1D )
2168       _regular1D = new Algo1D( _gen );
2169     _regular1D->SetSegmentLength( minSegLen );
2170
2171     vector<double> maParams;
2172     if ( ! divideMA( helper, ma, sinuFace, _regular1D, minSegLen, maParams ))
2173       return error(COMPERR_BAD_SHAPE);
2174
2175     _progress = 0.8;
2176     if ( _hyp2D )
2177       _regular1D->SetRadialDistribution( _hyp2D );
2178
2179     if ( !computeShortEdges( helper, sinuFace._shortSide[0], _regular1D, _hyp2D, 0 ) ||
2180          !computeShortEdges( helper, sinuFace._shortSide[1], _regular1D, _hyp2D, 1 ))
2181       return error("Failed to mesh short edges");
2182
2183     _progress = 0.85;
2184
2185     if ( !computeSinuEdges( helper, minSegLen, ma, maParams, sinuFace, _regular1D ))
2186       return error("Failed to mesh sinuous edges");
2187
2188     _progress = 0.9;
2189
2190     bool ok = computeQuads( helper, sinuFace._quad );
2191
2192     if ( ok )
2193       mergeNodes( helper, sinuFace );
2194
2195     _progress = 1.;
2196
2197     return ok;
2198   }
2199
2200   return error(COMPERR_BAD_SHAPE, "Not implemented so far");
2201 }
2202
2203 //================================================================================
2204 /*!
2205  * \brief Predict nb of elements
2206  */
2207 //================================================================================
2208
2209 bool StdMeshers_QuadFromMedialAxis_1D2D::Evaluate(SMESH_Mesh &         theMesh,
2210                                                   const TopoDS_Shape & theShape,
2211                                                   MapShapeNbElems&     theResMap)
2212 {
2213   return StdMeshers_Quadrangle_2D::Evaluate(theMesh,theShape,theResMap);
2214 }
2215
2216 //================================================================================
2217 /*!
2218  * \brief Return true if the algorithm can mesh this shape
2219  *  \param [in] aShape - shape to check
2220  *  \param [in] toCheckAll - if true, this check returns OK if all shapes are OK,
2221  *              else, returns OK if at least one shape is OK
2222  */
2223 //================================================================================
2224
2225 bool StdMeshers_QuadFromMedialAxis_1D2D::IsApplicable( const TopoDS_Shape & aShape,
2226                                                        bool                 toCheckAll )
2227 {
2228   TmpMesh tmpMesh;
2229   SMESH_MesherHelper helper( tmpMesh );
2230
2231   int nbFoundFaces = 0;
2232   for (TopExp_Explorer exp( aShape, TopAbs_FACE ); exp.More(); exp.Next(), ++nbFoundFaces )
2233   {
2234     const TopoDS_Face& face = TopoDS::Face( exp.Current() );
2235     SinuousFace sinuFace( face );
2236     bool isApplicable = getSinuousEdges( helper, sinuFace );
2237
2238     if ( toCheckAll  && !isApplicable ) return false;
2239     if ( !toCheckAll &&  isApplicable ) return true;
2240   }
2241   return ( toCheckAll && nbFoundFaces != 0 );
2242 }
2243