]> SALOME platform Git repositories - modules/geom.git/blob - src/OCC2VTK/GEOM_WireframeFace.cxx
Salome HOME
OCCT dev version porting (6.7.2)
[modules/geom.git] / src / OCC2VTK / GEOM_WireframeFace.cxx
1 // Copyright (C) 2007-2014  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 #include "OCC2VTK_internal.h"
22
23 #include <GEOMUtils_Hatcher.hxx>
24  
25 #include <vtkObjectFactory.h> 
26  
27 #include <vtkPoints.h> 
28 #include <vtkCellArray.h> 
29
30 #include <vtkPolyDataMapper.h>  
31 #include <vtkPolyData.h>  
32 #include <vtkInformation.h>
33 #include <vtkInformationVector.h>
34  
35 #include <Adaptor3d_HCurve.hxx>
36 #include <BRep_Tool.hxx>
37 #include <TColStd_Array1OfReal.hxx>
38
39 vtkStandardNewMacro(GEOM_WireframeFace);
40  
41 GEOM_WireframeFace::GEOM_WireframeFace(): 
42   Discret(15)
43
44   NbIso[0] = 1;
45   NbIso[1] = 1;
46
47   this->SetNumberOfInputPorts(0);
48
49  
50 GEOM_WireframeFace::~GEOM_WireframeFace() 
51
52
53  
54 int GEOM_WireframeFace::RequestData(vtkInformation *vtkNotUsed(request),
55                                     vtkInformationVector **vtkNotUsed(inputVector),
56                                     vtkInformationVector *outputVector)
57 {
58   vtkInformation *outInfo = outputVector->GetInformationObject(0);
59   vtkPolyData *aPolyData = vtkPolyData::SafeDownCast(
60     outInfo->Get(vtkDataObject::DATA_OBJECT()));
61
62   aPolyData->Allocate();
63   vtkPoints* aPts = vtkPoints::New();
64   aPolyData->SetPoints(aPts);
65   aPts->Delete();
66
67   TFaceSet::Iterator anIter(myData->myFaceSet);
68   for(; anIter.More(); anIter.Next()){
69     const TopoDS_Face& aFace = anIter.Value();
70     OCC2VTK(aFace,aPolyData,aPts,NbIso,Discret);
71   }
72   return 1;
73 }
74
75 void GEOM_WireframeFace::SetNbIso(const int theNb[2])
76 {
77   if ( theNb[0] == NbIso[0] && theNb[1] == NbIso[1])
78     return;
79
80   NbIso[0] = theNb[0];
81   NbIso[1] = theNb[1];
82
83   Modified();
84 }
85
86 void GEOM_WireframeFace::GetNbIso(int &theNbU,int &theNbV)
87 {
88   theNbU = NbIso[0];
89   theNbV = NbIso[1];
90 }
91
92 void  
93 GEOM_WireframeFace:: 
94 OCC2VTK(const TopoDS_Face& theFace,
95         vtkPolyData* thePolyData,
96                     vtkPoints* thePts,  
97         const int theNbIso[2], 
98         const int theDiscret) 
99
100   TopoDS_Face aFace = theFace; 
101   aFace.Orientation(TopAbs_FORWARD);
102   CreateIso(aFace,theNbIso,theDiscret,thePolyData,thePts);
103 }
104
105 void
106 GEOM_WireframeFace::
107 CreateIso(const TopoDS_Face& theFace,
108           const int theNbIso[2],
109           const int theDiscret,
110           vtkPolyData* thePolyData,
111           vtkPoints* thePts)
112 {
113   GEOMUtils_Hatcher aHatcher(theFace);
114
115   aHatcher.Init(theNbIso[0], theNbIso[1]);
116   aHatcher.Perform();
117
118   if (aHatcher.IsDone()) {
119     // Push iso lines in vtk kernel
120     CreateIso(aHatcher, Standard_True, theDiscret, thePolyData, thePts);
121     CreateIso(aHatcher, Standard_False, theDiscret, thePolyData, thePts);
122   }
123 }
124
125
126
127 void
128 GEOM_WireframeFace::
129 CreateIso(const GEOMUtils_Hatcher &theHatcher,
130           const Standard_Boolean   IsUIso,
131           const int                theDiscret,
132           vtkPolyData       *thePolyData,
133           vtkPoints         *thePts)
134 {
135   Handle(TColStd_HArray1OfInteger) anIndices;
136   Handle(TColStd_HArray1OfReal)    aParams;
137
138   if (IsUIso) {
139     // U-isolines
140     anIndices = theHatcher.GetUIndices();
141     aParams   = theHatcher.GetUParams();
142   } else {
143     // V-isolines
144     anIndices = theHatcher.GetVIndices();
145     aParams   = theHatcher.GetVParams();
146   }
147
148   if (anIndices.IsNull() == Standard_False &&
149       aParams.IsNull()   == Standard_False) {
150     const GeomAbs_IsoType aType    = (IsUIso ? GeomAbs_IsoU : GeomAbs_IsoV);
151     Standard_Integer      anIsoInd = anIndices->Lower();
152
153     for (; anIsoInd <= anIndices->Upper(); anIsoInd++) {
154       const Standard_Integer aHatchingIndex = anIndices->Value(anIsoInd);
155
156       if (aHatchingIndex != 0) {
157         const Standard_Real    aParam     = aParams->Value(anIsoInd);
158         const Standard_Integer aNbDomains =
159           theHatcher.GetNbDomains(aHatchingIndex);
160
161         if (aNbDomains >= 0) {
162           Standard_Integer anIDom = 1;
163           Standard_Real    aV1;
164           Standard_Real    aV2;
165
166           for (; anIDom <= aNbDomains; anIDom++) {
167             if (theHatcher.GetDomain(aHatchingIndex, anIDom, aV1, aV2)) {
168               CreateIso_(theHatcher.GetFace(), aType, aParam, aV1, aV2,
169                          theDiscret, thePolyData, thePts);
170             }
171           }
172         }
173       }
174     }
175   }
176 }
177
178
179
180 void 
181 GEOM_WireframeFace:: 
182 CreateIso_(const TopoDS_Face& theFace,
183            GeomAbs_IsoType theIsoType, 
184            Standard_Real Par, 
185            Standard_Real T1,
186            Standard_Real T2,
187            const int theDiscret, 
188            vtkPolyData* thePolyData,
189            vtkPoints* thePts)
190 {
191   Standard_Real U1, U2, V1, V2, stepU=0., stepV=0.;
192   Standard_Integer j;
193   gp_Pnt P;
194
195   TopLoc_Location aLoc;
196   const Handle(Geom_Surface)& S = BRep_Tool::Surface(theFace,aLoc);
197
198   if(!S.IsNull()){
199     BRepAdaptor_Surface S(theFace,Standard_False);
200       
201     GeomAbs_SurfaceType SurfType = S.GetType();
202
203     GeomAbs_CurveType CurvType = GeomAbs_OtherCurve;
204
205     Standard_Integer Intrv, nbIntv;
206     Standard_Integer nbUIntv = S.NbUIntervals(GeomAbs_CN);
207     Standard_Integer nbVIntv = S.NbVIntervals(GeomAbs_CN);
208     TColStd_Array1OfReal TI(1,Max(nbUIntv, nbVIntv)+1);
209
210     if(theIsoType == GeomAbs_IsoU){
211       S.VIntervals(TI, GeomAbs_CN);
212       V1 = Max(T1, TI(1));
213       V2 = Min(T2, TI(2));
214       U1 = Par;
215       U2 = Par;
216       stepU = 0;
217       nbIntv = nbVIntv;
218     }else{
219       S.UIntervals(TI, GeomAbs_CN);
220       U1 = Max(T1, TI(1));
221       U2 = Min(T2, TI(2));
222       V1 = Par;
223       V2 = Par;
224       stepV = 0;
225       nbIntv = nbUIntv;
226     }   
227         
228     S.D0(U1,V1,P);
229     MoveTo(P,thePts);
230
231     for(Intrv = 1; Intrv <= nbIntv; Intrv++){
232       if(TI(Intrv) <= T1 && TI(Intrv + 1) <= T1)
233         continue;
234       if(TI(Intrv) >= T2 && TI(Intrv + 1) >= T2)
235               continue;
236       if(theIsoType == GeomAbs_IsoU){
237               V1 = Max(T1, TI(Intrv));
238               V2 = Min(T2, TI(Intrv + 1));
239               stepV = (V2 - V1) / theDiscret;
240       }else{
241               U1 = Max(T1, TI(Intrv));
242               U2 = Min(T2, TI(Intrv + 1));
243               stepU = (U2 - U1) / theDiscret;
244       }
245
246       switch (SurfType) {
247       case GeomAbs_Plane :
248               break;
249       case GeomAbs_Cylinder :
250       case GeomAbs_Cone :
251         if(theIsoType == GeomAbs_IsoV){
252                 for(j = 1; j < theDiscret; j++){
253                   U1 += stepU;
254                   V1 += stepV;
255                   S.D0(U1,V1,P);
256                   DrawTo(P,thePolyData,thePts);
257                 }
258               }
259               break;
260       case GeomAbs_Sphere :
261       case GeomAbs_Torus :
262       case GeomAbs_OffsetSurface :
263       case GeomAbs_OtherSurface :
264         for(j = 1; j < theDiscret; j++){
265                 U1 += stepU;
266                 V1 += stepV;
267                 S.D0(U1,V1,P);
268                 DrawTo(P,thePolyData,thePts);
269               }
270               break;
271       case GeomAbs_BezierSurface :
272       case GeomAbs_BSplineSurface :
273         for(j = 1; j <= theDiscret/2; j++){
274           Standard_Real aStep = (theIsoType == GeomAbs_IsoV) ? stepU*2. : stepV*2.;
275                 CreateIso__(S, theIsoType, U1, V1, aStep, thePolyData, thePts);
276                 U1 += stepU*2.;
277                 V1 += stepV*2.;
278               }
279               break;
280       case GeomAbs_SurfaceOfExtrusion :
281       case GeomAbs_SurfaceOfRevolution :
282         if((theIsoType == GeomAbs_IsoV && SurfType == GeomAbs_SurfaceOfRevolution) ||
283                  (theIsoType == GeomAbs_IsoU && SurfType == GeomAbs_SurfaceOfExtrusion)) 
284         {
285                 if(SurfType == GeomAbs_SurfaceOfExtrusion) 
286             break;
287                 for(j = 1; j < theDiscret; j++){
288                   U1 += stepU;
289                   V1 += stepV;
290                   S.D0(U1,V1,P);
291                   DrawTo(P,thePolyData,thePts);
292                 }
293               }else{
294                 CurvType = (S.BasisCurve())->GetType();
295                 switch(CurvType){
296                 case GeomAbs_Line :
297                   break;
298                 case GeomAbs_Circle :
299                 case GeomAbs_Ellipse :
300                   for (j = 1; j < theDiscret; j++) {
301                     U1 += stepU;
302                     V1 += stepV;
303                     S.D0(U1,V1,P);
304                     DrawTo(P,thePolyData,thePts);
305                   }
306                   break;
307                 case GeomAbs_Parabola :
308                 case GeomAbs_Hyperbola :
309                 case GeomAbs_BezierCurve :
310                 case GeomAbs_BSplineCurve :
311                 case GeomAbs_OtherCurve :
312                   for(j = 1; j <= theDiscret/2; j++){
313               Standard_Real aStep = (theIsoType == GeomAbs_IsoV) ? stepU*2. : stepV*2.;
314                   CreateIso__(S, theIsoType, U1, V1, aStep, thePolyData, thePts);
315                     U1 += stepU*2.;
316                     V1 += stepV*2.;
317                   }
318                   break;
319                 }
320               }
321       }
322     }
323     S.D0(U2,V2,P);
324     DrawTo(P,thePolyData,thePts);
325   }  
326 }
327  
328  
329  
330  
331 void 
332 GEOM_WireframeFace:: 
333 CreateIso__(const BRepAdaptor_Surface& theSurface, 
334             GeomAbs_IsoType theIsoType,
335                                     Standard_Real& theU, 
336                                     Standard_Real& theV, 
337                                     Standard_Real theStep, 
338             vtkPolyData* thePolyData,
339             vtkPoints* thePts)
340 {
341   gp_Pnt Pl, Pr, Pm;
342   if (theIsoType == GeomAbs_IsoU) {
343     theSurface.D0(theU, theV, Pl);
344     theSurface.D0(theU, theV + theStep/2., Pm);
345     theSurface.D0(theU, theV + theStep, Pr);
346   } else {
347     theSurface.D0(theU, theV, Pl);
348     theSurface.D0(theU + theStep/2., theV, Pm);
349     theSurface.D0(theU + theStep, theV, Pr);
350   }
351
352   static Standard_Real ISO_RATIO = 1.001;
353   if (Pm.Distance(Pl) + Pm.Distance(Pr) <= ISO_RATIO*Pl.Distance(Pr)) {
354     DrawTo(Pr,thePolyData,thePts);
355   } else {
356     if (theIsoType == GeomAbs_IsoU) {
357       CreateIso__(theSurface, theIsoType, theU, theV, theStep/2, thePolyData, thePts);
358       Standard_Real aLocalV = theV + theStep/2 ;
359       CreateIso__(theSurface, theIsoType, theU, aLocalV , theStep/2, thePolyData, thePts);
360     } else {
361       CreateIso__(theSurface, theIsoType, theU, theV, theStep/2, thePolyData, thePts);
362       Standard_Real aLocalU = theU + theStep/2 ;
363       CreateIso__(theSurface, theIsoType, aLocalU , theV, theStep/2, thePolyData, thePts);
364     }
365   }
366 }