Salome HOME
Update copyrights 2014.
[modules/smesh.git] / src / StdMeshers / StdMeshers_Adaptive1D.cxx
1 // Copyright (C) 2007-2014  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_Adaptive1D.cxx
23 //  Module : SMESH
24 //
25 #include "StdMeshers_Adaptive1D.hxx"
26
27 #include "SMESH_Gen.hxx"
28 #include "SMESH_Mesh.hxx"
29 #include "SMESH_MesherHelper.hxx"
30 #include "SMESH_Octree.hxx"
31 #include "SMESH_subMesh.hxx"
32 #include "SMESH_HypoFilter.hxx"
33
34 #include <Utils_SALOME_Exception.hxx>
35
36 #include <BRepAdaptor_Curve.hxx>
37 #include <BRepAdaptor_Surface.hxx>
38 #include <BRepMesh_IncrementalMesh.hxx>
39 #include <BRep_Tool.hxx>
40 #include <Bnd_B3d.hxx>
41 #include <GCPnts_AbscissaPoint.hxx>
42 #include <GeomAdaptor_Curve.hxx>
43 #include <Geom_Curve.hxx>
44 #include <Poly_Array1OfTriangle.hxx>
45 #include <Poly_PolygonOnTriangulation.hxx>
46 #include <Poly_Triangulation.hxx>
47 #include <TColgp_Array1OfPnt.hxx>
48 #include <TopExp.hxx>
49 #include <TopExp_Explorer.hxx>
50 #include <TopLoc_Location.hxx>
51 #include <TopTools_IndexedMapOfShape.hxx>
52 #include <TopoDS.hxx>
53 #include <TopoDS_Edge.hxx>
54 #include <TopoDS_Face.hxx>
55 #include <TopoDS_Vertex.hxx>
56 #include <gp_Lin.hxx>
57 #include <gp_Pnt.hxx>
58
59 #include <limits>
60 #include <vector>
61 #include <set>
62
63 using namespace std;
64
65 namespace // internal utils
66 {
67   //================================================================================
68   /*!
69    * \brief Bnd_B3d with access to its center and half-size
70    */
71   struct BBox : public Bnd_B3d
72   {
73     gp_XYZ Center() const { return gp_XYZ( myCenter[0], myCenter[1], myCenter[2] ); }
74     gp_XYZ HSize()  const { return gp_XYZ( myHSize[0],  myHSize[1],  myHSize[2]  ); }
75     double Size()   const { return 2 * myHSize[0]; }
76   };
77   //================================================================================
78   /*!
79    * \brief Working data of an EDGE
80    */
81   struct EdgeData
82   {
83     struct ProbePnt
84     {
85       gp_Pnt myP;
86       double myU;
87       double mySegSize;
88       ProbePnt( gp_Pnt p, double u, double sz=1e100 ): myP( p ), myU( u ), mySegSize( sz ) {}
89     };
90     BRepAdaptor_Curve myC3d;
91     double            myLength;
92     list< ProbePnt >  myPoints;
93     BBox              myBBox;
94
95     typedef list< ProbePnt >::iterator TPntIter;
96     void AddPoint( TPntIter where, double u )
97     {
98       TPntIter it = myPoints.insert( where, ProbePnt( myC3d.Value( u ), u ));
99       myBBox.Add( it->myP.XYZ() );
100     }
101     const ProbePnt& First() const { return myPoints.front(); }
102     const ProbePnt& Last()  const { return myPoints.back(); }
103     const TopoDS_Edge& Edge() const { return myC3d.Edge(); }
104     bool IsTooDistant( const BBox& faceBox, double maxSegSize ) const
105     {
106       gp_XYZ hsize = myBBox.HSize() + gp_XYZ( maxSegSize, maxSegSize, maxSegSize );
107       return faceBox.IsOut ( Bnd_B3d( myBBox.Center(), hsize ));
108     }
109   };
110   //================================================================================
111   /*!
112    * \brief Octree of local segment size
113    */
114   class SegSizeTree : public SMESH_Octree
115   {
116     double mySegSize; // segment size
117
118     // structure holding some common parameters of SegSizeTree
119     struct _CommonData : public SMESH_TreeLimit
120     {
121       double myGrading, myMinSize, myMaxSize;
122     };
123     _CommonData* getData() const { return (_CommonData*) myLimit; }
124
125     SegSizeTree(double size): SMESH_Octree(), mySegSize(size)
126     {
127       allocateChildren();
128     }
129     void allocateChildren()
130     {
131       myChildren = new SMESH_Octree::TBaseTree*[nbChildren()];
132       for ( int i = 0; i < nbChildren(); ++i )
133         myChildren[i] = NULL;
134     }
135     virtual box_type* buildRootBox() { return 0; }
136     virtual SegSizeTree* newChild() const { return 0; }
137     virtual void buildChildrenData() {}
138
139   public:
140
141     SegSizeTree( Bnd_B3d & bb, double grading, double mixSize, double maxSize);
142     void   SetSize( const gp_Pnt& p, double size );
143     double SetSize( const gp_Pnt& p1, const gp_Pnt& p2 );
144     double GetSize( const gp_Pnt& p ) const;
145     const BBox* GetBox() const { return (BBox*) getBox(); }
146     double GetMinSize() { return getData()->myMinSize; }
147   };
148   //================================================================================
149   /*!
150    * \brief Adaptive wire discertizator.
151    */
152   class AdaptiveAlgo : public StdMeshers_Regular_1D
153   {
154   public:
155     AdaptiveAlgo(int hypId, int studyId, SMESH_Gen* gen);
156     virtual bool Compute(SMESH_Mesh & aMesh, const TopoDS_Shape & aShape );
157     virtual bool Evaluate(SMESH_Mesh &         theMesh,
158                           const TopoDS_Shape & theShape,
159                           MapShapeNbElems&     theResMap);
160     void SetHypothesis( const StdMeshers_Adaptive1D* hyp );
161   private:
162
163     bool makeSegments();
164
165     const StdMeshers_Adaptive1D* myHyp;
166     SMESH_Mesh*                  myMesh;
167     vector< EdgeData >           myEdges;
168     SegSizeTree*                 mySizeTree;
169   };
170
171   //================================================================================
172   /*!
173    * \brief Segment of Poly_PolygonOnTriangulation
174    */
175   struct Segment
176   {
177     gp_XYZ myPos, myDir;
178     double myLength;
179
180     void Init( const gp_Pnt& p1, const gp_Pnt& p2 )
181     {
182       myPos    = p1.XYZ();
183       myDir    = p2.XYZ() - p1.XYZ();
184       myLength = myDir.Modulus();
185       if ( myLength > std::numeric_limits<double>::min() )
186         myDir /= myLength;
187     }
188     bool Distance( const gp_Pnt& P, double& dist ) const // returns length of normal projection
189     {
190       gp_XYZ p = P.XYZ();
191       p.Subtract( myPos );
192       double proj = p.Dot( myDir );
193       if ( 0 < proj && proj < myLength )
194       {
195         p.Cross( myDir );
196         dist = p.Modulus();
197         return true;
198       }
199       return false;
200     }
201   };
202   //================================================================================
203   /*!
204    * \brief Data of triangle used to locate it in an octree and to find distance
205    *        to a point
206    */
207   struct Triangle
208   {
209     Bnd_B3d  myBox;
210     bool     myIsChecked; // to mark treated trias instead of using std::set
211     bool     myHasNodeOnVertex;
212     Segment* mySegments[3];
213     // data for DistToProjection()
214     gp_XYZ   myN0, myEdge1, myEdge2, myNorm, myPVec;
215     double   myInvDet, myMaxSize2;
216
217     void Init( const gp_Pnt& n1, const gp_Pnt& n2, const gp_Pnt& n3 ); 
218     bool DistToProjection( const gp_Pnt& p, double& dist ) const;
219     bool DistToSegment   ( const gp_Pnt& p, double& dist ) const;
220   };
221   //================================================================================
222   /*!
223    * \brief Element data held by ElementBndBoxTree + algorithm computing a distance
224    *        from a point to element
225    */
226   class ElementBndBoxTree;
227   struct ElemTreeData : public SMESH_TreeLimit
228   {
229     vector< int >                myWorkIDs[8];// to speed up filling ElementBndBoxTree::_elementIDs
230     virtual const Bnd_B3d* GetBox(int elemID) const = 0;
231   };
232   struct TriaTreeData : public ElemTreeData
233   {
234     vector< Triangle >           myTrias;
235     vector< Segment >            mySegments;
236     double                       myFaceTol;
237     double                       myTriasDeflection;
238     BBox                         myBBox;
239     BRepAdaptor_Surface          mySurface;
240     ElementBndBoxTree*           myTree;
241     const Poly_Array1OfTriangle* myPolyTrias;
242     const TColgp_Array1OfPnt*    myNodes;
243     bool                         myOwnNodes;
244
245     typedef vector<int> IntVec;
246     IntVec                       myFoundTriaIDs;
247
248     TriaTreeData( const TopoDS_Face& face, ElementBndBoxTree* triaTree );
249     ~TriaTreeData() { if ( myOwnNodes ) delete myNodes; myNodes = NULL; }
250     virtual const Bnd_B3d* GetBox(int elemID) const { return &myTrias[elemID].myBox; }
251     void PrepareToTriaSearch();
252     void SetSizeByTrias( SegSizeTree& sizeTree, double deflection ) const;
253     double GetMinDistInSphere(const gp_Pnt& p,
254                               const double  radius,
255                               const bool    projectedOnly,
256                               const gp_Pnt* avoidP=0) const;
257   };
258   //================================================================================
259   /*!
260    * \brief Octree of triangles or segments
261    */
262   class ElementBndBoxTree : public SMESH_Octree
263   {
264   public:
265     ElementBndBoxTree(const TopoDS_Face& face);
266     void GetElementsInSphere( const gp_XYZ& center,
267                               const double  radius, vector<int> & foundElemIDs) const;
268     void FillIn();
269     ElemTreeData* GetElemData() const { return (ElemTreeData*) myLimit; }
270     TriaTreeData* GetTriaData() const { return (TriaTreeData*) myLimit; }
271
272   protected:
273     ElementBndBoxTree() {}
274     SMESH_Octree* newChild() const { return new ElementBndBoxTree; }
275     void          buildChildrenData();
276     Bnd_B3d*      buildRootBox();
277   private:
278     vector< int > _elementIDs;
279   };
280   //================================================================================
281   /*!
282    * \brief BRepMesh_IncrementalMesh with access to its protected Bnd_Box
283    */
284   struct IncrementalMesh : public BRepMesh_IncrementalMesh
285   {
286     IncrementalMesh(const TopoDS_Shape& shape,
287                     const Standard_Real deflection,
288                     const bool          relative):
289       BRepMesh_IncrementalMesh( shape, deflection, relative )
290     {
291     }
292     Bnd_B3d GetBox() const
293     {
294       Standard_Real TXmin, TYmin, TZmin, TXmax, TYmax, TZmax;
295       myBox.Get(TXmin, TYmin, TZmin, TXmax, TYmax, TZmax);
296       Bnd_B3d bb;
297       bb.Add( gp_XYZ( TXmin, TYmin, TZmin ));
298       bb.Add( gp_XYZ( TXmax, TYmax, TZmax ));
299       return bb;
300     }
301   };
302   //================================================================================
303   /*!
304    * \brief Link of two nodes
305    */
306   struct NLink : public std::pair< int, int >
307   {
308     NLink( int n1, int n2 )
309     {
310       if ( n1 < n2 )
311       {
312         first  = n1;
313         second = n2;
314       }
315       else
316       {
317         first  = n2;
318         second = n1;
319       }
320     }
321     int N1() const { return first; }
322     int N2() const { return second; }
323   };
324
325   //================================================================================
326   /*!
327    * \brief Initialize TriaTreeData
328    */
329   //================================================================================
330
331   TriaTreeData::TriaTreeData( const TopoDS_Face& face, ElementBndBoxTree* triaTree )
332     : myTriasDeflection(0), mySurface( face ),
333       myTree(NULL), myPolyTrias(NULL), myNodes(NULL), myOwnNodes(false)
334   {
335     TopLoc_Location loc;
336     Handle(Poly_Triangulation) tr = BRep_Tool::Triangulation( face, loc );
337     if ( !tr.IsNull() )
338     {
339       myFaceTol         = SMESH_MesherHelper::MaxTolerance( face );
340       myTree            = triaTree;
341       myNodes           = & tr->Nodes();
342       myPolyTrias       = & tr->Triangles();
343       myTriasDeflection = tr->Deflection();
344       if ( !loc.IsIdentity() ) // transform nodes if necessary
345       {
346         TColgp_Array1OfPnt* trsfNodes = new TColgp_Array1OfPnt( myNodes->Lower(), myNodes->Upper() );
347         trsfNodes->Assign( *myNodes );
348         myNodes    = trsfNodes;
349         myOwnNodes = true;
350         const gp_Trsf& trsf = loc;
351         for ( int i = trsfNodes->Lower(); i <= trsfNodes->Upper(); ++i )
352           trsfNodes->ChangeValue(i).Transform( trsf );
353       }
354       for ( int i = myNodes->Lower(); i <= myNodes->Upper(); ++i )
355         myBBox.Add( myNodes->Value(i).XYZ() );
356     }
357   }
358   //================================================================================
359   /*!
360    * \brief Prepare data for search of trinagles in GetMinDistInSphere()
361    */
362   //================================================================================
363
364   void TriaTreeData::PrepareToTriaSearch()
365   {
366     if ( !myTrias.empty() ) return; // already done
367     if ( !myPolyTrias ) return;
368
369     // get all boundary links and nodes on VERTEXes
370     map< NLink, Segment* > linkToSegMap;
371     map< NLink, Segment* >::iterator l2s;
372     set< int > vertexNodes;
373     TopLoc_Location loc;
374     Handle(Poly_Triangulation) tr = BRep_Tool::Triangulation( mySurface.Face(), loc );
375     if ( !tr.IsNull() )
376     {
377       TopTools_IndexedMapOfShape edgeMap;
378       TopExp::MapShapes( mySurface.Face(), TopAbs_EDGE, edgeMap );
379       for ( int iE = 1; iE <= edgeMap.Extent(); ++iE )
380       {
381         const TopoDS_Edge& edge = TopoDS::Edge( edgeMap( iE ));
382         Handle(Poly_PolygonOnTriangulation) polygon =
383           BRep_Tool::PolygonOnTriangulation( edge, tr, loc );
384         if ( polygon.IsNull()  )
385           continue;
386         const TColStd_Array1OfInteger& nodes = polygon->Nodes();
387         for ( int i = nodes.Lower(); i < nodes.Upper(); ++i )
388           linkToSegMap.insert( make_pair( NLink( nodes(i), nodes(i+1)), (Segment*)0 ));
389         vertexNodes.insert( nodes( nodes.Lower()));
390         vertexNodes.insert( nodes( nodes.Upper()));
391       }
392       // fill mySegments by boundary links
393       mySegments.resize( linkToSegMap.size() );
394       int iS = 0;
395       for ( l2s = linkToSegMap.begin(); l2s != linkToSegMap.end(); ++l2s, ++iS )
396       {
397         const NLink& link = (*l2s).first;
398         (*l2s).second = & mySegments[ iS ];
399         mySegments[ iS ].Init( myNodes->Value( link.N1() ),
400                                myNodes->Value( link.N2() ));
401       }
402     }
403
404     // initialize myTrias
405     myTrias.resize( myPolyTrias->Length() );
406     Standard_Integer n1,n2,n3;
407     for ( int i = 1; i <= myPolyTrias->Upper(); ++i )
408     {
409       Triangle & t = myTrias[ i-1 ];
410       myPolyTrias->Value( i ).Get( n1,n2,n3 );
411       t.Init( myNodes->Value( n1 ),
412               myNodes->Value( n2 ),
413               myNodes->Value( n3 ));
414       int nbSeg = 0;
415       if (( l2s = linkToSegMap.find( NLink( n1, n2 ))) != linkToSegMap.end())
416         t.mySegments[ nbSeg++ ] = l2s->second;
417       if (( l2s = linkToSegMap.find( NLink( n2, n3 ))) != linkToSegMap.end())
418         t.mySegments[ nbSeg++ ] = l2s->second;
419       if (( l2s = linkToSegMap.find( NLink( n3, n1 ))) != linkToSegMap.end())
420         t.mySegments[ nbSeg++ ] = l2s->second;
421       while ( nbSeg < 3 )
422         t.mySegments[ nbSeg++ ] = NULL;
423
424       t.myIsChecked = false;
425       t.myHasNodeOnVertex = ( vertexNodes.count( n1 ) ||
426                               vertexNodes.count( n2 ) ||
427                               vertexNodes.count( n3 ));
428     }
429
430     // fill the tree of triangles
431     myTree->FillIn();
432   }
433
434   //================================================================================
435   /*!
436    * \brief Set size of segments by size of triangles
437    */
438   //================================================================================
439
440   void TriaTreeData::SetSizeByTrias( SegSizeTree& sizeTree, double hypDeflection ) const
441   {
442     if ( mySurface.GetType() == GeomAbs_Plane ||
443          myTriasDeflection   <= 1e-100 )
444       return;
445     const double factor = hypDeflection / myTriasDeflection;
446
447     bool isConstSize;
448     switch( mySurface.GetType() ) {
449     case GeomAbs_Cylinder:
450     case GeomAbs_Sphere:
451     case GeomAbs_Torus:
452       isConstSize = true; break;
453     default:
454       isConstSize = false;
455     }
456
457     map< NLink, double >           lenOfDoneLink;
458     map< NLink, double >::iterator link2len;
459
460     Standard_Integer n[4];
461     gp_Pnt p[4];
462     double a[3];
463     bool   isDone[3];
464     double size = -1., maxLinkLen;
465     int    jLongest;
466
467     //int nbLinks = 0;
468     for ( int i = 1; i <= myPolyTrias->Upper(); ++i )
469     {
470       // get corners of a triangle
471       myPolyTrias->Value( i ).Get( n[0],n[1],n[2] );
472       n[3] = n[0];
473       p[0] = myNodes->Value( n[0] );
474       p[1] = myNodes->Value( n[1] );
475       p[2] = myNodes->Value( n[2] );
476       p[3] = p[0];
477       // get length of links and find the longest one
478       maxLinkLen = 0;
479       for ( int j = 0; j < 3; ++j )
480       {
481         link2len  = lenOfDoneLink.insert( make_pair( NLink( n[j], n[j+1] ), -1. )).first;
482         isDone[j] = !((*link2len).second < 0 );
483         a[j]      = isDone[j] ? (*link2len).second : (*link2len).second = p[j].Distance( p[j+1] );
484         if ( isDone[j] )
485           lenOfDoneLink.erase( link2len );
486         if ( a[j] > maxLinkLen )
487         {
488           maxLinkLen = a[j];
489           jLongest   = j;
490         }
491       }
492       // compute minimal altitude of a triangle
493       if ( !isConstSize || size < 0. )
494       {
495         double s    = 0.5 * ( a[0] + a[1] + a[2] );
496         double area = sqrt( s * (s - a[0]) * (s - a[1]) * (s - a[2]));
497         size        = 2 * area / maxLinkLen; // minimal altitude
498       }
499       // set size to the size tree
500       if ( !isDone[ jLongest ] || !isConstSize )
501       {
502         //++nbLinks;
503         if ( size < numeric_limits<double>::min() )
504           continue;
505         int nb = Max( 1, int( maxLinkLen / size / 2 ));
506         for ( int k = 0; k <= nb; ++k )
507         {
508           double r = double( k ) / nb;
509           sizeTree.SetSize( r * p[ jLongest ].XYZ() + ( 1-r ) * p[ jLongest+1 ].XYZ(),
510                             size * factor );
511         }
512       }
513       //cout << "SetSizeByTrias, i="<< i << " " << sz * factor << endl;
514     }
515     // cout << "SetSizeByTrias, nn tria="<< myPolyTrias->Upper()
516     //      << " nb links" << nbLinks << " isConstSize="<<isConstSize
517     //      << " " << size * factor << endl;
518   }
519   //================================================================================
520   /*!
521    * \brief Return minimal distance from a given point to a trinangle but not more
522    *        distant than a given radius. Triangles with a node at avoidPnt are ignored.
523    *        If projectedOnly, 
524    */
525   //================================================================================
526
527   double TriaTreeData::GetMinDistInSphere(const gp_Pnt& p,
528                                           const double  radius,
529                                           const bool    projectedOnly,
530                                           const gp_Pnt* avoidPnt) const
531   {
532     double minDist2 = 1e100;
533     const double tol2 = myFaceTol * myFaceTol;
534     const double dMin2 = myTriasDeflection * myTriasDeflection;
535
536     TriaTreeData* me = const_cast<TriaTreeData*>( this );
537     me->myFoundTriaIDs.clear();
538     myTree->GetElementsInSphere( p.XYZ(), radius, me->myFoundTriaIDs );
539     if ( myFoundTriaIDs.empty() )
540       return minDist2;
541
542     Standard_Integer n[ 3 ];
543     for ( size_t i = 0; i < myFoundTriaIDs.size(); ++i )
544     {
545       Triangle& t = me->myTrias[ myFoundTriaIDs[i] ];
546       if ( t.myIsChecked )
547         continue;
548       t.myIsChecked = true;
549
550       double d, minD2 = minDist2;
551       myPolyTrias->Value( myFoundTriaIDs[i]+1 ).Get( n[0],n[1],n[2] );
552       if ( avoidPnt && t.myHasNodeOnVertex )
553       {
554         bool avoidTria = false;
555         for ( int i = 0; i < 3; ++i )
556         {
557           const gp_Pnt& pn = myNodes->Value(n[i]);
558           if ( avoidTria = ( pn.SquareDistance( *avoidPnt ) <= tol2 ))
559             break;
560           if ( !projectedOnly )
561             minD2 = Min( minD2, pn.SquareDistance( p ));
562         }
563         if ( avoidTria )
564           continue;
565         if (( projectedOnly || minD2 < t.myMaxSize2 ) &&
566             ( t.DistToProjection( p, d ) || t.DistToSegment( p, d )))
567           minD2 = Min( minD2, d*d );
568         minDist2 = Min( minDist2, minD2 );
569       }
570       else if ( projectedOnly )
571       {
572         if ( t.DistToProjection( p, d ) && d*d > dMin2 )
573           minDist2 = Min( minDist2, d*d );
574       }
575       else
576       {
577         for ( int i = 0; i < 3; ++i )
578           minD2 = Min( minD2, p.SquareDistance( myNodes->Value(n[i]) ));
579         if ( minD2 < t.myMaxSize2  && ( t.DistToProjection( p, d ) || t.DistToSegment( p, d )))
580           minD2 = Min( minD2, d*d );
581         minDist2 = Min( minDist2, minD2 );
582       }
583     }
584
585     for ( size_t i = 0; i < myFoundTriaIDs.size(); ++i )
586       me->myTrias[ myFoundTriaIDs[i] ].myIsChecked = false;
587
588     return sqrt( minDist2 );
589   }
590   //================================================================================
591   /*!
592    * \brief Prepare Triangle data
593    */
594   //================================================================================
595
596   void Triangle::Init( const gp_Pnt& p1, const gp_Pnt& p2, const gp_Pnt& p3 )
597   {
598     myBox.Add( p1 );
599     myBox.Add( p2 );
600     myBox.Add( p3 );
601     myN0 = p1.XYZ();
602     myEdge1 = p2.XYZ() - myN0;
603     myEdge2 = p3.XYZ() - myN0;
604     myNorm = myEdge1 ^ myEdge2;
605     double normSize = myNorm.Modulus();
606     if ( normSize > std::numeric_limits<double>::min() )
607     {
608       myNorm /= normSize;
609       myPVec = myNorm ^ myEdge2;
610       myInvDet = 1. / ( myEdge1 * myPVec );
611     }
612     else
613     {
614       myInvDet = 0.;
615     }
616     myMaxSize2 = Max( p2.SquareDistance( p3 ),
617                       Max( myEdge2.SquareModulus(), myEdge1.SquareModulus() ));
618   }
619   //================================================================================
620   /*!
621    * \brief Compute distance from a point to the triangle. Return false if the point
622    *        is not projected inside the triangle
623    */
624   //================================================================================
625
626   bool Triangle::DistToProjection( const gp_Pnt& p, double& dist ) const
627   {
628     if ( myInvDet == 0 )
629       return false; // degenerated triangle
630
631     /* distance from n0 to the point */
632     gp_XYZ tvec = p.XYZ() - myN0;
633
634     /* calculate U parameter and test bounds */
635     double u = ( tvec * myPVec ) * myInvDet;
636     if (u < 0.0 || u > 1.0)
637       return false; // projected outside the triangle
638
639     /* calculate V parameter and test bounds */
640     gp_XYZ qvec = tvec ^ myEdge1;
641     double v = ( myNorm * qvec) * myInvDet;
642     if ( v < 0.0 || u + v > 1.0 )
643       return false; // projected outside the triangle
644
645     dist = ( myEdge2 * qvec ) * myInvDet;
646     return true;
647   }
648
649   //================================================================================
650   /*!
651    * \brief Compute distance from a point to either of mySegments. Return false if the point
652    *        is not projected on a segment
653    */
654   //================================================================================
655
656   bool Triangle::DistToSegment( const gp_Pnt& p, double& dist ) const
657   {
658     dist = 1e100;
659     bool res = false;
660     double d;
661     for ( int i = 0; i < 3; ++i )
662     {
663       if ( !mySegments[ i ])
664         break;
665       if ( mySegments[ i ]->Distance( p, d ))
666       {
667         res = true;
668         dist = Min( dist, d );
669       }
670     }
671     return res;
672   }
673
674   //================================================================================
675   /*!
676    * \brief Consturct ElementBndBoxTree of Poly_Triangulation of a FACE
677    */
678   //================================================================================
679
680   ElementBndBoxTree::ElementBndBoxTree(const TopoDS_Face& face)
681     :SMESH_Octree()
682   {
683     TriaTreeData* data = new TriaTreeData( face, this );
684     data->myMaxLevel = 5;
685     myLimit = data;
686   }
687   //================================================================================
688   /*!
689    * \brief Fill all levels of octree of Poly_Triangulation of a FACE
690    */
691   //================================================================================
692
693   void ElementBndBoxTree::FillIn()
694   {
695     if ( myChildren ) return;
696     TriaTreeData* data = GetTriaData();
697     if ( !data->myTrias.empty() )
698     {
699       for ( size_t i = 0; i < data->myTrias.size(); ++i )
700         _elementIDs.push_back( i );
701
702       compute();
703     }
704   }
705   //================================================================================
706   /*!
707    * \brief Return the maximal box
708    */
709   //================================================================================
710
711   Bnd_B3d* ElementBndBoxTree::buildRootBox()
712   {
713     TriaTreeData* data = GetTriaData();
714     Bnd_B3d*       box = new Bnd_B3d( data->myBBox );
715     return box;
716   }
717   //================================================================================
718   /*!
719    * \brief Redistrubute element boxes among children
720    */
721   //================================================================================
722
723   void ElementBndBoxTree::buildChildrenData()
724   {
725     ElemTreeData* data = GetElemData();
726     for ( int i = 0; i < _elementIDs.size(); ++i )
727     {
728       const Bnd_B3d* elemBox = data->GetBox( _elementIDs[i] );
729       for (int j = 0; j < 8; j++)
730         if ( !elemBox->IsOut( *myChildren[ j ]->getBox() ))
731           data->myWorkIDs[ j ].push_back( _elementIDs[i] );
732     }
733     SMESHUtils::FreeVector( _elementIDs ); // = _elements.clear() + free memory
734
735     const int theMaxNbElemsInLeaf = 7;
736
737     for (int j = 0; j < 8; j++)
738     {
739       ElementBndBoxTree* child = static_cast<ElementBndBoxTree*>( myChildren[j] );
740       child->_elementIDs = data->myWorkIDs[ j ];
741       if ( child->_elementIDs.size() <= theMaxNbElemsInLeaf )
742         child->myIsLeaf = true;
743       data->myWorkIDs[ j ].clear();
744     }
745   }
746   //================================================================================
747   /*!
748    * \brief Return elements from leaves intersecting the sphere
749    */
750   //================================================================================
751
752   void ElementBndBoxTree::GetElementsInSphere( const gp_XYZ& center,
753                                                const double  radius,
754                                                vector<int> & foundElemIDs) const
755   {
756     if ( const box_type* box = getBox() )
757     {
758       if ( box->IsOut( center, radius ))
759         return;
760
761       if ( isLeaf() )
762       {
763         ElemTreeData* data = GetElemData();
764         for ( int i = 0; i < _elementIDs.size(); ++i )
765           if ( !data->GetBox( _elementIDs[i] )->IsOut( center, radius ))
766             foundElemIDs.push_back( _elementIDs[i] );
767       }
768       else
769       {
770         for (int i = 0; i < 8; i++)
771           ((ElementBndBoxTree*) myChildren[i])->GetElementsInSphere( center, radius, foundElemIDs );
772       }
773     }
774   }
775
776   //================================================================================
777   /*!
778    * \brief Constructor of SegSizeTree
779    *  \param [in,out] bb - bounding box enclosing all EDGEs to discretize 
780    *  \param [in] grading - factor to get max size of the neighbour segment by 
781    *         size of a current one.
782    */
783   //================================================================================
784
785   SegSizeTree::SegSizeTree( Bnd_B3d & bb, double grading, double minSize, double maxSize )
786     : SMESH_Octree( new _CommonData() )
787   {
788     // make cube myBox from the box bb
789     gp_XYZ pmin = bb.CornerMin(), pmax = bb.CornerMax();
790     double maxBoxHSize = 0.5 * Max( pmax.X()-pmin.X(), Max( pmax.Y()-pmin.Y(), pmax.Z()-pmin.Z() ));
791     maxBoxHSize *= 1.01;
792     bb.SetHSize( gp_XYZ( maxBoxHSize, maxBoxHSize, maxBoxHSize ));
793     myBox = new box_type( bb );
794
795     mySegSize = Min( 2 * maxBoxHSize, maxSize );
796
797     getData()->myGrading = grading;
798     getData()->myMinSize = Max( minSize, 2*maxBoxHSize / 1.e6 );
799     getData()->myMaxSize = maxSize;
800     allocateChildren();
801   }
802
803   //================================================================================
804   /*!
805    * \brief Set segment size at a given point
806    */
807   //================================================================================
808
809   void SegSizeTree::SetSize( const gp_Pnt& p, double size )
810   {
811     // check if the point is out of the largest cube
812     SegSizeTree* root = this;
813     while ( root->myFather )
814       root = (SegSizeTree*) root->myFather;
815     if ( root->getBox()->IsOut( p.XYZ() ))
816       return;
817
818     // keep size whthin the valid range
819     size = Max( size, getData()->myMinSize );
820     //size = Min( size, getData()->myMaxSize );
821
822     // find an existing leaf at the point
823     SegSizeTree* leaf = (SegSizeTree*) root;
824     int iChild;
825     while ( true )
826     {
827       iChild = SMESH_Octree::getChildIndex( p.X(), p.Y(), p.Z(), leaf->GetBox()->Center() );
828       if ( leaf->myChildren[ iChild ] )
829         leaf = (SegSizeTree*) leaf->myChildren[ iChild ];
830       else
831         break;
832     }
833     // don't increase the current size
834     if ( leaf->mySegSize <= 1.1 * size )
835       return;
836
837     // split the found leaf until its box size is less than the given size
838     const double rootSize = root->GetBox()->Size();
839     while ( leaf->GetBox()->Size() > size )
840     {
841       const BBox* bb = leaf->GetBox();
842       iChild   = SMESH_Octree::getChildIndex( p.X(), p.Y(), p.Z(), bb->Center() );
843       SegSizeTree* newLeaf = new SegSizeTree( bb->Size() / 2 );
844       leaf->myChildren[iChild] = newLeaf;
845       newLeaf->myFather = leaf;
846       newLeaf->myLimit  = leaf->myLimit;
847       newLeaf->myLevel  = leaf->myLevel + 1;
848       newLeaf->myBox    = leaf->newChildBox( iChild );
849       newLeaf->myBox->Enlarge( rootSize * 1e-10 );
850       //newLeaf->myIsLeaf = ( newLeaf->mySegSize <= size );
851       leaf = newLeaf;
852     }
853     leaf->mySegSize = size;
854
855     // propagate increased size out from the leaf
856     double boxSize = leaf->GetBox()->Size();
857     double sizeInc = size + boxSize * getData()->myGrading;
858     for ( int iDir = 1; iDir <= 3; ++iDir )
859     {
860       gp_Pnt outPnt = p;
861       outPnt.SetCoord( iDir, p.Coord( iDir ) + boxSize );
862       SetSize( outPnt, sizeInc );
863       outPnt.SetCoord( iDir, p.Coord( iDir ) - boxSize );
864       SetSize( outPnt, sizeInc );
865     }
866   }
867   //================================================================================
868   /*!
869    * \brief Set size of a segment given by two end points
870    */
871   //================================================================================
872
873   double SegSizeTree::SetSize( const gp_Pnt& p1, const gp_Pnt& p2 )
874   {
875     const double size = p1.Distance( p2 );
876     gp_XYZ p = 0.5 * ( p1.XYZ() + p2.XYZ() );
877     SetSize( p, size );
878     SetSize( p1, size );
879     SetSize( p2, size );
880     //cout << "SetSize " << p1.Distance( p2 ) << " at " << p.X() <<", "<< p.Y()<<", "<<p.Z()<< endl;
881     return GetSize( p );
882   }
883
884   //================================================================================
885   /*!
886    * \brief Return segment size at a point
887    */
888   //================================================================================
889
890   double SegSizeTree::GetSize( const gp_Pnt& p ) const
891   {
892     const SegSizeTree* leaf = this;
893     while ( true )
894     {
895       int iChild = SMESH_Octree::getChildIndex( p.X(), p.Y(), p.Z(), leaf->GetBox()->Center() );
896       if ( leaf->myChildren[ iChild ] )
897         leaf = (SegSizeTree*) leaf->myChildren[ iChild ];
898       else
899         return leaf->mySegSize;
900     }
901     return mySegSize; // just to return anything
902   }
903
904   //================================================================================
905   /*!
906    * \brief Evaluate curve deflection between two points
907    * \param theCurve - the curve
908    * \param theU1 - the parameter of the first point
909    * \param theU2 - the parameter of the second point
910    * \retval double - square deflection value
911    */
912   //================================================================================
913
914   double deflection2(const BRepAdaptor_Curve & theCurve,
915                      double                    theU1,
916                      double                    theU2)
917   {
918     // line between theU1 and theU2
919     gp_Pnt p1 = theCurve.Value( theU1 ), p2 = theCurve.Value( theU2 );
920     gp_Lin segment( p1, gp_Vec( p1, p2 ));
921
922     // evaluate square distance of theCurve from the segment
923     Standard_Real dist2 = 0;
924     const int nbPnt = 5;
925     const double step = ( theU2 - theU1 ) / nbPnt;
926     while (( theU1 += step ) < theU2 )
927       dist2 = Max( dist2, segment.SquareDistance( theCurve.Value( theU1 )));
928
929     return dist2;
930   }
931
932 } // namespace
933
934 //=======================================================================
935 //function : StdMeshers_Adaptive1D
936 //purpose  : Constructor
937 StdMeshers_Adaptive1D::StdMeshers_Adaptive1D(int         hypId,
938                                              int         studyId,
939                                              SMESH_Gen * gen)
940   :SMESH_Hypothesis(hypId, studyId, gen)
941 {
942   myMinSize       = 1e-10;
943   myMaxSize       = 1e+10;
944   myDeflection    = 1e-2;
945   myAlgo          = NULL;
946   _name           = "Adaptive1D";
947   _param_algo_dim = 1; // is used by SMESH_Regular_1D
948 }
949 //=======================================================================
950 //function : ~StdMeshers_Adaptive1D
951 //purpose  : Destructor
952 StdMeshers_Adaptive1D::~StdMeshers_Adaptive1D()
953 {
954   delete myAlgo; myAlgo = NULL;
955 }
956 //=======================================================================
957 //function : SetDeflection
958 //purpose  : 
959 void StdMeshers_Adaptive1D::SetDeflection(double value)
960   throw(SALOME_Exception)
961 {
962   if (value <= std::numeric_limits<double>::min() )
963     throw SALOME_Exception("Deflection must be greater that zero");
964   if (myDeflection != value)
965   {
966     myDeflection = value;
967     NotifySubMeshesHypothesisModification();
968   }
969 }
970 //=======================================================================
971 //function : SetMinSize
972 //purpose  : Sets minimal allowed segment length
973 void StdMeshers_Adaptive1D::SetMinSize(double minSize)
974   throw(SALOME_Exception)
975 {
976   if (minSize <= std::numeric_limits<double>::min() )
977     throw SALOME_Exception("Min size must be greater that zero");
978
979   if (myMinSize != minSize )
980   {
981     myMinSize = minSize;
982     NotifySubMeshesHypothesisModification();
983   }
984 }
985 //=======================================================================
986 //function : SetMaxSize
987 //purpose  : Sets maximal allowed segment length
988 void StdMeshers_Adaptive1D::SetMaxSize(double maxSize)
989   throw(SALOME_Exception)
990 {
991   if (maxSize <= std::numeric_limits<double>::min() )
992     throw SALOME_Exception("Max size must be greater that zero");
993
994   if (myMaxSize != maxSize )
995   {
996     myMaxSize = maxSize;
997     NotifySubMeshesHypothesisModification();
998   }
999 }
1000 //=======================================================================
1001 //function : SaveTo
1002 //purpose  : Persistence
1003 ostream & StdMeshers_Adaptive1D::SaveTo(ostream & save)
1004 {
1005   save << myMinSize << " " << myMaxSize << " " << myDeflection;
1006   save << " " << -1 << " " << -1; // preview addition of parameters
1007   return save;
1008 }
1009 //=======================================================================
1010 //function : LoadFrom
1011 //purpose  : Persistence
1012 istream & StdMeshers_Adaptive1D::LoadFrom(istream & load)
1013 {
1014   int dummyParam;
1015   bool isOK = (load >> myMinSize >> myMaxSize >> myDeflection >> dummyParam >> dummyParam);
1016   if (!isOK)
1017     load.clear(ios::badbit | load.rdstate());
1018   return load;
1019 }
1020 //=======================================================================
1021 //function : SetParametersByMesh
1022 //purpose  : Initialize parameters by the mesh built on the geometry
1023 //param theMesh - the built mesh
1024 //param theShape - the geometry of interest
1025 //retval bool - true if parameter values have been successfully defined
1026 bool StdMeshers_Adaptive1D::SetParametersByMesh(const SMESH_Mesh*   theMesh,
1027                                                 const TopoDS_Shape& theShape)
1028 {
1029   if ( !theMesh || theShape.IsNull() )
1030     return false;
1031
1032   int nbEdges = 0;
1033   TopTools_IndexedMapOfShape edgeMap;
1034   TopExp::MapShapes( theShape, TopAbs_EDGE, edgeMap );
1035
1036   SMESH_MesherHelper helper( (SMESH_Mesh&) *theMesh );
1037   double minSz2 = 1e100, maxSz2 = 0, sz2, maxDefl2 = 0;
1038   for ( int iE = 1; iE <= edgeMap.Extent(); ++iE )
1039   {
1040     const TopoDS_Edge& edge = TopoDS::Edge( edgeMap( iE ));
1041     SMESHDS_SubMesh* smDS = theMesh->GetMeshDS()->MeshElements( edge );
1042     if ( !smDS ) continue;
1043     ++nbEdges;
1044
1045     helper.SetSubShape( edge );
1046     BRepAdaptor_Curve curve( edge );
1047
1048     SMDS_ElemIteratorPtr segIt = smDS->GetElements();
1049     while ( segIt->more() )
1050     {
1051       const SMDS_MeshElement* seg = segIt->next();
1052       const SMDS_MeshNode* n1 = seg->GetNode(0);
1053       const SMDS_MeshNode* n2 = seg->GetNode(1);
1054       sz2 = SMESH_TNodeXYZ( n1 ).SquareDistance( n2 );
1055       minSz2 = Min( minSz2, sz2 );
1056       maxSz2 = Max( maxSz2, sz2 );
1057       if ( curve.GetType() != GeomAbs_Line )
1058       {
1059         double u1 = helper.GetNodeU( edge, n1, n2 );
1060         double u2 = helper.GetNodeU( edge, n2, n1 );
1061         maxDefl2 = Max( maxDefl2, deflection2( curve, u1, u2 ));
1062       }
1063     }
1064   }
1065   if ( nbEdges )
1066   {
1067     myMinSize = sqrt( minSz2 );
1068     myMaxSize = sqrt( maxSz2 );
1069     if ( maxDefl2 > 0 )
1070       myDeflection = maxDefl2;
1071   }
1072   return nbEdges;
1073 }
1074
1075 //=======================================================================
1076 //function : SetParametersByDefaults
1077 //purpose  : Initialize my parameter values by default parameters.
1078 //retval   : bool - true if parameter values have been successfully defined
1079 bool StdMeshers_Adaptive1D::SetParametersByDefaults(const TDefaults&  dflts,
1080                                                     const SMESH_Mesh* /*theMesh*/)
1081 {
1082   myMinSize = dflts._elemLength / 10;
1083   myMaxSize = dflts._elemLength * 2;
1084   myDeflection = myMinSize / 7;
1085   return true;
1086 }
1087
1088 //=======================================================================
1089 //function : GetAlgo
1090 //purpose  : Returns an algorithm that works using this hypothesis
1091 //=======================================================================
1092
1093 SMESH_Algo* StdMeshers_Adaptive1D::GetAlgo() const
1094 {
1095   if ( !myAlgo )
1096   {
1097     AdaptiveAlgo* newAlgo = 
1098       new AdaptiveAlgo( _gen->GetANewId(), _studyId, _gen );
1099     newAlgo->SetHypothesis( this );
1100
1101     ((StdMeshers_Adaptive1D*) this)->myAlgo = newAlgo;
1102   }
1103   return myAlgo;
1104 }
1105
1106 //================================================================================
1107 /*!
1108  * \brief Constructor
1109  */
1110 //================================================================================
1111
1112 AdaptiveAlgo::AdaptiveAlgo(int        hypId,
1113                            int        studyId,
1114                            SMESH_Gen* gen)
1115   : StdMeshers_Regular_1D( hypId, studyId, gen ),
1116     myHyp(NULL)
1117 {
1118   _name = "AdaptiveAlgo_1D";
1119 }
1120
1121 //================================================================================
1122 /*!
1123  * \brief Sets the hypothesis
1124  */
1125 //================================================================================
1126
1127 void AdaptiveAlgo::SetHypothesis( const StdMeshers_Adaptive1D* hyp )
1128 {
1129   myHyp = hyp;
1130 }
1131
1132 //================================================================================
1133 /*!
1134  * \brief Creates segments on all given EDGEs
1135  */
1136 //================================================================================
1137
1138 bool AdaptiveAlgo::Compute(SMESH_Mesh &         theMesh,
1139                            const TopoDS_Shape & theShape)
1140 {
1141   // *theProgress = 0.01;
1142
1143   if ( myHyp->GetMinSize() > myHyp->GetMaxSize() )
1144     return error( "Bad parameters: min size > max size" );
1145
1146   myMesh = &theMesh;
1147   SMESH_MesherHelper helper( theMesh );
1148   const double grading = 0.7;
1149
1150   TopTools_IndexedMapOfShape edgeMap, faceMap;
1151   TopExp::MapShapes( theShape,                 TopAbs_EDGE, edgeMap );
1152   TopExp::MapShapes( theMesh.GetShapeToMesh(), TopAbs_FACE, faceMap );
1153
1154   // Triangulate the shape with the given deflection ?????????
1155   Bnd_B3d box;
1156   {
1157     IncrementalMesh im( theMesh.GetShapeToMesh(), myHyp->GetDeflection(), /*Relatif=*/false);
1158     box = im.GetBox();
1159   }
1160   // *theProgress = 0.3;
1161
1162   // holder of segment size at each point
1163   SegSizeTree sizeTree( box, grading, myHyp->GetMinSize(), myHyp->GetMaxSize() );
1164   mySizeTree = & sizeTree;
1165
1166   // minimal segment size that sizeTree can store with reasonable tree height
1167   const double minSize = Max( myHyp->GetMinSize(), 1.1 * sizeTree.GetMinSize() );
1168
1169
1170   // fill myEdges - working data of EDGEs
1171   {
1172     // sort EDGEs by length
1173     multimap< double, TopoDS_Edge > edgeOfLength;
1174     for ( int iE = 1; iE <= edgeMap.Extent(); ++iE )
1175     {
1176       const TopoDS_Edge & edge = TopoDS::Edge( edgeMap( iE ));
1177       if ( !SMESH_Algo::isDegenerated( edge) )
1178         edgeOfLength.insert( make_pair( EdgeLength( edge ), edge ));
1179     }
1180     myEdges.clear();
1181     myEdges.resize( edgeOfLength.size() );
1182     multimap< double, TopoDS_Edge >::const_iterator len2edge = edgeOfLength.begin();
1183     for ( int iE = 0; len2edge != edgeOfLength.end(); ++len2edge, ++iE )
1184     {
1185       const TopoDS_Edge & edge = len2edge->second;
1186       EdgeData& eData = myEdges[ iE ];
1187       eData.myC3d.Initialize( edge );
1188       eData.myLength  = EdgeLength( edge );
1189       eData.AddPoint( eData.myPoints.end(), eData.myC3d.FirstParameter() );
1190       eData.AddPoint( eData.myPoints.end(), eData.myC3d.LastParameter() );
1191     }
1192   }
1193   if ( _computeCanceled ) return false;
1194
1195   // Take into account size of already existing segments
1196   SMDS_EdgeIteratorPtr segIterator = theMesh.GetMeshDS()->edgesIterator();
1197   while ( segIterator->more() )
1198   {
1199     const SMDS_MeshElement* seg = segIterator->next();
1200     sizeTree.SetSize( SMESH_TNodeXYZ( seg->GetNode( 0 )), SMESH_TNodeXYZ( seg->GetNode( 1 )));
1201   }
1202   if ( _computeCanceled ) return false;
1203
1204   // Set size of segments according to the deflection
1205
1206   StdMeshers_Regular_1D::_hypType = DEFLECTION;
1207   StdMeshers_Regular_1D::_value[ DEFLECTION_IND ] = myHyp->GetDeflection();
1208
1209   list< double > params;
1210   for ( int iE = 0; iE < myEdges.size(); ++iE )
1211   {
1212     EdgeData& eData = myEdges[ iE ];
1213     //cout << "E " << theMesh.GetMeshDS()->ShapeToIndex( eData.Edge() ) << endl;
1214
1215     double f = eData.First().myU, l = eData.Last().myU;
1216     if ( !computeInternalParameters( theMesh, eData.myC3d, eData.myLength, f,l, params, false, false ))
1217       continue;
1218     if ( params.size() <= 1 && helper.IsClosedEdge( eData.Edge() ) ) // 2 segments on a circle
1219     {
1220       params.clear();
1221       for ( int i = 1; i < 6; ++i )
1222         params.push_back(( l - f ) * i/6. );
1223     }
1224     EdgeData::TPntIter where = --eData.myPoints.end();
1225     list< double >::const_iterator param = params.begin();
1226     for ( ; param != params.end(); ++param )
1227       eData.AddPoint( where, *param );
1228
1229     EdgeData::TPntIter pIt2 = eData.myPoints.begin(), pIt1 = pIt2++;
1230     for ( ; pIt2 != eData.myPoints.end(); ++pIt1, ++pIt2 )
1231     {
1232       double sz = sizeTree.SetSize( (*pIt1).myP, (*pIt2).myP );
1233       sz = Min( sz, myHyp->GetMaxSize() );
1234       pIt1->mySegSize = Min( sz, pIt1->mySegSize );
1235       pIt2->mySegSize = Min( sz, pIt2->mySegSize );
1236     }
1237
1238     if ( _computeCanceled ) return false;
1239   }
1240
1241   // Limit size of segments according to distance to closest FACE
1242
1243   for ( int iF = 1; iF <= faceMap.Extent(); ++iF )
1244   {
1245     if ( _computeCanceled ) return false;
1246
1247     const TopoDS_Face & face = TopoDS::Face( faceMap( iF ));
1248     // cout << "FACE " << iF << "/" << faceMap.Extent()
1249     //      << " id-" << theMesh.GetMeshDS()->ShapeToIndex( face ) << endl;
1250
1251     ElementBndBoxTree triaTree( face ); // tree of FACE triangulation
1252     TriaTreeData*     triaSearcher = triaTree.GetTriaData();
1253
1254     triaSearcher->SetSizeByTrias( sizeTree, myHyp->GetDeflection() );
1255
1256     for ( int iE = 0; iE < myEdges.size(); ++iE )
1257     {
1258       EdgeData& eData = myEdges[ iE ];
1259
1260       // check if the face is in topological contact with the edge
1261       bool isAdjFace = ( helper.IsSubShape( helper.IthVertex( 0, eData.Edge()), face ) ||
1262                          helper.IsSubShape( helper.IthVertex( 1, eData.Edge()), face ));
1263
1264       if ( isAdjFace && triaSearcher->mySurface.GetType() == GeomAbs_Plane )
1265         continue;
1266
1267       bool sizeDecreased = true;
1268       for (int iLoop = 0; sizeDecreased; ++iLoop ) //repeat until segment size along the edge becomes stable
1269       {
1270         double maxSegSize = 0;
1271
1272         // get points to check distance to the face
1273         EdgeData::TPntIter pIt2 = eData.myPoints.begin(), pIt1 = pIt2++;
1274         maxSegSize = pIt1->mySegSize = Min( pIt1->mySegSize, sizeTree.GetSize( pIt1->myP ));
1275         for ( ; pIt2 != eData.myPoints.end(); )
1276         {
1277           pIt2->mySegSize = Min( pIt2->mySegSize, sizeTree.GetSize( pIt2->myP ));
1278           double curSize  = Min( pIt1->mySegSize, pIt2->mySegSize );
1279           maxSegSize      = Max( pIt2->mySegSize, maxSegSize );
1280           if ( pIt1->myP.Distance( pIt2->myP ) > curSize )
1281           {
1282             double midU  = 0.5*( pIt1->myU + pIt2->myU );
1283             gp_Pnt midP  = eData.myC3d.Value( midU );
1284             double midSz = sizeTree.GetSize( midP );
1285             pIt2 = eData.myPoints.insert( pIt2, EdgeData::ProbePnt( midP, midU, midSz ));
1286             eData.myBBox.Add( midP.XYZ() );
1287           }
1288           else
1289           {
1290             ++pIt1, ++pIt2;
1291           }
1292         }
1293         // check if the face is more distant than a half of the current segment size,
1294         // if not, segment size is decreased
1295
1296         if ( iLoop == 0 && eData.IsTooDistant( triaSearcher->myBBox, maxSegSize ))
1297           break;
1298         triaSearcher->PrepareToTriaSearch();
1299
1300         //cout << "E " << theMesh.GetMeshDS()->ShapeToIndex( eData.Edge() ) << endl;
1301         sizeDecreased = false;
1302         const gp_Pnt* avoidPnt = & eData.First().myP;
1303         for ( pIt1 = eData.myPoints.begin(); pIt1 != eData.myPoints.end();  )
1304         {
1305           double distToFace =
1306             triaSearcher->GetMinDistInSphere( pIt1->myP, pIt1->mySegSize, isAdjFace, avoidPnt );
1307           double allowedSize = Max( minSize, distToFace*( 1. + grading ));
1308           if ( allowedSize < pIt1->mySegSize  )
1309           {
1310             // cout << "E " << theMesh.GetMeshDS()->ShapeToIndex( eData.Edge() )
1311             //      << "\t closure detected " << endl;
1312             if ( 1.1 * allowedSize < pIt1->mySegSize  )
1313             {
1314               sizeDecreased = true;
1315               sizeTree.SetSize( pIt1->myP, allowedSize );
1316               // cout << "E " << theMesh.GetMeshDS()->ShapeToIndex( eData.Edge() )
1317               //      << "\t SetSize " << allowedSize << " at "
1318               //      << pIt1->myP.X() <<", "<< pIt1->myP.Y()<<", "<<pIt1->myP.Z() << endl;
1319               pIt2 = pIt1;
1320               if ( --pIt2 != eData.myPoints.end() && pIt2->mySegSize > allowedSize )
1321                 sizeTree.SetSize( eData.myC3d.Value( 0.6*pIt2->myU + 0.4*pIt1->myU ), allowedSize );
1322               pIt2 = pIt1;
1323               if ( ++pIt2 != eData.myPoints.end() && pIt2->mySegSize > allowedSize )
1324                 sizeTree.SetSize( eData.myC3d.Value( 0.6*pIt2->myU + 0.4*pIt1->myU ), allowedSize );
1325             }
1326             pIt1->mySegSize = allowedSize;
1327           }
1328           ++pIt1;
1329           if ( & (*pIt1) == & eData.Last() )
1330             avoidPnt = & eData.Last().myP;
1331           else
1332             avoidPnt = NULL;
1333
1334           if ( iLoop > 20 )
1335           {
1336 #ifdef _DEBUG_
1337             cout << "Infinite loop in AdaptiveAlgo::Compute()" << endl;
1338 #endif
1339             sizeDecreased = false;
1340             break;
1341           }
1342         }
1343       } // while ( sizeDecreased )
1344     } // loop on myEdges
1345
1346     // *theProgress = 0.3 + 0.3 * iF / double( faceMap.Extent() );
1347
1348   } // loop on faceMap
1349
1350   return makeSegments();
1351 }
1352
1353 //================================================================================
1354 /*!
1355  * \brief Create segments
1356  */
1357 //================================================================================
1358
1359 bool AdaptiveAlgo::makeSegments()
1360 {
1361   SMESH_HypoFilter quadHyp( SMESH_HypoFilter::HasName( "QuadraticMesh" ));
1362   _quadraticMesh = myMesh->GetHypothesis( myEdges[0].Edge(), quadHyp, /*andAncestors=*/true );
1363
1364   SMESH_MesherHelper helper( *myMesh );
1365   helper.SetIsQuadratic( _quadraticMesh );
1366
1367   vector< double > nbSegs, params;
1368
1369   for ( int iE = 0; iE < myEdges.size(); ++iE )
1370   {
1371     EdgeData& eData = myEdges[ iE ];
1372
1373     // estimate roughly min segment size on the EDGE
1374     double edgeMinSize = myHyp->GetMaxSize();
1375     EdgeData::TPntIter pIt1 = eData.myPoints.begin();
1376     for ( ; pIt1 != eData.myPoints.end(); ++pIt1 )
1377       edgeMinSize = Min( edgeMinSize,
1378                          Min( pIt1->mySegSize, mySizeTree->GetSize( pIt1->myP )));
1379
1380     const double f = eData.myC3d.FirstParameter(), l = eData.myC3d.LastParameter();
1381     const double parLen = l - f;
1382     const int  nbDivSeg = 5;
1383     int           nbDiv = Max( 1, int ( eData.myLength / edgeMinSize * nbDivSeg ));
1384
1385     // compute nb of segments
1386     bool toRecompute = true;
1387     double maxSegSize = 0;
1388     size_t i = 1, segCount;
1389     //cout << "E " << theMesh.GetMeshDS()->ShapeToIndex( eData.Edge() ) << endl;
1390     while ( toRecompute ) // recompute if segment size at some point is less than edgeMinSize/nbDivSeg
1391     {
1392       nbSegs.resize( nbDiv + 1 );
1393       nbSegs[0] = 0;
1394       toRecompute = false;
1395
1396       // fill nbSegs with segment size stored in EdgeData::ProbePnt::mySegSize which can
1397       // be less than size in mySizeTree
1398       pIt1 = eData.myPoints.begin();
1399       EdgeData::ProbePnt* pp1 = &(*pIt1), *pp2;
1400       for ( ++pIt1; pIt1 != eData.myPoints.end(); ++pIt1 )
1401       {
1402         pp2 = &(*pIt1);
1403         double size1 = Min( pp1->mySegSize, myHyp->GetMaxSize() );
1404         double size2 = Min( pp2->mySegSize, myHyp->GetMaxSize() );
1405         double r, u, du = pp2->myU - pp1->myU;
1406         while(( u = f + parLen * i / nbDiv ) < pp2->myU )
1407         {
1408           r = ( u - pp1->myU ) / du;
1409           nbSegs[i] = (1-r) * size1 + r * size2;
1410           ++i;
1411         }
1412         if ( i < nbSegs.size() )
1413           nbSegs[i] = size2;
1414         pp1 = pp2;
1415       }
1416       // fill nbSegs with local nb of segments
1417       gp_Pnt p1 = eData.First().myP, p2, pDiv = p1;
1418       for ( i = 1, segCount = 1; i < nbSegs.size(); ++i )
1419       {
1420         p2 = eData.myC3d.Value( f + parLen * i / nbDiv );
1421         double locSize = Min( mySizeTree->GetSize( p2 ), nbSegs[i] );
1422         double nb      = p1.Distance( p2 ) / locSize;
1423         // if ( nbSegs.size() < 30 )
1424         //   cout << "locSize " << locSize << " nb " << nb << endl;
1425         if ( nb > 1. )
1426         {
1427           toRecompute = true;
1428           edgeMinSize = locSize;
1429           nbDiv = int ( eData.myLength / edgeMinSize * nbDivSeg );
1430           break;
1431         }
1432         nbSegs[i] = nbSegs[i-1] + nb;
1433         p1 = p2;
1434         if ( nbSegs[i] >= segCount )
1435         {
1436           maxSegSize = Max( maxSegSize, pDiv.Distance( p2 ));
1437           pDiv = p2;
1438           ++segCount;
1439         }
1440       }
1441     }
1442
1443     // compute parameters of nodes
1444     int nbSegFinal = Max( 1, int(floor( nbSegs.back() + 0.5 )));
1445     double fact = nbSegFinal / nbSegs.back();
1446     if ( maxSegSize / fact > myHyp->GetMaxSize() )
1447       fact = ++nbSegFinal / nbSegs.back();
1448     //cout << "nbSegs.back() " << nbSegs.back() << " nbSegFinal " << nbSegFinal << endl;
1449     params.clear();
1450     for ( i = 0, segCount = 1; segCount < nbSegFinal; ++segCount )
1451     {
1452       while ( nbSegs[i] * fact < segCount )
1453         ++i;
1454       if ( i < nbDiv )
1455       {
1456         double d = i - ( nbSegs[i] - segCount/fact ) / ( nbSegs[i] - nbSegs[i-1] );
1457         params.push_back( f + parLen * d / nbDiv );
1458         //params.push_back( f + parLen * i / nbDiv );
1459       }
1460       else
1461         break;
1462     }
1463     // get nodes on VERTEXes
1464     TopoDS_Vertex vf = helper.IthVertex( 0, eData.Edge(), false );
1465     TopoDS_Vertex vl = helper.IthVertex( 1, eData.Edge(), false );
1466     myMesh->GetSubMesh( vf )->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1467     myMesh->GetSubMesh( vl )->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1468     const SMDS_MeshNode * nf = VertexNode( vf, myMesh->GetMeshDS() );
1469     const SMDS_MeshNode * nl = VertexNode( vl, myMesh->GetMeshDS() );
1470     if ( !nf || !nl )
1471       return error("No node on vertex");
1472
1473     // create segments
1474     helper.SetSubShape( eData.Edge() );
1475     helper.SetElementsOnShape( true );
1476     const int ID = 0;
1477     const SMDS_MeshNode *n1 = nf, *n2;
1478     for ( i = 0; i < params.size(); ++i, n1 = n2 )
1479     {
1480       gp_Pnt p2 = eData.myC3d.Value( params[i] );
1481       n2 = helper.AddNode( p2.X(), p2.Y(), p2.Z(), ID, params[i] );
1482       helper.AddEdge( n1, n2, ID, /*force3d=*/false );
1483     }
1484     helper.AddEdge( n1, nl, ID, /*force3d=*/false );
1485
1486     eData.myPoints.clear();
1487
1488     //*theProgress = 0.6 + 0.4 * iE / double( myEdges.size() );
1489     if ( _computeCanceled )
1490       return false;
1491
1492   } // loop on EDGEs
1493
1494   SMESHUtils::FreeVector( myEdges );
1495
1496   return true;
1497 }
1498
1499 //================================================================================
1500 /*!
1501  * \brief Predict number of segments on all given EDGEs
1502  */
1503 //================================================================================
1504
1505 bool AdaptiveAlgo::Evaluate(SMESH_Mesh &         theMesh,
1506                             const TopoDS_Shape & theShape,
1507                             MapShapeNbElems&     theResMap)
1508 {
1509   // initialize fields of inherited StdMeshers_Regular_1D
1510   StdMeshers_Regular_1D::_hypType = DEFLECTION;
1511   StdMeshers_Regular_1D::_value[ DEFLECTION_IND ] = myHyp->GetDeflection();
1512
1513   TopExp_Explorer edExp( theShape, TopAbs_EDGE );
1514
1515   for ( ; edExp.More(); edExp.Next() )
1516   {
1517     const TopoDS_Edge & edge = TopoDS::Edge( edExp.Current() );
1518     StdMeshers_Regular_1D::Evaluate( theMesh, theShape, theResMap );
1519   }
1520   return true;
1521 }
1522