Geant4  10.03
G4VSceneHandler.hh
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4VSceneHandler.hh 101714 2016-11-22 08:53:13Z gcosmo $
28 //
29 //
30 // John Allison 19th July 1996.
31 //
32 // Class description
33 //
34 // Abstract interface class for graphics scene handlers.
35 // Inherits from G4VGraphicsScene, in the intercoms category, which is
36 // a minimal abstract interface for the GEANT4 kernel.
37 
38 #ifndef G4VSCENEHANDLER_HH
39 #define G4VSCENEHANDLER_HH
40 
41 #include "globals.hh"
42 
43 #include "G4VGraphicsScene.hh"
44 #include "G4ViewerList.hh"
45 #include "G4ViewParameters.hh"
46 #include "G4THitsMap.hh"
47 
48 class G4Scene;
49 class G4VViewer;
50 class G4Colour;
51 class G4Visible;
53 class G4VModel;
54 class G4VGraphicsSystem;
55 class G4LogicalVolume;
56 class G4VPhysicalVolume;
57 class G4Material;
58 class G4Event;
59 class G4AttHolder;
60 
62 
63  friend class G4VViewer;
64  friend std::ostream& operator << (std::ostream& os, const G4VSceneHandler& s);
65 
66 public: // With description
67 
69 
71  G4int id,
72  const G4String& name = "");
73 
74  virtual ~G4VSceneHandler ();
75 
77  // Methods for adding raw GEANT4 objects to the scene handler. They
78  // must always be called in the triplet PreAddSolid, AddSolid and
79  // PostAddSolid. The transformation and visualization attributes
80  // must be set by the call to PreAddSolid. If your graphics system
81  // is sophisticated enough to handle a particular solid shape as a
82  // primitive, in your derived class write a function to override one
83  // or more of the following. See the implementation of
84  // G4VSceneHandler::AddSolid (const G4Box& box) for more
85  // suggestions. If not, please implement the base class invocation.
86 
87  virtual void PreAddSolid (const G4Transform3D& objectTransformation,
88  const G4VisAttributes&);
89  // objectTransformation is the transformation in the world
90  // coordinate system of the object about to be added, and visAttribs
91  // is its visualization attributes.
92  // IMPORTANT: invoke this from your polymorphic versions, e.g.:
93  // void MyXXXSceneHandler::PreAddSolid
94  // (const G4Transform3D& objectTransformation,
95  // const G4VisAttributes& visAttribs) {
96  // G4VSceneHandler::PreAddSolid (objectTransformation, visAttribs);
97  // ...
98  // }
99 
100  virtual void PostAddSolid ();
101  // IMPORTANT: invoke this from your polymorphic versions, e.g.:
102  // void MyXXXSceneHandler::PostAddSolid () {
103  // ...
104  // G4VSceneHandler::PostAddSolid ();
105  // }
106 
107  // From geometry/solids/CSG
108  virtual void AddSolid (const G4Box&);
109  virtual void AddSolid (const G4Cons&);
110  virtual void AddSolid (const G4Orb&);
111  virtual void AddSolid (const G4Para&);
112  virtual void AddSolid (const G4Sphere&);
113  virtual void AddSolid (const G4Torus&);
114  virtual void AddSolid (const G4Trap&);
115  virtual void AddSolid (const G4Trd&);
116  virtual void AddSolid (const G4Tubs&);
117 
118  // From geometry/solids/specific
119  virtual void AddSolid (const G4Ellipsoid&);
120  virtual void AddSolid (const G4Polycone&);
121  virtual void AddSolid (const G4Polyhedra&);
122 
123  // For solids not above.
124  virtual void AddSolid (const G4VSolid&);
125 
127  // Methods for adding "compound" GEANT4 objects to the scene
128  // handler. These methods may either (a) invoke "user code" that
129  // uses the "user interface", G4VVisManager (see, for example,
130  // G4VSceneHandler, which for trajectories uses
131  // G4VTrajectory::DrawTrajectory, via G4TrajectoriesModel in the
132  // Modeling Category) or (b) invoke AddPrimitives below (between
133  // calls to Begin/EndPrimitives) or (c) use graphics-system-specific
134  // code or (d) any combination of the above.
135 
136  virtual void AddCompound (const G4VTrajectory&);
137  virtual void AddCompound (const G4VHit&);
138  virtual void AddCompound (const G4VDigi&);
139  virtual void AddCompound (const G4THitsMap<G4double>&);
140  virtual void AddCompound (const G4THitsMap<G4StatDouble>&);
141 
143  // Functions for adding primitives.
144 
145  virtual void BeginModeling ();
146  // IMPORTANT: invoke this from your polymorphic versions, e.g.:
147  // void MyXXXSceneHandler::BeginModeling () {
148  // G4VSceneHandler::BeginModeling ();
149  // ...
150  // }
151 
152  virtual void EndModeling ();
153  // IMPORTANT: invoke this from your polymorphic versions, e.g.:
154  // void MyXXXSceneHandler::EndModeling () {
155  // ...
156  // G4VSceneHandler::EndModeling ();
157  // }
158 
159  virtual void BeginPrimitives
160  (const G4Transform3D& objectTransformation);
161  // IMPORTANT: invoke this from your polymorphic versions, e.g.:
162  // void MyXXXSceneHandler::BeginPrimitives
163  // (const G4Transform3D& objectTransformation) {
164  // G4VSceneHandler::BeginPrimitives (objectTransformation);
165  // ...
166  // }
167 
168  virtual void EndPrimitives ();
169  // IMPORTANT: invoke this from your polymorphic versions, e.g.:
170  // void MyXXXSceneHandler::EndPrimitives () {
171  // ...
172  // G4VSceneHandler::EndPrimitives ();
173  // }
174 
175  virtual void BeginPrimitives2D
176  (const G4Transform3D& objectTransformation);
177  // The x,y coordinates of the primitives passed to AddPrimitive are
178  // intrepreted as screen coordinates, -1 < x,y < 1. The
179  // z-coordinate is ignored.
180  // IMPORTANT: invoke this from your polymorphic versions, e.g.:
181  // void MyXXXSceneHandler::BeginPrimitives2D
182  // (const G4Transform3D& objectTransformation) {
183  // G4VSceneHandler::BeginPrimitives2D ();
184  // ...
185  // }
186 
187  virtual void EndPrimitives2D ();
188  // IMPORTANT: invoke this from your polymorphic versions, e.g.:
189  // void MyXXXSceneHandler::EndPrimitives2D () {
190  // ...
191  // G4VSceneHandler::EndPrimitives2D ();
192  // }
193 
194  virtual void AddPrimitive (const G4Polyline&) = 0;
195  virtual void AddPrimitive (const G4Scale&);
196  // Default implementation in this class but can be over-ridden.
197  virtual void AddPrimitive (const G4Text&) = 0;
198  virtual void AddPrimitive (const G4Circle&) = 0;
199  virtual void AddPrimitive (const G4Square&) = 0;
200  virtual void AddPrimitive (const G4Polymarker&);
201  // Default implementation in this class but can be over-ridden.
202  virtual void AddPrimitive (const G4Polyhedron&) = 0;
203 
204  virtual const G4VisExtent& GetExtent() const;
205 
207  // Access functions.
208  const G4String& GetName () const;
209  G4int GetSceneHandlerId () const;
210  G4int GetViewCount () const;
212  G4Scene* GetScene () const;
213  const G4ViewerList& GetViewerList () const;
214  G4VModel* GetModel () const;
215  G4VViewer* GetCurrentViewer () const;
217  G4bool IsReadyForTransients () const;
220  const G4Transform3D& GetObjectTransformation () const;
221  void SetName (const G4String&);
222  void SetCurrentViewer (G4VViewer*);
223  virtual void SetScene (G4Scene*);
224  G4ViewerList& SetViewerList (); // Non-const so you can change.
225  void SetModel (G4VModel*);
227  // Sets flag which will cause transient store to be cleared at the
228  // next call to BeginPrimitives(). Maintained by vis manager.
232  // Maintained by vis manager.
233 
235  // Public utility functions.
236 
237  const G4Colour& GetColour ();
238  const G4Colour& GetColor ();
239  // Returns colour - checks fpVisAttribs and gets applicable colour.
240 
241  const G4Colour& GetTextColour (const G4Text&);
242  const G4Colour& GetTextColor (const G4Text&);
243  // Returns colour of G4Text object, or default text colour.
244 
246  // Returns line width of G4VisAttributes multiplied by GlobalLineWidthScale.
247 
249  // Returns drawing style from current view parameters, unless the user
250  // has forced through the vis attributes, thereby over-riding the
251  // current view parameter.
252 
254  // Returns auxiliary edge visibility from current view parameters,
255  // unless the user has forced through the vis attributes, thereby
256  // over-riding the current view parameter.
257 
259  // Returns no. of sides (lines segments per circle) from current
260  // view parameters, unless the user has forced through the vis
261  // attributes, thereby over-riding the current view parameter.
262 
264  // Returns applicable marker size (diameter) and type (in second
265  // argument). Uses global default marker if marker sizes are not
266  // set. Multiplies by GlobalMarkerScale.
267 
269  // Alias for GetMarkerSize.
270 
272  // GetMarkerSize / 2.
273 
275  // Only the scene handler and view know what the Modeling Parameters should
276  // be. For historical reasons, the GEANT4 Visualization Environment
277  // maintains its own Scene Data and View Parameters, which must be
278  // converted, when needed, to Modeling Parameters.
279 
280  void DrawEvent(const G4Event*);
281  // Checks scene's end-of-event model list and draws trajectories,
282  // hits, etc.
283 
284  void DrawEndOfRunModels();
285  // Draws end-of-run models.
286 
288  // Administration functions.
289 
290  template <class T> void AddSolidT (const T& solid);
291  template <class T> void AddSolidWithAuxiliaryEdges (const T& solid);
292 
294 
295  virtual void ClearStore ();
296  // Clears graphics database (display lists) if any.
297 
298  virtual void ClearTransientStore ();
299  // Clears transient part of graphics database (display lists) if any.
300 
301  void AddViewerToList (G4VViewer* pView); // Add view to view List.
302  void RemoveViewerFromList (G4VViewer* pView); // Remove view from view List.
303 
304 protected:
305 
307  // Core routine for looping over models, redrawing stored events, etc.
308  // Overload with care (see, for example,
309  // G4OpenGLScenehandler::ProcessScene).
310  virtual void ProcessScene ();
311 
313  // Default routine used by default AddSolid ().
314  virtual void RequestPrimitives (const G4VSolid& solid);
315 
317  // Other internal routines...
318 
319  virtual G4VSolid* CreateSectionSolid ();
320  virtual G4VSolid* CreateCutawaySolid ();
321  // Generic clipping using the BooleanProcessor in graphics_reps is
322  // implemented in this class. Subclasses that implement their own
323  // clipping should provide an override that returns zero.
324 
325  void LoadAtts(const G4Visible&, G4AttHolder*);
326  // Load G4AttValues and G4AttDefs associated with the G4Visible
327  // object onto the G4AttHolder object. It checks fpModel, and also
328  // loads the G4AttValues and G4AttDefs from G4PhysicalVolumeModel,
329  // G4VTrajectory, G4VTrajectoryPoint, G4VHit or G4VDigi, as
330  // appropriate. The G4AttHolder object is an object of a class that
331  // publicly inherits G4AttHolder - see, e.g., SoG4Polyhedron in the
332  // Open Inventor driver. G4AttHolder deletes G4AttValues in its
333  // destructor to ensure proper clean-up of G4AttValues.
334 
336  // Data members
337 
338  G4VGraphicsSystem& fSystem; // Graphics system.
339  const G4int fSceneHandlerId; // Id of this instance.
341  G4int fViewCount; // To determine view ids.
343  G4VViewer* fpViewer; // Current viewer.
344  G4Scene* fpScene; // Scene for this scene handler.
346  G4bool fReadyForTransients; // I.e., not processing the
347  // run-duration part of scene.
348  G4bool fTransientsDrawnThisEvent; // Maintained by vis
350  G4bool fProcessingSolid; // True if within Pre/PostAddSolid.
351  G4bool fProcessing2D; // True for 2D.
352  G4VModel* fpModel; // Current model.
353  G4Transform3D fObjectTransformation; // Current accumulated
354  // object transformation.
355  G4int fNestingDepth; // For Begin/EndPrimitives.
356  const G4VisAttributes* fpVisAttribs; // Working vis attributes.
358 
359 private:
360 
363 };
364 
365 #include "G4VSceneHandler.icc"
366 
367 #endif
G4VModel * GetModel() const
virtual void ClearStore()
virtual ~G4VSceneHandler()
Definition: G4Para.hh:77
virtual void AddSolid(const G4Box &)
Definition: G4Text.hh:73
virtual G4VSolid * CreateSectionSolid()
G4double GetMarkerDiameter(const G4VMarker &, MarkerSizeType &)
static constexpr double s
Definition: G4SIunits.hh:169
virtual void BeginModeling()
G4ModelingParameters * CreateModelingParameters()
virtual void BeginPrimitives(const G4Transform3D &objectTransformation)
G4ViewerList & SetViewerList()
G4int GetViewCount() const
G4bool IsReadyForTransients() const
void LoadAtts(const G4Visible &, G4AttHolder *)
Definition: G4Box.hh:64
G4VViewer * fpViewer
G4VViewer * GetCurrentViewer() const
Definition: G4Tubs.hh:85
G4int IncrementViewCount()
void RemoveViewerFromList(G4VViewer *pView)
virtual void PostAddSolid()
const G4String & GetName() const
G4Transform3D fObjectTransformation
Definition: G4VHit.hh:48
const char * name(G4int ptype)
Definition: G4Trd.hh:72
G4int GetNoOfSides(const G4VisAttributes *)
G4bool GetMarkForClearingTransientStore() const
virtual const G4VisExtent & GetExtent() const
int G4int
Definition: G4Types.hh:78
virtual void AddPrimitive(const G4Polyline &)=0
virtual G4VSolid * CreateCutawaySolid()
G4VSceneHandler(G4VGraphicsSystem &system, G4int id, const G4String &name="")
const G4Transform3D & GetObjectTransformation() const
const G4ViewerList & GetViewerList() const
void SetTransientsDrawnThisEvent(G4bool)
G4double GetLineWidth(const G4VisAttributes *)
const G4int fSceneHandlerId
void AddViewerToList(G4VViewer *pView)
const G4Colour & GetColour()
bool G4bool
Definition: G4Types.hh:79
friend std::ostream & operator<<(std::ostream &os, const G4VSceneHandler &s)
const G4Transform3D fIdentityTransformation
Definition: G4Cons.hh:83
void SetName(const G4String &)
G4bool GetTransientsDrawnThisRun() const
G4double GetMarkerRadius(const G4VMarker &, MarkerSizeType &)
G4ViewerList fViewerList
virtual void EndModeling()
G4VGraphicsSystem * GetGraphicsSystem() const
virtual void EndPrimitives()
HepGeom::Transform3D G4Transform3D
G4bool GetTransientsDrawnThisEvent() const
void SetObjectTransformation(const G4Transform3D &)
Definition: G4Orb.hh:61
virtual void BeginPrimitives2D(const G4Transform3D &objectTransformation)
virtual void SetScene(G4Scene *)
const G4Colour & GetColor()
const G4VisAttributes * fpVisAttribs
G4Scene * GetScene() const
G4bool fMarkForClearingTransientStore
virtual void AddCompound(const G4VTrajectory &)
void DrawEvent(const G4Event *)
virtual void ProcessScene()
G4bool fTransientsDrawnThisEvent
const G4Colour & GetTextColor(const G4Text &)
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)
G4bool GetAuxEdgeVisible(const G4VisAttributes *)
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
virtual void RequestPrimitives(const G4VSolid &solid)
void AddSolidWithAuxiliaryEdges(const T &solid)
void SetModel(G4VModel *)
double G4double
Definition: G4Types.hh:76
virtual void ClearTransientStore()
G4bool fTransientsDrawnThisRun
void SetTransientsDrawnThisRun(G4bool)
G4int GetSceneHandlerId() const
G4VGraphicsSystem & fSystem
G4VSceneHandler & operator=(const G4VSceneHandler &)
void SetCurrentViewer(G4VViewer *)
void SetMarkForClearingTransientStore(G4bool)
virtual void EndPrimitives2D()
void AddSolidT(const T &solid)
const G4Colour & GetTextColour(const G4Text &)