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