Salome HOME
164a5a8ad41364ac7313beea1a77913d468aa590
[modules/geom.git] / src / OCC2VTK / GEOM_WireframeFace.cxx
1 // Copyright (C) 2007-2023  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "GEOM_WireframeFace.h"
21
22 #include <GEOMUtils_Hatcher.hxx>
23  
24 #include <vtkObjectFactory.h> 
25  
26 #include <vtkPoints.h> 
27 #include <vtkCellArray.h> 
28
29 #include <vtkPolyDataMapper.h>  
30 #include <vtkPolyData.h>  
31 #include <vtkInformation.h>
32 #include <vtkInformationVector.h>
33  
34 #include <BRep_Tool.hxx>
35 #include <TColStd_Array1OfReal.hxx>
36
37 #include <Basics_OCCTVersion.hxx>
38
39 #if OCC_VERSION_LARGE < 0x07070000
40 #include <Adaptor3d_HCurve.hxx>
41 #else
42 #include <Adaptor3d_Curve.hxx>
43 #endif
44
45 vtkStandardNewMacro(GEOM_WireframeFace)
46  
47 GEOM_WireframeFace::GEOM_WireframeFace(): 
48   Discret(15)
49
50   NbIso[0] = 1;
51   NbIso[1] = 1;
52
53   this->SetNumberOfInputPorts(0);
54
55  
56 GEOM_WireframeFace::~GEOM_WireframeFace() 
57
58
59  
60 int GEOM_WireframeFace::RequestData(vtkInformation *vtkNotUsed(request),
61                                     vtkInformationVector **vtkNotUsed(inputVector),
62                                     vtkInformationVector *outputVector)
63 {
64   vtkInformation *outInfo = outputVector->GetInformationObject(0);
65   vtkPolyData *aPolyData = vtkPolyData::SafeDownCast(
66     outInfo->Get(vtkDataObject::DATA_OBJECT()));
67
68   aPolyData->Allocate();
69   vtkPoints* aPts = vtkPoints::New();
70   aPolyData->SetPoints(aPts);
71   aPts->Delete();
72
73   TFaceSet::Iterator anIter(myFaceSet);
74   for(; anIter.More(); anIter.Next()){
75     const TopoDS_Face& aFace = anIter.Value();
76     OCC2VTK(aFace,aPolyData,aPts,NbIso,Discret);
77   }
78   return 1;
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   GEOMUtils::Hatcher aHatcher(theFace);
120
121   aHatcher.Init(theNbIso[0], theNbIso[1]);
122   aHatcher.Perform();
123
124   if (aHatcher.IsDone()) {
125     // Push iso lines in vtk kernel
126     CreateIso(aHatcher, Standard_True, theDiscret, thePolyData, thePts);
127     CreateIso(aHatcher, Standard_False, theDiscret, thePolyData, thePts);
128   }
129 }
130
131
132
133 void
134 GEOM_WireframeFace::
135 CreateIso(const GEOMUtils::Hatcher &theHatcher,
136           const Standard_Boolean   IsUIso,
137           const int                theDiscret,
138           vtkPolyData              *thePolyData,
139           vtkPoints                *thePts)
140 {
141   Handle(TColStd_HArray1OfInteger) anIndices;
142   Handle(TColStd_HArray1OfReal)    aParams;
143
144   if (IsUIso) {
145     // U-isolines
146     anIndices = theHatcher.GetUIndices();
147     aParams   = theHatcher.GetUParams();
148   } else {
149     // V-isolines
150     anIndices = theHatcher.GetVIndices();
151     aParams   = theHatcher.GetVParams();
152   }
153
154   if (anIndices.IsNull() == Standard_False &&
155       aParams.IsNull()   == Standard_False) {
156     const GeomAbs_IsoType aType    = (IsUIso ? GeomAbs_IsoU : GeomAbs_IsoV);
157     Standard_Integer      anIsoInd = anIndices->Lower();
158
159     for (; anIsoInd <= anIndices->Upper(); anIsoInd++) {
160       const Standard_Integer aHatchingIndex = anIndices->Value(anIsoInd);
161
162       if (aHatchingIndex != 0) {
163         const Standard_Real    aParam     = aParams->Value(anIsoInd);
164         const Standard_Integer aNbDomains =
165           theHatcher.GetNbDomains(aHatchingIndex);
166
167         if (aNbDomains >= 0) {
168           Standard_Integer anIDom = 1;
169           Standard_Real    aV1;
170           Standard_Real    aV2;
171
172           for (; anIDom <= aNbDomains; anIDom++) {
173             if (theHatcher.GetDomain(aHatchingIndex, anIDom, aV1, aV2)) {
174               CreateIso_(theHatcher.GetFace(), aType, aParam, aV1, aV2,
175                          theDiscret, thePolyData, thePts);
176             }
177           }
178         }
179       }
180     }
181   }
182 }
183
184
185
186 void 
187 GEOM_WireframeFace:: 
188 CreateIso_(const TopoDS_Face& theFace,
189            GeomAbs_IsoType theIsoType, 
190            Standard_Real Par, 
191            Standard_Real T1,
192            Standard_Real T2,
193            const int theDiscret, 
194            vtkPolyData* thePolyData,
195            vtkPoints* thePts)
196 {
197   Standard_Real U1, U2, V1, V2, stepU=0., stepV=0.;
198   Standard_Integer j;
199   gp_Pnt P;
200
201   TopLoc_Location aLoc;
202   const Handle(Geom_Surface)& S = BRep_Tool::Surface(theFace,aLoc);
203
204   if(!S.IsNull()){
205     BRepAdaptor_Surface S(theFace,Standard_False);
206       
207     GeomAbs_SurfaceType SurfType = S.GetType();
208
209     GeomAbs_CurveType CurvType = GeomAbs_OtherCurve;
210
211     Standard_Integer Intrv, nbIntv;
212     Standard_Integer nbUIntv = S.NbUIntervals(GeomAbs_CN);
213     Standard_Integer nbVIntv = S.NbVIntervals(GeomAbs_CN);
214     TColStd_Array1OfReal TI(1,Max(nbUIntv, nbVIntv)+1);
215
216     if(theIsoType == GeomAbs_IsoU){
217       S.VIntervals(TI, GeomAbs_CN);
218       V1 = Max(T1, TI(1));
219       V2 = Min(T2, TI(2));
220       U1 = Par;
221       U2 = Par;
222       stepU = 0;
223       nbIntv = nbVIntv;
224     }else{
225       S.UIntervals(TI, GeomAbs_CN);
226       U1 = Max(T1, TI(1));
227       U2 = Min(T2, TI(2));
228       V1 = Par;
229       V2 = Par;
230       stepV = 0;
231       nbIntv = nbUIntv;
232     }   
233         
234     S.D0(U1,V1,P);
235     MoveTo(P,thePts);
236
237     for(Intrv = 1; Intrv <= nbIntv; Intrv++){
238       if(TI(Intrv) <= T1 && TI(Intrv + 1) <= T1)
239         continue;
240       if(TI(Intrv) >= T2 && TI(Intrv + 1) >= T2)
241               continue;
242       if(theIsoType == GeomAbs_IsoU){
243               V1 = Max(T1, TI(Intrv));
244               V2 = Min(T2, TI(Intrv + 1));
245               stepV = (V2 - V1) / theDiscret;
246       }else{
247               U1 = Max(T1, TI(Intrv));
248               U2 = Min(T2, TI(Intrv + 1));
249               stepU = (U2 - U1) / theDiscret;
250       }
251
252       switch (SurfType) {
253       case GeomAbs_Plane :
254               break;
255       case GeomAbs_Cylinder :
256       case GeomAbs_Cone :
257         if(theIsoType == GeomAbs_IsoV){
258                 for(j = 1; j < theDiscret; j++){
259                   U1 += stepU;
260                   V1 += stepV;
261                   S.D0(U1,V1,P);
262                   DrawTo(P,thePolyData,thePts);
263                 }
264               }
265               break;
266       case GeomAbs_Sphere :
267       case GeomAbs_Torus :
268       case GeomAbs_OffsetSurface :
269       case GeomAbs_OtherSurface :
270         for(j = 1; j < theDiscret; j++){
271                 U1 += stepU;
272                 V1 += stepV;
273                 S.D0(U1,V1,P);
274                 DrawTo(P,thePolyData,thePts);
275               }
276               break;
277       case GeomAbs_BezierSurface :
278       case GeomAbs_BSplineSurface :
279         for(j = 1; j <= theDiscret/2; j++){
280           Standard_Real aStep = (theIsoType == GeomAbs_IsoV) ? stepU*2. : stepV*2.;
281                 CreateIso__(S, theIsoType, U1, V1, aStep, thePolyData, thePts);
282                 U1 += stepU*2.;
283                 V1 += stepV*2.;
284               }
285               break;
286       case GeomAbs_SurfaceOfExtrusion :
287       case GeomAbs_SurfaceOfRevolution :
288         if((theIsoType == GeomAbs_IsoV && SurfType == GeomAbs_SurfaceOfRevolution) ||
289                  (theIsoType == GeomAbs_IsoU && SurfType == GeomAbs_SurfaceOfExtrusion)) 
290         {
291                 if(SurfType == GeomAbs_SurfaceOfExtrusion) 
292             break;
293                 for(j = 1; j < theDiscret; j++){
294                   U1 += stepU;
295                   V1 += stepV;
296                   S.D0(U1,V1,P);
297                   DrawTo(P,thePolyData,thePts);
298                 }
299               }else{
300                 CurvType = (S.BasisCurve())->GetType();
301                 switch(CurvType){
302                 case GeomAbs_Line :
303                   break;
304                 case GeomAbs_Circle :
305                 case GeomAbs_Ellipse :
306                   for (j = 1; j < theDiscret; j++) {
307                     U1 += stepU;
308                     V1 += stepV;
309                     S.D0(U1,V1,P);
310                     DrawTo(P,thePolyData,thePts);
311                   }
312                   break;
313                 case GeomAbs_Parabola :
314                 case GeomAbs_Hyperbola :
315                 case GeomAbs_BezierCurve :
316                 case GeomAbs_BSplineCurve :
317                 case GeomAbs_OtherCurve :
318                   for(j = 1; j <= theDiscret/2; j++){
319                     Standard_Real aStep = (theIsoType == GeomAbs_IsoV) ? stepU*2. : stepV*2.;
320                     CreateIso__(S, theIsoType, U1, V1, aStep, thePolyData, thePts);
321                     U1 += stepU*2.;
322                     V1 += stepV*2.;
323                   }
324                   break;
325                 default:
326                   break;
327                 }
328               }
329       }
330     }
331     S.D0(U2,V2,P);
332     DrawTo(P,thePolyData,thePts);
333   }  
334 }
335  
336  
337  
338  
339 void 
340 GEOM_WireframeFace:: 
341 CreateIso__(const BRepAdaptor_Surface& theSurface, 
342             GeomAbs_IsoType theIsoType,
343                                     Standard_Real& theU, 
344                                     Standard_Real& theV, 
345                                     Standard_Real theStep, 
346             vtkPolyData* thePolyData,
347             vtkPoints* thePts)
348 {
349   gp_Pnt Pl, Pr, Pm;
350   if (theIsoType == GeomAbs_IsoU) {
351     theSurface.D0(theU, theV, Pl);
352     theSurface.D0(theU, theV + theStep/2., Pm);
353     theSurface.D0(theU, theV + theStep, Pr);
354   } else {
355     theSurface.D0(theU, theV, Pl);
356     theSurface.D0(theU + theStep/2., theV, Pm);
357     theSurface.D0(theU + theStep, theV, Pr);
358   }
359
360   static Standard_Real ISO_RATIO = 1.001;
361   if (Pm.Distance(Pl) + Pm.Distance(Pr) <= ISO_RATIO*Pl.Distance(Pr)) {
362     DrawTo(Pr,thePolyData,thePts);
363   } else {
364     if (theIsoType == GeomAbs_IsoU) {
365       CreateIso__(theSurface, theIsoType, theU, theV, theStep/2, thePolyData, thePts);
366       Standard_Real aLocalV = theV + theStep/2 ;
367       CreateIso__(theSurface, theIsoType, theU, aLocalV , theStep/2, thePolyData, thePts);
368     } else {
369       CreateIso__(theSurface, theIsoType, theU, theV, theStep/2, thePolyData, thePts);
370       Standard_Real aLocalU = theU + theStep/2 ;
371       CreateIso__(theSurface, theIsoType, aLocalU , theV, theStep/2, thePolyData, thePts);
372     }
373   }
374 }