BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
Simulation/BOOST/TofSim/include/TofSim/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 "G4Allocator.hh"
16#include "G4Material.hh"
17#include "G4THitsCollection.hh"
18#include "G4ThreeVector.hh"
19#include "G4VHit.hh"
20
21class BesTofHit : public G4VHit {
22public:
24 virtual ~BesTofHit();
25
27 const BesTofHit& operator=( const BesTofHit& );
28
29 void AddHit( const BesTofHit& );
30
31 virtual G4int operator==( const BesTofHit& ) const;
32 inline void* operator new( size_t );
33 inline void operator delete( void* );
34
35 virtual void Draw();
36 virtual void Print();
37
38public:
39 void SetEvent( G4double event ) { m_event = event; };
40 G4double GetEvent() { return m_event; }
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 // ETF(MRPC)
58 void SetModule( G4int module ) { m_scinNb = module; }
59 void SetStrip( G4int strip ) { m_strip = strip; }
60 void SetIons( G4int ions ) { m_number_of_ions = ions; }
61 void SetLocPos( G4ThreeVector locPos ) { m_locPos = locPos; }
62 void SetGapNb( G4int gapNb ) { m_gapNb = gapNb; }
63
64 G4int GetTrackIndex() { return m_trackIndex; }
65 G4int GetG4Index() { return m_g4Index; }
66 G4int GetPartId() { return m_partId; }
67 G4int GetScinNb() { return m_scinNb; }
68 G4double GetEdep() { return m_edep; }
69 G4double GetStepL() { return m_stepL; }
70 G4double GetTrackL() { return m_trackL; }
71 G4ThreeVector GetPos() { return m_pos; }
72 G4double GetTime() { return m_time; }
73 G4double GetDeltaT() { return m_deltaT; }
74 G4ThreeVector GetPDirection() { return m_pDirection; }
75 G4ThreeVector GetMomentum() { return m_momentum; }
76 // G4String GetPName() {return m_pName; }
77 // G4Material* GetMaterial() {return m_scinMaterial; }
78 G4double GetCharge() { return m_charge; }
79 G4int GetPDGcode() { return m_pdgcode; }
80 // ETF(MRPC)
81 G4int GetModule() { return m_scinNb; }
82 G4int GetStrip() { return m_strip; }
83 G4int GetIons() { return m_number_of_ions; }
84 G4ThreeVector GetLocPos() { return m_locPos; }
85 G4int GetGapNb() { return m_gapNb; }
86
87 inline void AddEdep( G4double de ) { m_edep += de; }
88
89private:
90 G4double m_event;
91 G4int m_trackIndex;
92 G4int m_g4Index;
93 G4int m_partId;
94 G4int m_scinNb;
95 G4int m_gapNb;
96 G4double m_edep;
97 G4double m_stepL;
98 G4double m_trackL;
99 G4ThreeVector m_pos;
100 G4double m_time;
101 G4double m_deltaT;
102 G4ThreeVector m_pDirection;
103 G4ThreeVector m_momentum;
104 // G4String m_pName;
105 // G4Material* m_scinMaterial;
106 G4double m_charge;
107 G4int m_number_of_ions;
108 G4int m_pdgcode;
109 G4int m_strip;
110 G4ThreeVector m_locPos;
111};
112
113typedef G4THitsCollection<BesTofHit> BesTofHitsCollection;
114
115extern G4Allocator<BesTofHit> BesTofHitAllocator;
116
117inline void* BesTofHit::operator new( size_t ) {
118 void* aHit;
119 aHit = (void*)BesTofHitAllocator.MallocSingle();
120 return aHit;
121}
122
123inline void BesTofHit::operator delete( void* aHit ) {
124 BesTofHitAllocator.FreeSingle( (BesTofHit*)aHit );
125}
126
127#endif
**********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
Double_t time
G4THitsCollection< BesTofHit > BesTofHitsCollection
G4Allocator< BesTofHit > BesTofHitAllocator
Definition BesTofHit.cc:20
virtual ~BesTofHit()
virtual G4int operator==(const BesTofHit &) const
BesTofHit(const BesTofHit &)
virtual void Print()
virtual void Draw()
void SetPDirection(G4ThreeVector pDirection)
const BesTofHit & operator=(const BesTofHit &)
void SetMomentum(G4ThreeVector momentum)
void AddHit(const BesTofHit &)