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