Geant4-11
G4RichTrajectory.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// G4RichTrajectory
27//
28// Class description:
29//
30// This class extends G4Trajectory, which includes the following:
31// 1) List of trajectory points which compose the trajectory,
32// 2) static information of particle which generated the
33// trajectory,
34// 3) trackID and parent particle ID of the trajectory.
35// The extended information, only publicly accessible through AttValues,
36// includes:
37// 4) physical volume and next physical volume;
38// 5) creator process;
39// ...and much more.
40
41// Contact:
42// Questions and comments on G4Trajectory should be sent to
43// Katsuya Amako (e-mail: Katsuya.Amako@kek.jp)
44// Makoto Asai (e-mail: asai@slac.stanford.edu)
45// Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
46// and on the extended code to:
47// John Allison (e-mail: John.Allison@manchester.ac.uk)
48// Joseph Perl (e-mail: perl@slac.stanford.edu)
49// --------------------------------------------------------------------
50#ifndef G4RICHTRAJECTORY_HH
51#define G4RICHTRAJECTORY_HH
52
53#include "trkgdefs.hh"
54#include "G4Trajectory.hh"
55#include "G4TouchableHandle.hh"
56
58{
59 using RichTrajectoryPointsContainer = std::vector<G4VTrajectoryPoint*>;
60
61 public:
62
63 // Constructors/destructor
64 //
66 G4RichTrajectory(const G4Track* aTrack);
68 virtual ~G4RichTrajectory();
69
70 // Operators
71 //
73 G4int operator == (const G4RichTrajectory& r) const { return (this==&r); }
74
75 inline void* operator new(size_t);
76 inline void operator delete(void*);
77
78 // Other (virtual) member functions
79 //
80 void ShowTrajectory(std::ostream& os=G4cout) const;
81 void DrawTrajectory() const;
82 void AppendStep(const G4Step* aStep);
83 void MergeTrajectory(G4VTrajectory* secondTrajectory);
84 inline G4int GetPointEntries() const;
85 inline G4VTrajectoryPoint* GetPoint(G4int i) const;
86
87 // Get methods for HepRep style attributes
88 //
89 virtual const std::map<G4String, G4AttDef>* GetAttDefs() const;
90 virtual std::vector<G4AttValue>* CreateAttValues() const;
91
92 private:
93
94 // Extended information (only publicly accessible through AttValues)...
95 //
99 const G4VProcess* fpCreatorProcess = nullptr;
103 const G4VProcess* fpEndingProcess = nullptr;
105};
106
108
109inline void* G4RichTrajectory::operator new(size_t)
110{
111 if (aRichTrajectoryAllocator() == nullptr)
112 {
114 }
115 return (void*)aRichTrajectoryAllocator()->MallocSingle();
116}
117
118inline void G4RichTrajectory::operator delete(void* aRichTrajectory)
119{
120 aRichTrajectoryAllocator()->FreeSingle((G4RichTrajectory*)aRichTrajectory);
121}
122
124{
125 return G4int(fpRichPointsContainer->size());
126}
127
129{
130 return (*fpRichPointsContainer)[i];
131}
132
133#endif
G4TRACKING_DLL G4Allocator< G4RichTrajectory > *& aRichTrajectoryAllocator()
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
virtual ~G4RichTrajectory()
G4TouchableHandle fpFinalVolume
void ShowTrajectory(std::ostream &os=G4cout) const
G4TouchableHandle fpFinalNextVolume
G4TouchableHandle fpInitialNextVolume
void AppendStep(const G4Step *aStep)
G4double fFinalKineticEnergy
G4VTrajectoryPoint * GetPoint(G4int i) const
std::vector< G4VTrajectoryPoint * > RichTrajectoryPointsContainer
virtual std::vector< G4AttValue > * CreateAttValues() const
RichTrajectoryPointsContainer * fpRichPointsContainer
const G4VProcess * fpCreatorProcess
const G4VProcess * fpEndingProcess
void DrawTrajectory() const
void MergeTrajectory(G4VTrajectory *secondTrajectory)
G4TouchableHandle fpInitialVolume
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
G4int GetPointEntries() const
G4RichTrajectory & operator=(const G4RichTrajectory &)=delete
G4int operator==(const G4RichTrajectory &r) const
Definition: G4Step.hh:62
#define G4TRACKING_DLL
Definition: trkgdefs.hh:45