Geant4 11.1.1
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#include "G4PseudoScene.hh"
47
48class G4Scene;
50class G4AttHolder;
51
53
54 friend class G4VViewer;
55 friend std::ostream& operator << (std::ostream& os, const G4VSceneHandler& s);
56
57public: // With description
58
60
62 G4int id,
63 const G4String& name = "");
64
65 virtual ~G4VSceneHandler ();
66
67 ///////////////////////////////////////////////////////////////////
68 // Methods for adding raw GEANT4 objects to the scene handler. They
69 // must always be called in the triplet PreAddSolid, AddSolid and
70 // PostAddSolid. The transformation and visualization attributes
71 // must be set by the call to PreAddSolid. If your graphics system
72 // is sophisticated enough to handle a particular solid shape as a
73 // primitive, in your derived class write a function to override one
74 // or more of the following. See the implementation of
75 // G4VSceneHandler::AddSolid (const G4Box& box) for more
76 // suggestions. If not, please implement the base class invocation.
77
78 virtual void PreAddSolid (const G4Transform3D& objectTransformation,
79 const G4VisAttributes&);
80 // objectTransformation is the transformation in the world
81 // coordinate system of the object about to be added, and visAttribs
82 // is its visualization attributes.
83 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
84 // void MyXXXSceneHandler::PreAddSolid
85 // (const G4Transform3D& objectTransformation,
86 // const G4VisAttributes& visAttribs) {
87 // G4VSceneHandler::PreAddSolid (objectTransformation, visAttribs);
88 // ...
89 // }
90
91 virtual void PostAddSolid ();
92 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
93 // void MyXXXSceneHandler::PostAddSolid () {
94 // ...
95 // G4VSceneHandler::PostAddSolid ();
96 // }
97
98 // From geometry/solids/CSG
99 virtual void AddSolid (const G4Box&);
100 virtual void AddSolid (const G4Cons&);
101 virtual void AddSolid (const G4Orb&);
102 virtual void AddSolid (const G4Para&);
103 virtual void AddSolid (const G4Sphere&);
104 virtual void AddSolid (const G4Torus&);
105 virtual void AddSolid (const G4Trap&);
106 virtual void AddSolid (const G4Trd&);
107 virtual void AddSolid (const G4Tubs&);
108
109 // From geometry/solids/specific
110 virtual void AddSolid (const G4Ellipsoid&);
111 virtual void AddSolid (const G4Polycone&);
112 virtual void AddSolid (const G4Polyhedra&);
113 virtual void AddSolid (const G4TessellatedSolid&);
114
115 // For solids not above.
116 virtual void AddSolid (const G4VSolid&);
117
118 ///////////////////////////////////////////////////////////////////
119 // Methods for adding "compound" GEANT4 objects to the scene
120 // handler. These methods may either (a) invoke "user code" that
121 // uses the "user interface", G4VVisManager (see, for example,
122 // G4VSceneHandler, which for trajectories uses
123 // G4VTrajectory::DrawTrajectory, via G4TrajectoriesModel in the
124 // Modeling Category) or (b) invoke AddPrimitives below (between
125 // calls to Begin/EndPrimitives) or (c) use graphics-system-specific
126 // code or (d) any combination of the above.
127
128 virtual void AddCompound (const G4VTrajectory&);
129 virtual void AddCompound (const G4VHit&);
130 virtual void AddCompound (const G4VDigi&);
131 virtual void AddCompound (const G4THitsMap<G4double>&);
132 virtual void AddCompound (const G4THitsMap<G4StatDouble>&);
133 virtual void AddCompound (const G4Mesh&);
134
135 //////////////////////////////////////////////////////////////
136 // Functions for adding primitives.
137
138 virtual void BeginModeling ();
139 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
140 // void MyXXXSceneHandler::BeginModeling () {
141 // G4VSceneHandler::BeginModeling ();
142 // ...
143 // }
144
145 virtual void EndModeling ();
146 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
147 // void MyXXXSceneHandler::EndModeling () {
148 // ...
149 // G4VSceneHandler::EndModeling ();
150 // }
151
152 virtual void BeginPrimitives
153 (const G4Transform3D& objectTransformation = G4Transform3D());
154 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
155 // void MyXXXSceneHandler::BeginPrimitives
156 // (const G4Transform3D& objectTransformation) {
157 // G4VSceneHandler::BeginPrimitives (objectTransformation);
158 // ...
159 // }
160
161 virtual void EndPrimitives ();
162 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
163 // void MyXXXSceneHandler::EndPrimitives () {
164 // ...
165 // G4VSceneHandler::EndPrimitives ();
166 // }
167
168 virtual void BeginPrimitives2D
169 (const G4Transform3D& objectTransformation = G4Transform3D());
170 // The x,y coordinates of the primitives passed to AddPrimitive are
171 // intrepreted as screen coordinates, -1 < x,y < 1. The
172 // z-coordinate is ignored.
173 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
174 // void MyXXXSceneHandler::BeginPrimitives2D
175 // (const G4Transform3D& objectTransformation) {
176 // G4VSceneHandler::BeginPrimitives2D (objectTransformation);
177 // ...
178 // }
179
180 virtual void EndPrimitives2D ();
181 // IMPORTANT: invoke this from your polymorphic versions, e.g.:
182 // void MyXXXSceneHandler::EndPrimitives2D () {
183 // ...
184 // G4VSceneHandler::EndPrimitives2D ();
185 // }
186
187 virtual void AddPrimitive (const G4Polyline&) = 0;
188 virtual void AddPrimitive (const G4Text&) = 0;
189 virtual void AddPrimitive (const G4Circle&) = 0;
190 virtual void AddPrimitive (const G4Square&) = 0;
191 virtual void AddPrimitive (const G4Polymarker&);
192 virtual void AddPrimitive (const G4Polyhedron&) = 0;
193 virtual void AddPrimitive (const G4Plotter&);
194
195 // Other virtual functions
196 virtual const G4VisExtent& GetExtent() const;
197
198 //////////////////////////////////////////////////////////////
199 // Access functions.
200 const G4String& GetName () const;
204 G4Scene* GetScene () const;
213 void SetName (const G4String&);
215 virtual void SetScene (G4Scene*);
216 G4ViewerList& SetViewerList (); // Non-const so you can change.
219 // Sets flag which will cause transient store to be cleared at the
220 // next call to BeginPrimitives(). Maintained by vis manager.
224 // Maintained by vis manager.
225
226 //////////////////////////////////////////////////////////////
227 // Public utility functions.
228
229 const G4Colour& GetColour (); // To be deprecated?
231 // Returns colour - checks fpVisAttribs and gets applicable colour.
232 // Assumes fpVisAttribs point to the G4VisAttributes of the current object.
233 // If the pointer is null, the colour is obtained from the default view
234 // parameters of the current viewer.
235
236 const G4Colour& GetColour (const G4Visible&);
237 const G4Colour& GetColor (const G4Visible&);
238 // Returns colour, or viewer default colour.
239 // Makes no assumptions about the validity of data member fpVisAttribs.
240 // If the G4Visible has no vis attributes, i.e., the pointer is null,
241 // the colour is obtained from the default view parameters of the
242 // current viewer.
243
244 const G4Colour& GetTextColour (const G4Text&);
245 const G4Colour& GetTextColor (const G4Text&);
246 // Returns colour of G4Text object, or default text colour.
247
249 // Returns line width of G4VisAttributes multiplied by GlobalLineWidthScale.
250
252 // Returns drawing style from current view parameters, unless the user
253 // has forced through the vis attributes, thereby over-riding the
254 // current view parameter.
255
257 // Returns no of cloud points from current view parameters, unless the user
258 // has forced through the vis attributes, thereby over-riding the
259 // current view parameter.
260
262 // Returns auxiliary edge visibility from current view parameters,
263 // unless the user has forced through the vis attributes, thereby
264 // over-riding the current view parameter.
265
267 // Returns no. of sides (lines segments per circle) from current
268 // view parameters, unless the user has forced through the vis
269 // attributes, thereby over-riding the current view parameter.
270
272 // Returns applicable marker size (diameter) and type (in second
273 // argument). Uses global default marker if marker sizes are not
274 // set. Multiplies by GlobalMarkerScale.
275
277 // Alias for GetMarkerSize.
278
280 // GetMarkerSize / 2.
281
283 // Only the scene handler and view know what the Modeling Parameters should
284 // be. For historical reasons, the GEANT4 Visualization Environment
285 // maintains its own Scene Data and View Parameters, which must be
286 // converted, when needed, to Modeling Parameters.
287
288 void DrawEvent(const G4Event*);
289 // Checks scene's end-of-event model list and draws trajectories,
290 // hits, etc.
291
292 void DrawEndOfRunModels();
293 // Draws end-of-run models.
294
295 //////////////////////////////////////////////////////////////
296 // Administration functions.
297
298 template <class T> void AddSolidT (const T& solid);
299 template <class T> void AddSolidWithAuxiliaryEdges (const T& solid);
300
302
303 virtual void ClearStore ();
304 // Clears graphics database (display lists) if any.
305
306 virtual void ClearTransientStore ();
307 // Clears transient part of graphics database (display lists) if any.
308
309 void AddViewerToList (G4VViewer* pView); // Add view to view List.
310 void RemoveViewerFromList (G4VViewer* pView); // Remove view from view List.
311
312protected:
313
314 //////////////////////////////////////////////////////////////
315 // Core routine for looping over models, redrawing stored events, etc.
316 // Overload with care (see, for example,
317 // G4OpenGLScenehandler::ProcessScene).
318 virtual void ProcessScene ();
319
320 //////////////////////////////////////////////////////////////
321 // Default routine used by default AddSolid ().
322 virtual void RequestPrimitives (const G4VSolid& solid);
323
324 //////////////////////////////////////////////////////////////
325 // Other internal routines...
326
329 // Generic clipping using the BooleanProcessor in graphics_reps is
330 // implemented in this class. Subclasses that implement their own
331 // clipping should provide an override that returns zero.
332
333 void LoadAtts(const G4Visible&, G4AttHolder*);
334 // Load G4AttValues and G4AttDefs associated with the G4Visible
335 // object onto the G4AttHolder object. It checks fpModel, and also
336 // loads the G4AttValues and G4AttDefs from G4PhysicalVolumeModel,
337 // G4VTrajectory, G4VTrajectoryPoint, G4VHit or G4VDigi, as
338 // appropriate. The G4AttHolder object is an object of a class that
339 // publicly inherits G4AttHolder - see, e.g., SoG4Polyhedron in the
340 // Open Inventor driver. G4AttHolder deletes G4AttValues in its
341 // destructor to ensure proper clean-up of G4AttValues.
342
343 //////////////////////////////////////////////////////////////
344 // Special mesh rendering utilities...
345
347 NameAndVisAtts(const G4String& name = "",const G4VisAttributes& visAtts = G4VisAttributes())
348 : fName(name),fVisAtts(visAtts) {}
351 };
352
354 public:
356 (G4PhysicalVolumeModel* pvModel // input
357 , G4int depth // input...the following are outputs by reference
358 , std::multimap<const G4Material*,const G4ThreeVector>& positionByMaterial
359 , std::map<const G4Material*,NameAndVisAtts>& nameAndVisAttsByMaterial)
360 : fpPVModel(pvModel)
361 , fDepth(depth)
362 , fPositionByMaterial(positionByMaterial)
363 , fNameAndVisAttsByMaterial(nameAndVisAttsByMaterial)
364 {}
365 private:
366 using G4PseudoScene::AddSolid; // except for...
367 void AddSolid(const G4Box&) override;
368 void ProcessVolume(const G4VSolid&) override {
369 // Do nothing if uninteresting solids found, e.g., the container if not marked invisible.
370 }
371 G4PhysicalVolumeModel* fpPVModel;
372 G4int fDepth;
373 std::multimap<const G4Material*,const G4ThreeVector>& fPositionByMaterial;
374 std::map<const G4Material*,NameAndVisAtts>& fNameAndVisAttsByMaterial;
375 };
376
378 public:
380 (G4PhysicalVolumeModel* pvModel // input
381 , G4int depth // input...the following are outputs by reference
382 , std::multimap<const G4Material*,std::vector<G4ThreeVector>>& verticesByMaterial
383 , std::map<const G4Material*,NameAndVisAtts>& nameAndVisAttsByMaterial)
384 : fpPVModel(pvModel)
385 , fDepth(depth)
386 , fVerticesByMaterial(verticesByMaterial)
387 , fNameAndVisAttsByMaterial(nameAndVisAttsByMaterial)
388 {}
389 private:
390 using G4PseudoScene::AddSolid; // except for...
391 void AddSolid(const G4VSolid& solid) override;
392 void ProcessVolume(const G4VSolid&) override {
393 // Do nothing if uninteresting solids found, e.g., the container if not marked invisible.
394 }
395 G4PhysicalVolumeModel* fpPVModel;
396 G4int fDepth;
397 std::multimap<const G4Material*,std::vector<G4ThreeVector>>& fVerticesByMaterial;
398 std::map<const G4Material*,NameAndVisAtts>& fNameAndVisAttsByMaterial;
399 };
400
402 // Standard way of special mesh rendering.
403 // MySceneHandler::AddCompound(const G4Mesh& mesh) may use this if
404 // appropriate or implement its own special mesh rendereing.
405
406 void Draw3DRectMeshAsDots(const G4Mesh&);
407 // For a rectangular 3-D mesh, draw as coloured dots by colour and material,
408 // one dot randomly placed in each visible mesh cell.
409
410 void Draw3DRectMeshAsSurfaces(const G4Mesh&);
411 // For a rectangular 3-D mesh, draw as surfaces by colour and material
412 // with inner shared faces removed.
413
414 void DrawTetMeshAsDots(const G4Mesh&);
415 // For a tetrahedron mesh, draw as coloured dots by colour and material,
416 // one dot randomly placed in each visible mesh cell.
417
418 void DrawTetMeshAsSurfaces(const G4Mesh&);
419 // For a tetrahedron mesh, draw as surfaces by colour and material
420 // with inner shared faces removed.
421
423 G4double halfX,
424 G4double halfY,
425 G4double halfZ) const;
426 // Sample a random point inside the box
427
428 G4ThreeVector GetPointInTet(const std::vector<G4ThreeVector>& vertices) const;
429 // Sample a random point inside the tetrahedron
430
431 //////////////////////////////////////////////////////////////
432 // Data members
433
434 G4VGraphicsSystem& fSystem; // Graphics system.
435 const G4int fSceneHandlerId; // Id of this instance.
437 G4int fViewCount; // To determine view ids.
439 G4VViewer* fpViewer; // Current viewer.
440 G4Scene* fpScene; // Scene for this scene handler.
442 G4bool fReadyForTransients; // I.e., not processing the
443 // run-duration part of scene.
444 G4bool fTransientsDrawnThisEvent; // Maintained by vis
446 G4bool fProcessingSolid; // True if within Pre/PostAddSolid.
447 G4bool fProcessing2D; // True for 2D.
448 G4VModel* fpModel; // Current model.
449 G4Transform3D fObjectTransformation; // Current accumulated
450 // object transformation.
451 G4int fNestingDepth; // For Begin/EndPrimitives.
452 const G4VisAttributes* fpVisAttribs; // Working vis attributes.
454
455private:
456
458 G4VSceneHandler& operator = (const G4VSceneHandler&);
459};
460
461#include "G4VSceneHandler.icc"
462
463#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: G4Mesh.hh:48
Definition: G4Orb.hh:56
Definition: G4Para.hh:79
void AddSolid(const G4Box &solid)
Definition: G4Text.hh:72
Definition: G4Trd.hh:63
Definition: G4Tubs.hh:75
Definition: G4VHit.hh:48
PseudoSceneFor3DRectMeshPositions(G4PhysicalVolumeModel *pvModel, G4int depth, std::multimap< const G4Material *, const G4ThreeVector > &positionByMaterial, std::map< const G4Material *, NameAndVisAtts > &nameAndVisAttsByMaterial)
PseudoSceneForTetVertices(G4PhysicalVolumeModel *pvModel, G4int depth, std::multimap< const G4Material *, std::vector< G4ThreeVector > > &verticesByMaterial, std::map< const G4Material *, NameAndVisAtts > &nameAndVisAttsByMaterial)
virtual void BeginModeling()
G4int GetNumberOfCloudPoints(const G4VisAttributes *) const
const G4Colour & GetColor(const G4Visible &)
friend std::ostream & operator<<(std::ostream &os, const G4VSceneHandler &s)
G4int GetNoOfSides(const G4VisAttributes *)
virtual void AddPrimitive(const G4Polyhedron &)=0
void DrawTetMeshAsSurfaces(const G4Mesh &)
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
void Draw3DRectMeshAsDots(const G4Mesh &)
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())
G4ThreeVector GetPointInBox(const G4ThreeVector &pos, G4double halfX, G4double halfY, G4double halfZ) const
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()
void DrawTetMeshAsDots(const G4Mesh &)
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())
G4double GetLineWidth(const G4VisAttributes *)
G4ThreeVector GetPointInTet(const std::vector< G4ThreeVector > &vertices) const
G4int GetViewCount() const
G4VGraphicsSystem & fSystem
virtual void AddSolid(const G4Box &)
virtual void ClearStore()
void Draw3DRectMeshAsSurfaces(const G4Mesh &)
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 *)
void StandardSpecialMeshRendering(const G4Mesh &)
NameAndVisAtts(const G4String &name="", const G4VisAttributes &visAtts=G4VisAttributes())