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