Salome HOME
23304: [EDF 10304] Radial Quadrangle on ellipse
[modules/smesh.git] / src / StdMeshers / StdMeshers_RadialQuadrangle_1D2D.cxx
1 // Copyright (C) 2007-2016  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 // File  : StdMeshers_RadialQuadrangle_1D2D.cxx
21 // Module: SMESH
22
23 #include "StdMeshers_RadialQuadrangle_1D2D.hxx"
24
25 #include "StdMeshers_NumberOfLayers.hxx"
26 #include "StdMeshers_LayerDistribution.hxx"
27 #include "StdMeshers_Regular_1D.hxx"
28 #include "StdMeshers_NumberOfSegments.hxx"
29 #include "StdMeshers_FaceSide.hxx"
30
31 #include "SMDS_MeshNode.hxx"
32 #include "SMESHDS_Mesh.hxx"
33 #include "SMESHDS_SubMesh.hxx"
34 #include "SMESH_Gen.hxx"
35 #include "SMESH_HypoFilter.hxx"
36 #include "SMESH_Mesh.hxx"
37 #include "SMESH_MesherHelper.hxx"
38 #include "SMESH_subMesh.hxx"
39 #include "SMESH_subMeshEventListener.hxx"
40 #include "SMESH_Block.hxx"
41
42 #include "utilities.h"
43
44 #include <BRepAdaptor_CompCurve.hxx>
45 #include <BRepAdaptor_Curve.hxx>
46 #include <BRepBuilderAPI_MakeEdge.hxx>
47 #include <BRep_Builder.hxx>
48 #include <BRep_Tool.hxx>
49 #include <GCPnts_AbscissaPoint.hxx>
50 #include <Geom2d_Line.hxx>
51 #include <Geom2d_TrimmedCurve.hxx>
52 #include <GeomAPI_ProjectPointOnSurf.hxx>
53 #include <Geom_Circle.hxx>
54 #include <Geom_Line.hxx>
55 #include <Geom_TrimmedCurve.hxx>
56 #include <ShapeFix_Edge.hxx>
57 #include <TColgp_SequenceOfPnt.hxx>
58 #include <TColgp_SequenceOfPnt2d.hxx>
59 #include <TopExp.hxx>
60 #include <TopExp_Explorer.hxx>
61 #include <TopTools_ListIteratorOfListOfShape.hxx>
62 #include <TopoDS.hxx>
63
64
65 using namespace std;
66
67 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
68 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
69
70
71 //=======================================================================
72 //function : StdMeshers_RadialQuadrangle_1D2D
73 //purpose  : 
74 //=======================================================================
75
76 StdMeshers_RadialQuadrangle_1D2D::StdMeshers_RadialQuadrangle_1D2D(int        hypId,
77                                                                    int        studyId,
78                                                                    SMESH_Gen* gen)
79   :StdMeshers_Quadrangle_2D( hypId, studyId, gen )
80 {
81   _name = "RadialQuadrangle_1D2D";
82   _shapeType = (1 << TopAbs_FACE);        // 1 bit per shape type
83
84   _compatibleHypothesis.push_back("LayerDistribution2D");
85   _compatibleHypothesis.push_back("NumberOfLayers2D");
86   _requireDiscreteBoundary = false;
87   _supportSubmeshes = true;
88   _neededLowerHyps[ 1 ] = true;  // suppress warning on hiding a global 1D algo
89
90   myNbLayerHypo      = 0;
91   myDistributionHypo = 0;
92 }
93
94
95 //================================================================================
96 /*!
97  * \brief Destructor
98  */
99 //================================================================================
100
101 StdMeshers_RadialQuadrangle_1D2D::~StdMeshers_RadialQuadrangle_1D2D()
102 {}
103
104
105 //=======================================================================
106 //function : CheckHypothesis
107 //purpose  : 
108 //=======================================================================
109
110 bool StdMeshers_RadialQuadrangle_1D2D::CheckHypothesis
111                            (SMESH_Mesh&                          aMesh,
112                             const TopoDS_Shape&                  aShape,
113                             SMESH_Hypothesis::Hypothesis_Status& aStatus)
114 {
115   // check aShape 
116   myNbLayerHypo = 0;
117   myDistributionHypo = 0;
118
119   list <const SMESHDS_Hypothesis * >::const_iterator itl;
120
121   const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
122   if ( hyps.size() == 0 ) {
123     aStatus = SMESH_Hypothesis::HYP_OK;
124     return true;  // can work with no hypothesis
125   }
126
127   if ( hyps.size() > 1 ) {
128     aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
129     return false;
130   }
131
132   const SMESHDS_Hypothesis *theHyp = hyps.front();
133
134   string hypName = theHyp->GetName();
135
136   if (hypName == "NumberOfLayers2D") {
137     myNbLayerHypo = static_cast<const StdMeshers_NumberOfLayers *>(theHyp);
138     aStatus = SMESH_Hypothesis::HYP_OK;
139     return true;
140   }
141   if (hypName == "LayerDistribution2D") {
142     myDistributionHypo = static_cast<const StdMeshers_LayerDistribution *>(theHyp);
143     aStatus = SMESH_Hypothesis::HYP_OK;
144     return true;
145   }
146   aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
147   return true;
148 }
149
150 namespace
151 {
152   // ------------------------------------------------------------------------------
153   /*!
154    * \brief Listener used to mark edges meshed by StdMeshers_RadialQuadrangle_1D2D
155    */
156   class TEdgeMarker : public SMESH_subMeshEventListener
157   {
158     TEdgeMarker(): SMESH_subMeshEventListener(/*isDeletable=*/false,
159                                               "StdMeshers_RadialQuadrangle_1D2D::TEdgeMarker") {}
160   public:
161     //!<  Return static listener
162     static SMESH_subMeshEventListener* getListener()
163     {
164       static TEdgeMarker theEdgeMarker;
165       return &theEdgeMarker;
166     }
167     //! Clear face sumbesh if something happens on edges
168     void ProcessEvent(const int          event,
169                       const int          eventType,
170                       SMESH_subMesh*     edgeSubMesh,
171                       EventListenerData* data,
172                       const SMESH_Hypothesis*  /*hyp*/)
173     {
174       if ( data && !data->mySubMeshes.empty() && eventType == SMESH_subMesh::ALGO_EVENT)
175       {
176         ASSERT( data->mySubMeshes.front() != edgeSubMesh );
177         SMESH_subMesh* faceSubMesh = data->mySubMeshes.front();
178         faceSubMesh->ComputeStateEngine( SMESH_subMesh::CLEAN );
179       }
180     }
181   };
182
183   //================================================================================
184   /*!
185    * \brief Mark an edge as computed by StdMeshers_RadialQuadrangle_1D2D
186    */
187   //================================================================================
188
189   void markEdgeAsComputedByMe(const TopoDS_Edge& edge, SMESH_subMesh* faceSubMesh)
190   {
191     if ( SMESH_subMesh* edgeSM = faceSubMesh->GetFather()->GetSubMeshContaining( edge ))
192     {
193       if ( !edgeSM->GetEventListenerData( TEdgeMarker::getListener() ))
194         faceSubMesh->SetEventListener( TEdgeMarker::getListener(),
195                                        SMESH_subMeshEventListenerData::MakeData(faceSubMesh),
196                                        edgeSM);
197     }
198   }
199
200   //================================================================================
201   /*!
202    * \brief Return sides of the face connected in the order: aCircEdge, aLinEdge1, aLinEdge2
203    *  \retval int - nb of sides
204    */
205   //================================================================================
206
207   int analyseFace(const TopoDS_Shape&     aShape,
208                   SMESH_Mesh*             aMesh,
209                   StdMeshers_FaceSidePtr& aCircSide,
210                   StdMeshers_FaceSidePtr& aLinSide1,
211                   StdMeshers_FaceSidePtr& aLinSide2)
212   {
213     const TopoDS_Face& face = TopoDS::Face( aShape );
214     aCircSide.reset(); aLinSide1.reset(); aLinSide2.reset();
215
216     list< TopoDS_Edge > edges;
217     list< int > nbEdgesInWire;
218     int nbWire = SMESH_Block::GetOrderedEdges ( face, edges, nbEdgesInWire );
219     if ( nbWire > 2 || nbEdgesInWire.front() < 1 ) return 0;
220
221     // remove degenerated EDGEs
222     list<TopoDS_Edge>::iterator edge = edges.begin();
223     while ( edge != edges.end() )
224       if ( SMESH_Algo::isDegenerated( *edge ))
225         edge = edges.erase( edge );
226       else
227         ++edge;
228     int nbEdges = edges.size();
229
230     // find VERTEXes between continues EDGEs
231     TopTools_MapOfShape contVV;
232     if ( nbEdges > 1 )
233     {
234       TopoDS_Edge ePrev = edges.back();
235       for ( edge = edges.begin(); edge != edges.end(); ++edge )
236       {
237         if ( SMESH_Algo::IsContinuous( ePrev, *edge ))
238           contVV.Add( SMESH_MesherHelper::IthVertex( 0, *edge ));
239         ePrev = *edge;
240       }
241     }
242     // make edges start from a non-continues VERTEX
243     if ( 1 < contVV.Extent() && contVV.Extent() < nbEdges )
244     {
245       while ( contVV.Contains( SMESH_MesherHelper::IthVertex( 0, edges.front() )))
246         edges.splice( edges.end(), edges, edges.begin() );
247     }
248
249     // make face sides
250     TSideVector sides;
251     while ( !edges.empty() )
252     {
253       list< TopoDS_Edge > sideEdges;
254       sideEdges.splice( sideEdges.end(), edges, edges.begin() );
255       while ( !edges.empty() &&
256               contVV.Contains( SMESH_MesherHelper::IthVertex( 0, edges.front() )))
257         sideEdges.splice( sideEdges.end(), edges, edges.begin() );
258
259       StdMeshers_FaceSidePtr side;
260       if ( aMesh ) 
261         side = StdMeshers_FaceSide::New( face, sideEdges, aMesh,
262                                          /*isFwd=*/true, /*skipMedium=*/ true );
263       sides.push_back( side );
264     }
265
266     if ( !aMesh ) // call from IsApplicable()
267       return sides.size();
268
269     if ( sides.size() > 3 )
270       return sides.size();
271
272     if ( nbWire == 2 && (( sides.size() != 2 ) ||
273                          ( sides[0]->IsClosed() && sides[1]->IsClosed() ) ||
274                          ( !sides[0]->IsClosed() && !sides[1]->IsClosed() )))
275       return -1;
276
277     // detect an elliptic side
278
279     if ( sides.size() == 1 )
280     {
281       aCircSide = sides[0];
282       return sides.size();
283     }
284
285     // sort sides by deviation from a straight line
286     multimap< double, int > deviation2sideInd;
287     const double nbSamples = 7;
288     for ( size_t iS = 0; iS < sides.size(); ++iS )
289     {
290       gp_Pnt pf = BRep_Tool::Pnt( sides[iS]->FirstVertex() );
291       gp_Pnt pl = BRep_Tool::Pnt( sides[iS]->LastVertex() );
292       gp_Vec v1( pf, pl );
293       double v1Len = v1.Magnitude();
294       if ( v1Len < std::numeric_limits< double >::min() )
295       {
296         deviation2sideInd.insert( make_pair( sides[iS]->Length(), iS )); // the side seems closed
297         continue;
298       }
299       double devia = 0;
300       for ( int i = 0; i < nbSamples; ++i )
301       {
302         gp_Pnt pi( sides[iS]->Value3d(( i + 1 ) / nbSamples ));
303         gp_Vec vi( pf, pi );
304         double h = 0.5 * v1.Crossed( vi ).Magnitude() / v1Len;
305         devia = Max( devia, h );
306       }
307       deviation2sideInd.insert( make_pair( devia, iS ));
308     }
309
310     int iCirc = deviation2sideInd.rbegin()->second; 
311     aCircSide = sides[ iCirc ];
312     aLinSide1 = sides[( iCirc + 1 ) % sides.size() ];
313     if ( sides.size() > 2 )
314     {
315       aLinSide2 = sides[( iCirc + 2 ) % sides.size() ];
316       aLinSide2->Reverse(); // to be "parallel" to aLinSide1
317     }
318
319     if (( nbWire == 2 && aLinSide1 ) &&
320         ( aLinSide1->Edge(0).Orientation() == TopAbs_INTERNAL ) &&
321         ( aCircSide->IsClosed() ))
322     {
323       // assure that aCircSide starts at aLinSide1
324       TopoDS_Vertex v0 = aLinSide1->FirstVertex();
325       TopoDS_Vertex v1 = aLinSide1->LastVertex();
326       if ( ! aCircSide->FirstVertex().IsSame( v0 ) &&
327            ! aCircSide->FirstVertex().IsSame( v1 ))
328       {
329         int iE = 0;
330         for ( ; iE < aCircSide->NbEdges(); ++iE )
331           if ( aCircSide->FirstVertex(iE).IsSame( v0 ) ||
332                aCircSide->FirstVertex(iE).IsSame( v1 ))
333             break;
334         if ( iE == aCircSide->NbEdges() )
335           return -2;
336
337         edges.clear();
338         for ( int i = 0; i < aCircSide->NbEdges(); ++i, ++iE )
339           edges.push_back( aCircSide->Edge( iE % aCircSide->NbEdges() ));
340
341         aCircSide = StdMeshers_FaceSide::New( face, edges, aMesh,
342                                               /*isFwd=*/true, /*skipMedium=*/ true );
343       }
344     }
345
346     return sides.size();
347   }
348
349   //================================================================================
350   /*!
351    * \brief Checks if the common vertex between LinSide's lies inside the circle
352    *  and not outside
353    *  \return bool - false if there are 3 EDGEs and the corner is outside
354    */
355   //================================================================================
356
357   bool isCornerInsideCircle(const StdMeshers_FaceSidePtr& CircSide,
358                             const StdMeshers_FaceSidePtr& LinSide1,
359                             const StdMeshers_FaceSidePtr& LinSide2)
360   {
361     // if ( CircSide && LinSide1 && LinSide2 )
362     // {
363     //   Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircSide ));
364     //   TopoDS_Vertex aCommonV;
365     //   if ( !aCirc.IsNull() &&
366     //        TopExp::CommonVertex( LinSide1, LinSide2, aCommonV ))
367     //   {
368     //     gp_Pnt aCommonP = BRep_Tool::Pnt( aCommonV );
369     //     gp_Pnt  aCenter = aCirc->Location();
370     //     double     dist = aCenter.Distance( aCommonP );
371     //     return dist < 0.1 * aCirc->Radius();
372     //   }
373     // }
374     return true;
375   }
376
377   //================================================================================
378   /*!
379    * \brief Create an EDGE connecting the ellipse center with the most distant point
380    *        of the ellipse.
381    *  \param [in] circSide - the elliptic side
382   *  \param [in] face - the FACE
383   *  \param [out] circNode - a node on circSide most distant from the center
384   *  \return TopoDS_Edge - the create EDGE
385   */
386   //================================================================================
387
388   TopoDS_Edge makeEdgeToCenter( StdMeshers_FaceSidePtr& circSide,
389                                 const TopoDS_Face&      face,
390                                 const SMDS_MeshNode*&   circNode)
391   {
392     // find the center and a point most distant from it
393
394     double maxDist = 0, normPar;
395     gp_XY uv1, uv2;
396     for ( int i = 0; i < 32; ++i )
397     {
398       double    u = 0.5 * i / 32.;
399       gp_Pnt2d p1 = circSide->Value2d( u );
400       gp_Pnt2d p2 = circSide->Value2d( u + 0.5 );
401       double dist = p1.SquareDistance( p2 );
402       if ( dist > maxDist )
403       {
404         maxDist = dist;
405         uv1 = p1.XY();
406         uv2 = p2.XY();
407         normPar = u;
408       }
409     }
410     gp_XY center = 0.5 * ( uv1 + uv2 );
411     double   len = 0.5 * Sqrt( maxDist );
412     bool  isCirc = ( Abs( len - circSide->Value2d( 0 ).Distance( center )) < 1e-3 * len );
413
414     // find a node closest to the most distant point
415
416     size_t iDist = 0;
417     const UVPtStructVec& circNodes = circSide->GetUVPtStruct();
418     if ( !isCirc )
419     {
420       double minDist = 1e100;
421       for ( size_t i = 0; i <= circNodes.size(); ++i )
422       {
423         double dist = Abs( circNodes[i].normParam - normPar );
424         if ( dist < minDist )
425         {
426           iDist = i;
427           minDist = dist;
428         }
429       }
430     }
431     circNode = circNodes[iDist].node;
432     uv1      = circNodes[iDist].UV();
433     len      = ( uv1 - center ).Modulus();
434
435     // make the EDGE between the most distant point and the center
436
437     Handle(Geom2d_Line) line = new Geom2d_Line( uv1, gp_Dir2d( center - uv1 ) );
438     Handle(Geom2d_Curve) pcu = new Geom2d_TrimmedCurve( line, 0, len );
439
440     Handle(Geom_Surface) surface = BRep_Tool::Surface( face );
441     TopoDS_Edge edge = BRepBuilderAPI_MakeEdge( pcu, surface, 0, len );
442
443     BRep_Builder().UpdateEdge( edge, pcu, face, 1e-7 );
444     ShapeFix_Edge().FixAddCurve3d( edge );
445
446
447     // assure that circSide starts at circNode
448     if ( iDist != 0  &&  iDist != circNodes.size()-1 )
449     {
450       // create a new circSide
451       UVPtStructVec nodesNew;
452       nodesNew.reserve( circNodes.size() );
453       nodesNew.insert( nodesNew.end(), circNodes.begin() + iDist, circNodes.end() );
454       nodesNew.insert( nodesNew.end(), circNodes.begin() + 1, circNodes.begin() + iDist + 1 );
455       circSide = StdMeshers_FaceSide::New( nodesNew );
456     }
457
458     return edge;
459   }
460
461   //================================================================================
462   /*!
463    * \brief Set nodes existing on a linSide to UVPtStructVec and create missing nodes
464    *        corresponding to layerPositions
465    */
466   //================================================================================
467
468   void makeMissingMesh( StdMeshers_FaceSidePtr& linSide,
469                         UVPtStructVec&          nodes,
470                         const vector< double >& layerPositions,
471                         SMESH_MesherHelper*     helper )
472   {
473     // tolerance to compare normParam
474     double tol = 1e100;
475     for ( size_t i = 1; i < layerPositions.size(); ++i )
476       tol = Min( tol, layerPositions[i] - layerPositions[i-1] );
477     tol *= 0.05;
478
479     // merge existing nodes with layerPositions into UVPtStructVec
480     // ------------------------------------------------------------
481
482     const UVPtStructVec& exiNodes = linSide->GetUVPtStruct();
483     nodes.clear();
484     nodes.reserve( layerPositions.size() + exiNodes.size() );
485     vector< double >::const_iterator pos = layerPositions.begin(), posEnd = layerPositions.end();
486     for ( size_t i = 0; i < exiNodes.size(); ++i )
487     {
488       switch ( exiNodes[i].node->GetPosition()->GetTypeOfPosition() )
489       {
490       case SMDS_TOP_VERTEX:
491       {
492         // allocate UVPtStruct's for non-existing nodes
493         while ( pos != posEnd && *pos < exiNodes[i].normParam - tol )
494         {
495           UVPtStruct uvPS;
496           uvPS.normParam = *pos++;
497           nodes.push_back( uvPS );
498         }
499         // save existing node on a VERTEX
500         nodes.push_back( exiNodes[i] );
501         break;
502       }
503       case SMDS_TOP_EDGE:
504       {
505         // save existing nodes on an EDGE
506         while ( i < exiNodes.size() && exiNodes[i].node->GetPosition()->GetDim() == 1 )
507         {
508           nodes.push_back( exiNodes[i++] );
509         }
510         // save existing node on a VERTEX
511         if ( i < exiNodes.size() && exiNodes[i].node->GetPosition()->GetDim() == 0 )
512         {
513           nodes.push_back( exiNodes[i] );
514         }
515         break;
516       }
517       default:;
518       }
519
520       // skip layer positions covered by saved nodes
521       while ( pos != posEnd && *pos < nodes.back().normParam + tol )
522       {
523         ++pos;
524       }
525     }
526     // allocate UVPtStruct's for the rest non-existing nodes
527     while ( pos != posEnd )
528     {
529       UVPtStruct uvPS;
530       uvPS.normParam = *pos++;
531       nodes.push_back( uvPS );
532     }
533
534     // create missing nodes
535     // ---------------------
536
537     SMESHDS_Mesh * meshDS = helper->GetMeshDS();
538     const TopoDS_Face& F  = TopoDS::Face( helper->GetSubShape() );
539     Handle(ShapeAnalysis_Surface) surface = helper->GetSurface( F );
540
541     helper->SetElementsOnShape( false ); // we create nodes on EDGEs, not on the FACE
542
543     for ( size_t i = 0; i < nodes.size(); ++i )
544     {
545       if ( nodes[ i ].node ) continue;
546
547       gp_Pnt2d uv = linSide->Value2d( nodes[i].normParam );
548       gp_Pnt  xyz = surface->Value( uv.X(), uv.Y() );
549
550       nodes[ i ].SetUV( uv.XY() );
551       nodes[ i ].node = helper->AddNode( xyz.X(), xyz.Y(), xyz.Z() );
552     }
553
554     // set nodes on VERTEXes
555     for ( int iE = 0; iE < linSide->NbEdges(); ++iE )
556     {
557       TopoDS_Vertex v = linSide->LastVertex( iE );
558       if ( SMESH_Algo::VertexNode( v, meshDS ))
559         continue;
560
561       double normPar = linSide->LastParameter( iE );
562       size_t i = 0;
563       while ( nodes[ i ].normParam < normPar )
564         ++i;
565       if (( nodes[ i ].normParam - normPar ) > ( normPar - nodes[ i-1 ].normParam ))
566         --i;
567       meshDS->SetNodeOnVertex( nodes[ i ].node, v );
568     }
569
570     // set nodes on EDGEs
571     int edgeID;
572     for ( size_t i = 0; i < nodes.size(); ++i )
573     {
574       if ( nodes[ i ].node->getshapeId() > 0 ) continue;
575
576       double u = linSide->Parameter( nodes[i].normParam, edgeID );
577       meshDS->SetNodeOnEdge( nodes[ i ].node, edgeID, u );
578     }
579
580     // create segments
581     for ( size_t i = 1; i < nodes.size(); ++i )
582     {
583       if ( meshDS->FindEdge( nodes[i].node, nodes[i-1].node )) continue;
584
585       const SMDS_MeshElement* seg = helper->AddEdge( nodes[i].node, nodes[i-1].node );
586
587       double normParam = 0.5 * ( nodes[i].normParam + nodes[i-1].normParam );
588       edgeID = linSide->EdgeID( linSide->EdgeIndex( normParam ));
589       meshDS->SetMeshElementOnShape( seg, edgeID );
590     }
591
592     helper->SetElementsOnShape( true );
593
594   } // end makeMissingMesh()
595
596 //================================================================================
597 //================================================================================
598 /*!
599  * \brief Class computing layers distribution using data of
600  *        StdMeshers_LayerDistribution hypothesis
601  */
602 //================================================================================
603 //================================================================================
604
605 class TNodeDistributor: public StdMeshers_Regular_1D
606 {
607   list <const SMESHDS_Hypothesis *> myUsedHyps;
608 public:
609   // -----------------------------------------------------------------------------
610   static TNodeDistributor* GetDistributor(SMESH_Mesh& aMesh)
611   {
612     const int myID = -1001;
613     TNodeDistributor* myHyp = dynamic_cast<TNodeDistributor*>( aMesh.GetHypothesis( myID ));
614     if ( !myHyp )
615       myHyp = new TNodeDistributor( myID, 0, aMesh.GetGen() );
616     return myHyp;
617   }
618   // -----------------------------------------------------------------------------
619   //! Computes distribution of nodes on a straight line ending at pIn and pOut
620   bool Compute( vector< double > &      positions,
621                 const TopoDS_Edge&      edge,
622                 Adaptor3d_Curve&        curve,
623                 double                  f,
624                 double                  l,
625                 SMESH_Mesh&             mesh,
626                 const SMESH_Hypothesis* hyp1d)
627   {
628     if ( !hyp1d ) return error( "Invalid LayerDistribution hypothesis");
629
630     myUsedHyps.clear();
631     myUsedHyps.push_back( hyp1d );
632
633     SMESH_Hypothesis::Hypothesis_Status aStatus;
634     if ( !StdMeshers_Regular_1D::CheckHypothesis( mesh, edge, aStatus ))
635       return error( "StdMeshers_Regular_1D::CheckHypothesis() failed "
636                     "with LayerDistribution hypothesis");
637
638     double len = GCPnts_AbscissaPoint::Length( curve, f, l );
639
640     list< double > params;
641     if ( !StdMeshers_Regular_1D::computeInternalParameters( mesh, curve, len, f, l, params, false ))
642       return error("StdMeshers_Regular_1D failed to compute layers distribution");
643
644     params.push_front( f );
645     params.push_back ( l );
646     positions.clear();
647     positions.reserve( params.size() );
648     for (list<double>::iterator itU = params.begin(); itU != params.end(); itU++)
649       positions.push_back(( *itU - f ) / ( l - f ));
650     return true;
651   }
652   // -----------------------------------------------------------------------------
653   //! Make mesh on an adge using assigned 1d hyp or defaut nb of segments
654   bool ComputeCircularEdge( SMESH_Mesh&                   aMesh,
655                             const StdMeshers_FaceSidePtr& aSide )
656   {
657     bool ok = true;
658     for ( int i = 0; i < aSide->NbEdges(); ++i )
659     {
660       const TopoDS_Edge& edge = aSide->Edge( i );
661       _gen->Compute( aMesh, edge );
662       SMESH_subMesh *sm = aMesh.GetSubMesh( edge );
663       if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK)
664       {
665         // find any 1d hyp assigned (there can be a hyp w/o algo)
666         myUsedHyps = SMESH_Algo::GetUsedHypothesis( aMesh, edge, /*ignoreAux=*/true );
667         Hypothesis_Status aStatus;
668         if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, edge, aStatus ))
669         {
670           // no valid 1d hyp assigned, use default nb of segments
671           _hypType                    = NB_SEGMENTS;
672           _ivalue[ DISTR_TYPE_IND ]   = StdMeshers_NumberOfSegments::DT_Regular;
673           _ivalue[ NB_SEGMENTS_IND  ] = _gen->GetDefaultNbSegments();
674         }
675         ok &= StdMeshers_Regular_1D::Compute( aMesh, edge );
676       }
677     }
678     return ok;
679   }
680   // -----------------------------------------------------------------------------
681   //! Make mesh on an adge using assigned 1d hyp or defaut nb of segments
682   bool EvaluateCircularEdge(SMESH_Mesh&                  aMesh,
683                             const StdMeshers_FaceSidePtr aSide,
684                             MapShapeNbElems&             aResMap)
685   {
686     bool ok = true;
687     for ( int i = 0; i < aSide->NbEdges(); ++i )
688     {
689       const TopoDS_Edge& anEdge = aSide->Edge( i );
690       _gen->Evaluate( aMesh, anEdge, aResMap );
691       if ( aResMap.count( aMesh.GetSubMesh( anEdge )))
692         continue;
693
694       // find any 1d hyp assigned
695       myUsedHyps = SMESH_Algo::GetUsedHypothesis(aMesh, anEdge, /*ignoreAux=*/true);
696       Hypothesis_Status aStatus;
697       if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, anEdge, aStatus ))
698       {
699         // no valid 1d hyp assigned, use default nb of segments
700         _hypType                    = NB_SEGMENTS;
701         _ivalue[ DISTR_TYPE_IND ]   = StdMeshers_NumberOfSegments::DT_Regular;
702         _ivalue[ NB_SEGMENTS_IND  ] = _gen->GetDefaultNbSegments();
703       }
704       ok &= StdMeshers_Regular_1D::Evaluate( aMesh, anEdge, aResMap );
705     }
706     return ok;
707   }
708 protected:
709   // -----------------------------------------------------------------------------
710   TNodeDistributor( int hypId, int studyId, SMESH_Gen* gen)
711     : StdMeshers_Regular_1D( hypId, studyId, gen)
712   {
713   }
714   // -----------------------------------------------------------------------------
715   virtual const list <const SMESHDS_Hypothesis *> &
716     GetUsedHypothesis(SMESH_Mesh &, const TopoDS_Shape &, const bool)
717   {
718     return myUsedHyps;
719   }
720   // -----------------------------------------------------------------------------
721 };
722 }
723
724 //=======================================================================
725 /*!
726  * \brief Allow algo to do something after persistent restoration
727  * \param subMesh - restored submesh
728  *
729  * call markEdgeAsComputedByMe()
730  */
731 //=======================================================================
732
733 void StdMeshers_RadialQuadrangle_1D2D::SubmeshRestored(SMESH_subMesh* faceSubMesh)
734 {
735   if ( !faceSubMesh->IsEmpty() )
736   {
737     for ( TopExp_Explorer e( faceSubMesh->GetSubShape(), TopAbs_EDGE ); e.More(); e.Next() )
738     {
739       markEdgeAsComputedByMe( TopoDS::Edge( e.Current() ), faceSubMesh );
740     }
741   }
742 }
743
744 //=======================================================================
745 //function : Compute
746 //purpose  :
747 //=======================================================================
748
749 bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh&         aMesh,
750                                                const TopoDS_Shape& aShape)
751 {
752   SMESH_MesherHelper helper( aMesh );
753   StdMeshers_Quadrangle_2D::myHelper     = & helper;
754   StdMeshers_Quadrangle_2D::myNeedSmooth = false;
755   StdMeshers_Quadrangle_2D::myCheckOri   = false;
756   StdMeshers_Quadrangle_2D::myQuadList.clear();
757
758   StdMeshers_FaceSidePtr circSide, linSide1, linSide2;
759   int nbSides = analyseFace( aShape, &aMesh, circSide, linSide1, linSide2 );
760   if( nbSides > 3 || nbSides < 1 )
761     return error("The face must be a full ellipse or a part of ellipse (i.e. the number "
762                  "of edges is less or equal to 3 and one of them is an ellipse curve)");
763
764   // get not yet computed EDGEs
765   list< TopoDS_Edge > emptyEdges;
766   for ( TopExp_Explorer e( aShape, TopAbs_EDGE ); e.More(); e.Next() )
767   {
768     if ( aMesh.GetSubMesh( e.Current() )->IsEmpty() )
769       emptyEdges.push_back( TopoDS::Edge( e.Current() ));
770   }
771
772   TNodeDistributor* algo1d = TNodeDistributor::GetDistributor(aMesh);
773
774   if ( !algo1d->ComputeCircularEdge( aMesh, circSide ))
775     return error( algo1d->GetComputeError() );
776
777
778   TopoDS_Face F = TopoDS::Face(aShape);
779   Handle(Geom_Surface) S = BRep_Tool::Surface(F);
780
781   myHelper->IsQuadraticSubMesh( aShape );
782   myHelper->SetElementsOnShape( true );
783
784   vector< double > layerPositions; // [0,1]
785
786   const SMDS_MeshNode* centerNode = 0;
787   gp_Pnt2d             centerUV(0,0);
788
789   // ------------------------------------------------------------------------------------------
790   if ( nbSides == 1 ) // full ellipse
791   {
792     const SMDS_MeshNode* circNode;
793     TopoDS_Edge linEdge = makeEdgeToCenter( circSide, F, circNode );
794
795     StdMeshers_FaceSidePtr tmpSide =
796       StdMeshers_FaceSide::New( F, linEdge, &aMesh, /*isFrw=*/true, /*skipMedium=*/true );
797
798     if ( !computeLayerPositions( tmpSide, layerPositions ))
799       return false;
800
801     UVPtStructVec nodes( layerPositions.size() );
802     nodes[0].node = circNode;
803     for ( size_t i = 0; i < layerPositions.size(); ++i )
804     {
805       gp_Pnt2d uv = tmpSide->Value2d( layerPositions[i] );
806       gp_Pnt  xyz = S->Value( uv.X(), uv.Y() );
807
808       nodes[ i ].SetUV( uv.XY() );
809       nodes[ i ].normParam = layerPositions[i];
810       if ( i )
811         nodes[ i ].node = myHelper->AddNode( xyz.X(), xyz.Y(), xyz.Z(), 0, uv.X(), uv.Y() );
812     }
813
814     linSide1 = StdMeshers_FaceSide::New( nodes );
815     linSide2 = StdMeshers_FaceSide::New( nodes );
816
817     centerNode = nodes.back().node;
818     centerUV   = nodes.back().UV();
819   }
820   // ------------------------------------------------------------------------------------------
821   else if ( nbSides == 2 && linSide1->Edge(0).Orientation() == TopAbs_INTERNAL )
822   {
823     // full ellipse with an internal radial side
824
825     // eliminate INTERNAL orientation
826     list< TopoDS_Edge > edges;
827     for ( int iE = 0; iE < linSide1->NbEdges(); ++iE )
828     {
829       edges.push_back( linSide1->Edge( iE ));
830       edges.back().Orientation( TopAbs_FORWARD );
831     }
832
833     // orient the internal side
834     bool isVIn0Shared = false;
835     TopoDS_Vertex vIn0 = myHelper->IthVertex( 0, edges.front() );
836     for ( int iE = 0; iE < circSide->NbEdges() && !isVIn0Shared; ++iE )
837       isVIn0Shared = vIn0.IsSame( circSide->FirstVertex( iE ));
838
839     linSide1 = StdMeshers_FaceSide::New( F, edges, &aMesh,
840                                          /*isFrw=*/isVIn0Shared, /*skipMedium=*/true );
841
842     int nbMeshedEdges;
843     if ( !computeLayerPositions( linSide1, layerPositions, &nbMeshedEdges ))
844       return false;
845
846     // merge existing nodes with new nodes at layerPositions into a UVPtStructVec
847     UVPtStructVec nodes;
848     if ( nbMeshedEdges != linSide1->NbEdges() )
849       makeMissingMesh( linSide1, nodes, layerPositions, myHelper );
850     else
851       nodes = linSide1->GetUVPtStruct();
852
853     linSide1 = StdMeshers_FaceSide::New( nodes );
854     linSide2 = StdMeshers_FaceSide::New( nodes );
855
856     centerNode = nodes.back().node;
857     centerUV   = nodes.back().UV();
858   }
859   // ------------------------------------------------------------------------------------------
860   else if ( nbSides == 2 )
861   {
862     // find positions of layers for the first half of linSide1
863     int nbMeshedEdges;
864     if ( !computeLayerPositions( linSide1, layerPositions, &nbMeshedEdges, /*useHalf=*/true ))
865       return false;
866
867     // make positions for the whole linSide1
868     for ( size_t i = 0; i < layerPositions.size(); ++i )
869     {
870       layerPositions[i] *= 0.5;
871     }
872     layerPositions.reserve( layerPositions.size() * 2 );
873     for ( int nb = layerPositions.size()-1; nb > 0; --nb )
874       layerPositions.push_back( layerPositions.back() + layerPositions[nb] - layerPositions[nb-1] );
875
876     // merge existing nodes with new nodes at layerPositions into a UVPtStructVec
877     UVPtStructVec nodes;
878     if ( nbMeshedEdges != linSide1->NbEdges() )
879       makeMissingMesh( linSide1, nodes, layerPositions, myHelper );
880     else
881       nodes = linSide1->GetUVPtStruct();
882
883     // find a central node
884     size_t i = 0;
885     while ( nodes[ i ].normParam < 0.5 ) ++i;
886     if (( nodes[ i ].normParam - 0.5 ) > ( 0.5 - nodes[ i-1 ].normParam )) --i;
887
888     // distribute nodes between two linear sides
889     UVPtStructVec nodes2( nodes.rbegin(), nodes.rbegin() + nodes.size() - i );
890     nodes.resize( i + 1 );
891
892     linSide1 = StdMeshers_FaceSide::New( nodes );
893     linSide2 = StdMeshers_FaceSide::New( nodes2 );
894
895     centerNode = nodes.back().node;
896     centerUV   = nodes.back().UV();
897   }
898   // ------------------------------------------------------------------------------------------
899   else // nbSides == 3 
900   {
901     // one curve must be a part of ellipse and 2 other curves must be segments of line
902
903     int nbMeshedEdges1, nbMeshedEdges2;
904     vector< double > layerPositions2;
905     bool ok1 = computeLayerPositions( linSide1, layerPositions,  &nbMeshedEdges1 );
906     bool ok2 = computeLayerPositions( linSide2, layerPositions2, &nbMeshedEdges2 );
907     if ( !ok1 && !ok2 )
908       return false;
909
910     bool linSide1Computed = ( nbMeshedEdges1 == linSide1->NbEdges() );
911     bool linSide2Computed = ( nbMeshedEdges2 == linSide2->NbEdges() );
912
913     UVPtStructVec nodes;
914
915     if ( linSide1Computed && !linSide2Computed )
916     {
917       // use layer positions of linSide1 to mesh linSide2
918       makeMissingMesh( linSide2, nodes, layerPositions, myHelper );
919       linSide2 = StdMeshers_FaceSide::New( nodes );
920     }
921     else if ( linSide2Computed && !linSide1Computed )
922     {
923       // use layer positions of linSide2 to mesh linSide1
924       makeMissingMesh( linSide1, nodes, layerPositions2, myHelper );
925       linSide1 = StdMeshers_FaceSide::New( nodes );
926     }
927     else if ( !linSide2Computed && !linSide1Computed )
928     {
929       // use layer positions of a longer side to mesh the shorter side
930       vector< double >& lp =
931         ( linSide1->Length() > linSide2->Length() ) ? layerPositions : layerPositions2;
932
933       makeMissingMesh( linSide1, nodes, lp, myHelper );
934       linSide1 = StdMeshers_FaceSide::New( nodes );
935       makeMissingMesh( linSide2, nodes, lp, myHelper );
936       linSide2 = StdMeshers_FaceSide::New( nodes );
937     }
938
939     const UVPtStructVec& nodes2 = linSide2->GetUVPtStruct();
940     centerNode = nodes2.back().node;
941     centerUV   = nodes2.back().UV();
942   }
943
944   list< TopoDS_Edge >::iterator ee = emptyEdges.begin();
945   for ( ; ee != emptyEdges.end(); ++ee )
946     markEdgeAsComputedByMe( *ee, aMesh.GetSubMesh( F ));
947
948   circSide->GetUVPtStruct(); // let sides take into account just computed nodes
949   linSide1->GetUVPtStruct();
950   linSide2->GetUVPtStruct();
951
952   FaceQuadStruct::Ptr quad( new FaceQuadStruct );
953   quad->side.resize( 4 );
954   quad->side[0] = circSide;
955   quad->side[1] = linSide1;
956   quad->side[2] = StdMeshers_FaceSide::New( circSide.get(), centerNode, &centerUV );
957   quad->side[3] = linSide2;
958
959   myQuadList.push_back( quad );
960
961   // create quadrangles
962   bool ok;
963   if ( linSide1->NbPoints() == linSide2->NbPoints() )
964     ok = StdMeshers_Quadrangle_2D::computeQuadDominant( aMesh, F, quad );
965   else
966     ok = StdMeshers_Quadrangle_2D::computeTriangles( aMesh, F, quad );
967
968   StdMeshers_Quadrangle_2D::myHelper = 0;
969
970   return ok;
971 }
972
973 //================================================================================
974 /*!
975  * \brief Compute nodes on the radial edge
976  * \retval int - nb of segments on the linSide
977  */
978 //================================================================================
979
980 int StdMeshers_RadialQuadrangle_1D2D::computeLayerPositions(StdMeshers_FaceSidePtr linSide,
981                                                             vector< double > &     positions,
982                                                             int*                   nbMeshedEdges,
983                                                             bool                   useHalf)
984 {
985   // First, try to compute positions of layers
986
987   positions.clear();
988
989   SMESH_Mesh * mesh = myHelper->GetMesh();
990
991   const SMESH_Hypothesis* hyp1D = myDistributionHypo ? myDistributionHypo->GetLayerDistribution() : 0;
992   int                  nbLayers = myNbLayerHypo ? myNbLayerHypo->GetNumberOfLayers() : 0;
993
994   if ( !hyp1D && !nbLayers )
995   {
996     // No own algo hypotheses assigned, so first try to find any 1D hypothesis.
997     // find a hyp usable by TNodeDistributor
998     TopoDS_Shape edge = linSide->Edge(0);
999     const SMESH_HypoFilter* hypKind =
1000       TNodeDistributor::GetDistributor(*mesh)->GetCompatibleHypoFilter(/*ignoreAux=*/true);
1001     hyp1D = mesh->GetHypothesis( edge, *hypKind, /*fromAncestors=*/true);
1002   }
1003   if ( hyp1D ) // try to compute with hyp1D
1004   {
1005     BRepAdaptor_CompCurve* curve = linSide->GetCurve3d();
1006     SMESHUtils::Deleter< BRepAdaptor_CompCurve > delCurve( curve );
1007     double f = curve->FirstParameter();
1008     double l = curve->LastParameter();
1009
1010     if ( useHalf )
1011       l = 0.5 * ( f + l ); // first part of linSide is used
1012
1013     if ( !TNodeDistributor::GetDistributor(*mesh)->Compute( positions, linSide->Edge(0),
1014                                                             *curve, f, l, *mesh, hyp1D ))
1015     {
1016       if ( myDistributionHypo ) { // bad hyp assigned 
1017         return error( TNodeDistributor::GetDistributor(*mesh)->GetComputeError() );
1018       }
1019       else {
1020         // bad hyp found, its Ok, lets try with default nb of segments
1021       }
1022     }
1023   }
1024   
1025   if ( positions.empty() ) // try to use nb of layers
1026   {
1027     if ( !nbLayers )
1028       nbLayers = _gen->GetDefaultNbSegments();
1029
1030     if ( nbLayers )
1031     {
1032       positions.resize( nbLayers + 1 );
1033       for ( int z = 0; z < nbLayers; ++z )
1034         positions[ z ] = double( z )/ double( nbLayers );
1035       positions.back() = 1;
1036     }
1037   }
1038
1039   // Second, check presence of a mesh built by other algo on linSide
1040
1041   int nbEdgesComputed = 0;
1042   for ( int i = 0; i < linSide->NbEdges(); ++i )
1043   {
1044     nbEdgesComputed += ( !mesh->GetSubMesh( linSide->Edge(i))->IsEmpty() );
1045   }
1046
1047   if ( nbEdgesComputed == linSide->NbEdges() )
1048   {
1049     const UVPtStructVec& points = linSide->GetUVPtStruct();
1050     if ( points.size() >= 2 )
1051     {
1052       positions.resize( points.size() );
1053       for ( size_t i = 0; i < points.size(); ++i )
1054         positions[ i ] = points[i].normParam;
1055     }
1056   }
1057
1058   if ( nbMeshedEdges ) *nbMeshedEdges = nbEdgesComputed;
1059
1060   return positions.size();
1061 }
1062
1063
1064 //=======================================================================
1065 //function : Evaluate
1066 //purpose  : 
1067 //=======================================================================
1068
1069 bool StdMeshers_RadialQuadrangle_1D2D::Evaluate(SMESH_Mesh&         aMesh,
1070                                                 const TopoDS_Shape& aShape,
1071                                                 MapShapeNbElems&    aResMap)
1072 {
1073   if( aShape.ShapeType() != TopAbs_FACE ) {
1074     return false;
1075   }
1076   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
1077   if( aResMap.count(sm) )
1078     return false;
1079
1080   vector<int>& aResVec =
1081     aResMap.insert( make_pair(sm, vector<int>(SMDSEntity_Last,0))).first->second;
1082
1083   myHelper = new SMESH_MesherHelper( aMesh );
1084   myHelper->SetSubShape( aShape );
1085   SMESHUtils::Deleter<SMESH_MesherHelper> helperDeleter( myHelper );
1086
1087   TNodeDistributor* algo1d = TNodeDistributor::GetDistributor(aMesh);
1088
1089   StdMeshers_FaceSidePtr circSide, linSide1, linSide2;
1090   int nbSides = analyseFace( aShape, &aMesh, circSide, linSide1, linSide2 );
1091   if( nbSides > 3 || nbSides < 1 )
1092     return false;
1093
1094   if ( algo1d->EvaluateCircularEdge( aMesh, circSide, aResMap ))
1095     return false;
1096
1097   vector< double > layerPositions; // [0,1]
1098
1099   // ------------------------------------------------------------------------------------------
1100   if ( nbSides == 1 )
1101   {
1102     const TopoDS_Face& F = TopoDS::Face( aShape );
1103
1104     const SMDS_MeshNode* circNode;
1105     TopoDS_Edge linEdge = makeEdgeToCenter( circSide, F, circNode );
1106
1107     StdMeshers_FaceSidePtr tmpSide =
1108       StdMeshers_FaceSide::New( F, linEdge, &aMesh, /*isFrw=*/true, /*skipMedium=*/true );
1109
1110     if ( !computeLayerPositions( tmpSide, layerPositions ))
1111       return false;
1112   }
1113   // ------------------------------------------------------------------------------------------
1114   else if ( nbSides == 2 && linSide1->Edge(0).Orientation() == TopAbs_INTERNAL )
1115   {
1116     if ( !computeLayerPositions( linSide1, layerPositions ))
1117       return false;
1118   }
1119   // ------------------------------------------------------------------------------------------
1120   else if ( nbSides == 2 )
1121   {
1122     // find positions of layers for the first half of linSide1
1123     if ( !computeLayerPositions( linSide1, layerPositions, 0, /*useHalf=*/true ))
1124       return false;
1125   }
1126   // ------------------------------------------------------------------------------------------
1127   else // nbSides == 3 
1128   {
1129     if ( !computeLayerPositions(( linSide1->Length() > linSide2->Length() ) ? linSide1 : linSide2,
1130                                  layerPositions ))
1131       return false;
1132   }
1133
1134   bool isQuadratic = false;
1135   for ( TopExp_Explorer edge( aShape, TopAbs_EDGE ); edge.More() &&  !isQuadratic ; edge.Next() )
1136   {
1137     sm = aMesh.GetSubMesh( edge.Current() );
1138     vector<int>& nbElems = aResMap[ sm ];
1139     if ( SMDSEntity_Quad_Edge < (int) nbElems.size() )
1140       isQuadratic = nbElems[ SMDSEntity_Quad_Edge ];
1141   }
1142
1143   int nbCircSegments = 0;
1144   for ( int iE = 0; iE < circSide->NbEdges(); ++iE )
1145   {
1146     sm = aMesh.GetSubMesh( circSide->Edge( iE ));
1147     vector<int>& nbElems = aResMap[ sm ];
1148     if ( SMDSEntity_Quad_Edge < (int) nbElems.size() )
1149       nbCircSegments += ( nbElems[ SMDSEntity_Edge ] + nbElems[ SMDSEntity_Quad_Edge ]);
1150   }
1151   
1152   int nbQuads = nbCircSegments * ( layerPositions.size() - 1 );
1153   int nbTria  = nbCircSegments;
1154   int nbNodes = ( nbCircSegments - 1 ) * ( layerPositions.size() - 2 );
1155   if ( isQuadratic )
1156   {
1157     nbNodes += (( nbCircSegments - 1 ) * ( layerPositions.size() - 1 ) + // radial
1158                 ( nbCircSegments     ) * ( layerPositions.size() - 2 )); // circular
1159     aResVec[SMDSEntity_Quad_Triangle  ] = nbTria;
1160     aResVec[SMDSEntity_Quad_Quadrangle] = nbQuads;
1161   }
1162   else
1163   {
1164     aResVec[SMDSEntity_Triangle  ] = nbTria;
1165     aResVec[SMDSEntity_Quadrangle] = nbQuads;
1166   }
1167   aResVec[SMDSEntity_Node] = nbNodes;
1168
1169   if ( linSide1 )
1170   {
1171     // evaluation for linSides
1172     vector<int> aResVec(SMDSEntity_Last, 0);
1173     if ( isQuadratic ) {
1174       aResVec[SMDSEntity_Node     ] = 2 * ( layerPositions.size() - 1 ) + 1;
1175       aResVec[SMDSEntity_Quad_Edge] = layerPositions.size() - 1;
1176     }
1177     else {
1178       aResVec[SMDSEntity_Node] = layerPositions.size() - 2;
1179       aResVec[SMDSEntity_Edge] = layerPositions.size() - 1;
1180     }
1181     sm = aMesh.GetSubMesh( linSide1->Edge(0) );
1182     aResMap[sm] = aResVec;
1183     if ( linSide2 )
1184     {
1185       sm = aMesh.GetSubMesh( linSide2->Edge(0) );
1186       aResMap[sm] = aResVec;
1187     }
1188   }
1189
1190   return true;
1191 }
1192
1193 //================================================================================
1194 /*!
1195  * \brief Return true if the algorithm can compute mesh on this shape
1196  */
1197 //================================================================================
1198
1199 bool StdMeshers_RadialQuadrangle_1D2D::IsApplicable( const TopoDS_Shape & aShape, bool toCheckAll )
1200 {
1201   int nbFoundFaces = 0;
1202   for (TopExp_Explorer exp( aShape, TopAbs_FACE ); exp.More(); exp.Next(), ++nbFoundFaces )
1203   {
1204     StdMeshers_FaceSidePtr circSide, linSide1, linSide2;
1205     int nbSides = analyseFace( aShape, NULL, circSide, linSide1, linSide2 );
1206     bool ok = ( 0 < nbSides && nbSides <= 3 &&
1207                 isCornerInsideCircle( circSide, linSide1, linSide2 ));
1208     if( toCheckAll  && !ok ) return false;
1209     if( !toCheckAll && ok  ) return true;
1210   }
1211   if( toCheckAll && nbFoundFaces != 0 ) return true;
1212   return false;
1213 };