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