]> SALOME platform Git repositories - modules/visu.git/blob - src/VISU_I/VISU_Prs3d_i.hh
Salome HOME
Join modifications from BR_Dev_For_4_0 tag V4_1_1.
[modules/visu.git] / src / VISU_I / VISU_Prs3d_i.hh
1 //  VISU OBJECT : interactive object for VISU entities implementation
2 //
3 //  Copyright (C) 2003  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 //
23 //  File   : VISU_PrsObject_i.hxx
24 //  Author : Alexey PETROV
25 //  Module : VISU
26
27 #ifndef VISU_Prs3d_i_HeaderFile
28 #define VISU_Prs3d_i_HeaderFile
29
30 #include "VISU_PrsObject_i.hh"
31
32 #include "VISU_ActorFactory.h"
33 #include "VISU_ConvertorDef.hxx"
34
35 #include "SALOME_GenericObj_i.hh"
36 #include "SALOME_GenericObjPointer.hh"
37 #include "SALOME_InteractiveObject.hxx"
38
39 #include "VTKViewer.h"
40
41 #include <vtkSmartPointer.h>
42 #include <vtkTimeStamp.h>
43
44 class VISU_PipeLine;
45 class VISU_Actor;
46
47 class vtkPlane;
48 class vtkActorCollection;
49 class vtkDataSet;
50
51 namespace VISU
52 {
53   //----------------------------------------------------------------------------
54   class Result_i;
55
56
57   //----------------------------------------------------------------------------
58   struct TResultObserver: public virtual boost::signalslib::trackable
59   {
60     virtual
61     void
62     UpdateFromResult(Result_i* theResult) = 0;
63   };
64
65
66   //----------------------------------------------------------------------------
67   //! Base class for all VTK 3D presentations.
68   /*!
69     It is a root class for a middle level of VISU functionality.
70     Almost all functionality of the the class implemented through redirection 
71     external requestes to its VISU_PipeLine.
72     It define an interface and implement the following topics:
73     - provide persistent mechanism;
74     - implement basic actor management (CreateActor, UpdateActor, UpdateActors, RemoveActor and RemoveActors);
75     - implement common 3D functionality like "clipping planes" and offset.
76   */
77   class VISU_I_EXPORT Prs3d_i : public virtual POA_VISU::Prs3d,
78                                 public virtual SALOME::GenericObj_i,
79                                 public virtual TActorFactory,
80                                 public virtual PrsObject_i
81
82   {
83     Prs3d_i(const Prs3d_i&);
84
85   public:
86     //----------------------------------------------------------------------------
87     typedef PrsObject_i TSuperClass;
88     typedef VISU::Prs3d TInterface;
89
90     //----------------------------------------------------------------------------
91     //! A constructor to create a fresh instance of the class
92     Prs3d_i();
93
94     //! To create a deep copy from another instance of the class
95     virtual
96     void
97     SameAs(const Prs3d_i* theOrigin);
98
99     virtual
100     ~Prs3d_i();
101
102     //----------------------------------------------------------------------------
103     virtual
104     CORBA::Boolean 
105     Apply(bool theReInit);
106
107     //----------------------------------------------------------------------------
108     void 
109     SetCResult(Result_i* theResult);
110
111     Result_i* 
112     GetCResult() const;
113
114     virtual
115     void 
116     SetResultObject(VISU::Result_ptr theResult);
117
118     virtual
119     VISU::Result_ptr
120     GetResultObject();
121
122     //----------------------------------------------------------------------------
123     virtual
124     void 
125     SetMeshName(const char* theMeshName);
126
127     virtual
128     char*
129     GetMeshName();
130
131     std::string
132     GetCMeshName() const;
133
134     //----------------------------------------------------------------------------
135     //! To generate an unique type name for the class (used into persistent functionality)
136     virtual
137     const char* 
138     GetComment() const = 0;
139
140     //! To generate an unique name for the instance of the class
141     virtual
142     QString
143     GenerateName() = 0;
144
145     //! To save paramters of the instance to std::ostringstream
146     virtual
147     void
148     ToStream(std::ostringstream& theStr);
149
150     //! To restore paramters of the instance from Storable::TRestoringMap
151     virtual
152     Storable* 
153     Restore(SALOMEDS::SObject_ptr theSObject,
154             const Storable::TRestoringMap& theMap);
155
156     //----------------------------------------------------------------------------
157     //! Get corresponding SALOMEDS::SObject
158     virtual
159     SALOMEDS::SObject_var 
160     GetSObject();
161
162     //----------------------------------------------------------------------------
163     //! To update is internal state
164     virtual 
165     void
166     Update();
167
168     //! To remove the instance from study
169     virtual
170     void
171     RemoveFromStudy();
172
173     //----------------------------------------------------------------------------
174     //! Get corresponding VISU_PipeLine
175     VISU_PipeLine* 
176     GetPipeLine() const;
177
178     bool
179     IsPipeLineExists();
180
181     //! Get input of the VISU_PipeLine
182     vtkDataSet* 
183     GetInput();
184
185     //----------------------------------------------------------------------------
186     //! To define a way to create VTK representation of the instance
187     virtual 
188     VISU_Actor* 
189     CreateActor() = 0;
190
191     //! To unregister the pointed actor
192     virtual 
193     void
194     RemoveActor(VISU_Actor* theActor);
195
196     //! To unregister all actors of the instance
197     virtual 
198     void
199     RemoveActors();
200
201     //! To update the pointed actor
202     virtual 
203     void
204     UpdateActor(VISU_Actor* theActor);
205
206     //! To update all actors of the instance
207     virtual 
208     void
209     UpdateActors();
210
211     //----------------------------------------------------------------------------
212     // Clipping planes
213     void
214     RemoveAllClippingPlanes();
215
216     bool
217     AddClippingPlane(vtkPlane* thePlane);
218
219     vtkIdType
220     GetNumberOfClippingPlanes() const;
221
222     vtkPlane* 
223     GetClippingPlane(vtkIdType theID) const;
224
225     void
226     SetPlaneParam(vtkFloatingPointType theDir[3], 
227                   vtkFloatingPointType theDist, 
228                   vtkPlane* thePlane);
229
230     //----------------------------------------------------------------------------
231     void
232     GetBounds(vtkFloatingPointType aBounds[6]);
233
234     int
235     GetNumberOfActors ();
236
237     //! Move the 3D presentation according to the given offset parameters
238     void
239     SetOffset(const CORBA::Float* theOffsets);
240
241     //! Move the 3D presentation according to the given offset parameters
242     virtual
243     void
244     SetOffset(CORBA::Float theDx, 
245               CORBA::Float theDy, 
246               CORBA::Float theDz);
247
248     //! Gets offset parameters for the 3D presentation
249     void
250     GetOffset(CORBA::Float* theOffsets);
251
252     //! Gets offset parameters for the 3D presentation
253     virtual
254     void
255     GetOffset(CORBA::Float& theDx, 
256               CORBA::Float& theDy, 
257               CORBA::Float& theDz);
258
259     //----------------------------------------------------------------------------
260     //! Gets memory size actually used by the presentation (Mb).
261     virtual
262     CORBA::Float
263     GetMemorySize();
264
265     //----------------------------------------------------------------------------
266     //! Gets know whether the factory instance can be used for actor management or not
267     virtual
268     bool 
269     GetActiveState();
270
271     //----------------------------------------------------------------------------
272     //! Return modified time of the presentation
273     virtual
274     unsigned long int 
275     GetMTime();
276
277     //----------------------------------------------------------------------------
278     //! Create and return the interactive object
279     virtual
280     Handle(SALOME_InteractiveObject)
281     GetIO();
282
283     //! Used in derived classes to initilize the IO for actors
284     virtual
285     std::string
286     GetActorEntry();
287
288   protected:
289     /*! 
290       Used in Apply method to get know whether it is possible to create presentation
291       with the input parameters or not. The derived classes can use this method 
292       to customize Apply behaviour.
293     */
294     virtual 
295     bool 
296     SetInput(bool theReInit);
297
298     //! Restore input parameters if Apply function fails
299     virtual 
300     void 
301     OnRestoreInput();
302
303     //! Used in derived classes to initilize the myPipeLine member
304     void
305     SetPipeLine(VISU_PipeLine* thePipeLine);
306
307     //! Used in derived classes to initilize the myPipeLine member
308     void
309     CreateActor(VISU_Actor* theActor);
310
311     //! Gets or creates VISU_PipeLine instance for actor initilization
312     virtual 
313     VISU_PipeLine* 
314     GetActorPipeLine();
315
316     //! To check dataset validity, throws std::exception if not valid
317     virtual
318     void
319     CheckDataSet();
320
321   protected:
322     vtkTimeStamp myUpdateTime;
323     vtkTimeStamp myParamsTime;
324
325   private:
326     void
327     SetResultEntry(const std::string& theResultEntry);
328
329     std::string
330     GetResultEntry();
331
332     typedef SALOME::GenericObjPtr<VISU::Result_i> TResultPtr;
333     TResultPtr myResult;
334     TResultPtr myPreviousResult;
335
336     std::string myMeshName;
337     std::string myPreviousMeshName;
338
339     CORBA::Float myOffset[3];
340
341     boost::signal0<void> myUpdateActorsSignal;
342     boost::signal0<void> myRemoveActorsFromRendererSignal;
343     vtkSmartPointer<vtkActorCollection> myActorCollection;
344
345     vtkSmartPointer<VISU_PipeLine> myPipeLine;
346
347     Handle(SALOME_InteractiveObject) myIO;
348
349   private:
350     friend class ColoredPrs3dCache_i;
351
352     //! Sets activity flag for the factory instance
353     void
354     SetActiveState(bool theState);
355     
356     bool myIsActiveSatate;
357   };
358
359   //----------------------------------------------------------------------------
360 }
361
362 #endif