CGEM BOSS 6.6.5.g
BESIII Offline Software System
Loading...
Searching...
No Matches
BesTofHit.hh
Go to the documentation of this file.
1//---------------------------------------------------------------------------//
2// BOOST --- BESIII Object_Oriented Simulation Tool //
3//---------------------------------------------------------------------------//
4//Description:
5//Author: Dengzy
6//Created: Mar, 2004
7//Modified:
8//Comment:
9//---------------------------------------------------------------------------//
10// $Id: BesTofHit.hh
11
12#ifndef BesTofHit_h
13#define BesTofHit_h 1
14
15#include "G4VHit.hh"
16#include "G4THitsCollection.hh"
17#include "G4Allocator.hh"
18#include "G4ThreeVector.hh"
19#include "G4Material.hh"
20
21
22class BesTofHit : public G4VHit
23{
24 public:
25 BesTofHit();
26 virtual ~BesTofHit();
27
28 BesTofHit(const BesTofHit&);
29 const BesTofHit& operator=(const BesTofHit&);
30
31 void AddHit(const BesTofHit&);
32
33 virtual G4int operator==(const BesTofHit&) const;
34 inline void* operator new(size_t);
35 inline void operator delete(void*);
36
37 virtual void Draw();
38 virtual void Print();
39
40 public:
41 void SetTrackIndex(G4int trackIndex) { m_trackIndex = trackIndex; };
42 void SetG4Index(G4int index) {m_g4Index = index;}
43 void SetPartId(G4int partId) {m_partId = partId;}
44 void SetScinNb(G4int scinNb) { m_scinNb = scinNb; }
45 void SetEdep(G4double edep) { m_edep = edep; }
46 void SetStepL(G4double stepL) { m_stepL = stepL;}
47 void SetTrackL(G4double length) { m_trackL = length; }
48 void SetPos(G4ThreeVector pos) { m_pos = pos; }
49 void SetTime(G4double time) { m_time=time;}
50 void SetDeltaT(G4double deltaT) { m_deltaT = deltaT;}
51 void SetPDirection(G4ThreeVector pDirection) { m_pDirection = pDirection; }
52 void SetMomentum(G4ThreeVector momentum) { m_momentum = momentum; }
53 //void SetPName(G4String pName) { m_pName = pName; }
54 //void SetMaterial(G4Material* myMaterial){m_scinMaterial = myMaterial;}
55 void SetCharge(G4double charge) {m_charge = charge;}
56 void SetPDGcode(G4int pdgcode){m_pdgcode=pdgcode;}
57
58 G4int GetTrackIndex() { return m_trackIndex; }
59 G4int GetG4Index() {return m_g4Index;}
60 G4int GetPartId() { return m_partId; }
61 G4int GetScinNb() { return m_scinNb; }
62 G4double GetEdep() { return m_edep; }
63 G4double GetStepL() { return m_stepL; }
64 G4double GetTrackL() {return m_trackL; }
65 G4ThreeVector GetPos() { return m_pos; }
66 G4double GetTime() { return m_time;}
67 G4double GetDeltaT() {return m_deltaT;}
68 G4ThreeVector GetPDirection() {return m_pDirection;}
69 G4ThreeVector GetMomentum() {return m_momentum;}
70 //G4String GetPName() {return m_pName; }
71 //G4Material* GetMaterial() {return m_scinMaterial; }
72 G4double GetCharge() {return m_charge; }
73 G4int GetPDGcode(){return m_pdgcode;}
74
75
76
77 inline void AddEdep(G4double de) { m_edep += de; }
78
79 void SetIons(G4int ions){m_number_of_ions=ions;}
80 G4int GetIons(){return m_number_of_ions;}
81
82 //#Matthias
83 G4int GetModule_mrpc(){return m_scinNb;}
84 void SetModule_mrpc(G4int module_mrpc){m_scinNb=module_mrpc;}
85
86 private:
87 G4int m_trackIndex;
88 G4int m_g4Index;
89 G4int m_partId;
90 G4int m_scinNb;
91 G4double m_edep;
92 G4double m_stepL;
93 G4double m_trackL;
94 G4ThreeVector m_pos;
95 G4double m_time;
96 G4double m_deltaT;
97 G4ThreeVector m_pDirection;
98 G4ThreeVector m_momentum;
99 //G4String m_pName;
100 //G4Material* m_scinMaterial;
101 G4double m_charge;
102 G4int m_number_of_ions;
103 G4int m_pdgcode;
104
105};
106
107
108typedef G4THitsCollection<BesTofHit> BesTofHitsCollection; //BesTofHitsCollection ist nur neuer Name für die Klasse BesTofHit!
109
110extern G4Allocator<BesTofHit> BesTofHitAllocator;
111
112
113inline void* BesTofHit::operator new(size_t)
114{
115 void *aHit;
116 aHit = (void *) BesTofHitAllocator.MallocSingle();
117 return aHit;
118}
119
120
121inline void BesTofHit::operator delete(void *aHit)
122{
123 BesTofHitAllocator.FreeSingle((BesTofHit*) aHit);
124}
125
126
127#endif
128
129
**********INTEGER nmxhep !maximum number of particles DOUBLE PRECISION vhep INTEGER jdahep COMMON hepevt $ !serial number $ !number of particles $ !status code $ !particle ident KF $ !parent particles $ !childreen particles $ !four momentum
G4THitsCollection< BesTofHit > BesTofHitsCollection
Definition: BesTofHit.hh:108
G4Allocator< BesTofHit > BesTofHitAllocator
Definition: BesTofHit.cc:20
double length
Double_t time
G4double GetDeltaT()
Definition: BesTofHit.hh:67
void AddEdep(G4double de)
Definition: BesTofHit.hh:77
void SetPos(G4ThreeVector pos)
Definition: BesTofHit.hh:48
G4ThreeVector GetPDirection()
Definition: BesTofHit.hh:68
void SetDeltaT(G4double deltaT)
Definition: BesTofHit.hh:50
void SetCharge(G4double charge)
Definition: BesTofHit.hh:55
void SetPDGcode(G4int pdgcode)
Definition: BesTofHit.hh:56
G4int GetIons()
Definition: BesTofHit.hh:80
virtual G4int operator==(const BesTofHit &) const
Definition: BesTofHit.cc:70
void SetTrackIndex(G4int trackIndex)
Definition: BesTofHit.hh:41
G4int GetPDGcode()
Definition: BesTofHit.hh:73
void SetPDirection(G4ThreeVector pDirection)
Definition: BesTofHit.hh:51
G4double GetStepL()
Definition: BesTofHit.hh:63
G4double GetCharge()
Definition: BesTofHit.hh:72
virtual void Print()
Definition: BesTofHit.cc:101
G4double GetTime()
Definition: BesTofHit.hh:66
G4int GetModule_mrpc()
Definition: BesTofHit.hh:83
void SetIons(G4int ions)
Definition: BesTofHit.hh:79
G4double GetEdep()
Definition: BesTofHit.hh:62
void SetPartId(G4int partId)
Definition: BesTofHit.hh:43
void SetModule_mrpc(G4int module_mrpc)
Definition: BesTofHit.hh:84
virtual ~BesTofHit()
Definition: BesTofHit.cc:24
G4ThreeVector GetPos()
Definition: BesTofHit.hh:65
const BesTofHit & operator=(const BesTofHit &)
Definition: BesTofHit.cc:47
void SetScinNb(G4int scinNb)
Definition: BesTofHit.hh:44
G4int GetScinNb()
Definition: BesTofHit.hh:61
G4int GetPartId()
Definition: BesTofHit.hh:60
void SetStepL(G4double stepL)
Definition: BesTofHit.hh:46
G4double GetTrackL()
Definition: BesTofHit.hh:64
void SetTrackL(G4double length)
Definition: BesTofHit.hh:47
G4int GetG4Index()
Definition: BesTofHit.hh:59
virtual void Draw()
Definition: BesTofHit.cc:85
void SetTime(G4double time)
Definition: BesTofHit.hh:49
void SetEdep(G4double edep)
Definition: BesTofHit.hh:45
G4ThreeVector GetMomentum()
Definition: BesTofHit.hh:69
void SetMomentum(G4ThreeVector momentum)
Definition: BesTofHit.hh:52
void SetG4Index(G4int index)
Definition: BesTofHit.hh:42
G4int GetTrackIndex()
Definition: BesTofHit.hh:58
void AddHit(const BesTofHit &)
Definition: BesTofHit.cc:75