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