BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
BesEmcHit.cc
Go to the documentation of this file.
1//---------------------------------------------------------------------------//
2// BOOST --- BESIII Object_Oreiented Simulation Tool //
3//---------------------------------------------------------------------------//
4// Descpirtion: EMC detector
5// Author: Fu Chengdong
6// Created: Sep 4, 2003
7// Comment:
8//---------------------------------------------------------------------------//
9//
10#include "EmcSim/BesEmcHit.hh"
11
12#include "G4UnitsTable.hh"
13#include "G4ios.hh"
14#include <iomanip>
15
16G4Allocator<BesEmcHit> BesEmcHitAllocator;
17G4Allocator<BesEmcTruthHit> BesEmcTruthHitAllocator;
18
19using namespace std;
20
21//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
22
24 TotalEdepCrystal = 0.;
25 TotalTrackLengthCrystal = 0.;
26 EdepCrystal = 0.;
27 TrackLengthCrystal = 0.;
28 EdepCasing = 0.;
29 PositionCrystal = G4ThreeVector( 0, 0, 0 );
30 TimeCrystal = 0.;
31 PartId = 0;
32 NumTheta = 0;
33 NumPhi = 0;
34 trackIndex = 0;
35 g4Index = 0;
36 momentum = G4ThreeVector( 0, 0, 0 );
37}
38
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
40
42
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
44
45BesEmcHit::BesEmcHit( const BesEmcHit& right ) : G4VHit() {
46 EdepCrystal = right.EdepCrystal;
47 TrackLengthCrystal = right.TrackLengthCrystal;
48 EdepCasing = right.EdepCasing;
49 PositionCrystal = right.PositionCrystal;
50 TimeCrystal = right.TimeCrystal;
51 PartId = right.PartId;
52 NumTheta = right.NumTheta;
53 NumPhi = right.NumPhi;
54 trackIndex = right.trackIndex;
55 g4Index = right.g4Index;
56 momentum = right.momentum;
57}
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
62 EdepCrystal = right.EdepCrystal;
63 TrackLengthCrystal = right.TrackLengthCrystal;
64 EdepCasing = right.EdepCasing;
65 PositionCrystal = right.PositionCrystal;
66 TimeCrystal = right.TimeCrystal;
67 PartId = right.PartId;
68 NumTheta = right.NumTheta;
69 NumPhi = right.NumPhi;
70 trackIndex = right.trackIndex;
71 g4Index = right.g4Index;
72 momentum = right.momentum;
73 return *this;
74}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
77
78int BesEmcHit::operator==( const BesEmcHit& right ) const {
79 return ( this == &right ) ? 1 : 0;
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87
88void BesEmcHit::Print( G4int verboseLevel ) {
89 if ( verboseLevel > 0 )
90 G4cout << "Hit in crystal:" << NumTheta << "," << NumPhi << G4endl
91 << "Energy deposited:" << G4BestUnit( EdepCrystal, "Energy" ) << G4endl;
92 if ( verboseLevel > 1 )
93 G4cout << "Hit time :" << G4BestUnit( TimeCrystal, "Time" ) << G4endl
94 << " position :" << G4BestUnit( PositionCrystal, "Length" ) << G4endl;
95 if ( verboseLevel > 2 )
96 G4cout << "Track length :" << G4BestUnit( TrackLengthCrystal, "Length" ) << G4endl;
97}
98
100 G4cout << "time: " << TimeCrystal << " edep: " << EdepCrystal << G4endl;
101}
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103
105 m_trackIndex = -1;
106 m_g4TrackId = -1;
107 m_hitEmc = -1;
108 m_PDGCode = 0;
109 m_PDGCharge = 0;
110 m_particleName = G4String();
111 m_edep = 0;
112 m_time = 0;
113 m_momentum = G4ThreeVector( 0, 0, 0 );
114 m_position = G4ThreeVector( 0, 0, 0 );
115 m_hitMap.clear();
116}
117
119
120BesEmcTruthHit::BesEmcTruthHit( const BesEmcTruthHit& right ) : G4VHit() { *this = right; }
121
123 m_hitMap.clear();
124 std::map<Identifier, G4double>::const_iterator iHitMap;
125 if ( this != &right )
126 {
127 for ( iHitMap = right.Begin(); iHitMap != right.End(); iHitMap++ )
128 { Insert( iHitMap->first, iHitMap->second ); }
129
130 m_identify = right.m_identify;
131 m_trackIndex = right.m_trackIndex;
132 m_g4TrackId = right.m_g4TrackId;
133 m_hitEmc = right.m_hitEmc;
134 m_PDGCode = right.m_PDGCode;
135 m_PDGCharge = right.m_PDGCharge;
136 m_particleName = right.m_particleName;
137 m_edep = right.m_edep;
138 m_time = right.m_time;
139 m_momentum = right.m_momentum;
140 m_position = right.m_position;
141 }
142
143 return *this;
144}
145
147 G4cout << "Id: " << m_identify << "\tTrack Index: " << m_trackIndex
148 << "\tG4 Track Id: " << m_g4TrackId << "\tHit Emc: " << m_hitEmc
149 << "\tTotal Energy: " << m_edep << "\nPDGCode: " << m_PDGCode
150 << "\tCharge: " << m_PDGCharge << "\tParticle Name: " << m_particleName
151 << "\nGloble Time: " << m_time << "\tMomentum: " << m_momentum.mag()
152 << "\tPosition: " << m_position << G4endl;
153
154 std::map<Identifier, G4double>::iterator iHitMap;
155 for ( iHitMap = m_hitMap.begin(); iHitMap != m_hitMap.end(); iHitMap++ )
156 { G4cout << iHitMap->first << "\t" << iHitMap->second << G4endl; }
157}
158
159std::map<Identifier, G4double>::const_iterator BesEmcTruthHit::Begin() const {
160 return m_hitMap.begin();
161}
162
163std::map<Identifier, G4double>::const_iterator BesEmcTruthHit::End() const {
164 return m_hitMap.end();
165}
166
167std::map<Identifier, G4double>::const_iterator BesEmcTruthHit::Find( Identifier id ) const {
168 return m_hitMap.find( id );
169}
170
171G4double BesEmcTruthHit::GetEHit( Identifier id ) { return m_hitMap[id]; }
172
174 if ( energy > 0 ) m_hitMap[id] += energy;
175}
176
178 if ( energy > 0 ) m_hitMap[id] = energy;
179}
180
181G4int BesEmcTruthHit::Size() const { return m_hitMap.size(); }
G4Allocator< BesEmcTruthHit > BesEmcTruthHitAllocator
Definition BesEmcHit.cc:17
G4Allocator< BesEmcHit > BesEmcHitAllocator
Definition BesEmcHit.cc:16
************Class m_ypar INTEGER m_KeyWgt INTEGER m_KeyIHVP INTEGER m_KeyGPS INTEGER m_IsBeamPolarized INTEGER m_EvtGenInterface DOUBLE PRECISION m_Emin DOUBLE PRECISION m_sphot DOUBLE PRECISION m_Xenph DOUBLE PRECISION m_q2 DOUBLE PRECISION m_PolBeam2 DOUBLE PRECISION m_xErrPb *COMMON c_KK2f $ !CMS energy average $ !Spin Polarization vector first beam $ !Spin Polarization vector second beam $ !Beam energy spread[GeV] $ !minimum hadronization energy[GeV] $ !input READ never touch them !$ !debug facility $ !maximum weight $ !inverse alfaQED $ !minimum real photon energy
Definition KK2f.h:50
void Draw()
Definition BesEmcHit.cc:84
int operator==(const BesEmcHit &) const
Definition BesEmcHit.cc:78
void Print()
Definition BesEmcHit.cc:99
const BesEmcHit & operator=(const BesEmcHit &)
Definition BesEmcHit.cc:61
void AddEHit(Identifier, G4double)
Definition BesEmcHit.cc:173
std::map< Identifier, G4double >::const_iterator End() const
Definition BesEmcHit.cc:163
std::map< Identifier, G4double >::const_iterator Find(Identifier) const
Definition BesEmcHit.cc:167
const BesEmcTruthHit & operator=(const BesEmcTruthHit &)
Definition BesEmcHit.cc:122
void Insert(Identifier, G4double)
Definition BesEmcHit.cc:177
G4double GetEHit(Identifier)
Definition BesEmcHit.cc:171
std::map< Identifier, G4double >::const_iterator Begin() const
Definition BesEmcHit.cc:159
virtual ~BesEmcTruthHit()
Definition BesEmcHit.cc:118
G4int Size() const
Definition BesEmcHit.cc:181