Salome HOME
23072: [CEA 1500] Split biquadratic elements into linear elements
[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 "StdMeshers_FaceSide.hxx"
37 #include "StdMeshers_Regular_1D.hxx"
38 #include "StdMeshers_ViscousLayers2D.hxx"
39
40 #include <BRepBuilderAPI_MakeEdge.hxx>
41 #include <BRepTools.hxx>
42 #include <BRep_Tool.hxx>
43 #include <GeomAPI_Interpolate.hxx>
44 #include <Geom_Surface.hxx>
45 #include <Precision.hxx>
46 #include <TColgp_HArray1OfPnt.hxx>
47 #include <TopExp.hxx>
48 #include <TopLoc_Location.hxx>
49 #include <TopTools_MapOfShape.hxx>
50 #include <TopoDS.hxx>
51 #include <TopoDS_Edge.hxx>
52 #include <TopoDS_Face.hxx>
53 #include <TopoDS_Vertex.hxx>
54 #include <gp_Pnt.hxx>
55
56 #include <list>
57 #include <vector>
58
59 //================================================================================
60 /*!
61  * \brief 1D algo
62  */
63 class StdMeshers_QuadFromMedialAxis_1D2D::Algo1D : public StdMeshers_Regular_1D
64 {
65 public:
66   Algo1D(int studyId, SMESH_Gen* gen):
67     StdMeshers_Regular_1D( gen->GetANewId(), studyId, gen )
68   {
69   }
70   void SetSegmentLength( double len )
71   {
72     _value[ BEG_LENGTH_IND ] = len;
73     _value[ PRECISION_IND  ] = 1e-7;
74     _hypType = LOCAL_LENGTH;
75   }
76 };
77  
78 //================================================================================
79 /*!
80  * \brief Constructor sets algo features
81  */
82 //================================================================================
83
84 StdMeshers_QuadFromMedialAxis_1D2D::StdMeshers_QuadFromMedialAxis_1D2D(int        hypId,
85                                                                        int        studyId,
86                                                                        SMESH_Gen* gen)
87   : StdMeshers_Quadrangle_2D(hypId, studyId, gen),
88     _regular1D( 0 )
89 {
90   _name = "QuadFromMedialAxis_1D2D";
91   _shapeType = (1 << TopAbs_FACE);
92   _onlyUnaryInput          = true;  // FACE by FACE so far
93   _requireDiscreteBoundary = false; // make 1D by myself
94   _supportSubmeshes        = true; // make 1D by myself
95   _neededLowerHyps[ 1 ]    = true;  // suppress warning on hiding a global 1D algo
96   _neededLowerHyps[ 2 ]    = true;  // suppress warning on hiding a global 2D algo
97   _compatibleHypothesis.clear();
98   _compatibleHypothesis.push_back("ViscousLayers2D");
99 }
100
101 //================================================================================
102 /*!
103  * \brief Destructor
104  */
105 //================================================================================
106
107 StdMeshers_QuadFromMedialAxis_1D2D::~StdMeshers_QuadFromMedialAxis_1D2D()
108 {
109   delete _regular1D;
110   _regular1D = 0;
111 }
112
113 //================================================================================
114 /*!
115  * \brief Check if needed hypotheses are present
116  */
117 //================================================================================
118
119 bool StdMeshers_QuadFromMedialAxis_1D2D::CheckHypothesis(SMESH_Mesh&         aMesh,
120                                                          const TopoDS_Shape& aShape,
121                                                          Hypothesis_Status&  aStatus)
122 {
123   aStatus = HYP_OK;
124   return true; // does not require hypothesis
125 }
126
127 namespace
128 {
129   typedef map< const SMDS_MeshNode*, list< const SMDS_MeshNode* > > TMergeMap;
130   
131   //================================================================================
132   /*!
133    * \brief Sinuous face
134    */
135   struct SinuousFace
136   {
137     FaceQuadStruct::Ptr          _quad;
138     vector< TopoDS_Edge >        _edges;
139     vector< TopoDS_Edge >        _sinuSide[2], _shortSide[2];
140     vector< TopoDS_Edge >        _sinuEdges;
141     vector< Handle(Geom_Curve) > _sinuCurves;
142     int                          _nbWires;
143     list< int >                  _nbEdgesInWire;
144     TMergeMap                    _nodesToMerge;
145
146     SinuousFace( const TopoDS_Face& f ): _quad( new FaceQuadStruct )
147     {
148       list< TopoDS_Edge > edges;
149       _nbWires = SMESH_Block::GetOrderedEdges (f, edges, _nbEdgesInWire);
150       _edges.assign( edges.begin(), edges.end() );
151
152       _quad->side.resize( 4 );
153       _quad->face = f;
154     }
155     const TopoDS_Face& Face() const { return _quad->face; }
156   };
157
158   //================================================================================
159   /*!
160    * \brief Temporary mesh
161    */
162   struct TmpMesh : public SMESH_Mesh
163   {
164     TmpMesh()
165     {
166       _myMeshDS = new SMESHDS_Mesh(/*id=*/0, /*isEmbeddedMode=*/true);
167     }
168   };
169
170   //================================================================================
171   /*!
172    * \brief Return a member of a std::pair
173    */
174   //================================================================================
175
176   template< typename T >
177   T& get( std::pair< T, T >& thePair, bool is2nd )
178   {
179     return is2nd ? thePair.second : thePair.first;
180   }
181
182   //================================================================================
183   /*!
184    * \brief Select two EDGEs from a map, either mapped to least values or to max values
185    */
186   //================================================================================
187
188   // template< class TVal2EdgesMap >
189   // void getTwo( bool                 least,
190   //              TVal2EdgesMap&       map,
191   //              vector<TopoDS_Edge>& twoEdges,
192   //              vector<TopoDS_Edge>& otherEdges)
193   // {
194   //   twoEdges.clear();
195   //   otherEdges.clear();
196   //   if ( least )
197   //   {
198   //     TVal2EdgesMap::iterator i = map.begin();
199   //     twoEdges.push_back( i->second );
200   //     twoEdges.push_back( ++i->second );
201   //     for ( ; i != map.end(); ++i )
202   //       otherEdges.push_back( i->second );
203   //   }
204   //   else
205   //   {
206   //     TVal2EdgesMap::reverse_iterator i = map.rbegin();
207   //     twoEdges.push_back( i->second );
208   //     twoEdges.push_back( ++i->second );
209   //     for ( ; i != map.rend(); ++i )
210   //       otherEdges.push_back( i->second );
211   //   }
212   //   TopoDS_Vertex v;
213   //   if ( TopExp::CommonVertex( twoEdges[0], twoEdges[1], v ))
214   //   {
215   //     twoEdges.clear(); // two EDGEs must not be connected
216   //     otherEdges.clear();
217   //   }
218   // }
219
220   //================================================================================
221   /*!
222    * \brief Finds out a minimal segment length given EDGEs will be divided into.
223    *        This length is further used to discretize the Medial Axis
224    */
225   //================================================================================
226
227   double getMinSegLen(SMESH_MesherHelper&        theHelper,
228                       const vector<TopoDS_Edge>& theEdges)
229   {
230     TmpMesh tmpMesh;
231     SMESH_Mesh* mesh = theHelper.GetMesh();
232
233     vector< SMESH_Algo* > algos( theEdges.size() );
234     for ( size_t i = 0; i < theEdges.size(); ++i )
235     {
236       SMESH_subMesh* sm = mesh->GetSubMesh( theEdges[i] );
237       algos[i] = sm->GetAlgo();
238     }
239
240     const int nbSegDflt = mesh->GetGen()->GetDefaultNbSegments();
241     double minSegLen    = Precision::Infinite();
242
243     for ( size_t i = 0; i < theEdges.size(); ++i )
244     {
245       SMESH_subMesh* sm = mesh->GetSubMesh( theEdges[i] );
246       if ( SMESH_Algo::IsStraight( theEdges[i], /*degenResult=*/true ))
247         continue;
248       // get algo
249       size_t iOpp = ( theEdges.size() == 4 ? (i+2)%4 : i );
250       SMESH_Algo*  algo = sm->GetAlgo();
251       if ( !algo ) algo = algos[ iOpp ];
252       // get hypo
253       SMESH_Hypothesis::Hypothesis_Status status = SMESH_Hypothesis::HYP_MISSING;
254       if ( algo )
255       {
256         if ( !algo->CheckHypothesis( *mesh, theEdges[i], status ))
257           algo->CheckHypothesis( *mesh, theEdges[iOpp], status );
258       }
259       // compute
260       if ( status != SMESH_Hypothesis::HYP_OK )
261       {
262         minSegLen = Min( minSegLen, SMESH_Algo::EdgeLength( theEdges[i] ) / nbSegDflt );
263       }
264       else
265       {
266         tmpMesh.Clear();
267         tmpMesh.ShapeToMesh( TopoDS_Shape());
268         tmpMesh.ShapeToMesh( theEdges[i] );
269         try {
270           mesh->GetGen()->Compute( tmpMesh, theEdges[i], true, true ); // make nodes on VERTEXes
271           if ( !algo->Compute( tmpMesh, theEdges[i] ))
272             continue;
273         }
274         catch (...) {
275           continue;
276         }
277         SMDS_EdgeIteratorPtr segIt = tmpMesh.GetMeshDS()->edgesIterator();
278         while ( segIt->more() )
279         {
280           const SMDS_MeshElement* seg = segIt->next();
281           double len = SMESH_TNodeXYZ( seg->GetNode(0) ).Distance( seg->GetNode(1) );
282           minSegLen = Min( minSegLen, len );
283         }
284       }
285     }
286     if ( Precision::IsInfinite( minSegLen ))
287       minSegLen = mesh->GetShapeDiagonalSize() / nbSegDflt;
288
289     return minSegLen;
290   }
291
292   //================================================================================
293   /*!
294    * \brief Returns EDGEs located between two VERTEXes at which given MA branches end
295    *  \param [in] br1 - one MA branch
296    *  \param [in] br2 - one more MA branch
297    *  \param [in] allEdges - all EDGEs of a FACE
298    *  \param [out] shortEdges - the found EDGEs
299    *  \return bool - is OK or not
300    */
301   //================================================================================
302
303   bool getConnectedEdges( const SMESH_MAT2d::Branch* br1,
304                           const SMESH_MAT2d::Branch* br2,
305                           const vector<TopoDS_Edge>& allEdges,
306                           vector<TopoDS_Edge>&       shortEdges)
307   {
308     vector< size_t > edgeIDs[4];
309     br1->getGeomEdges( edgeIDs[0], edgeIDs[1] );
310     br2->getGeomEdges( edgeIDs[2], edgeIDs[3] );
311
312     // EDGEs returned by a Branch form a connected chain with a VERTEX where
313     // the Branch ends at the chain middle. One of end EDGEs of the chain is common
314     // with either end EDGE of the chain of the other Branch, or the chains are connected
315     // at a common VERTEX;
316
317     // Get indices of end EDGEs of the branches
318     bool vAtStart1 = ( br1->getEnd(0)->_type == SMESH_MAT2d::BE_ON_VERTEX );
319     bool vAtStart2 = ( br2->getEnd(0)->_type == SMESH_MAT2d::BE_ON_VERTEX );
320     size_t iEnd[4] = {
321       vAtStart1 ? edgeIDs[0].back() : edgeIDs[0][0],
322       vAtStart1 ? edgeIDs[1].back() : edgeIDs[1][0],
323       vAtStart2 ? edgeIDs[2].back() : edgeIDs[2][0],
324       vAtStart2 ? edgeIDs[3].back() : edgeIDs[3][0]
325     };
326
327     set< size_t > connectedIDs;
328     TopoDS_Vertex vCommon;
329     // look for the same EDGEs
330     for ( int i = 0; i < 2; ++i )
331       for ( int j = 2; j < 4; ++j )
332         if ( iEnd[i] == iEnd[j] )
333         {
334           connectedIDs.insert( edgeIDs[i].begin(), edgeIDs[i].end() );
335           connectedIDs.insert( edgeIDs[j].begin(), edgeIDs[j].end() );
336           i = j = 4;
337         }
338     if ( connectedIDs.empty() )
339       // look for connected EDGEs
340       for ( int i = 0; i < 2; ++i )
341         for ( int j = 2; j < 4; ++j )
342           if ( TopExp::CommonVertex( allEdges[ iEnd[i]], allEdges[ iEnd[j]], vCommon ))
343           {
344             connectedIDs.insert( edgeIDs[i].begin(), edgeIDs[i].end() );
345             connectedIDs.insert( edgeIDs[j].begin(), edgeIDs[j].end() );
346             i = j = 4;
347           }
348     if ( connectedIDs.empty() ||                     // nothing
349          allEdges.size() - connectedIDs.size() < 2 ) // too many
350       return false;
351
352     // set shortEdges in the order as in allEdges
353     if ( connectedIDs.count( 0 ) &&
354          connectedIDs.count( allEdges.size()-1 ))
355     {
356       size_t iE = allEdges.size()-1;
357       while ( connectedIDs.count( iE-1 ))
358         --iE;
359       for ( size_t i = 0; i < connectedIDs.size(); ++i )
360       {
361         shortEdges.push_back( allEdges[ iE ]);
362         iE = ( iE + 1 ) % allEdges.size();
363       }
364     }
365     else
366     {
367       set< size_t >::iterator i = connectedIDs.begin();
368       for ( ; i != connectedIDs.end(); ++i )
369         shortEdges.push_back( allEdges[ *i ]);
370     }
371     return true;
372   }
373
374   //================================================================================
375   /*!
376    * \brief Find EDGEs to discretize using projection from MA
377    *  \param [in,out] theSinuFace - the FACE to be meshed
378    *  \return bool - OK or not
379    *
380    * It separates all EDGEs into four sides of a quadrangle connected in the order:
381    * theSinuEdges[0], theShortEdges[0], theSinuEdges[1], theShortEdges[1]
382    */
383   //================================================================================
384
385   bool getSinuousEdges( SMESH_MesherHelper& theHelper,
386                         SinuousFace&        theSinuFace)
387   {
388     vector<TopoDS_Edge> * theSinuEdges  = & theSinuFace._sinuSide [0];
389     vector<TopoDS_Edge> * theShortEdges = & theSinuFace._shortSide[0];
390     theSinuEdges[0].clear();
391     theSinuEdges[1].clear();
392     theShortEdges[0].clear();
393     theShortEdges[1].clear();
394    
395     vector<TopoDS_Edge> & allEdges = theSinuFace._edges;
396     const size_t nbEdges = allEdges.size();
397     if ( nbEdges < 4 && theSinuFace._nbWires == 1 )
398       return false;
399
400     if ( theSinuFace._nbWires == 2 ) // ring
401     {
402       size_t nbOutEdges = theSinuFace._nbEdgesInWire.front();
403       theSinuEdges[0].assign ( allEdges.begin(), allEdges.begin() + nbOutEdges );
404       theSinuEdges[1].assign ( allEdges.begin() + nbOutEdges, allEdges.end() );
405       return true;
406     }
407     if ( theSinuFace._nbWires > 2 )
408       return false;
409
410     // create MedialAxis to find short edges by analyzing MA branches
411     double minSegLen = getMinSegLen( theHelper, allEdges );
412     SMESH_MAT2d::MedialAxis ma( theSinuFace.Face(), allEdges, minSegLen * 3 );
413
414     // in an initial request case, theFace represents a part of a river with almost parallel banks
415     // so there should be two branch points
416     using SMESH_MAT2d::BranchEnd;
417     using SMESH_MAT2d::Branch;
418     const vector< const BranchEnd* >& braPoints = ma.getBranchPoints();
419     if ( braPoints.size() < 2 )
420       return false;
421     TopTools_MapOfShape shortMap;
422     size_t nbBranchPoints = 0;
423     for ( size_t i = 0; i < braPoints.size(); ++i )
424     {
425       vector< const Branch* > vertBranches; // branches with an end on VERTEX
426       for ( size_t ib = 0; ib < braPoints[i]->_branches.size(); ++ib )
427       {
428         const Branch* branch = braPoints[i]->_branches[ ib ];
429         if ( branch->hasEndOfType( SMESH_MAT2d::BE_ON_VERTEX ))
430           vertBranches.push_back( branch );
431       }
432       if ( vertBranches.size() != 2 || braPoints[i]->_branches.size() != 3)
433         continue;
434
435       // get common EDGEs of two branches
436       if ( !getConnectedEdges( vertBranches[0], vertBranches[1],
437                                allEdges, theShortEdges[ nbBranchPoints > 0 ] ))
438         return false;
439
440       for ( size_t iS = 0; iS < theShortEdges[ nbBranchPoints ].size(); ++iS )
441         shortMap.Add( theShortEdges[ nbBranchPoints ][ iS ]);
442
443       ++nbBranchPoints;
444     }
445
446     if ( nbBranchPoints != 2 )
447       return false;
448
449     // add to theSinuEdges all edges that are not theShortEdges
450     vector< vector<TopoDS_Edge> > sinuEdges(1);
451     TopoDS_Vertex vCommon;
452     for ( size_t i = 0; i < allEdges.size(); ++i )
453     {
454       if ( !shortMap.Contains( allEdges[i] ))
455       {
456         if ( !sinuEdges.back().empty() )
457           if ( !TopExp::CommonVertex( sinuEdges.back().back(), allEdges[ i ], vCommon ))
458             sinuEdges.resize( sinuEdges.size() + 1 );
459
460         sinuEdges.back().push_back( allEdges[i] );
461       }
462     }
463     if ( sinuEdges.size() == 3 )
464     {
465       if ( !TopExp::CommonVertex( sinuEdges.back().back(), sinuEdges[0][0], vCommon ))
466         return false;
467       vector<TopoDS_Edge>& last = sinuEdges.back();
468       last.insert( last.end(), sinuEdges[0].begin(), sinuEdges[0].end() );
469       sinuEdges[0].swap( last );
470       sinuEdges.resize( 2 );
471     }
472     if ( sinuEdges.size() != 2 )
473       return false;
474
475     theSinuEdges[0].swap( sinuEdges[0] );
476     theSinuEdges[1].swap( sinuEdges[1] );
477
478     if ( !TopExp::CommonVertex( theSinuEdges[0].back(), theShortEdges[0][0], vCommon ) ||
479          !vCommon.IsSame( theHelper.IthVertex( 1, theSinuEdges[0].back() )))
480       theShortEdges[0].swap( theShortEdges[1] );
481
482     theSinuFace._sinuEdges = theSinuEdges[0];
483     theSinuFace._sinuEdges.insert( theSinuFace._sinuEdges.end(),
484                                    theSinuEdges[1].begin(), theSinuEdges[1].end() );
485
486     return ( theShortEdges[0].size() > 0 && theShortEdges[1].size() > 0 &&
487              theSinuEdges [0].size() > 0 && theSinuEdges [1].size() > 0 );
488
489     // the sinuous EDGEs can be composite and C0 continuous,
490     // therefor we use a complex criterion to find TWO short non-sinuous EDGEs
491     // and the rest EDGEs will be treated as sinuous.
492     // A short edge should have the following features:
493     // a) straight
494     // b) short
495     // c) with convex corners at ends
496     // d) far from the other short EDGE
497
498     // vector< double > isStraightEdge( nbEdges, 0 ); // criterion value
499
500     // // a0) evaluate continuity
501     // const double contiWgt = 0.5; // weight of continuity in the criterion
502     // multimap< int, TopoDS_Edge > continuity;
503     // for ( size_t i = 0; i < nbEdges; ++I )
504     // {
505     //   BRepAdaptor_Curve curve( allEdges[i] );
506     //   GeomAbs_Shape C = GeomAbs_CN;
507     //   try:
508     //     C = curve.Continuity(); // C0, G1, C1, G2, C2, C3, CN
509     //   catch ( Standard_Failure ) {}
510     //   continuity.insert( make_pair( C, allEdges[i] ));
511     //   isStraight[i] += double( C ) / double( CN ) * contiWgt;
512     // }
513
514     // // try to choose by continuity
515     // int mostStraight = (int) continuity.rbegin()->first;
516     // int lessStraight = (int) continuity.begin()->first;
517     // if ( mostStraight != lessStraight )
518     // {
519     //   int nbStraight = continuity.count( mostStraight );
520     //   if ( nbStraight == 2 )
521     //   {
522     //     getTwo( /*least=*/false, continuity, theShortEdges, theSinuEdges );
523     //   }
524     //   else if ( nbStraight == 3 && nbEdges == 4 )
525     //   {
526     //     theSinuEdges.push_back( continuity.begin()->second );
527     //     vector<TopoDS_Edge>::iterator it =
528     //       std::find( allEdges.begin(), allEdges.end(), theSinuEdges[0] );
529     //     int i = std::distance( allEdges.begin(), it );
530     //     theSinuEdges .push_back( allEdges[( i+2 )%4 ]);
531     //     theShortEdges.push_back( allEdges[( i+1 )%4 ]);
532     //     theShortEdges.push_back( allEdges[( i+3 )%4 ]);
533     //   }
534     //   if ( theShortEdges.size() == 2 )
535     //     return true;
536     // }
537
538     // // a) curvature; evaluate aspect ratio
539     // {
540     //   const double curvWgt = 0.5;
541     //   for ( size_t i = 0; i < nbEdges; ++I )
542     //   {
543     //     BRepAdaptor_Curve curve( allEdges[i] );
544     //     double curvature = 1;
545     //     if ( !curve.IsClosed() )
546     //     {
547     //       const double f = curve.FirstParameter(), l = curve.LastParameter();
548     //       gp_Pnt pf = curve.Value( f ), pl = curve.Value( l );
549     //       gp_Lin line( pf, pl.XYZ() - pf.XYZ() );
550     //       double distMax = 0;
551     //       for ( double u = f; u < l; u += (l-f)/30. )
552     //         distMax = Max( distMax, line.SquareDistance( curve.Value( u )));
553     //       curvature = Sqrt( distMax ) / ( pf.Distance( pl ));
554     //     }
555     //     isStraight[i] += curvWgt / (              curvature + 1e-20 );
556     //   }
557     // }
558     // // b) length
559     // {
560     //   const double lenWgt = 0.5;
561     //   for ( size_t i = 0; i < nbEdges; ++I )
562     //   {
563     //     double length = SMESH_Algo::Length( allEdges[i] );
564     //     if ( length > 0 )
565     //       isStraight[i] += lenWgt / length;
566     //   }
567     // }
568     // // c) with convex corners at ends
569     // {
570     //   const double cornerWgt = 0.25;
571     //   for ( size_t i = 0; i < nbEdges; ++I )
572     //   {
573     //     double convex = 0;
574     //     int iPrev = SMESH_MesherHelper::WrapIndex( int(i)-1, nbEdges );
575     //     int iNext = SMESH_MesherHelper::WrapIndex( int(i)+1, nbEdges );
576     //     TopoDS_Vertex v = helper.IthVertex( 0, allEdges[i] );
577     //     double angle = SMESH_MesherHelper::GetAngle( allEdges[iPrev], allEdges[i], theFace, v );
578     //     if ( angle < M_PI ) // [-PI; PI]
579     //       convex += ( angle + M_PI ) / M_PI / M_PI;
580     //     v = helper.IthVertex( 1, allEdges[i] );
581     //     angle = SMESH_MesherHelper::GetAngle( allEdges[iNext], allEdges[i], theFace, v );
582     //     if ( angle < M_PI ) // [-PI; PI]
583     //       convex += ( angle + M_PI ) / M_PI / M_PI;
584     //     isStraight[i] += cornerWgt * convex;
585     //   }
586     // }
587   }
588
589   //================================================================================
590   /*!
591    * \brief Creates an EDGE from a sole branch of MA
592    */
593   //================================================================================
594
595   TopoDS_Edge makeEdgeFromMA( SMESH_MesherHelper&            theHelper,
596                               const SMESH_MAT2d::MedialAxis& theMA,
597                               const double                   theMinSegLen)
598   {
599     if ( theMA.nbBranches() != 1 )
600       return TopoDS_Edge();
601
602     vector< gp_XY > uv;
603     theMA.getPoints( theMA.getBranch(0), uv );
604     if ( uv.size() < 2 )
605       return TopoDS_Edge();
606
607     TopoDS_Face face = TopoDS::Face( theHelper.GetSubShape() );
608     Handle(Geom_Surface) surface = BRep_Tool::Surface( face );
609
610     vector< gp_Pnt > pnt;
611     pnt.reserve( uv.size() * 2 );
612     pnt.push_back( surface->Value( uv[0].X(), uv[0].Y() ));
613     for ( size_t i = 1; i < uv.size(); ++i )
614     {
615       gp_Pnt p = surface->Value( uv[i].X(), uv[i].Y() );
616       int nbDiv = int( p.Distance( pnt.back() ) / theMinSegLen );
617       for ( int iD = 1; iD < nbDiv; ++iD )
618       {
619         double  R = iD / double( nbDiv );
620         gp_XY uvR = uv[i-1] * (1 - R) + uv[i] * R;
621         pnt.push_back( surface->Value( uvR.X(), uvR.Y() ));
622       }
623       pnt.push_back( p );
624     }
625
626     // cout << "from salome.geom import geomBuilder" << endl;
627     // cout << "geompy = geomBuilder.New(salome.myStudy)" << endl;
628     Handle(TColgp_HArray1OfPnt) points = new TColgp_HArray1OfPnt(1, pnt.size());
629     for ( size_t i = 0; i < pnt.size(); ++i )
630     {
631       gp_Pnt& p = pnt[i];
632       points->SetValue( i+1, p );
633       // cout << "geompy.MakeVertex( "<< p.X()<<", " << p.Y()<<", " << p.Z()
634       //      <<" theName = 'p_" << i << "')" << endl;
635     }
636
637     GeomAPI_Interpolate interpol( points, /*isClosed=*/false, gp::Resolution());
638     interpol.Perform();
639     if ( !interpol.IsDone())
640       return TopoDS_Edge();
641
642     TopoDS_Edge branchEdge = BRepBuilderAPI_MakeEdge(interpol.Curve());
643     return branchEdge;
644   }
645
646   //================================================================================
647   /*!
648    * \brief Returns a type of shape, to which a hypothesis used to mesh a given edge is assigned
649    */
650   //================================================================================
651
652   TopAbs_ShapeEnum getHypShape( SMESH_Mesh* mesh, const TopoDS_Shape& edge )
653   {
654     TopAbs_ShapeEnum shapeType = TopAbs_SHAPE;
655
656     SMESH_subMesh* sm = mesh->GetSubMesh( edge );
657     SMESH_Algo*  algo = sm->GetAlgo();
658     if ( !algo ) return shapeType;
659
660     const list <const SMESHDS_Hypothesis *> & hyps =
661       algo->GetUsedHypothesis( *mesh, edge, /*ignoreAuxiliary=*/true );
662     if ( hyps.empty() ) return shapeType;
663
664     TopoDS_Shape shapeOfHyp =
665       SMESH_MesherHelper::GetShapeOfHypothesis( hyps.front(), edge, mesh);
666
667     return SMESH_MesherHelper::GetGroupType( shapeOfHyp, /*woCompound=*/true);
668   }
669
670   //================================================================================
671   /*!
672    * \brief Discretize a sole branch of MA an returns parameters of divisions on MA
673    */
674   //================================================================================
675
676   bool divideMA( SMESH_MesherHelper&            theHelper,
677                  const SMESH_MAT2d::MedialAxis& theMA,
678                  const SinuousFace&             theSinuFace,
679                  SMESH_Algo*                    the1dAlgo,
680                  const double                   theMinSegLen,
681                  vector<double>&                theMAParams )
682   {
683     // check if all EDGEs of one size are meshed, then MA discretization is not needed
684     SMESH_Mesh* mesh = theHelper.GetMesh();
685     size_t nbComputedEdges[2] = { 0, 0 };
686     for ( size_t iS = 0; iS < 2; ++iS )
687       for ( size_t i = 0; i < theSinuFace._sinuSide[iS].size(); ++i )
688       {
689         bool isComputed = ( ! mesh->GetSubMesh( theSinuFace._sinuSide[iS][i] )->IsEmpty() );
690         nbComputedEdges[ iS ] += isComputed;
691       }
692     if ( nbComputedEdges[0] == theSinuFace._sinuSide[0].size() ||
693          nbComputedEdges[1] == theSinuFace._sinuSide[1].size() )
694       return true; // discretization is not needed
695
696
697     TopoDS_Edge branchEdge = makeEdgeFromMA( theHelper, theMA, theMinSegLen );
698     if ( branchEdge.IsNull() )
699       return false;
700
701     // const char* file = "/misc/dn25/salome/eap/salome/misc/tmp/MAedge.brep";
702     // BRepTools::Write( branchEdge, file);
703     // cout << "Write " << file << endl;
704
705     // look for a most local hyps assigned to theSinuEdges
706     TopoDS_Edge edge = theSinuFace._sinuEdges[0];
707     int mostSimpleShape = (int) getHypShape( mesh, edge );
708     for ( size_t i = 1; i < theSinuFace._sinuEdges.size(); ++i )
709     {
710       int shapeType = (int) getHypShape( mesh, theSinuFace._sinuEdges[i] );
711       if ( shapeType > mostSimpleShape )
712         edge = theSinuFace._sinuEdges[i];
713     }
714
715     SMESH_Algo* algo = the1dAlgo;
716     if ( mostSimpleShape != TopAbs_SHAPE )
717     {
718       algo = mesh->GetSubMesh( edge )->GetAlgo();
719       SMESH_Hypothesis::Hypothesis_Status status;
720       if ( !algo->CheckHypothesis( *mesh, edge, status ))
721         algo = the1dAlgo;
722     }
723
724     TmpMesh tmpMesh;
725     tmpMesh.ShapeToMesh( branchEdge );
726     try {
727       mesh->GetGen()->Compute( tmpMesh, branchEdge, true, true ); // make nodes on VERTEXes
728       if ( !algo->Compute( tmpMesh, branchEdge ))
729         return false;
730     }
731     catch (...) {
732       return false;
733     }
734     return SMESH_Algo::GetNodeParamOnEdge( tmpMesh.GetMeshDS(), branchEdge, theMAParams );
735   }
736
737   //================================================================================
738   /*!
739    * \brief Select division parameters on MA and make them coincide at ends with
740    *        projections of VERTEXes to MA for a given pair of opposite EDGEs
741    *  \param [in] theEdgePairInd - index of the EDGE pair
742    *  \param [in] theDivPoints - the BranchPoint's dividing MA into parts each
743    *         corresponding to a unique pair of opposite EDGEs
744    *  \param [in] theMAParams - the MA division parameters
745    *  \param [out] theSelectedMAParams - the selected MA parameters
746    *  \return bool - is OK
747    */
748   //================================================================================
749
750   bool getParamsForEdgePair( const size_t                              theEdgePairInd,
751                              const vector< SMESH_MAT2d::BranchPoint >& theDivPoints,
752                              const vector<double>&                     theMAParams,
753                              vector<double>&                           theSelectedMAParams)
754   {
755     if ( theDivPoints.empty() )
756     {
757       theSelectedMAParams = theMAParams;
758       return true;
759     }
760     if ( theEdgePairInd > theDivPoints.size() || theMAParams.empty() )
761       return false;
762
763     // find a range of params to copy
764
765     double par1 = 0;
766     size_t iPar1 = 0;
767     if ( theEdgePairInd > 0 )
768     {
769       const SMESH_MAT2d::BranchPoint& bp = theDivPoints[ theEdgePairInd-1 ];
770       bp._branch->getParameter( bp, par1 );
771       while ( theMAParams[ iPar1 ] < par1 ) ++iPar1;
772       if ( par1 - theMAParams[ iPar1-1 ] < theMAParams[ iPar1 ] - par1 )
773         --iPar1;
774     }
775
776     double par2 = 1;
777     size_t iPar2 = theMAParams.size() - 1;
778     if ( theEdgePairInd < theDivPoints.size() )
779     {
780       const SMESH_MAT2d::BranchPoint& bp = theDivPoints[ theEdgePairInd ];
781       bp._branch->getParameter( bp, par2 );
782       iPar2 = iPar1;
783       while ( theMAParams[ iPar2 ] < par2 ) ++iPar2;
784       if ( par2 - theMAParams[ iPar2-1 ] < theMAParams[ iPar2 ] - par2 )
785         --iPar2;
786     }
787
788     theSelectedMAParams.assign( theMAParams.begin() + iPar1,
789                                 theMAParams.begin() + iPar2 + 1 );
790
791     // adjust theSelectedMAParams to fit between par1 and par2
792
793     double d = par1 - theSelectedMAParams[0];
794     double f = ( par2 - par1 ) / ( theSelectedMAParams.back() - theSelectedMAParams[0] );
795
796     for ( size_t i = 0; i < theSelectedMAParams.size(); ++i )
797     {
798       theSelectedMAParams[i] += d;
799       theSelectedMAParams[i] = par1 + ( theSelectedMAParams[i] - par1 ) * f;
800     }
801
802     return true;
803   }
804
805   //--------------------------------------------------------------------------------
806   // node or node parameter on EDGE
807   struct NodePoint
808   {
809     const SMDS_MeshNode* _node;
810     double               _u;
811     int                  _edgeInd; // index in theSinuEdges vector
812
813     NodePoint(): _node(0), _u(0), _edgeInd(-1) {}
814     NodePoint(const SMDS_MeshNode* n, double u, size_t iEdge ): _node(n), _u(u), _edgeInd(iEdge) {}
815     NodePoint(double u, size_t iEdge) : _node(0), _u(u), _edgeInd(iEdge) {}
816     NodePoint(const SMESH_MAT2d::BoundaryPoint& p) : _node(0), _u(p._param), _edgeInd(p._edgeIndex) {}
817     gp_Pnt Point(const vector< Handle(Geom_Curve) >& curves) const
818     {
819       return curves[ _edgeInd ]->Value( _u );
820     }
821   };
822
823   //================================================================================
824   /*!
825    * \brief Finds a VERTEX corresponding to a point on EDGE, which is also filled
826    *        with a node on the VERTEX, present or created
827    *  \param [in,out] theNodePnt - the node position on the EDGE
828    *  \param [in] theSinuEdges - the sinuous EDGEs
829    *  \param [in] theMeshDS - the mesh
830    *  \return bool - true if the \a theBndPnt is on VERTEX
831    */
832   //================================================================================
833
834   bool findVertexAndNode( NodePoint&                  theNodePnt,
835                           const vector<TopoDS_Edge>&  theSinuEdges,
836                           SMESHDS_Mesh*               theMeshDS = 0,
837                           size_t                      theEdgeIndPrev = 0,
838                           size_t                      theEdgeIndNext = 0)
839   {
840     if ( theNodePnt._edgeInd >= theSinuEdges.size() )
841       return false;
842
843     double f,l;
844     BRep_Tool::Range( theSinuEdges[ theNodePnt._edgeInd ], f,l );
845     const double tol = 1e-3 * ( l - f );
846
847     TopoDS_Vertex V;
848     if      ( Abs( f - theNodePnt._u ) < tol )
849       V = SMESH_MesherHelper::IthVertex( 0, theSinuEdges[ theNodePnt._edgeInd ], /*CumOri=*/false);
850     else if ( Abs( l - theNodePnt._u ) < tol )
851       V = SMESH_MesherHelper::IthVertex( 1, theSinuEdges[ theNodePnt._edgeInd ], /*CumOri=*/false);
852     else if ( theEdgeIndPrev != theEdgeIndNext )
853       TopExp::CommonVertex( theSinuEdges[theEdgeIndPrev], theSinuEdges[theEdgeIndNext], V );
854
855     if ( !V.IsNull() && theMeshDS )
856     {
857       theNodePnt._node = SMESH_Algo::VertexNode( V, theMeshDS );
858       if ( !theNodePnt._node )
859       {
860         gp_Pnt p = BRep_Tool::Pnt( V );
861         theNodePnt._node = theMeshDS->AddNode( p.X(), p.Y(), p.Z() );
862         theMeshDS->SetNodeOnVertex( theNodePnt._node, V );
863       }
864     }
865     return !V.IsNull();
866   }
867
868   //================================================================================
869   /*!
870    * \brief Add to the map of NodePoint's those on VERTEXes
871    *  \param [in,out] theHelper - the helper
872    *  \param [in] theMA - Medial Axis
873    *  \param [in] theMinSegLen - minimal segment length
874    *  \param [in] theDivPoints - projections of VERTEXes to MA
875    *  \param [in] theSinuEdges - the sinuous EDGEs
876    *  \param [in] theSideEdgeIDs - indices of sinuous EDGEs per side
877    *  \param [in] theIsEdgeComputed - is sinuous EGDE is meshed
878    *  \param [in,out] thePointsOnE - the map to fill
879    *  \param [out] theNodes2Merge - the map of nodes to merge
880    */
881   //================================================================================
882
883   bool projectVertices( SMESH_MesherHelper&                           theHelper,
884                         const SMESH_MAT2d::MedialAxis&                theMA,
885                         const vector< SMESH_MAT2d::BranchPoint >&     theDivPoints,
886                         const vector< std::size_t > &                 theEdgeIDs1,
887                         const vector< std::size_t > &                 theEdgeIDs2,
888                         const vector< bool >&                         theIsEdgeComputed,
889                         map< double, pair< NodePoint, NodePoint > > & thePointsOnE,
890                         SinuousFace&                                  theSinuFace)
891   {
892     SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
893     const vector<TopoDS_Edge>&       theSinuEdges = theSinuFace._sinuEdges;
894     const vector< Handle(Geom_Curve) >& theCurves = theSinuFace._sinuCurves;
895
896     double uMA;
897     SMESH_MAT2d::BoundaryPoint bp[2];
898     const SMESH_MAT2d::Branch& branch = *theMA.getBranch(0);
899
900     // add to thePointsOnE NodePoint's of ends of theSinuEdges
901     if ( !branch.getBoundaryPoints( 0., bp[0], bp[1] ) ||
902          !theMA.getBoundary().moveToClosestEdgeEnd( bp[0] ) ||
903          !theMA.getBoundary().moveToClosestEdgeEnd( bp[1] )) return false;
904     NodePoint np0( bp[0]), np1( bp[1] );
905     findVertexAndNode( np0, theSinuEdges, meshDS );
906     findVertexAndNode( np1, theSinuEdges, meshDS );
907     thePointsOnE.insert( make_pair( -0.1, make_pair( np0, np1 )));
908
909     if ( !branch.getBoundaryPoints( 1., bp[0], bp[1] ) ||
910          !theMA.getBoundary().moveToClosestEdgeEnd( bp[0] ) ||
911          !theMA.getBoundary().moveToClosestEdgeEnd( bp[1] )) return false;
912     np0 = bp[0]; np1 = bp[1];
913     findVertexAndNode( np0, theSinuEdges, meshDS );
914     findVertexAndNode( np1, theSinuEdges, meshDS );
915     thePointsOnE.insert( make_pair( 1.1, make_pair( np0, np1)));
916
917     // project theDivPoints
918
919     if ( theDivPoints.empty() )
920       return true;
921
922     for ( size_t i = 0; i < theDivPoints.size(); ++i )
923     {
924       if ( !branch.getParameter( theDivPoints[i], uMA ))
925         return false;
926       if ( !branch.getBoundaryPoints( theDivPoints[i], bp[0], bp[1] ))
927         return false;
928
929       NodePoint  np[2] = {
930         NodePoint( bp[0] ),
931         NodePoint( bp[1] )
932       };
933       bool isVertex[2] = {
934         findVertexAndNode( np[0], theSinuEdges, meshDS, theEdgeIDs1[i], theEdgeIDs1[i+1] ),
935         findVertexAndNode( np[1], theSinuEdges, meshDS, theEdgeIDs2[i], theEdgeIDs2[i+1] )
936       };
937
938       map< double, pair< NodePoint, NodePoint > >::iterator u2NP =
939         thePointsOnE.insert( make_pair( uMA, make_pair( np[0], np[1]))).first;
940
941       if ( !isVertex[0] && !isVertex[1] ) return false; // error
942       if ( isVertex[0] && isVertex[1] )
943         continue;
944       const size_t iVert = isVertex[0] ? 0 : 1;
945       const size_t iNode   = 1 - iVert;
946
947       bool isOppComputed = theIsEdgeComputed[ np[ iNode ]._edgeInd ];
948       if ( !isOppComputed )
949         continue;
950
951       // a VERTEX is projected on a meshed EDGE; there are two options:
952       // 1) a projected point is joined with a closet node if a strip between this and neighbor
953       // projection is WIDE enough; joining is done by creating a node coincident with the
954       //  existing node which will be merged together after all;
955       // 2) a neighbor projection is merged with this one if it is TOO CLOSE; a node of deleted
956       // projection is set to the BoundaryPoint of this projection
957
958       // evaluate distance to neighbor projections
959       const double rShort = 0.2;
960       bool isShortPrev[2], isShortNext[2];
961       map< double, pair< NodePoint, NodePoint > >::iterator u2NPPrev = u2NP, u2NPNext = u2NP;
962       --u2NPPrev; ++u2NPNext;
963       // bool hasPrev = ( u2NP     != thePointsOnE.begin() );
964       // bool hasNext = ( u2NPNext != thePointsOnE.end() );
965       // if ( !hasPrev ) u2NPPrev = u2NP0;
966       // if ( !hasNext ) u2NPNext = u2NP1;
967       for ( int iS = 0; iS < 2; ++iS ) // side with Vertex and side with Nodes
968       {
969         NodePoint np     = get( u2NP->second,     iS );
970         NodePoint npPrev = get( u2NPPrev->second, iS );
971         NodePoint npNext = get( u2NPNext->second, iS );
972         gp_Pnt     p     = np    .Point( theCurves );
973         gp_Pnt     pPrev = npPrev.Point( theCurves );
974         gp_Pnt     pNext = npNext.Point( theCurves );
975         double  distPrev = p.Distance( pPrev );
976         double  distNext = p.Distance( pNext );
977         double         r = distPrev / ( distPrev + distNext );
978         isShortPrev[iS] = ( r < rShort );
979         isShortNext[iS] = (( 1 - r ) > ( 1 - rShort ));
980       }
981       // if ( !hasPrev ) isShortPrev[0] = isShortPrev[1] = false;
982       // if ( !hasNext ) isShortNext[0] = isShortNext[1] = false;
983
984       map< double, pair< NodePoint, NodePoint > >::iterator u2NPClose;
985
986       if (( isShortPrev[0] && isShortPrev[1] ) || // option 2) -> remove a too close projection
987           ( isShortNext[0] && isShortNext[1] ))
988       {
989         u2NPClose = isShortPrev[0] ? u2NPPrev : u2NPNext;
990         NodePoint& npProj  = get( u2NP->second,      iNode ); // NP of VERTEX projection
991         NodePoint npCloseN = get( u2NPClose->second, iNode ); // NP close to npProj
992         NodePoint npCloseV = get( u2NPClose->second, iVert ); // NP close to VERTEX
993         if ( !npCloseV._node )
994         {
995           npProj = npCloseN;
996           thePointsOnE.erase( isShortPrev[0] ? u2NPPrev : u2NPNext );
997           continue;
998         }
999         else
1000         {
1001           // can't remove the neighbor projection as it is also from VERTEX, -> option 1)
1002         }
1003       }
1004       // else: option 1) - wide enough -> "duplicate" existing node
1005       {
1006         u2NPClose = isShortPrev[ iNode ] ? u2NPPrev : u2NPNext;
1007         NodePoint& npProj   = get( u2NP->second,      iNode ); // NP of VERTEX projection
1008         NodePoint& npCloseN = get( u2NPClose->second, iNode ); // NP close to npProj
1009         npProj = npCloseN;
1010         npProj._node = 0;
1011         //npProj._edgeInd = npCloseN._edgeInd;
1012         // npProj._u       = npCloseN._u + 1e-3 * Abs( get( u2NPPrev->second, iNode )._u -
1013         //                                             get( u2NPNext->second, iNode )._u );
1014         // gp_Pnt        p = npProj.Point( theCurves );
1015         // npProj._node    = meshDS->AddNode( p.X(), p.Y(), p.Z() );
1016         // meshDS->SetNodeOnEdge( npProj._node, theSinuEdges[ npProj._edgeInd ], npProj._u  );
1017
1018         //theNodes2Merge[ npCloseN._node ].push_back( npProj._node );
1019       }
1020     }
1021     return true;
1022   }
1023
1024   //================================================================================
1025   /*!
1026    * \brief Move coincident nodes to make node params on EDGE unique
1027    *  \param [in] theHelper - the helper
1028    *  \param [in] thePointsOnE - nodes on two opposite river sides
1029    *  \param [in] theSinuFace - the sinuous FACE
1030    *  \param [out] theNodes2Merge - the map of nodes to merge
1031    */
1032   //================================================================================
1033
1034   void separateNodes( SMESH_MesherHelper&                           theHelper,
1035                       map< double, pair< NodePoint, NodePoint > > & thePointsOnE,
1036                       SinuousFace&                                  theSinuFace )
1037   {
1038     if ( thePointsOnE.size() < 2 )
1039       return;
1040
1041     SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
1042     const vector<TopoDS_Edge>& theSinuEdges = theSinuFace._sinuEdges;
1043
1044     typedef map< double, pair< NodePoint, NodePoint > >::iterator TIterator;
1045
1046     for ( int iSide = 0; iSide < 2; ++iSide )
1047     {
1048       TIterator u2NP0, u2NP1, u2NP = thePointsOnE.begin();
1049       while ( u2NP != thePointsOnE.end() )
1050       {
1051         while ( u2NP != thePointsOnE.end() &&
1052                 get( u2NP->second, iSide )._node )
1053           ++u2NP; // skip NP with an existing node (VERTEXes must be meshed)
1054         if ( u2NP == thePointsOnE.end() )
1055           break;
1056
1057         // find a range of not meshed NP on one EDGE
1058         u2NP0 = u2NP;
1059         if ( !findVertexAndNode( get( u2NP0->second, iSide ), theSinuEdges ))
1060           --u2NP0;
1061         int iCurEdge = get( u2NP->second, iSide )._edgeInd;
1062         int nbNP = 1;
1063         while ( get( u2NP->second, iSide )._edgeInd == iCurEdge &&
1064                 get( u2NP->second, iSide )._node == 0 )
1065           ++u2NP, ++nbNP;
1066         u2NP1 = u2NP; // end of not meshed NP on iCurEdge
1067
1068         // fix parameters of extremity NP of the range
1069         NodePoint* np0 = & get( u2NP0->second, iSide );
1070         NodePoint* np1 = & get( u2NP1->second, iSide );
1071         const TopoDS_Edge& edge = TopoDS::Edge( theSinuFace._sinuEdges[ iCurEdge ]);
1072         if ( np0->_node && np0->_edgeInd != iCurEdge )
1073         {
1074           np0->_u       = theHelper.GetNodeU( edge, np0->_node );
1075           np0->_edgeInd = iCurEdge;
1076         }
1077         if ( np1->_node && np1->_edgeInd != iCurEdge )
1078         {
1079           np1->_u       = theHelper.GetNodeU( edge, np1->_node );
1080           np1->_edgeInd = iCurEdge;
1081         }
1082
1083         // find coincident NPs
1084         double f,l;
1085         BRep_Tool::Range( edge, f,l );
1086         double tol = 1e-2* (l-f) / nbNP;
1087         TIterator u2NPEq = thePointsOnE.end();
1088         u2NP = u2NP0;
1089         for ( ++u2NP; u2NP0 != u2NP1; ++u2NP, ++u2NP0 )
1090         {
1091           np0 = & get( u2NP0->second, iSide );
1092           np1 = & get( u2NP->second,  iSide );
1093           bool coincides = ( Abs( np0->_u - np1->_u ) < tol );
1094           if ( coincides && u2NPEq == thePointsOnE.end() )
1095             u2NPEq = u2NP0;
1096
1097           if (( u2NPEq != thePointsOnE.end() ) &&
1098               ( u2NP == u2NP1 || !coincides ))
1099           {
1100             if ( !get( u2NPEq->second, iSide )._node )
1101               --u2NPEq;
1102             if ( coincides && !get( u2NP->second, iSide )._node && u2NP0 != u2NP1 )
1103               ++u2NP;
1104
1105             // distribute nodes between u2NPEq and u2NP
1106             size_t nbSeg = std::distance( u2NPEq, u2NP );
1107             double    du = 1. / nbSeg * ( get( u2NP->second,   iSide )._u -
1108                                           get( u2NPEq->second, iSide )._u );
1109             double     u = get( u2NPEq->second, iSide )._u + du;
1110
1111             const SMDS_MeshNode* closeNode =
1112               get(( coincides ? u2NP : u2NPEq )->second, iSide )._node;
1113             list< const SMDS_MeshNode* >& eqNodes = theSinuFace._nodesToMerge[ closeNode ];
1114
1115             for ( ++u2NPEq; u2NPEq != u2NP; ++u2NPEq, u += du )
1116             {
1117               np0       = & get( u2NPEq->second, iSide );
1118               np0->_u   = u;
1119               gp_Pnt p  = np0->Point( theSinuFace._sinuCurves );
1120               np0->_node = meshDS->AddNode( p.X(), p.Y(), p.Z() );
1121               meshDS->SetNodeOnEdge( np0->_node, theSinuEdges[ np0->_edgeInd ], np0->_u  );
1122               if ( !closeNode )
1123                 eqNodes = theSinuFace._nodesToMerge[ closeNode = np0->_node ];
1124               else
1125                 eqNodes.push_back( np0->_node );
1126             }
1127           }
1128         }
1129         u2NP = u2NP1;
1130         while ( get( u2NP->second, iSide )._edgeInd != iCurEdge )
1131           --u2NP;
1132         u2NP++;
1133       }
1134     }
1135   }
1136
1137   //================================================================================
1138   /*!
1139    * \brief Divide the sinuous EDGEs by projecting the division point of Medial
1140    *        Axis to the EGDEs
1141    *  \param [in] theHelper - the helper
1142    *  \param [in] theMinSegLen - minimal segment length
1143    *  \param [in] theMA - the Medial Axis
1144    *  \param [in] theMAParams - parameters of division points of \a theMA
1145    *  \param [in] theSinuEdges - the EDGEs to make nodes on
1146    *  \param [in] theSinuSide0Size - the number of EDGEs in the 1st sinuous side
1147    *  \return bool - is OK or not
1148    */
1149   //================================================================================
1150
1151   bool computeSinuEdges( SMESH_MesherHelper&        theHelper,
1152                          double                     /*theMinSegLen*/,
1153                          SMESH_MAT2d::MedialAxis&   theMA,
1154                          vector<double>&            theMAParams,
1155                          SinuousFace&               theSinuFace)
1156   {
1157     if ( theMA.nbBranches() != 1 )
1158       return false;
1159
1160     // normalize theMAParams
1161     for ( size_t i = 0; i < theMAParams.size(); ++i )
1162       theMAParams[i] /= theMAParams.back();
1163
1164
1165     SMESH_Mesh*     mesh = theHelper.GetMesh();
1166     SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
1167     double f,l;
1168
1169     // get data of sinuous EDGEs and remove unnecessary nodes
1170     const vector< TopoDS_Edge >& theSinuEdges = theSinuFace._sinuEdges;
1171     vector< Handle(Geom_Curve) >& curves      = theSinuFace._sinuCurves;
1172     vector< int >                edgeIDs( theSinuEdges.size() );
1173     vector< bool >            isComputed( theSinuEdges.size() );
1174     curves.resize( theSinuEdges.size(), 0 );
1175     for ( size_t i = 0; i < theSinuEdges.size(); ++i )
1176     {
1177       curves[i] = BRep_Tool::Curve( theSinuEdges[i], f,l );
1178       if ( !curves[i] )
1179         return false;
1180       SMESH_subMesh* sm = mesh->GetSubMesh( theSinuEdges[i] );
1181       edgeIDs   [i] = sm->GetId();
1182       isComputed[i] = ( !sm->IsEmpty() );
1183       if ( isComputed[i] )
1184       {
1185         TopAbs_ShapeEnum shape = getHypShape( mesh, theSinuEdges[i] );
1186         if ( shape == TopAbs_SHAPE || shape <= TopAbs_FACE )
1187         {
1188           // EDGE computed using global hypothesis -> clear it
1189           bool hasComputedFace = false;
1190           PShapeIteratorPtr faceIt = theHelper.GetAncestors( theSinuEdges[i], *mesh, TopAbs_FACE );
1191           while ( const TopoDS_Shape* face = faceIt->next() )
1192             if (( !face->IsSame( theSinuFace.Face())) &&
1193                 ( hasComputedFace = !mesh->GetSubMesh( *face )->IsEmpty() ))
1194               break;
1195           if ( !hasComputedFace )
1196             sm->ComputeStateEngine( SMESH_subMesh::CLEAN );
1197           isComputed[i] = false;
1198         }
1199       }
1200     }
1201
1202     const SMESH_MAT2d::Branch& branch = *theMA.getBranch(0);
1203     SMESH_MAT2d::BoundaryPoint bp[2];
1204
1205     vector< std::size_t > edgeIDs1, edgeIDs2;
1206     vector< SMESH_MAT2d::BranchPoint > divPoints;
1207     branch.getOppositeGeomEdges( edgeIDs1, edgeIDs2, divPoints );
1208     for ( size_t i = 0; i < edgeIDs1.size(); ++i )
1209       if ( isComputed[ edgeIDs1[i]] &&
1210            isComputed[ edgeIDs2[i]])
1211         return false;
1212
1213     // map param on MA to parameters of nodes on a pair of theSinuEdges
1214     typedef map< double, pair< NodePoint, NodePoint > > TMAPar2NPoints;
1215     TMAPar2NPoints pointsOnE;
1216     vector<double> maParams;
1217
1218     // compute params of nodes on EDGEs by projecting division points from MA
1219     //const double tol = 1e-5 * theMAParams.back();
1220     size_t iEdgePair = 0;
1221     while ( iEdgePair < edgeIDs1.size() )
1222     {
1223       if ( isComputed[ edgeIDs1[ iEdgePair ]] ||
1224            isComputed[ edgeIDs2[ iEdgePair ]])
1225       {
1226         // "projection" from one side to the other
1227
1228         size_t iEdgeComputed = edgeIDs1[iEdgePair], iSideComputed = 0;
1229         if ( !isComputed[ iEdgeComputed ])
1230           ++iSideComputed, iEdgeComputed = edgeIDs2[iEdgePair];
1231
1232         map< double, const SMDS_MeshNode* > nodeParams; // params of existing nodes
1233         if ( !SMESH_Algo::GetSortedNodesOnEdge( meshDS, theSinuEdges[ iEdgeComputed ], /*skipMedium=*/true, nodeParams ))
1234           return false;
1235
1236         SMESH_MAT2d::BoundaryPoint& bndPnt = bp[ 1-iSideComputed ];
1237         SMESH_MAT2d::BranchPoint brp;
1238         NodePoint npN, npB;
1239         NodePoint& np0 = iSideComputed ? npB : npN;
1240         NodePoint& np1 = iSideComputed ? npN : npB;
1241
1242         double maParam1st, maParamLast, maParam;
1243         if ( !theMA.getBoundary().getBranchPoint( iEdgeComputed, nodeParams.begin()->first, brp ))
1244             return false;
1245         branch.getParameter( brp, maParam1st );
1246         if ( !theMA.getBoundary().getBranchPoint( iEdgeComputed, nodeParams.rbegin()->first, brp ))
1247             return false;
1248         branch.getParameter( brp, maParamLast );
1249
1250         map< double, const SMDS_MeshNode* >::iterator u2n = nodeParams.begin(), u2nEnd = --nodeParams.end();
1251         TMAPar2NPoints::iterator end = pointsOnE.end(), pos = end;
1252         TMAPar2NPoints::iterator & hint = (maParamLast > maParam1st) ? end : pos;
1253         for ( ++u2n; u2n != u2nEnd; ++u2n )
1254         {
1255           if ( !theMA.getBoundary().getBranchPoint( iEdgeComputed, u2n->first, brp ))
1256             return false;
1257           if ( !branch.getBoundaryPoints( brp, bp[0], bp[1] ))
1258             return false;
1259           if ( !branch.getParameter( brp, maParam ))
1260             return false;
1261
1262           npN = NodePoint( u2n->second, u2n->first, iEdgeComputed );
1263           npB = NodePoint( bndPnt );
1264           pos = pointsOnE.insert( hint, make_pair( maParam, make_pair( np0, np1 )));
1265         }
1266
1267         // move iEdgePair forward
1268         while ( iEdgePair < edgeIDs1.size() )
1269           if ( edgeIDs1[ iEdgePair ] == bp[0]._edgeIndex &&
1270                edgeIDs2[ iEdgePair ] == bp[1]._edgeIndex )
1271             break;
1272           else
1273             ++iEdgePair;
1274       }
1275       else
1276       {
1277         // projection from MA
1278         maParams.clear();
1279         if ( !getParamsForEdgePair( iEdgePair, divPoints, theMAParams, maParams ))
1280           return false;
1281
1282         for ( size_t i = 1; i < maParams.size()-1; ++i )
1283         {
1284           if ( !branch.getBoundaryPoints( maParams[i], bp[0], bp[1] ))
1285             return false;
1286
1287           pointsOnE.insert( pointsOnE.end(), make_pair( maParams[i], make_pair( NodePoint(bp[0]),
1288                                                                                 NodePoint(bp[1]))));
1289         }
1290       }
1291       ++iEdgePair;
1292     }
1293
1294     if ( !projectVertices( theHelper, theMA, divPoints, edgeIDs1, edgeIDs2,
1295                            isComputed, pointsOnE, theSinuFace ))
1296       return false;
1297
1298     separateNodes( theHelper, pointsOnE, theSinuFace );
1299
1300     // create nodes
1301     TMAPar2NPoints::iterator u2np = pointsOnE.begin();
1302     for ( ; u2np != pointsOnE.end(); ++u2np )
1303     {
1304       NodePoint* np[2] = { & u2np->second.first, & u2np->second.second };
1305       for ( int iSide = 0; iSide < 2; ++iSide )
1306       {
1307         if ( np[ iSide ]->_node ) continue;
1308         size_t       iEdge = np[ iSide ]->_edgeInd;
1309         double           u = np[ iSide ]->_u;
1310         gp_Pnt           p = curves[ iEdge ]->Value( u );
1311         np[ iSide ]->_node = meshDS->AddNode( p.X(), p.Y(), p.Z() );
1312         meshDS->SetNodeOnEdge( np[ iSide ]->_node, edgeIDs[ iEdge ], u );
1313       }
1314     }
1315
1316     // create mesh segments on EDGEs
1317     theHelper.SetElementsOnShape( false );
1318     TopoDS_Face face = TopoDS::Face( theHelper.GetSubShape() );
1319     for ( size_t i = 0; i < theSinuEdges.size(); ++i )
1320     {
1321       SMESH_subMesh* sm = mesh->GetSubMesh( theSinuEdges[i] );
1322       if ( sm->GetSubMeshDS() && sm->GetSubMeshDS()->NbElements() > 0 )
1323         continue;
1324
1325       StdMeshers_FaceSide side( face, theSinuEdges[i], mesh,
1326                                 /*isFwd=*/true, /*skipMediumNodes=*/true );
1327       vector<const SMDS_MeshNode*> nodes = side.GetOrderedNodes();
1328       for ( size_t in = 1; in < nodes.size(); ++in )
1329       {
1330         const SMDS_MeshElement* seg = theHelper.AddEdge( nodes[in-1], nodes[in], 0, false );
1331         meshDS->SetMeshElementOnShape( seg, edgeIDs[ i ] );
1332       }
1333     }
1334
1335     // update sub-meshes on VERTEXes
1336     for ( size_t i = 0; i < theSinuEdges.size(); ++i )
1337     {
1338       mesh->GetSubMesh( theHelper.IthVertex( 0, theSinuEdges[i] ))
1339         ->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1340       mesh->GetSubMesh( theHelper.IthVertex( 1, theSinuEdges[i] ))
1341         ->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1342     }
1343
1344     return true;
1345   }
1346
1347   //================================================================================
1348   /*!
1349    * \brief Mesh short EDGEs
1350    */
1351   //================================================================================
1352
1353   bool computeShortEdges( SMESH_MesherHelper&        theHelper,
1354                           const vector<TopoDS_Edge>& theShortEdges,
1355                           SMESH_Algo*                the1dAlgo )
1356   {
1357     for ( size_t i = 0; i < theShortEdges.size(); ++i )
1358     {
1359       theHelper.GetGen()->Compute( *theHelper.GetMesh(), theShortEdges[i], true, true );
1360
1361       SMESH_subMesh* sm = theHelper.GetMesh()->GetSubMesh(theShortEdges[i] );
1362       if ( sm->IsEmpty() )
1363       {
1364         try {
1365           if ( !the1dAlgo->Compute( *theHelper.GetMesh(), theShortEdges[i] ))
1366             return false;
1367         }
1368         catch (...) {
1369           return false;
1370         }
1371         sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1372         if ( sm->IsEmpty() )
1373           return false;
1374       }
1375     }
1376     return true;
1377   }
1378
1379   inline double area( const UVPtStruct& p1, const UVPtStruct& p2, const UVPtStruct& p3 )
1380   {
1381     gp_XY v1 = p2.UV() - p1.UV();
1382     gp_XY v2 = p3.UV() - p1.UV();
1383     return v2 ^ v1;
1384   }
1385
1386   bool ellipticSmooth( FaceQuadStruct::Ptr quad, int nbLoops )
1387   {
1388     //nbLoops = 10;
1389     if ( quad->uv_grid.empty() )
1390       return true;
1391
1392     int nbhoriz  = quad->iSize;
1393     int nbvertic = quad->jSize;
1394
1395     const double dksi = 0.5, deta = 0.5;
1396     const double  dksi2 = dksi*dksi, deta2 = deta*deta;
1397     double err = 0., g11, g22, g12;
1398     int nbErr = 0;
1399
1400     FaceQuadStruct& q = *quad;
1401     UVPtStruct pNew;
1402
1403     double refArea = area( q.UVPt(0,0), q.UVPt(1,0), q.UVPt(1,1) );
1404
1405     for ( int iLoop = 0; iLoop < nbLoops; ++iLoop )
1406     {
1407       err = 0;
1408       for ( int i = 1; i < nbhoriz - 1; i++ )
1409         for ( int j = 1; j < nbvertic - 1; j++ )
1410         {
1411           g11 = ( (q.U(i,j+1) - q.U(i,j-1))*(q.U(i,j+1) - q.U(i,j-1))/dksi2 +
1412                   (q.V(i,j+1) - q.V(i,j-1))*(q.V(i,j+1) - q.V(i,j-1))/deta2 )/4;
1413
1414           g22 = ( (q.U(i+1,j) - q.U(i-1,j))*(q.U(i+1,j) - q.U(i-1,j))/dksi2 +
1415                   (q.V(i+1,j) - q.V(i-1,j))*(q.V(i+1,j) - q.V(i-1,j))/deta2 )/4;
1416
1417           g12 = ( (q.U(i+1,j) - q.U(i-1,j))*(q.U(i,j+1) - q.U(i,j-1))/dksi2 +
1418                   (q.V(i+1,j) - q.V(i-1,j))*(q.V(i,j+1) - q.V(i,j-1))/deta2 )/(4*dksi*deta);
1419
1420           pNew.u = dksi2/(2*(g11+g22)) * (g11*(q.U(i+1,j) + q.U(i-1,j))/dksi2 +
1421                                           g22*(q.U(i,j+1) + q.U(i,j-1))/dksi2
1422                                           - 0.5*g12*q.U(i+1,j+1) + 0.5*g12*q.U(i-1,j+1) +
1423                                           - 0.5*g12*q.U(i-1,j-1) + 0.5*g12*q.U(i+1,j-1));
1424
1425           pNew.v = deta2/(2*(g11+g22)) * (g11*(q.V(i+1,j) + q.V(i-1,j))/deta2 +
1426                                           g22*(q.V(i,j+1) + q.V(i,j-1))/deta2
1427                                           - 0.5*g12*q.V(i+1,j+1) + 0.5*g12*q.V(i-1,j+1) +
1428                                           - 0.5*g12*q.V(i-1,j-1) + 0.5*g12*q.V(i+1,j-1));
1429
1430           // if (( refArea * area( q.UVPt(i-1,j-1), q.UVPt(i,j-1), pNew ) > 0 ) &&
1431           //     ( refArea * area( q.UVPt(i+1,j-1), q.UVPt(i+1,j), pNew ) > 0 ) &&
1432           //     ( refArea * area( q.UVPt(i+1,j+1), q.UVPt(i,j+1), pNew ) > 0 ) &&
1433           //     ( refArea * area( q.UVPt(i-1,j), q.UVPt(i-1,j-1), pNew ) > 0 ))
1434           {
1435             err += sqrt(( q.U(i,j) - pNew.u ) * ( q.U(i,j) - pNew.u ) +
1436                         ( q.V(i,j) - pNew.v ) * ( q.V(i,j) - pNew.v ));
1437             q.U(i,j) = pNew.u;
1438             q.V(i,j) = pNew.v;
1439           }
1440           // else if ( ++nbErr < 10 )
1441           // {
1442           //   cout << i << ", " << j << endl;
1443           //   cout << "x = ["
1444           //        << "[ " << q.U(i-1,j-1) << ", " <<q.U(i,j-1) << ", " << q.U(i+1,j-1) << " ],"
1445           //        << "[ " << q.U(i-1,j-0) << ", " <<q.U(i,j-0) << ", " << q.U(i+1,j-0) << " ],"
1446           //        << "[ " << q.U(i-1,j+1) << ", " <<q.U(i,j+1) << ", " << q.U(i+1,j+1) << " ]]" << endl;
1447           //   cout << "y = ["
1448           //        << "[ " << q.V(i-1,j-1) << ", " <<q.V(i,j-1) << ", " << q.V(i+1,j-1) << " ],"
1449           //        << "[ " << q.V(i-1,j-0) << ", " <<q.V(i,j-0) << ", " << q.V(i+1,j-0) << " ],"
1450           //        << "[ " << q.V(i-1,j+1) << ", " <<q.V(i,j+1) << ", " << q.V(i+1,j+1) << " ]]" << endl<<endl;
1451           // }
1452         }
1453
1454       if ( err / ( nbhoriz - 2 ) / ( nbvertic - 2 ) < 1e-6 )
1455         break;
1456     }
1457     //cout << " ERR " << err / ( nbhoriz - 2 ) / ( nbvertic - 2 ) << endl;
1458
1459     return true;
1460   }
1461
1462   //================================================================================
1463   /*!
1464    * \brief Remove temporary node
1465    */
1466   //================================================================================
1467
1468   void mergeNodes( SMESH_MesherHelper& theHelper,
1469                    SinuousFace&        theSinuFace )
1470   {
1471     SMESH_MeshEditor editor( theHelper.GetMesh() );
1472     SMESH_MeshEditor::TListOfListOfNodes nodesGroups;
1473
1474     TMergeMap::iterator n2nn = theSinuFace._nodesToMerge.begin();
1475     for ( ; n2nn != theSinuFace._nodesToMerge.end(); ++n2nn )
1476     {
1477       if ( !n2nn->first ) continue;
1478       nodesGroups.push_back( list< const SMDS_MeshNode* >() );
1479       list< const SMDS_MeshNode* > & group = nodesGroups.back();
1480
1481       group.push_back( n2nn->first );
1482       group.splice( group.end(), n2nn->second );
1483     }
1484     editor.MergeNodes( nodesGroups );
1485   }
1486
1487 } // namespace
1488
1489 //================================================================================
1490 /*!
1491  * \brief Create quadrangle elements
1492  *  \param [in] theHelper - the helper
1493  *  \param [in] theFace - the face to mesh
1494  *  \param [in] theSinuEdges - the sinuous EDGEs
1495  *  \param [in] theShortEdges - the short EDGEs
1496  *  \return bool - is OK or not
1497  */
1498 //================================================================================
1499
1500 bool StdMeshers_QuadFromMedialAxis_1D2D::computeQuads( SMESH_MesherHelper&       theHelper,
1501                                                        const TopoDS_Face&        theFace,
1502                                                        const vector<TopoDS_Edge> theSinuEdges[2],
1503                                                        const vector<TopoDS_Edge> theShortEdges[2])
1504 {
1505   SMESH_Mesh* mesh = theHelper.GetMesh();
1506   SMESH_ProxyMesh::Ptr proxyMesh = StdMeshers_ViscousLayers2D::Compute( *mesh, theFace );
1507   if ( !proxyMesh )
1508     return false;
1509
1510   StdMeshers_Quadrangle_2D::myProxyMesh  = proxyMesh;
1511   StdMeshers_Quadrangle_2D::myHelper     = &theHelper;
1512   StdMeshers_Quadrangle_2D::myNeedSmooth = false;
1513   StdMeshers_Quadrangle_2D::myCheckOri   = false;
1514   StdMeshers_Quadrangle_2D::myQuadList.clear();
1515
1516   // fill FaceQuadStruct
1517
1518   list< TopoDS_Edge > side[4];
1519   side[0].insert( side[0].end(), theShortEdges[0].begin(), theShortEdges[0].end() );
1520   side[1].insert( side[1].end(), theSinuEdges[1].begin(),  theSinuEdges[1].end() );
1521   side[2].insert( side[2].end(), theShortEdges[1].begin(), theShortEdges[1].end() );
1522   side[3].insert( side[3].end(), theSinuEdges[0].begin(),  theSinuEdges[0].end() );
1523
1524   FaceQuadStruct::Ptr quad( new FaceQuadStruct );
1525   quad->side.resize( 4 );
1526   quad->face = theFace;
1527   for ( int i = 0; i < 4; ++i )
1528   {
1529     quad->side[i] = StdMeshers_FaceSide::New( theFace, side[i], mesh, i < QUAD_TOP_SIDE,
1530                                               /*skipMediumNodes=*/true, proxyMesh );
1531   }
1532   int nbNodesShort0 = quad->side[0].NbPoints();
1533   int nbNodesShort1 = quad->side[2].NbPoints();
1534
1535   // compute UV of internal points
1536   myQuadList.push_back( quad );
1537   if ( !StdMeshers_Quadrangle_2D::setNormalizedGrid( quad ))
1538     return false;
1539
1540   // elliptic smooth of internal points to get boundary cell normal to the boundary
1541   ellipticSmooth( quad, 1 );
1542
1543   // create quadrangles
1544   bool ok;
1545   if ( nbNodesShort0 == nbNodesShort1 )
1546     ok = StdMeshers_Quadrangle_2D::computeQuadDominant( *mesh, theFace, quad );
1547   else
1548     ok = StdMeshers_Quadrangle_2D::computeTriangles( *mesh, theFace, quad );
1549
1550   StdMeshers_Quadrangle_2D::myProxyMesh.reset();
1551   StdMeshers_Quadrangle_2D::myHelper = 0;
1552
1553   return ok;
1554 }
1555
1556 //================================================================================
1557 /*!
1558  * \brief Generate quadrangle mesh
1559  */
1560 //================================================================================
1561
1562 bool StdMeshers_QuadFromMedialAxis_1D2D::Compute(SMESH_Mesh&         theMesh,
1563                                                  const TopoDS_Shape& theShape)
1564 {
1565   SMESH_MesherHelper helper( theMesh );
1566   helper.SetSubShape( theShape );
1567
1568   TopoDS_Face F = TopoDS::Face( theShape );
1569   if ( F.Orientation() >= TopAbs_INTERNAL ) F.Orientation( TopAbs_FORWARD );
1570
1571   SinuousFace sinuFace( F );
1572
1573   _progress = 0.01;
1574
1575   if ( getSinuousEdges( helper, sinuFace ))
1576   {
1577     _progress = 0.2;
1578
1579     // if ( sinuFace._sinuEdges.size() > 2 )
1580     //   return error(COMPERR_BAD_SHAPE, "Not yet supported case" );
1581
1582     double minSegLen = getMinSegLen( helper, sinuFace._sinuEdges );
1583     SMESH_MAT2d::MedialAxis ma( F, sinuFace._sinuEdges, minSegLen, /*ignoreCorners=*/true );
1584
1585     if ( !_regular1D )
1586       _regular1D = new Algo1D( _studyId, _gen );
1587     _regular1D->SetSegmentLength( minSegLen );
1588
1589     vector<double> maParams;
1590     if ( ! divideMA( helper, ma, sinuFace, _regular1D, minSegLen, maParams ))
1591       return error(COMPERR_BAD_SHAPE);
1592
1593     _progress = 0.4;
1594
1595     if ( !computeShortEdges( helper, sinuFace._shortSide[0], _regular1D ) ||
1596          !computeShortEdges( helper, sinuFace._shortSide[1], _regular1D ))
1597       return error("Failed to mesh short edges");
1598
1599     _progress = 0.6;
1600
1601     if ( !computeSinuEdges( helper, minSegLen, ma, maParams, sinuFace ))
1602       return error("Failed to mesh sinuous edges");
1603
1604     _progress = 0.8;
1605
1606     bool ok = computeQuads( helper, F, sinuFace._sinuSide, sinuFace._shortSide );
1607
1608     if ( ok )
1609       mergeNodes( helper, sinuFace );
1610
1611     _progress = 1.;
1612
1613     return ok;
1614   }
1615
1616   return error(COMPERR_BAD_SHAPE, "Not implemented so far");
1617 }
1618
1619 //================================================================================
1620 /*!
1621  * \brief Predict nb of elements
1622  */
1623 //================================================================================
1624
1625 bool StdMeshers_QuadFromMedialAxis_1D2D::Evaluate(SMESH_Mesh &         theMesh,
1626                                                   const TopoDS_Shape & theShape,
1627                                                   MapShapeNbElems&     theResMap)
1628 {
1629   return StdMeshers_Quadrangle_2D::Evaluate(theMesh,theShape,theResMap);
1630 }
1631