BOSS 7.0.6
BESIII Offline Software System
Loading...
Searching...
No Matches
BesSensitiveManager.hh
Go to the documentation of this file.
1//---------------------------------------------------------------------------//
2//////// BOOST --- BESIII Object_Oriented Simulation Tool //
3////////---------------------------------------------------------------------------//
4////////Description:
5////////Author : Dengzy
6////
7//// ////Created: Aug, 2004
8//// ////Modified:
9//// ////Comment:
10//// ////---------------------------------------------------------------------------//
11//// //// $Id: BesSensitiveManager.hh
12
13#ifndef BesSensitiveManager_hh
14#define BesSensitiveManager_hh
15
17#include "G4TrackingManager.hh"
18#include "globals.hh"
19#include "HepMC/GenEvent.h"
20
21class G4Svc;
22
23class G4VProcess;
24class G4Track;
25class G4Event;
26
27
28class BesTruthVertex;
29class BesTruthTrack;
30// Since G4 lacks a "permanent stack" feature like GEANT3,
31// we need to keep track of the complete parentage of particles
32
33class BesTStats {
34
36
37public:
39 G4index(-1),
40 p4(0),
41 time(0),
42 trackIndex(-1),
43 originVertex(-1),
44 savedByDefault(false),
45 vertices(0) {;}
46 BesTStats( const G4int index, const CLHEP::HepLorentzVector theP4, G4double theTime )
47 : G4index(index),
48 p4(theP4),
49 time(theTime),
50 trackIndex(-1),
51 originVertex(-1),
52 savedByDefault(false),
53 vertices(0) {;}
55
56 bool operator==( const BesTStats &other ) const { return G4index == other.G4index; }
57
58private:
59 G4int G4index; // Index as in G4Track::fTrackID
60 CLHEP::HepLorentzVector p4; // 4-Momentum at origin
61 G4double time; // Global time at origin
62 G4int trackIndex; // Index into trackList, or -1 if no such
63 G4int originVertex; // Index into vertexList of origin vertex
64 bool savedByDefault; // True if saved to BesTruthTrack by default
65
66 // A list of vertices (index into vertexList)
67 // associated with an interaction
68 // of this track.
69 std::vector<int> vertices;
70
71};
72
74public:
76 {
77 m_hepmcevt = 0;
79 G4Exception("BesSensitiveManager constructed twice.");
80 m_sensitiveManager = this ;
81 }
82
84
86
87 // Add a new client
89 {
90 clients.push_back( detector );
91 }
92
93 G4int GetCurrentTrackIndex() const { return m_trackIndex; }
94 std::vector<BesTruthTrack*>* GetTrackList() {return m_trackList;}
95 std::vector<BesTruthVertex*>* GetVertexList() {return m_vertexList;}
96
97 // Take care of things at the beginning and end
98 // of event, including dumping all required stuff
99 // at the end of the event
100 void BeginOfTruthEvent(const G4Event*);
101 void EndOfTruthEvent(const G4Event*);
102
103 G4int CheckType(const HepMC::GenEvent* hepmcevt);
105 void SetVertex0(const G4Event*);
106 void UpdatePrimaryTrack(const G4Track*);
107 G4bool CheckDecayTrack(const G4Track*);
108 void UpdateVertex(BesTStats, const G4Track*);
109 G4bool MatchDaughterTrack(const G4Track*);
110 void GetDaughterVertexes(BesTruthTrack* pTrack, std::vector<int>* vDau);
111 G4bool MatchVertex(G4int vIndex, std::vector<int>* vDau);
112
113 // Take care of things at the beginning and end
114 // of tracks
115 void BeginOfTrack( const G4Track *track );
116 void EndOfTrack( const G4Track *track , G4TrackingManager* );
117
118 // Find out how many tracks and vertices we have
119 G4int GetNumberTracks() const { return m_trackList ? m_trackList->size() : 0; }
120 G4int GetNumberVertices() const { return m_vertexList ? m_vertexList->size() : 0; }
121
122 void ClearEvent();
123
124 void SetLogLevel(G4int level) {m_logLevel = level;}
125
126protected:
127
129
130 std::vector<BesSensitiveDetector*>::iterator iter;
131 std::vector<BesSensitiveDetector*>::iterator iter_end;
132
133 // Current track indices
136
137 // The following are the BesTruthTrack, BesTruthVertex
138 // lists. Don't forget to delete them at the end of the event.
139 std::vector<BesTruthTrack*> *m_trackList;
140 std::vector<BesTruthVertex*> *m_vertexList;
141
142 // Here is a list of clients. Note that
143 // these objects are decendents of G4VSensitiveDetector
144 // and are thus owned by GEANT4. We must not delete
145 // them ourselves!
146 std::vector<BesSensitiveDetector*> clients;
147
148 std::vector<BesTStats> chain; // Stats on the decay chain
149
150 // Make a new TruthTrack (and TruthVertex if necessary)
151 void MakeNewTrack( BesTStats &stat, const G4Track *track );
152
153 // Update trackList, vertexList, for a new track, and return the new track's stats
154 BesTStats FollowTrack( const G4Track *track );
155
156
157 // this data member is set for recording trackIndex
158 G4int m_count;
159
160 G4ThreeVector m_pos0;
161 G4double m_t0;
162
164 HepMC::GenEvent* m_hepmcevt;
165};
166
167#endif
void EndOfTruthEvent(const G4Event *)
BesTStats FollowTrack(const G4Track *track)
std::vector< BesSensitiveDetector * >::iterator iter_end
void UpdatePrimaryTrack(const G4Track *)
HepMC::GenEvent * m_hepmcevt
std::vector< BesTruthVertex * > * GetVertexList()
void SetVertex0(const G4Event *)
static BesSensitiveManager * m_sensitiveManager
G4bool MatchDaughterTrack(const G4Track *)
std::vector< BesSensitiveDetector * > clients
G4int GetCurrentTrackIndex() const
void SetLogLevel(G4int level)
void MakeNewTrack(BesTStats &stat, const G4Track *track)
void BeginOfTrack(const G4Track *track)
G4int GetNumberVertices() const
void AddSensitiveDetector(BesSensitiveDetector *detector)
void UpdateVertex(BesTStats, const G4Track *)
std::vector< BesTStats > chain
std::vector< BesTruthTrack * > * m_trackList
std::vector< BesTruthTrack * > * GetTrackList()
void EndOfTrack(const G4Track *track, G4TrackingManager *)
void BeginOfTruthEvent(const G4Event *)
G4bool MatchVertex(G4int vIndex, std::vector< int > *vDau)
std::vector< BesTruthVertex * > * m_vertexList
void GetDaughterVertexes(BesTruthTrack *pTrack, std::vector< int > *vDau)
std::vector< BesSensitiveDetector * >::iterator iter
static BesSensitiveManager * GetSensitiveManager()
G4bool CheckDecayTrack(const G4Track *)
G4int CheckType(const HepMC::GenEvent *hepmcevt)
bool operator==(const BesTStats &other) const
BesTStats(const G4int index, const CLHEP::HepLorentzVector theP4, G4double theTime)
Definition: G4Svc.h:32