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