]> SALOME platform Git repositories - modules/geom.git/blob - src/OBJECT/GEOM_OCCReader.cxx
Salome HOME
IPAL20715 impossible to select sub-shape after apply group
[modules/geom.git] / src / OBJECT / GEOM_OCCReader.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 //  GEOM OBJECT : interactive object for Geometry entities visualization
23 //  File   : GEOM_OCCReader.h
24 //  Author : Christophe ATTANASIO
25 //  Module : GEOM
26 //  $Header$
27 //
28 #include "GEOM_OCCReader.h"
29
30 // VTK Includes
31 #include <vtkPoints.h>
32 #include <vtkCellArray.h>
33
34 #include <vtkObjectFactory.h>
35 #include <vtkPolyData.h>
36 #include <vtkPolyDataMapper.h>
37 #include <vtkMergePoints.h>
38
39 #include <vtkTransform.h>
40 #include <vtkMatrix4x4.h>
41
42 // OpenCASCADE Includes
43 #include <BRepAdaptor_Surface.hxx>
44 #include <TopExp_Explorer.hxx>
45 #include <BRepMesh_IncrementalMesh.hxx>
46 #include <Poly_Triangulation.hxx>
47 #include <Poly_Polygon3D.hxx>
48 #include <BRep_Tool.hxx>
49 #include <TopoDS_Face.hxx>
50 #include <TopoDS_Edge.hxx>
51 #include <TopoDS_Wire.hxx>
52 #include <BRepBndLib.hxx>
53 #include <TopoDS.hxx>
54 #include <TopAbs.hxx>
55 #include <Precision.hxx>
56 #include <BRepTools.hxx>
57 #include <BRep_Tool.hxx>
58 #include <Geom2dAdaptor_Curve.hxx>
59 #include <Geom2dHatch_Intersector.hxx>
60 #include <Geom2dHatch_Hatcher.hxx>
61 #include <Geom2d_Curve.hxx>
62 #include <Geom2d_Line.hxx>
63 #include <Geom2d_TrimmedCurve.hxx>
64 #include <HatchGen_Domain.hxx>
65 #include <GeomAbs_IsoType.hxx>
66 #include <Precision.hxx>
67 #include <TopAbs_ShapeEnum.hxx>
68 #include <TopExp_Explorer.hxx>
69 #include <TopoDS.hxx>
70 #include <TopoDS_Edge.hxx>
71 #include <gp_Dir2d.hxx>
72 #include <gp_Pnt2d.hxx>
73 #include <TColStd_Array1OfInteger.hxx>
74 #include <TColStd_Array1OfReal.hxx>
75 #include <Adaptor3d_HCurve.hxx>
76
77 #include "utilities.h"
78
79 using namespace std;
80
81 #define MAX2(X, Y)      (  Abs(X) > Abs(Y)? Abs(X) : Abs(Y) )
82 #define MAX3(X, Y, Z)   ( MAX2 ( MAX2(X,Y) , Z) )
83
84 // Constante for iso building
85 static Standard_Real IntersectorConfusion = 1.e-10 ; // -8 ;
86 static Standard_Real IntersectorTangency  = 1.e-10 ; // -8 ;
87 static Standard_Real HatcherConfusion2d   = 1.e-8 ;
88 static Standard_Real HatcherConfusion3d   = 1.e-8 ;
89
90 static Standard_Integer lastVTKpoint = 0;
91 static Standard_Integer PlotCount = 0;
92 static Standard_Real IsoRatio = 1.001;
93 static Standard_Integer MaxPlotCount = 5; 
94
95 //=======================================================================
96 // Function : New
97 // Purpose  : 
98 //=======================================================================
99
100 GEOM_OCCReader* GEOM_OCCReader::New()
101 {
102   vtkObject* ret = vtkObjectFactory::CreateInstance("GEOM_OCCReader");
103   if(ret) {
104     return (GEOM_OCCReader*)ret;
105   }
106   return new GEOM_OCCReader;
107 }
108
109 //=======================================================================
110 // Function : GEOM_OCCReader
111 // Purpose  : 
112 //=======================================================================
113
114 GEOM_OCCReader::GEOM_OCCReader()
115 {
116   //this->myShape = NULL;
117   this->amode = 0;
118   this->forced = Standard_False;
119   this->discretiso = 15;
120   this->nbisos = 1;
121 }
122 //=======================================================================
123 // Function : ~GEOM_OCCReader
124 // Purpose  : 
125 //=======================================================================
126
127 GEOM_OCCReader::~GEOM_OCCReader()
128 {
129 }
130
131
132 //=======================================================================
133 // Function : Execute
134 // Purpose  : 
135 //=======================================================================
136
137
138 void GEOM_OCCReader::Execute() {
139
140   vtkPolyData* output = this->GetOutput();
141   vtkPoints* Pts = NULL;
142   vtkCellArray* Cells = NULL;
143   TopLoc_Location aLoc;
144
145   // Allocation
146   Pts = vtkPoints::New();
147   Cells = vtkCellArray::New();
148         
149   //Compute number of triangles and points
150   Standard_Integer nbpoly=0,nbpts=0;
151
152   if(amode==1) {
153     //for shading
154     
155     if(myShape.ShapeType() == TopAbs_FACE) {
156       // whole FACE 
157       const TopoDS_Face& aFace = TopoDS::Face(myShape);
158       Handle(Poly_Triangulation) aPoly = BRep_Tool::Triangulation(aFace,aLoc);
159       if(aPoly.IsNull()) {
160         Pts->Delete();
161         Cells->Delete();
162         return;
163       }
164
165       nbpts = aPoly->NbNodes();
166       nbpoly = aPoly->NbTriangles();
167
168       Pts->SetNumberOfPoints(nbpts);
169       Cells->Allocate(Cells->EstimateSize(nbpoly,3));
170     }
171     else { 
172         Cells->Delete();
173         Pts->Delete();
174         return; 
175     }
176   }
177
178   // Start computation
179   if(amode == 0) {
180     ComputeWireframe(Pts,Cells);
181     output->SetPoints(Pts);
182     output->SetLines(Cells);
183     output->Squeeze();
184   }
185   else {
186     if(myShape.ShapeType() == TopAbs_FACE) {
187       ComputeShading(Pts,Cells);
188
189       output->SetPoints(Pts);
190       output->SetPolys(Cells);
191       output->Squeeze();
192     }
193   }
194   Pts->Delete();
195   Cells->Delete();
196   
197 }
198
199 //=======================================================================
200 // Function : ComputeWireframe
201 // Purpose  : Compute the shape in CAD wireframe mode
202 //=======================================================================
203
204 void GEOM_OCCReader::ComputeWireframe(vtkPoints* Pts,vtkCellArray* Cells){
205
206   // Check the type of the shape:
207   if(myShape.ShapeType() == TopAbs_FACE) {
208     // Face
209     TransferFaceWData(TopoDS::Face(myShape),Pts,Cells);
210   } else if(myShape.ShapeType() == TopAbs_EDGE) {
211     // Edge
212     TransferEdgeWData(TopoDS::Edge(myShape),Pts,Cells);
213   } else {
214     if(myShape.ShapeType() == TopAbs_VERTEX) {
215       // Vertex
216       TransferVertexWData(TopoDS::Vertex(myShape),Pts,Cells);
217     }
218   }
219 }
220
221 //=======================================================================
222 // Function : TransferFaceWData
223 // Purpose  : Transfert wireframe data for FACE
224 //=======================================================================
225
226 void GEOM_OCCReader::TransferFaceWData(const TopoDS_Face& aFace,
227                                          vtkPoints* Pts,
228                                          vtkCellArray* Cells) 
229 {
230   TopoDS_Face aCopyFace = aFace; 
231   aCopyFace.Orientation (TopAbs_FORWARD);
232   createISO(aCopyFace,Precision::Infinite(),1,Pts,Cells);
233 }
234
235 //=======================================================================
236 // Function : createISO
237 // Purpose  : Create ISO for Face Wireframe representation 
238 //=======================================================================
239
240 void GEOM_OCCReader::createISO (const TopoDS_Face&     TopologicalFace,
241                                   const Standard_Real    Infinite,
242                                   const Standard_Integer NbIsos,
243                                   vtkPoints* Pts,
244                                   vtkCellArray* Cell)
245 {
246   Geom2dHatch_Hatcher aHatcher (Geom2dHatch_Intersector (IntersectorConfusion,
247                                                          IntersectorTangency),
248                                 HatcherConfusion2d,
249                                 HatcherConfusion3d,
250                                 Standard_True,
251                                 Standard_False);
252   
253   Standard_Real myInfinite,myUMin,myUMax,myVMin,myVMax;
254   //myInfinite = Precision::Infinite();
255   myInfinite = 1e38; // VTK uses float numbers - Precision::Infinite() is double and can not be accepted.
256
257   Standard_Integer myNbDom;
258   TColStd_Array1OfReal myUPrm(1, NbIsos),myVPrm(1, NbIsos);
259   TColStd_Array1OfInteger myUInd(1, NbIsos),myVInd(1, NbIsos);
260
261   myUInd.Init(0);
262   myVInd.Init(0);
263
264   //-----------------------------------------------------------------------
265   // If the Min Max bounds are infinite, there are bounded to Infinite
266   // value.
267   //-----------------------------------------------------------------------
268
269   BRepTools::UVBounds (TopologicalFace, myUMin, myUMax, myVMin, myVMax) ;
270   Standard_Boolean InfiniteUMin = Precision::IsNegativeInfinite (myUMin) ;
271   Standard_Boolean InfiniteUMax = Precision::IsPositiveInfinite (myUMax) ;
272   Standard_Boolean InfiniteVMin = Precision::IsNegativeInfinite (myVMin) ;
273   Standard_Boolean InfiniteVMax = Precision::IsPositiveInfinite (myVMax) ;
274   if (InfiniteUMin && InfiniteUMax) {
275     myUMin = - myInfinite ;
276     myUMax =   myInfinite ;
277   } else if (InfiniteUMin) {
278     myUMin = myUMax - myInfinite ;
279   } else if (InfiniteUMax) {
280     myUMax = myUMin + myInfinite ;
281   }
282   if (InfiniteVMin && InfiniteVMax) {
283     myVMin = - myInfinite ;
284     myVMax =   myInfinite ;
285   } else if (InfiniteVMin) {
286     myVMin = myVMax - myInfinite ;
287   } else if (InfiniteVMax) {
288     myVMax = myVMin + myInfinite ;
289   }
290
291   //-----------------------------------------------------------------------
292   // Retreiving the edges and loading them into the hatcher.
293   //-----------------------------------------------------------------------
294
295   TopExp_Explorer ExpEdges ;
296   for (ExpEdges.Init (TopologicalFace, TopAbs_EDGE) ; ExpEdges.More() ; ExpEdges.Next()) {
297     const TopoDS_Edge& TopologicalEdge = TopoDS::Edge (ExpEdges.Current()) ;
298     Standard_Real U1, U2 ;
299     const Handle(Geom2d_Curve) PCurve = BRep_Tool::CurveOnSurface (TopologicalEdge, TopologicalFace, U1, U2) ;
300
301     if ( PCurve.IsNull() ) {
302       return;
303     }
304
305     if ( U1==U2) {
306       return;
307     }
308
309     //-- Test if a TrimmedCurve is necessary
310     if(   Abs(PCurve->FirstParameter()-U1)<= Precision::PConfusion() 
311           && Abs(PCurve->LastParameter()-U2)<= Precision::PConfusion()) { 
312       aHatcher.AddElement (PCurve, TopologicalEdge.Orientation()) ;      
313     }
314     else { 
315       if (!PCurve->IsPeriodic()) {
316         Handle (Geom2d_TrimmedCurve) TrimPCurve =Handle(Geom2d_TrimmedCurve)::DownCast(PCurve);
317         if (!TrimPCurve.IsNull()) {
318           if (TrimPCurve->BasisCurve()->FirstParameter()-U1 > Precision::PConfusion() ||
319               U2-TrimPCurve->BasisCurve()->LastParameter()  > Precision::PConfusion()) {
320             aHatcher.AddElement (PCurve, TopologicalEdge.Orientation()) ;      
321             return;
322           }
323         }
324         else {
325           if (PCurve->FirstParameter()-U1 > Precision::PConfusion()){
326             U1=PCurve->FirstParameter();
327           }
328           if (U2-PCurve->LastParameter()  > Precision::PConfusion()){
329             U2=PCurve->LastParameter();
330           }
331         }
332       }
333       Handle (Geom2d_TrimmedCurve) TrimPCurve = new Geom2d_TrimmedCurve (PCurve, U1, U2) ;
334       aHatcher.AddElement (TrimPCurve, TopologicalEdge.Orientation()) ;
335     }
336   }
337
338
339   //-----------------------------------------------------------------------
340   // Loading and trimming the hatchings.
341   //-----------------------------------------------------------------------
342
343   Standard_Integer IIso ;
344   Standard_Real DeltaU = Abs (myUMax - myUMin) ;
345   Standard_Real DeltaV = Abs (myVMax - myVMin) ;
346   Standard_Real confusion = Min (DeltaU, DeltaV) * HatcherConfusion3d ;
347   aHatcher.Confusion3d (confusion) ;
348
349   Standard_Real StepU = DeltaU / (Standard_Real) NbIsos ;
350   if (StepU > confusion) {
351     Standard_Real UPrm = myUMin + StepU / 2. ;
352     gp_Dir2d Dir (0., 1.) ;
353     for (IIso = 1 ; IIso <= NbIsos ; IIso++) {
354       myUPrm(IIso) = UPrm ;
355       gp_Pnt2d Ori (UPrm, 0.) ;
356       Geom2dAdaptor_Curve HCur (new Geom2d_Line (Ori, Dir)) ;
357       myUInd(IIso) = aHatcher.AddHatching (HCur) ;
358       UPrm += StepU ;
359     }
360   }
361
362   Standard_Real StepV = DeltaV / (Standard_Real) NbIsos ;
363   if (StepV > confusion) {
364     Standard_Real VPrm = myVMin + StepV / 2. ;
365     gp_Dir2d Dir (1., 0.) ;
366     for (IIso = 1 ; IIso <= NbIsos ; IIso++) {
367       myVPrm(IIso) = VPrm ;
368       gp_Pnt2d Ori (0., VPrm) ;
369       Geom2dAdaptor_Curve HCur (new Geom2d_Line (Ori, Dir)) ;
370       myVInd(IIso) = aHatcher.AddHatching (HCur) ;
371       VPrm += StepV ;
372     }
373   }
374
375   //-----------------------------------------------------------------------
376   // Computation.
377   //-----------------------------------------------------------------------
378
379   aHatcher.Trim() ;
380
381   myNbDom = 0 ;
382   for (IIso = 1 ; IIso <= NbIsos ; IIso++) {
383     Standard_Integer Index ;
384
385     Index = myUInd(IIso) ;
386     if (Index != 0) {
387       if (aHatcher.TrimDone (Index) && !aHatcher.TrimFailed (Index)) {
388         aHatcher.ComputeDomains (Index);
389         if (aHatcher.IsDone (Index)) myNbDom = myNbDom + aHatcher.NbDomains (Index) ;
390       }
391     }
392
393     Index = myVInd(IIso) ;
394     if (Index != 0) {
395       if (aHatcher.TrimDone (Index) && !aHatcher.TrimFailed (Index)) {
396         aHatcher.ComputeDomains (Index);
397         if (aHatcher.IsDone (Index)) myNbDom = myNbDom + aHatcher.NbDomains (Index) ;
398       }
399     }
400   }
401
402   //-----------------------------------------------------------------------
403   // Push iso lines in vtk kernel
404   //-----------------------------------------------------------------------
405
406
407   Standard_Integer pt_start_idx = 0;
408
409   for (Standard_Integer UIso = myUPrm.Lower() ; UIso <= myUPrm.Upper() ; UIso++) {
410     Standard_Integer UInd = myUInd.Value (UIso) ;
411     if (UInd != 0) {
412       Standard_Real UPrm = myUPrm.Value (UIso) ;
413       if (!aHatcher.IsDone (UInd)) {
414         MESSAGE("DBRep_IsoBuilder:: U iso of parameter: "<<UPrm)
415         switch (aHatcher.Status (UInd)) {
416         case HatchGen_NoProblem          : MESSAGE("No Problem")          ; break ;
417         case HatchGen_TrimFailure        : MESSAGE("Trim Failure")        ; break ;
418         case HatchGen_TransitionFailure  : MESSAGE("Transition Failure")  ; break ;
419         case HatchGen_IncoherentParity   : MESSAGE("Incoherent Parity")   ; break ;
420         case HatchGen_IncompatibleStates : MESSAGE("Incompatible States") ; break ;
421         }
422       } else {
423         Standard_Integer NbDom = aHatcher.NbDomains (UInd) ;
424         for (Standard_Integer IDom = 1 ; IDom <= NbDom ; IDom++) {
425           const HatchGen_Domain& Dom = aHatcher.Domain (UInd, IDom) ;
426           Standard_Real V1 = Dom.HasFirstPoint()  ? Dom.FirstPoint().Parameter()  : myVMin - myInfinite ;
427           Standard_Real V2 = Dom.HasSecondPoint() ? Dom.SecondPoint().Parameter() : myVMax + myInfinite ;
428           DrawIso(GeomAbs_IsoU, UPrm, V1, V2, Pts, Cell,pt_start_idx);
429         }
430       }
431     }
432   }
433
434   for (Standard_Integer VIso = myVPrm.Lower() ; VIso <= myVPrm.Upper() ; VIso++) {
435     Standard_Integer VInd = myVInd.Value (VIso) ;
436     if (VInd != 0) {
437       Standard_Real VPrm = myVPrm.Value (VIso) ;
438       if (!aHatcher.IsDone (VInd)) {
439         MESSAGE("DBRep_IsoBuilder:: V iso of parameter: "<<VPrm)
440         switch (aHatcher.Status (VInd)) {
441         case HatchGen_NoProblem          : MESSAGE("No Problem")          ; break ;
442         case HatchGen_TrimFailure        : MESSAGE("Trim Failure")        ; break ;
443         case HatchGen_TransitionFailure  : MESSAGE("Transition Failure")  ; break ;
444         case HatchGen_IncoherentParity   : MESSAGE("Incoherent Parity")   ; break ;
445         case HatchGen_IncompatibleStates : MESSAGE("Incompatible States") ; break ;
446         }
447       } else {
448         Standard_Integer NbDom = aHatcher.NbDomains (VInd) ;
449         for (Standard_Integer IDom = 1 ; IDom <= NbDom ; IDom++) {
450           const HatchGen_Domain& Dom = aHatcher.Domain (VInd, IDom) ;
451           Standard_Real U1 = Dom.HasFirstPoint()  ? Dom.FirstPoint().Parameter()  : myVMin - myInfinite ;
452           Standard_Real U2 = Dom.HasSecondPoint() ? Dom.SecondPoint().Parameter() : myVMax + myInfinite ;
453           DrawIso(GeomAbs_IsoV, VPrm, U1, U2, Pts, Cell,pt_start_idx) ;
454         }
455       }
456     }
457   }
458
459 }
460
461 //=======================================================================
462 // Function : MoveTo
463 // Purpose  : Init VTK ISO PLOT
464 //=======================================================================
465 void GEOM_OCCReader::MoveTo(gp_Pnt P,
466                               vtkPoints* Pts)
467 {    
468   float coord[3];
469
470   coord[0] = P.X(); coord[1] = P.Y(); coord[2] = P.Z();
471   lastVTKpoint = Pts->InsertNextPoint(coord);
472     
473
474
475 //=======================================================================
476 // Function : DrawTo
477 // Purpose  : Plot point in VTK
478 //=======================================================================
479 void GEOM_OCCReader::DrawTo(gp_Pnt P,
480                               vtkPoints* Pts,
481                               vtkCellArray* Cells)
482 {
483   float coord[3];
484   coord[0] = P.X(); coord[1] = P.Y(); coord[2] = P.Z();
485   Standard_Integer NewVTKpoint =  Pts->InsertNextPoint(coord);
486
487   int pts[2];
488   pts[0] = lastVTKpoint;
489   pts[1] = NewVTKpoint;
490
491   Cells->InsertNextCell(2,pts);
492     
493   lastVTKpoint = NewVTKpoint;
494 }
495
496
497 //=======================================================================
498 // Function : DrawIso
499 // Purpose  : Draw an iso on vtk
500 //=======================================================================
501 void GEOM_OCCReader::DrawIso(GeomAbs_IsoType T, 
502                                Standard_Real Par, 
503                                Standard_Real T1,
504                                Standard_Real T2,
505                                vtkPoints* Pts,
506                                vtkCellArray* Cells,
507                                Standard_Integer& startidx)
508 {
509
510   Standard_Boolean halt = Standard_False;
511   Standard_Integer j,myDiscret = discretiso;
512   Standard_Real U1,U2,V1,V2,stepU=0.,stepV=0.;
513   gp_Pnt P;
514   TopLoc_Location l;
515
516   const Handle(Geom_Surface)& S = BRep_Tool::Surface(TopoDS::Face(myShape),l);
517   if (!S.IsNull()) {
518     BRepAdaptor_Surface S(TopoDS::Face(myShape),Standard_False);
519       
520     GeomAbs_SurfaceType SurfType = S.GetType();
521
522     GeomAbs_CurveType CurvType = GeomAbs_OtherCurve;
523
524     Standard_Integer Intrv, nbIntv;
525     Standard_Integer nbUIntv = S.NbUIntervals(GeomAbs_CN);
526     Standard_Integer nbVIntv = S.NbVIntervals(GeomAbs_CN);
527     TColStd_Array1OfReal TI(1,Max(nbUIntv, nbVIntv)+1);
528
529
530     if (T == GeomAbs_IsoU) {
531       S.VIntervals(TI, GeomAbs_CN);
532       V1 = Max(T1, TI(1));
533       V2 = Min(T2, TI(2));
534       U1 = Par;
535       U2 = Par;
536       stepU = 0;
537       nbIntv = nbVIntv;
538     }
539     else {
540       S.UIntervals(TI, GeomAbs_CN);
541       U1 = Max(T1, TI(1));
542       U2 = Min(T2, TI(2));
543       V1 = Par;
544       V2 = Par;
545       stepV = 0;
546       nbIntv = nbUIntv;
547     }   
548         
549     S.D0(U1,V1,P);
550     MoveTo(P,Pts);
551
552     for (Intrv = 1; Intrv <= nbIntv; Intrv++) {
553
554       if (TI(Intrv) <= T1 && TI(Intrv + 1) <= T1)
555         continue;
556       if (TI(Intrv) >= T2 && TI(Intrv + 1) >= T2)
557         continue;
558       if (T == GeomAbs_IsoU) {
559         V1 = Max(T1, TI(Intrv));
560         V2 = Min(T2, TI(Intrv + 1));
561         stepV = (V2 - V1) / myDiscret;
562       }
563       else {
564         U1 = Max(T1, TI(Intrv));
565         U2 = Min(T2, TI(Intrv + 1));
566         stepU = (U2 - U1) / myDiscret;
567       }
568
569       switch (SurfType) {
570         //-------------GeomAbs_Plane---------------
571       case GeomAbs_Plane :
572         break;
573         //----GeomAbs_Cylinder   GeomAbs_Cone------
574       case GeomAbs_Cylinder :
575       case GeomAbs_Cone :
576         if (T == GeomAbs_IsoV) {
577           for (j = 1; j < myDiscret; j++) {
578             U1 += stepU;
579             V1 += stepV;
580             S.D0(U1,V1,P);
581             DrawTo(P,Pts,Cells);
582           }
583         }
584         break;
585         //---GeomAbs_Sphere   GeomAbs_Torus--------
586         //GeomAbs_BezierSurface GeomAbs_BezierSurface
587       case GeomAbs_Sphere :
588       case GeomAbs_Torus :
589       case GeomAbs_OffsetSurface :
590       case GeomAbs_OtherSurface :
591         for (j = 1; j < myDiscret; j++) {
592           U1 += stepU;
593           V1 += stepV;
594           S.D0(U1,V1,P);
595           DrawTo(P,Pts,Cells);
596         }
597         break;
598         //-------------GeomAbs_BSplineSurface------
599       case GeomAbs_BezierSurface :
600       case GeomAbs_BSplineSurface :
601         for (j = 1; j <= myDiscret/2; j++) {
602
603           PlotCount = 0;
604
605           PlotIso ( S, T, U1, V1, (T == GeomAbs_IsoV) ? stepU*2. : stepV*2., halt, Pts, Cells);
606           U1 += stepU*2.;
607           V1 += stepV*2.;
608         }
609         break;
610         //-------------GeomAbs_SurfaceOfExtrusion--
611         //-------------GeomAbs_SurfaceOfRevolution-
612       case GeomAbs_SurfaceOfExtrusion :
613       case GeomAbs_SurfaceOfRevolution :
614         if ((T == GeomAbs_IsoV && SurfType == GeomAbs_SurfaceOfRevolution) ||
615             (T == GeomAbs_IsoU && SurfType == GeomAbs_SurfaceOfExtrusion)) {
616           if (SurfType == GeomAbs_SurfaceOfExtrusion) break;
617           for (j = 1; j < myDiscret; j++) {
618             U1 += stepU;
619             V1 += stepV;
620             S.D0(U1,V1,P);
621             DrawTo(P,Pts,Cells);
622           }
623         } else {
624           CurvType = (S.BasisCurve())->GetType();
625           switch (CurvType) {
626           case GeomAbs_Line :
627             break;
628           case GeomAbs_Circle :
629           case GeomAbs_Ellipse :
630             for (j = 1; j < myDiscret; j++) {
631               U1 += stepU;
632               V1 += stepV;
633               S.D0(U1,V1,P);
634               DrawTo(P,Pts,Cells);
635             }
636             break;
637           case GeomAbs_Parabola :
638           case GeomAbs_Hyperbola :
639           case GeomAbs_BezierCurve :
640           case GeomAbs_BSplineCurve :
641           case GeomAbs_OtherCurve :
642             for (j = 1; j <= myDiscret/2; j++) {
643
644               PlotCount = 0;
645
646               PlotIso ( S, T, U1, V1,(T == GeomAbs_IsoV) ? stepU*2. : stepV*2., halt, Pts, Cells);
647               U1 += stepU*2.;
648               V1 += stepV*2.;
649             }
650             break;
651           }
652         }
653       }
654     }
655     S.D0(U2,V2,P);
656     DrawTo(P,Pts,Cells);
657   }  
658 }
659
660 //=======================================================================
661 // Function : PlotIso
662 // Purpose  : Plot iso for other surface
663 //=======================================================================
664
665 void GEOM_OCCReader::PlotIso (BRepAdaptor_Surface& S, 
666                                 GeomAbs_IsoType T,
667                                 Standard_Real& U, 
668                                 Standard_Real& V, 
669                                 Standard_Real Step, 
670                                 Standard_Boolean& halt,
671                                 vtkPoints* Pts,
672                                 vtkCellArray* Cells)
673 {
674
675   ++PlotCount; 
676
677   gp_Pnt Pl, Pr, Pm;
678
679   if (T == GeomAbs_IsoU) {
680     S.D0(U, V, Pl);
681     S.D0(U, V + Step/2., Pm);
682     S.D0(U, V + Step, Pr);
683   } else {
684     S.D0(U, V, Pl);
685     S.D0(U + Step/2., V, Pm);
686     S.D0(U + Step, V, Pr);
687   }
688
689   if (PlotCount > MaxPlotCount) {
690     DrawTo(Pr,Pts,Cells);
691     return;
692   }
693
694   if (Pm.Distance(Pl) + Pm.Distance(Pr) <= IsoRatio*Pl.Distance(Pr)) {
695     DrawTo(Pr,Pts,Cells);
696   } else 
697     if (T == GeomAbs_IsoU) {
698       PlotIso ( S, T, U, V, Step/2, halt, Pts, Cells);
699       Standard_Real aLocalV = V + Step/2 ;
700       PlotIso ( S, T, U, aLocalV , Step/2, halt, Pts, Cells);
701     } else {
702       PlotIso ( S, T, U, V, Step/2, halt, Pts, Cells);
703       Standard_Real aLocalU = U + Step/2 ;
704       PlotIso ( S, T, aLocalU , V, Step/2, halt, Pts, Cells);
705     }
706 }
707
708 //=======================================================================
709 // Function : TransferEdgeWData
710 // Purpose  : Transfert wireframe data for EDGE
711 //=======================================================================
712
713 void GEOM_OCCReader::TransferEdgeWData(const TopoDS_Edge& aEdge,
714                                          vtkPoints* Pts,
715                                          vtkCellArray* Cells) {
716   
717   
718   Handle(Poly_PolygonOnTriangulation) aEdgePoly;
719   Standard_Integer i = 1;
720   Handle(Poly_Triangulation) T;
721   TopLoc_Location aEdgeLoc;
722   BRep_Tool::PolygonOnTriangulation(aEdge, aEdgePoly, T, aEdgeLoc, i);
723   
724   Handle(Poly_Polygon3D) P;
725   if(aEdgePoly.IsNull()) {
726     P = BRep_Tool::Polygon3D(aEdge, aEdgeLoc);
727   }
728
729   if(P.IsNull() && aEdgePoly.IsNull())
730     return;
731   
732   // Location edges
733   //---------------
734   
735   gp_Trsf edgeTransf;
736   Standard_Boolean isidtrsf = true;
737   if(!aEdgeLoc.IsIdentity())  {
738     isidtrsf = false;
739     edgeTransf = aEdgeLoc.Transformation();
740   }
741
742   gp_Pnt aP1, aP2;
743
744   Standard_Integer nbnodes;
745   if (aEdgePoly.IsNull()) {
746     nbnodes = P->NbNodes();
747     const TColgp_Array1OfPnt& theNodesP = P->Nodes();
748
749     aP1 = theNodesP(1);
750     aP2 = theNodesP(nbnodes);
751
752     float coord[3];
753     int pts[2];
754
755     for(int j=1;j<nbnodes;j++) {
756       gp_Pnt pt1 = theNodesP(j);
757       gp_Pnt pt2 = theNodesP(j+1);
758     
759       if(!isidtrsf) {
760         // apply edge transformation
761         pt1.Transform(edgeTransf);
762         pt2.Transform(edgeTransf);
763       }
764       
765       // insert pt1
766       coord[0] = pt1.X(); coord[1] = pt1.Y(); coord[2] = pt1.Z();
767       pts[0] = Pts->InsertNextPoint(coord);
768       
769       // insert pt2
770       coord[0] = pt2.X(); coord[1] = pt2.Y(); coord[2] = pt2.Z();
771       pts[1] = Pts->InsertNextPoint(coord);
772       
773       // insert line (pt1,pt2)
774       Cells->InsertNextCell(2,pts);
775     }
776   } else {
777     nbnodes = aEdgePoly->NbNodes();
778     const TColStd_Array1OfInteger& Nodesidx = aEdgePoly->Nodes();
779     const TColgp_Array1OfPnt& theNodesPoly = T->Nodes();
780
781     aP1 = theNodesPoly(1);
782     aP2 = theNodesPoly(nbnodes);
783
784     float coord[3];
785     int pts[2];
786     
787     for(int j=1;j<nbnodes;j++) {
788       Standard_Integer id1 = Nodesidx(j);
789       Standard_Integer id2 = Nodesidx(j+1);
790       
791       gp_Pnt pt1 = theNodesPoly(id1);
792       gp_Pnt pt2 = theNodesPoly(id2);
793           
794       if(!isidtrsf) {
795         // apply edge transformation
796         pt1.Transform(edgeTransf);
797         pt2.Transform(edgeTransf);
798       }
799       
800       // insert pt1
801       coord[0] = pt1.X(); coord[1] = pt1.Y(); coord[2] = pt1.Z();
802       pts[0] = Pts->InsertNextPoint(coord);
803       
804       // insert pt2
805       coord[0] = pt2.X(); coord[1] = pt2.Y(); coord[2] = pt2.Z();
806       pts[1] = Pts->InsertNextPoint(coord);
807       
808       // insert line (pt1,pt2)
809       Cells->InsertNextCell(2,pts);
810     }
811   }
812
813   // vector representation has an arrow on its end
814   if (myIsVector)
815   {
816     if (!isidtrsf) {
817       // apply edge transformation
818       aP1.Transform(edgeTransf);
819       aP2.Transform(edgeTransf);
820     }
821
822     // draw an arrow
823     gp_Vec aDirVec (aP1, aP2);
824     Standard_Real aDist = aDirVec.Magnitude();
825     if (aDist < gp::Resolution()) return;
826     gp_Dir aDirection (aDirVec);
827
828     Standard_Real anAngle = PI/180.*5.;
829     Standard_Real aLength = aDist/10.;
830
831     Standard_Real dx,dy,dz;
832     aDirection.Coord(dx,dy,dz);
833
834     // Pointe de la fleche
835     Standard_Real xo,yo,zo;
836     aP2.Coord(xo,yo,zo);
837
838     // Centre du cercle base de la fleche
839     gp_XYZ aPc = aP2.XYZ() - aDirection.XYZ() * aLength;
840
841     // Construction d'un repere i,j pour le cercle
842     gp_Dir aDirN;
843     if      (Abs(dx) <= Abs(dy) && Abs(dx) <= Abs(dz)) aDirN = gp::DX();
844     else if (Abs(dy) <= Abs(dz) && Abs(dy) <= Abs(dx)) aDirN = gp::DY();
845     else aDirN = gp::DZ();
846
847     gp_Dir aDirI = aDirection ^ aDirN;
848     gp_Dir aDirJ = aDirection ^ aDirI;
849
850     // Add points and segments, composing the arrow
851     Standard_Real cosinus, sinus, Tg = tan(anAngle);
852
853     float coord[3];
854     coord[0] = xo; coord[1] = yo; coord[2] = zo;
855
856     int ptLoc = Pts->InsertNextPoint(coord);
857     int ptFirst = 0;
858     int ptPrev = 0;
859     int ptCur = 0;
860
861     int pts[2];
862
863     int NbPoints = 15;
864     for (int i = 1; i <= NbPoints; i++, ptPrev = ptCur)
865     {
866       cosinus = cos(2. * PI / NbPoints * (i-1));   
867       sinus   = sin(2. * PI / NbPoints * (i-1));
868
869       gp_XYZ aP = aPc + (aDirI.XYZ() * cosinus + aDirJ.XYZ() * sinus) * aLength * Tg;
870       coord[0] = aP.X();
871       coord[1] = aP.Y();
872       coord[2] = aP.Z();
873
874       // insert pts
875       ptCur = Pts->InsertNextPoint(coord);
876       pts[0] = ptCur;
877
878       if (i == 1) {
879         ptFirst = ptCur;
880       }
881       else {
882         // insert line (ptCur,ptPrev)
883         pts[1] = ptPrev;
884         Cells->InsertNextCell(2,pts);
885       }
886
887       // insert line (ptCur,ptLoc)
888       pts[1] = ptLoc;
889       Cells->InsertNextCell(2,pts);
890     }
891
892     // insert line (ptCur,ptFirst)
893     pts[0] = ptCur;
894     pts[1] = ptFirst;
895     Cells->InsertNextCell(2,pts);
896   }
897 }
898
899 /*  Standard_Integer nbnodes = aEdgePoly->NbNodes();
900   const TColStd_Array1OfInteger& Nodesidx = aEdgePoly->Nodes();
901   const TColgp_Array1OfPnt& theNodes = T->Nodes();
902     
903   float coord[3];
904   int pts[2];
905     
906
907   // PUSH NODES
908   for(i=1;i<=nbnodes;i++) {
909     Standard_Integer id = Nodesidx(i);
910     gp_Pnt pt = theNodes(id);
911  
912     float coord[3];
913     if(!isidtrsf) pt.Transform(edgeTransf);
914
915     coord[0] = pt.X(); coord[1] = pt.Y(); coord[2] = pt.Z();
916
917     Pts->SetPoint(id-1,coord);
918
919   }
920
921   // PUSH EDGES
922   for(i=1;i<nbnodes;i++) {
923       
924     Standard_Integer id1 = Nodesidx(i);
925     Standard_Integer id2 = Nodesidx(i+1);
926     
927     int pts[2];
928     pts[0] = id1-1; pts[1] = id2-1;
929
930     // insert line (pt1,pt2)
931     Cells->InsertNextCell(2,pts);
932   }
933  
934   
935   }*/
936
937 //=======================================================================
938 // Function : TransferVertexWData
939 // Purpose  : Transfert wireframe data for VERTEX
940 //=======================================================================
941
942 void GEOM_OCCReader::TransferVertexWData(const TopoDS_Vertex& aVertex,
943                                          vtkPoints* Pts,
944                                          vtkCellArray* Cells) {
945 #define ZERO_COORD coord[0] = 0.0; coord[1] = 0.0; coord[2] = 0.0
946   
947   gp_Pnt P = BRep_Tool::Pnt( aVertex );
948   float delta = 1, coord[3];
949   int pts[2];
950   // insert pt
951   ZERO_COORD; coord[0] = +delta;
952   pts[0] = Pts->InsertNextPoint(coord);
953   coord[0] = -delta;
954   pts[1] = Pts->InsertNextPoint(coord);
955   // insert line (pt1,pt2)
956   Cells->InsertNextCell(2,pts);
957
958   ZERO_COORD; coord[1] = +delta;
959   pts[0] = Pts->InsertNextPoint(coord);
960   coord[1] = -delta;
961   pts[1] = Pts->InsertNextPoint(coord);
962   // insert line (pt1,pt2)
963   Cells->InsertNextCell(2,pts);
964
965   ZERO_COORD; coord[2] = +delta;
966   pts[0] = Pts->InsertNextPoint(coord);
967   coord[2] = -delta;
968   pts[1] = Pts->InsertNextPoint(coord);
969   // insert line (pt1,pt2)
970   Cells->InsertNextCell(2,pts);
971
972 #undef ZERO_COORD
973 }       
974
975 //=======================================================================
976 // Function : TransferEdgeSData(
977 // Purpose  : Transfert shading data for EDGE
978 //=======================================================================
979
980 void GEOM_OCCReader::TransferEdgeSData(const TopoDS_Edge& aFace,
981                                          vtkPoints* Pts,
982                                          vtkCellArray* Cells) 
983 {
984 }
985
986
987 //=======================================================================
988 // Function : TransferFaceSData
989 // Purpose  : Transfert shading data for FACE
990 //=======================================================================
991 void GEOM_OCCReader::TransferFaceSData(const TopoDS_Face& aFace,
992                                          vtkPoints* Pts,
993                                          vtkCellArray* Cells) {
994
995   TopLoc_Location aLoc;
996   Handle(Poly_Triangulation) aPoly = BRep_Tool::Triangulation(aFace,aLoc);
997   if(aPoly.IsNull()) return;
998   else {
999     
1000     gp_Trsf myTransf;
1001     Standard_Boolean identity = true;
1002     if(!aLoc.IsIdentity())  {
1003       identity = false;
1004       myTransf = aLoc.Transformation();
1005     }
1006
1007     Standard_Integer nbNodesInFace = aPoly->NbNodes();
1008     Standard_Integer nbTriInFace = aPoly->NbTriangles();
1009                 
1010     const Poly_Array1OfTriangle& Triangles = aPoly->Triangles();
1011     const TColgp_Array1OfPnt& Nodes = aPoly->Nodes();
1012       
1013     Standard_Integer i;
1014     for(i=1;i<=nbNodesInFace;i++) {
1015       gp_Pnt P = Nodes(i);
1016       float coord[3];
1017       if(!identity) P.Transform(myTransf);
1018       coord[0] = P.X(); coord[1] = P.Y(); coord[2] = P.Z();
1019       Pts->SetPoint(i-1,coord);
1020     }
1021
1022     for(i=1;i<=nbTriInFace;i++) {
1023       // Get the triangle
1024         
1025       Standard_Integer N1,N2,N3;
1026       Triangles(i).Get(N1,N2,N3);
1027         
1028       int pts[3];
1029       pts[0] = N1-1; pts[1] = N2-1; pts[2] = N3-1;
1030       Cells->InsertNextCell(3,pts);
1031
1032     }
1033   } 
1034 }
1035
1036 //=======================================================================
1037 // Function : ComputeShading
1038 // Purpose  : Compute the shape in shading mode
1039 //=======================================================================
1040 void GEOM_OCCReader::ComputeShading(vtkPoints* Pts,vtkCellArray* Cells){
1041
1042   // Check the type of the shape:
1043   if(myShape.ShapeType() == TopAbs_FACE) {
1044     // Face
1045     TransferFaceSData(TopoDS::Face(myShape),Pts,Cells);
1046   }
1047   else {
1048     if(myShape.ShapeType() == TopAbs_EDGE) {
1049       // Edge
1050       TransferEdgeSData(TopoDS::Edge(myShape),Pts,Cells);
1051     }
1052     else {
1053     }
1054
1055   } 
1056
1057 }
1058
1059 //=======================================================================
1060 // Function : 
1061 // Purpose  : Set parameters
1062 //=======================================================================
1063 void GEOM_OCCReader::setDisplayMode(int thenewmode) {
1064   amode = thenewmode;
1065 }
1066
1067 void GEOM_OCCReader::setTopo(const TopoDS_Shape& aShape, bool isVector) {
1068   myShape = aShape;
1069   myIsVector = isVector;
1070 }
1071
1072 void GEOM_OCCReader::setForceUpdate(Standard_Boolean bol) {
1073   forced = bol;
1074 }
1075
1076 //=======================================================================
1077 // Function : 
1078 // Purpose  : Get parameters
1079 //=======================================================================
1080 const TopoDS_Shape& GEOM_OCCReader::getTopo() {
1081   return myShape;
1082 }
1083
1084 int GEOM_OCCReader::getDisplayMode() {
1085   return amode;
1086 }
1087
1088