Salome HOME
0020044: EDF 866 GEOM: Extrusion along a path : impossible to select a Wire
[modules/geom.git] / src / OBJECT / GEOM_WireframeFace.cxx
1 #include "GEOM_WireframeFace.h" 
2  
3 #include <vtkObjectFactory.h> 
4  
5 #include <vtkPoints.h> 
6 #include <vtkCellArray.h> 
7
8 #include <vtkPolyDataMapper.h>  
9 #include <vtkPolyData.h>  
10  
11 #include <Precision.hxx>
12 #include <BRepTools.hxx>
13 #include <TopExp_Explorer.hxx>
14 #include <Geom2dHatch_Hatcher.hxx>
15 #include <Geom2dHatch_Intersector.hxx>
16 #include <TColStd_Array1OfReal.hxx>
17 #include <TColStd_Array1OfInteger.hxx>
18  
19 #include <TopoDS.hxx> 
20 #include <TopoDS_Edge.hxx> 
21 #include <BRep_Tool.hxx>
22 #include <Geom2d_TrimmedCurve.hxx>
23 #include <Geom2d_Line.hxx>
24 #include <gp_Dir2d.hxx>
25 #include <gp_Pnt2d.hxx>
26  
27 #include <Geom2dHatch_Hatcher.hxx>
28 #include <HatchGen_Domain.hxx>
29
30 #include <Adaptor3d_HCurve.hxx>
31
32 vtkStandardNewMacro(GEOM_WireframeFace);
33  
34 GEOM_WireframeFace::GEOM_WireframeFace(): 
35   Discret(15)
36
37   NbIso[0] = 1;
38   NbIso[1] = 1;
39
40  
41 GEOM_WireframeFace::~GEOM_WireframeFace() 
42
43
44  
45 void
46 GEOM_WireframeFace:: 
47 Execute()
48 {
49   vtkPolyData* aPolyData = GetOutput();
50   aPolyData->Allocate();
51   vtkPoints* aPts = vtkPoints::New();
52   aPolyData->SetPoints(aPts);
53   aPts->Delete();
54
55   TFaceSet::Iterator anIter(myFaceSet);
56   for(; anIter.More(); anIter.Next()){
57     const TopoDS_Face& aFace = anIter.Value();
58     OCC2VTK(aFace,aPolyData,aPts,NbIso,Discret);
59   }
60 }
61
62 void GEOM_WireframeFace::SetNbIso(const int theNb[2])
63 {
64   if ( theNb[0] == NbIso[0] && theNb[1] == NbIso[1])
65     return;
66
67   NbIso[0] = theNb[0];
68   NbIso[1] = theNb[1];
69
70   Modified();
71 }
72
73 void GEOM_WireframeFace::GetNbIso(int &theNbU,int &theNbV)
74 {
75   theNbU = NbIso[0];
76   theNbV = NbIso[1];
77 }
78
79 void  
80 GEOM_WireframeFace:: 
81 OCC2VTK(const TopoDS_Face& theFace,
82         vtkPolyData* thePolyData,
83                     vtkPoints* thePts,  
84         const int theNbIso[2], 
85         const int theDiscret) 
86
87   TopoDS_Face aFace = theFace; 
88   aFace.Orientation(TopAbs_FORWARD);
89   CreateIso(aFace,theNbIso,theDiscret,thePolyData,thePts);
90 }
91
92 void 
93 GEOM_WireframeFace:: 
94 CreateIso(const TopoDS_Face& theFace,
95           const int theNbIso[2], 
96           const int theDiscret, 
97           vtkPolyData* thePolyData,
98           vtkPoints* thePts)
99 {
100   // Constants for iso building
101   static Standard_Real INTERSECTOR_CONFUSION = 1.e-10 ; // -8 ;
102   static Standard_Real INTERSECTOR_TANGENCY  = 1.e-10 ; // -8 ;
103
104   static Standard_Real HATHCER_CONFUSION_2D = 1.e-8 ;
105   static Standard_Real HATHCER_CONFUSION_3D = 1.e-8 ;
106
107   Geom2dHatch_Hatcher 
108     aHatcher(Geom2dHatch_Intersector(INTERSECTOR_CONFUSION,
109                                      INTERSECTOR_TANGENCY),
110                          HATHCER_CONFUSION_2D,
111                          HATHCER_CONFUSION_3D,
112                                      Standard_True,
113                                      Standard_False);
114   
115   Standard_Real anUMin, anUMax, aVMin, aVMax;
116   TColStd_Array1OfReal anUPrm(1, theNbIso[0]), aVPrm(1, theNbIso[1]);
117   TColStd_Array1OfInteger anUInd(1, theNbIso[0]), aVInd(1, theNbIso[1]);
118
119   anUInd.Init(0);
120   aVInd.Init(0);
121
122   //-----------------------------------------------------------------------
123   // If the Min Max bounds are infinite, there are bounded to Infinite
124   // value.
125   //-----------------------------------------------------------------------
126   BRepTools::UVBounds(theFace, anUMin, anUMax, aVMin, aVMax) ;
127   Standard_Boolean InfiniteUMin = Precision::IsNegativeInfinite (anUMin) ;
128   Standard_Boolean InfiniteUMax = Precision::IsPositiveInfinite (anUMax) ;
129   Standard_Boolean InfiniteVMin = Precision::IsNegativeInfinite (aVMin) ;
130   Standard_Boolean InfiniteVMax = Precision::IsPositiveInfinite (aVMax) ;
131
132   static float VTKINFINITE = 1.0E38;
133   if(InfiniteUMin && InfiniteUMax){
134     anUMin = - VTKINFINITE ;
135     anUMax =   VTKINFINITE ;
136   }else if(InfiniteUMin){
137     anUMin = anUMax - VTKINFINITE ;
138   }else if(InfiniteUMax){
139     anUMax = anUMin + VTKINFINITE ;
140   }
141
142   if(InfiniteVMin && InfiniteVMax){
143     aVMin = - VTKINFINITE ;
144     aVMax =   VTKINFINITE ;
145   }else if(InfiniteVMin){
146     aVMin = aVMax - VTKINFINITE ;
147   }else if(InfiniteVMax){
148     aVMax = aVMin + VTKINFINITE ;
149   }
150
151   //-----------------------------------------------------------------------
152   // Retreiving the edges and loading them into the hatcher.
153   //-----------------------------------------------------------------------
154   TopExp_Explorer ExpEdges(theFace, TopAbs_EDGE);
155   for(; ExpEdges.More(); ExpEdges.Next()){
156     const TopoDS_Edge& anEdge = TopoDS::Edge(ExpEdges.Current());
157     Standard_Real U1, U2 ;
158     const Handle(Geom2d_Curve) PCurve = 
159       BRep_Tool::CurveOnSurface(anEdge, theFace, U1, U2) ;
160
161     if(PCurve.IsNull() || U1 == U2)
162       return;
163
164     //-- Test if a TrimmedCurve is necessary
165     if(Abs(PCurve->FirstParameter()-U1) <= Precision::PConfusion() &&
166              Abs(PCurve->LastParameter()-U2) <= Precision::PConfusion())
167     { 
168       aHatcher.AddElement(PCurve, anEdge.Orientation()) ;      
169     }else{ 
170       if(!PCurve->IsPeriodic()){
171               Handle(Geom2d_TrimmedCurve) TrimPCurve =
172           Handle(Geom2d_TrimmedCurve)::DownCast(PCurve);
173               if(!TrimPCurve.IsNull()){
174           Handle_Geom2d_Curve aBasisCurve = TrimPCurve->BasisCurve();
175                 if(aBasisCurve->FirstParameter()-U1 > Precision::PConfusion() ||
176                    U2-aBasisCurve->LastParameter() > Precision::PConfusion()) 
177           {
178                   aHatcher.AddElement(PCurve, anEdge.Orientation()) ;      
179                   return;
180                 }
181               }else{
182                 if(PCurve->FirstParameter()-U1 > Precision::PConfusion()){
183                   U1=PCurve->FirstParameter();
184                 }
185                 if(U2-PCurve->LastParameter()  > Precision::PConfusion()){
186                   U2=PCurve->LastParameter();
187                 }
188               }
189       }
190       Handle(Geom2d_TrimmedCurve) TrimPCurve = 
191         new Geom2d_TrimmedCurve(PCurve, U1, U2);
192       aHatcher.AddElement(TrimPCurve, anEdge.Orientation());
193     }
194   }
195
196
197   //-----------------------------------------------------------------------
198   // Loading and trimming the hatchings.
199   //-----------------------------------------------------------------------
200   Standard_Integer IIso;
201   Standard_Real DeltaU = Abs(anUMax - anUMin) ;
202   Standard_Real DeltaV = Abs(aVMax - aVMin) ;
203   Standard_Real confusion = Min(DeltaU, DeltaV) * HATHCER_CONFUSION_3D ;
204   aHatcher.Confusion3d (confusion) ;
205
206   Standard_Real StepU = DeltaU / (Standard_Real)theNbIso[0];
207   if(StepU > confusion){
208     Standard_Real UPrm = anUMin + StepU / 2.;
209     gp_Dir2d Dir(0., 1.) ;
210     for(IIso = 1 ; IIso <= theNbIso[0] ; IIso++) {
211       anUPrm(IIso) = UPrm ;
212       gp_Pnt2d Ori (UPrm, 0.) ;
213       Geom2dAdaptor_Curve HCur (new Geom2d_Line (Ori, Dir)) ;
214       anUInd(IIso) = aHatcher.AddHatching (HCur) ;
215       UPrm += StepU ;
216     }
217   }
218
219   Standard_Real StepV = DeltaV / (Standard_Real) theNbIso[1] ;
220   if(StepV > confusion){
221     Standard_Real VPrm = aVMin + StepV / 2.;
222     gp_Dir2d Dir(1., 0.);
223     for(IIso = 1 ; IIso <= theNbIso[1] ; IIso++){
224       aVPrm(IIso) = VPrm;
225       gp_Pnt2d Ori (0., VPrm);
226       Geom2dAdaptor_Curve HCur(new Geom2d_Line (Ori, Dir));
227       aVInd(IIso) = aHatcher.AddHatching (HCur) ;
228       VPrm += StepV ;
229     }
230   }
231
232   //-----------------------------------------------------------------------
233   // Computation.
234   //-----------------------------------------------------------------------
235   aHatcher.Trim() ;
236
237   Standard_Integer aNbDom = 0 ; // for debug purpose
238   Standard_Integer Index ;
239
240   for(IIso = 1 ; IIso <= theNbIso[0] ; IIso++){
241     Index = anUInd(IIso) ;
242     if(Index != 0){
243       if(aHatcher.TrimDone(Index) && !aHatcher.TrimFailed(Index)){
244               aHatcher.ComputeDomains(Index);
245               if(aHatcher.IsDone (Index)) 
246           aNbDom = aHatcher.NbDomains (Index);
247       }
248     }
249   }
250
251   for(IIso = 1 ; IIso <= theNbIso[1] ; IIso++){
252     Index = aVInd(IIso);
253     if(Index != 0){
254       if(aHatcher.TrimDone (Index) && !aHatcher.TrimFailed(Index)){
255               aHatcher.ComputeDomains (Index);
256               if(aHatcher.IsDone (Index)) 
257           aNbDom = aHatcher.NbDomains (Index);
258       }
259     }
260   }
261
262   //-----------------------------------------------------------------------
263   // Push iso lines in vtk kernel
264   //-----------------------------------------------------------------------
265   for(Standard_Integer UIso = anUPrm.Lower() ; UIso <= anUPrm.Upper(); UIso++){
266     Standard_Integer UInd = anUInd.Value(UIso);
267     if(UInd != 0){
268       Standard_Real UPrm = anUPrm.Value(UIso);
269       if(aHatcher.IsDone(UInd)){
270               Standard_Integer NbDom = aHatcher.NbDomains(UInd);
271               for(Standard_Integer IDom = 1 ; IDom <= NbDom ; IDom++){
272                 const HatchGen_Domain& Dom = aHatcher.Domain(UInd, IDom) ;
273                 Standard_Real V1 = Dom.HasFirstPoint()? Dom.FirstPoint().Parameter(): aVMin - VTKINFINITE;
274                 Standard_Real V2 = Dom.HasSecondPoint()? Dom.SecondPoint().Parameter(): aVMax + VTKINFINITE;
275                 CreateIso_(theFace, GeomAbs_IsoU, UPrm, V1, V2, theDiscret, thePolyData, thePts);
276         }
277             }
278     }
279   }
280
281   for(Standard_Integer VIso = aVPrm.Lower() ; VIso <= aVPrm.Upper(); VIso++){
282     Standard_Integer VInd = aVInd.Value(VIso);
283     if(VInd != 0){
284       Standard_Real VPrm = aVPrm.Value(VIso);
285       if(aHatcher.IsDone (VInd)){
286               Standard_Integer NbDom = aHatcher.NbDomains(VInd);
287               for (Standard_Integer IDom = 1 ; IDom <= NbDom ; IDom++){
288                 const HatchGen_Domain& Dom = aHatcher.Domain(VInd, IDom);
289                 Standard_Real U1 = Dom.HasFirstPoint()? Dom.FirstPoint().Parameter(): aVMin - VTKINFINITE;
290                 Standard_Real U2 = Dom.HasSecondPoint()? Dom.SecondPoint().Parameter(): aVMax + VTKINFINITE;
291             CreateIso_(theFace, GeomAbs_IsoV, VPrm, U1, U2, theDiscret, thePolyData, thePts);
292               }
293       }
294     }
295   }
296 }
297
298  
299
300 void 
301 GEOM_WireframeFace:: 
302 CreateIso_(const TopoDS_Face& theFace,
303            GeomAbs_IsoType theIsoType, 
304            Standard_Real Par, 
305            Standard_Real T1,
306            Standard_Real T2,
307            const int theDiscret, 
308            vtkPolyData* thePolyData,
309            vtkPoints* thePts)
310 {
311   Standard_Real U1, U2, V1, V2, stepU=0., stepV=0.;
312   Standard_Integer j;
313   gp_Pnt P;
314
315   TopLoc_Location aLoc;
316   const Handle(Geom_Surface)& S = BRep_Tool::Surface(theFace,aLoc);
317
318   if(!S.IsNull()){
319     BRepAdaptor_Surface S(theFace,Standard_False);
320       
321     GeomAbs_SurfaceType SurfType = S.GetType();
322
323     GeomAbs_CurveType CurvType = GeomAbs_OtherCurve;
324
325     Standard_Integer Intrv, nbIntv;
326     Standard_Integer nbUIntv = S.NbUIntervals(GeomAbs_CN);
327     Standard_Integer nbVIntv = S.NbVIntervals(GeomAbs_CN);
328     TColStd_Array1OfReal TI(1,Max(nbUIntv, nbVIntv)+1);
329
330     if(theIsoType == GeomAbs_IsoU){
331       S.VIntervals(TI, GeomAbs_CN);
332       V1 = Max(T1, TI(1));
333       V2 = Min(T2, TI(2));
334       U1 = Par;
335       U2 = Par;
336       stepU = 0;
337       nbIntv = nbVIntv;
338     }else{
339       S.UIntervals(TI, GeomAbs_CN);
340       U1 = Max(T1, TI(1));
341       U2 = Min(T2, TI(2));
342       V1 = Par;
343       V2 = Par;
344       stepV = 0;
345       nbIntv = nbUIntv;
346     }   
347         
348     S.D0(U1,V1,P);
349     MoveTo(P,thePts);
350
351     for(Intrv = 1; Intrv <= nbIntv; Intrv++){
352       if(TI(Intrv) <= T1 && TI(Intrv + 1) <= T1)
353         continue;
354       if(TI(Intrv) >= T2 && TI(Intrv + 1) >= T2)
355               continue;
356       if(theIsoType == GeomAbs_IsoU){
357               V1 = Max(T1, TI(Intrv));
358               V2 = Min(T2, TI(Intrv + 1));
359               stepV = (V2 - V1) / theDiscret;
360       }else{
361               U1 = Max(T1, TI(Intrv));
362               U2 = Min(T2, TI(Intrv + 1));
363               stepU = (U2 - U1) / theDiscret;
364       }
365
366       switch (SurfType) {
367       case GeomAbs_Plane :
368               break;
369       case GeomAbs_Cylinder :
370       case GeomAbs_Cone :
371         if(theIsoType == GeomAbs_IsoV){
372                 for(j = 1; j < theDiscret; j++){
373                   U1 += stepU;
374                   V1 += stepV;
375                   S.D0(U1,V1,P);
376                   DrawTo(P,thePolyData,thePts);
377                 }
378               }
379               break;
380       case GeomAbs_Sphere :
381       case GeomAbs_Torus :
382       case GeomAbs_OffsetSurface :
383       case GeomAbs_OtherSurface :
384         for(j = 1; j < theDiscret; j++){
385                 U1 += stepU;
386                 V1 += stepV;
387                 S.D0(U1,V1,P);
388                 DrawTo(P,thePolyData,thePts);
389               }
390               break;
391       case GeomAbs_BezierSurface :
392       case GeomAbs_BSplineSurface :
393         for(j = 1; j <= theDiscret/2; j++){
394           Standard_Real aStep = (theIsoType == GeomAbs_IsoV) ? stepU*2. : stepV*2.;
395                 CreateIso__(S, theIsoType, U1, V1, aStep, thePolyData, thePts);
396                 U1 += stepU*2.;
397                 V1 += stepV*2.;
398               }
399               break;
400       case GeomAbs_SurfaceOfExtrusion :
401       case GeomAbs_SurfaceOfRevolution :
402         if((theIsoType == GeomAbs_IsoV && SurfType == GeomAbs_SurfaceOfRevolution) ||
403                  (theIsoType == GeomAbs_IsoU && SurfType == GeomAbs_SurfaceOfExtrusion)) 
404         {
405                 if(SurfType == GeomAbs_SurfaceOfExtrusion) 
406             break;
407                 for(j = 1; j < theDiscret; j++){
408                   U1 += stepU;
409                   V1 += stepV;
410                   S.D0(U1,V1,P);
411                   DrawTo(P,thePolyData,thePts);
412                 }
413               }else{
414                 CurvType = (S.BasisCurve())->GetType();
415                 switch(CurvType){
416                 case GeomAbs_Line :
417                   break;
418                 case GeomAbs_Circle :
419                 case GeomAbs_Ellipse :
420                   for (j = 1; j < theDiscret; j++) {
421                     U1 += stepU;
422                     V1 += stepV;
423                     S.D0(U1,V1,P);
424                     DrawTo(P,thePolyData,thePts);
425                   }
426                   break;
427                 case GeomAbs_Parabola :
428                 case GeomAbs_Hyperbola :
429                 case GeomAbs_BezierCurve :
430                 case GeomAbs_BSplineCurve :
431                 case GeomAbs_OtherCurve :
432                   for(j = 1; j <= theDiscret/2; j++){
433               Standard_Real aStep = (theIsoType == GeomAbs_IsoV) ? stepU*2. : stepV*2.;
434                   CreateIso__(S, theIsoType, U1, V1, aStep, thePolyData, thePts);
435                     U1 += stepU*2.;
436                     V1 += stepV*2.;
437                   }
438                   break;
439                 }
440               }
441       }
442     }
443     S.D0(U2,V2,P);
444     DrawTo(P,thePolyData,thePts);
445   }  
446 }
447  
448  
449  
450  
451 void 
452 GEOM_WireframeFace:: 
453 CreateIso__(const BRepAdaptor_Surface& theSurface, 
454             GeomAbs_IsoType theIsoType,
455                                     Standard_Real& theU, 
456                                     Standard_Real& theV, 
457                                     Standard_Real theStep, 
458             vtkPolyData* thePolyData,
459             vtkPoints* thePts)
460 {
461   gp_Pnt Pl, Pr, Pm;
462   if (theIsoType == GeomAbs_IsoU) {
463     theSurface.D0(theU, theV, Pl);
464     theSurface.D0(theU, theV + theStep/2., Pm);
465     theSurface.D0(theU, theV + theStep, Pr);
466   } else {
467     theSurface.D0(theU, theV, Pl);
468     theSurface.D0(theU + theStep/2., theV, Pm);
469     theSurface.D0(theU + theStep, theV, Pr);
470   }
471
472   static Standard_Real ISO_RATIO = 1.001;
473   if (Pm.Distance(Pl) + Pm.Distance(Pr) <= ISO_RATIO*Pl.Distance(Pr)) {
474     DrawTo(Pr,thePolyData,thePts);
475   } else {
476     if (theIsoType == GeomAbs_IsoU) {
477       CreateIso__(theSurface, theIsoType, theU, theV, theStep/2, thePolyData, thePts);
478       Standard_Real aLocalV = theV + theStep/2 ;
479       CreateIso__(theSurface, theIsoType, theU, aLocalV , theStep/2, thePolyData, thePts);
480     } else {
481       CreateIso__(theSurface, theIsoType, theU, theV, theStep/2, thePolyData, thePts);
482       Standard_Real aLocalU = theU + theStep/2 ;
483       CreateIso__(theSurface, theIsoType, aLocalU , theV, theStep/2, thePolyData, thePts);
484     }
485   }
486 }