Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
28//
29// John Allison 19th July 1996.
30//
31// Class description
32//
33// Abstract interface class for graphics scene handlers.
34// Inherits from G4VGraphicsScene, in the intercoms category, which is
35// a minimal abstract interface for the GEANT4 kernel.
36
37#ifndef G4VSCENEHANDLER_HH
38#define G4VSCENEHANDLER_HH
39
40#include "globals.hh"
41
42#include "G4VGraphicsScene.hh"
43#include "G4ViewerList.hh"
44#include "G4ViewParameters.hh"
45#include "G4THitsMap.hh"
46
47class G4Scene;
48class G4VViewer;
49class G4Colour;
50class G4Visible;
52class G4VModel;
54class G4LogicalVolume;
56class G4Material;
57class G4Event;
58class G4AttHolder;
60
62
63 friend class G4VViewer;
64 friend std::ostream& operator << (std::ostream& os, const G4VSceneHandler& s);
65
66public: // With description
67
69
71 G4int id,
72 const G4String& name = "");
73
74 virtual ~G4VSceneHandler ();
75
76 ///////////////////////////////////////////////////////////////////
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 virtual void AddSolid (const G4TessellatedSolid&);
123
124 // For solids not above.
125 virtual void AddSolid (const G4VSolid&);
126
127 ///////////////////////////////////////////////////////////////////
128 // Methods for adding "compound" GEANT4 objects to the scene
129 // handler. These methods may either (a) invoke "user code" that
130 // uses the "user interface", G4VVisManager (see, for example,
131 // G4VSceneHandler, which for trajectories uses
132 // G4VTrajectory::DrawTrajectory, via G4TrajectoriesModel in the
133 // Modeling Category) or (b) invoke AddPrimitives below (between
134 // calls to Begin/EndPrimitives) or (c) use graphics-system-specific
135 // code or (d) any combination of the above.
136
137 virtual void AddCompound (const G4VTrajectory&);
138 virtual void AddCompound (const G4VHit&);
139 virtual void AddCompound (const G4VDigi&);
140 virtual void AddCompound (const G4THitsMap<G4double>&);
141 virtual void AddCompound (const G4THitsMap<G4StatDouble>&);
142
143 //////////////////////////////////////////////////////////////
144 // Functions for adding primitives.
145
146 virtual void BeginModeling ();
147 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
148 // void MyXXXSceneHandler::BeginModeling () {
149 // G4VSceneHandler::BeginModeling ();
150 // ...
151 // }
152
153 virtual void EndModeling ();
154 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
155 // void MyXXXSceneHandler::EndModeling () {
156 // ...
157 // G4VSceneHandler::EndModeling ();
158 // }
159
160 virtual void BeginPrimitives
161 (const G4Transform3D& objectTransformation = G4Transform3D());
162 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
163 // void MyXXXSceneHandler::BeginPrimitives
164 // (const G4Transform3D& objectTransformation) {
165 // G4VSceneHandler::BeginPrimitives (objectTransformation);
166 // ...
167 // }
168
169 virtual void EndPrimitives ();
170 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
171 // void MyXXXSceneHandler::EndPrimitives () {
172 // ...
173 // G4VSceneHandler::EndPrimitives ();
174 // }
175
176 virtual void BeginPrimitives2D
177 (const G4Transform3D& objectTransformation = G4Transform3D());
178 // The x,y coordinates of the primitives passed to AddPrimitive are
179 // intrepreted as screen coordinates, -1 < x,y < 1. The
180 // z-coordinate is ignored.
181 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
182 // void MyXXXSceneHandler::BeginPrimitives2D
183 // (const G4Transform3D& objectTransformation) {
184 // G4VSceneHandler::BeginPrimitives2D (objectTransformation);
185 // ...
186 // }
187
188 virtual void EndPrimitives2D ();
189 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
190 // void MyXXXSceneHandler::EndPrimitives2D () {
191 // ...
192 // G4VSceneHandler::EndPrimitives2D ();
193 // }
194
195 virtual void AddPrimitive (const G4Polyline&) = 0;
196 virtual void AddPrimitive (const G4Scale&);
197 // Default implementation in this class but can be over-ridden.
198 virtual void AddPrimitive (const G4Text&) = 0;
199 virtual void AddPrimitive (const G4Circle&) = 0;
200 virtual void AddPrimitive (const G4Square&) = 0;
201 virtual void AddPrimitive (const G4Polymarker&);
202 // Default implementation in this class but can be over-ridden.
203 virtual void AddPrimitive (const G4Polyhedron&) = 0;
204
205 virtual const G4VisExtent& GetExtent() const;
206
207 //////////////////////////////////////////////////////////////
208 // Access functions.
209 const G4String& GetName () const;
213 G4Scene* GetScene () const;
222 void SetName (const G4String&);
224 virtual void SetScene (G4Scene*);
225 G4ViewerList& SetViewerList (); // Non-const so you can change.
228 // Sets flag which will cause transient store to be cleared at the
229 // next call to BeginPrimitives(). Maintained by vis manager.
233 // Maintained by vis manager.
234
235 //////////////////////////////////////////////////////////////
236 // Public utility functions.
237
240 // Returns colour - checks fpVisAttribs and gets applicable colour.
241
242 const G4Colour& GetTextColour (const G4Text&);
243 const G4Colour& GetTextColor (const G4Text&);
244 // Returns colour of G4Text object, or default text colour.
245
247 // Returns line width of G4VisAttributes multiplied by GlobalLineWidthScale.
248
250 // Returns drawing style from current view parameters, unless the user
251 // has forced through the vis attributes, thereby over-riding the
252 // current view parameter.
253
255 // Returns no of cloud points from current view parameters, unless the user
256 // has forced through the vis attributes, thereby over-riding the
257 // current view parameter.
258
260 // Returns auxiliary edge visibility from current view parameters,
261 // unless the user has forced through the vis attributes, thereby
262 // over-riding the current view parameter.
263
265 // Returns no. of sides (lines segments per circle) from current
266 // view parameters, unless the user has forced through the vis
267 // attributes, thereby over-riding the current view parameter.
268
270 // Returns applicable marker size (diameter) and type (in second
271 // argument). Uses global default marker if marker sizes are not
272 // set. Multiplies by GlobalMarkerScale.
273
275 // Alias for GetMarkerSize.
276
278 // GetMarkerSize / 2.
279
281 // Only the scene handler and view know what the Modeling Parameters should
282 // be. For historical reasons, the GEANT4 Visualization Environment
283 // maintains its own Scene Data and View Parameters, which must be
284 // converted, when needed, to Modeling Parameters.
285
286 void DrawEvent(const G4Event*);
287 // Checks scene's end-of-event model list and draws trajectories,
288 // hits, etc.
289
290 void DrawEndOfRunModels();
291 // Draws end-of-run models.
292
293 //////////////////////////////////////////////////////////////
294 // Administration functions.
295
296 template <class T> void AddSolidT (const T& solid);
297 template <class T> void AddSolidWithAuxiliaryEdges (const T& solid);
298
300
301 virtual void ClearStore ();
302 // Clears graphics database (display lists) if any.
303
304 virtual void ClearTransientStore ();
305 // Clears transient part of graphics database (display lists) if any.
306
307 void AddViewerToList (G4VViewer* pView); // Add view to view List.
308 void RemoveViewerFromList (G4VViewer* pView); // Remove view from view List.
309
310protected:
311
312 //////////////////////////////////////////////////////////////
313 // Core routine for looping over models, redrawing stored events, etc.
314 // Overload with care (see, for example,
315 // G4OpenGLScenehandler::ProcessScene).
316 virtual void ProcessScene ();
317
318 //////////////////////////////////////////////////////////////
319 // Default routine used by default AddSolid ().
320 virtual void RequestPrimitives (const G4VSolid& solid);
321
322 //////////////////////////////////////////////////////////////
323 // Other internal routines...
324
327 // Generic clipping using the BooleanProcessor in graphics_reps is
328 // implemented in this class. Subclasses that implement their own
329 // clipping should provide an override that returns zero.
330
331 void LoadAtts(const G4Visible&, G4AttHolder*);
332 // Load G4AttValues and G4AttDefs associated with the G4Visible
333 // object onto the G4AttHolder object. It checks fpModel, and also
334 // loads the G4AttValues and G4AttDefs from G4PhysicalVolumeModel,
335 // G4VTrajectory, G4VTrajectoryPoint, G4VHit or G4VDigi, as
336 // appropriate. The G4AttHolder object is an object of a class that
337 // publicly inherits G4AttHolder - see, e.g., SoG4Polyhedron in the
338 // Open Inventor driver. G4AttHolder deletes G4AttValues in its
339 // destructor to ensure proper clean-up of G4AttValues.
340
341 //////////////////////////////////////////////////////////////
342 // Data members
343
344 G4VGraphicsSystem& fSystem; // Graphics system.
345 const G4int fSceneHandlerId; // Id of this instance.
347 G4int fViewCount; // To determine view ids.
349 G4VViewer* fpViewer; // Current viewer.
350 G4Scene* fpScene; // Scene for this scene handler.
352 G4bool fReadyForTransients; // I.e., not processing the
353 // run-duration part of scene.
354 G4bool fTransientsDrawnThisEvent; // Maintained by vis
356 G4bool fProcessingSolid; // True if within Pre/PostAddSolid.
357 G4bool fProcessing2D; // True for 2D.
358 G4VModel* fpModel; // Current model.
359 G4Transform3D fObjectTransformation; // Current accumulated
360 // object transformation.
361 G4int fNestingDepth; // For Begin/EndPrimitives.
362 const G4VisAttributes* fpVisAttribs; // Working vis attributes.
364
365private:
366
368 G4VSceneHandler& operator = (const G4VSceneHandler&);
369};
370
371#include "G4VSceneHandler.icc"
372
373#endif
HepGeom::Transform3D G4Transform3D
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
Definition: G4Box.hh:56
Definition: G4Cons.hh:78
Definition: G4Orb.hh:56
Definition: G4Para.hh:79
Definition: G4Text.hh:72
Definition: G4Trd.hh:63
Definition: G4Tubs.hh:75
Definition: G4VHit.hh:48
virtual void BeginModeling()
G4int GetNumberOfCloudPoints(const G4VisAttributes *) const
friend std::ostream & operator<<(std::ostream &os, const G4VSceneHandler &s)
G4int GetNoOfSides(const G4VisAttributes *)
virtual void AddPrimitive(const G4Polyhedron &)=0
virtual void ClearTransientStore()
G4bool GetTransientsDrawnThisEvent() const
void SetTransientsDrawnThisRun(G4bool)
void LoadAtts(const G4Visible &, G4AttHolder *)
void DrawEvent(const G4Event *)
G4ModelingParameters * CreateModelingParameters()
const G4Colour & GetTextColour(const G4Text &)
void SetMarkForClearingTransientStore(G4bool)
const G4Colour & GetTextColor(const G4Text &)
G4VModel * GetModel() const
const G4Colour & GetColour()
G4VGraphicsSystem * GetGraphicsSystem() const
virtual void AddPrimitive(const G4Circle &)=0
G4double GetMarkerRadius(const G4VMarker &, MarkerSizeType &)
const G4ViewerList & GetViewerList() const
void AddSolidT(const T &solid)
G4bool IsReadyForTransients() const
G4Scene * GetScene() const
G4Transform3D fObjectTransformation
virtual void EndPrimitives()
G4bool fTransientsDrawnThisEvent
G4VViewer * GetCurrentViewer() const
void AddSolidWithAuxiliaryEdges(const T &solid)
G4int IncrementViewCount()
virtual G4DisplacedSolid * CreateSectionSolid()
virtual void AddPrimitive(const G4Text &)=0
G4ViewerList & SetViewerList()
virtual void EndModeling()
const G4int fSceneHandlerId
void SetName(const G4String &)
virtual const G4VisExtent & GetExtent() const
G4int GetSceneHandlerId() const
G4ViewerList fViewerList
virtual void ProcessScene()
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)
void SetObjectTransformation(const G4Transform3D &)
void SetModel(G4VModel *)
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
G4bool fTransientsDrawnThisRun
G4VViewer * fpViewer
virtual void PostAddSolid()
const G4String & GetName() const
void AddViewerToList(G4VViewer *pView)
void SetTransientsDrawnThisEvent(G4bool)
virtual void EndPrimitives2D()
virtual void SetScene(G4Scene *)
G4bool fMarkForClearingTransientStore
const G4Transform3D & GetObjectTransformation() const
const G4VisAttributes * fpVisAttribs
virtual void BeginPrimitives2D(const G4Transform3D &objectTransformation=G4Transform3D())
virtual void RequestPrimitives(const G4VSolid &solid)
G4double GetMarkerDiameter(const G4VMarker &, MarkerSizeType &)
virtual void AddPrimitive(const G4Square &)=0
const G4Colour & GetColor()
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
G4bool GetTransientsDrawnThisRun() const
void RemoveViewerFromList(G4VViewer *pView)
virtual G4DisplacedSolid * CreateCutawaySolid()
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())
G4double GetLineWidth(const G4VisAttributes *)
G4int GetViewCount() const
G4VGraphicsSystem & fSystem
virtual void AddSolid(const G4Box &)
virtual void ClearStore()
void SetCurrentViewer(G4VViewer *)
virtual void AddCompound(const G4VTrajectory &)
virtual ~G4VSceneHandler()
virtual void AddPrimitive(const G4Polyline &)=0
const G4Transform3D fIdentityTransformation
G4bool GetMarkForClearingTransientStore() const
G4bool GetAuxEdgeVisible(const G4VisAttributes *)