BOSS 7.0.7
BESIII Offline Software System
Loading...
Searching...
No Matches
BesTofSD.cc
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: BesTofSD.cc
11
12#include "GaudiKernel/MsgStream.h"
13#include "GaudiKernel/Bootstrap.h"
14#include "GaudiKernel/PropertyMgr.h"
15#include "GaudiKernel/IJobOptionsSvc.h"
16#include "GaudiKernel/ISvcLocator.h"
17#include "GaudiKernel/IDataProviderSvc.h"
18
19#include "BesTofSD.hh"
20#include "ReadBoostRoot.hh"
21#include "G4HCofThisEvent.hh"
22#include "G4Step.hh"
23#include "G4ThreeVector.hh"
24#include "G4SDManager.hh"
25#include "G4ios.hh"
26#include "G4UnitsTable.hh"
27#include "G4Event.hh"
28#include "GaudiKernel/ISvcLocator.h"
29#include "GaudiKernel/Bootstrap.h"
30#include "GaudiKernel/MsgStream.h"
31#include "G4LossTableManager.hh"
32#include "G4ElectronIonPair.hh"
33
34//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
35
36BesTofSD::BesTofSD( G4String name ) : BesSensitiveDetector(name), m_besTofList(0) {
37 collectionName.insert("BesTofHitsCollection");
38 collectionName.insert("BesTofHitsList");
39 elIonPair = new G4ElectronIonPair();
40
41 PropertyMgr m_propMgr;
42 m_propMgr.declareProperty("FanoFactor", m_fanoFactor = 0.2);
43 m_propMgr.declareProperty("nionOff", m_nionOff = 0.5);
44
45 IJobOptionsSvc* jobSvc;
46 Gaudi::svcLocator()->service("JobOptionsSvc", jobSvc);
47 jobSvc->setMyProperties("BesTofSD", &m_propMgr);
48
49 cout<<"BesTofSD Property:"<<endl
50 <<" FanoFactor= "<<m_fanoFactor
51 <<" nionOff= "<<m_nionOff
52 <<endl;
53}
54
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56
58}
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61
62void BesTofSD::Initialize( G4HCofThisEvent* ) {
63 m_besTofCollection = new BesTofHitsCollection( SensitiveDetectorName, collectionName[0] );
64 //elIonPair = G4LossTableManager::Instance()->ElectronIonPair();
65}
66
67//for MC Truth
68void BesTofSD::BeginOfTruthEvent( const G4Event* evt) {
69 m_event = evt->GetEventID();
70 m_besTofList = new BesTofHitsCollection( SensitiveDetectorName, collectionName[1] );
71 m_trackIndexes.clear();
72 m_trackIndex = -99;
73}
74
75void BesTofSD::EndOfTruthEvent( const G4Event* evt ) {
76 static G4int HLID=-1;
77 if( HLID<0 ) {
78 HLID = G4SDManager::GetSDMpointer()->GetCollectionID( collectionName[1] );
79 }
80 G4HCofThisEvent* HCE = evt->GetHCofThisEvent();
81 HCE->AddHitsCollection( HLID, m_besTofList );
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85
86G4bool BesTofSD::ProcessHits( G4Step* aStep, G4TouchableHistory* ) {
87
88 G4double chg=aStep->GetTrack()->GetDefinition()->GetPDGCharge();
89 G4double edep = aStep->GetTotalEnergyDeposit();
90 G4double stepL=aStep->GetStepLength();
91 G4double deltaT=aStep->GetDeltaTime();
92 G4StepPoint* preStep = aStep->GetPreStepPoint();
93 G4ThreeVector pDirection=preStep->GetMomentumDirection();
94 G4String particleName = aStep->GetTrack()->GetDefinition()->GetParticleName();
95 G4Material* scinMaterial = aStep->GetTrack()->GetMaterial();
96 G4double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
97 G4int pdgcode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
98
99 if( chg==0 && edep==0 && stepL==0 ) { return false; }
100
101 BesTofHit *newHit = new BesTofHit();
102 G4int trackId = aStep->GetTrack()->GetTrackID();
103
104 newHit->SetEvent(m_event);
105 newHit->SetTrackIndex(-99);
106 newHit->SetG4Index(aStep->GetTrack()->GetTrackID());
107 newHit->SetEdep(edep);
108 newHit->SetStepL(stepL);
109 //newHit->SetPName(particleName);
110 newHit->SetTrackL(aStep->GetTrack()->GetTrackLength());
111 G4ThreeVector pos=preStep->GetPosition();
112 newHit->SetPos(pos);
113 G4double globalTime=preStep->GetGlobalTime();
114 newHit->SetTime(globalTime);
115 newHit->SetDeltaT(deltaT);
116 newHit->SetPDirection(pDirection);
117 newHit->SetMomentum(preStep->GetMomentum());
118 //newHit->SetMaterial(scinMaterial);
119 newHit->SetCharge(charge);
120 newHit->SetPDGcode(pdgcode);
121
122 //Get local position in sensitive detector, add by anff
123 G4ThreeVector locPos(0,0,0);
124 G4TouchableHistory* theTouchable = (G4TouchableHistory*)( preStep->GetTouchable() );
125 //partId: barrel(1), east endcap(0), west endcap(2);
126 //scinNb: barrel(0-175),east endcap(0-47), west endcap(0-47)
127 //partId: mrpc: east(3), west(4),
128
129 G4String name;
130 //for(G4int i=0; i<theTouchable->GetHistoryDepth(); i++)
131 //{
132 // name = theTouchable->GetVolume(i)->GetName();
133 // G4cout<<"********* name of the "<<i<<" volume: "<<name<<" ***************"<<G4endl;
134 //}
135 name = theTouchable->GetVolume(0)->GetName();
136
137 G4int partId=-1, scinNb=-1, gapNb=-1, number=-1;
138 G4int strip = -1;
139 gapNb = theTouchable->GetReplicaNumber(0);
140 number = theTouchable->GetReplicaNumber(2); //This is valid only for scintillator
141
142 //MRPC ENDCAPS construct by geant4 code
143 if(name.contains("physical_gasLayer"))
144 {
145 locPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(pos);
146 number = theTouchable->GetReplicaNumber(3);
147 scinNb = number;
148
149 G4String name1 = theTouchable->GetVolume(4)->GetName();
150 if(name1 == "physicalEcTofEast") partId=3;
151 else if(name1 == "physicalEcTofWest") partId=4;
152 //G4cout<<"scinNb= "<<scinNb<<" gapNb= "<<gapNb<<" partId="<<partId<<G4endl;
153 }
154 //MRPC ENDCAPS construct by GDML
155 else if(name=="logical_gasLayer")
156 {
157 locPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(pos);
158 number = theTouchable->GetReplicaNumber(3);
159 scinNb = 35-number;
160
161 G4String name1 = theTouchable->GetVolume(4)->GetName();
162 if(name1 == "logicalEcTofEast") partId=3;
163 else if(name1 == "logicalEcTofWest") partId=4;
164 //G4cout<<"scinNb= "<<scinNb<<" gapNb= "<<gapNb<<" partId="<<partId<<G4endl;
165 }
166
167
168 // Scintillator
169 else if( name=="physicalScinBr1" ) {
170 partId = 1;
171 scinNb = number;
172 //G4cout<<"scinNb= "<<scinNb<<" partId="<<partId<<G4endl;
173 }
174 else if( name=="physicalScinBr2" ) {
175 partId = 1;
176 scinNb = number+88;
177 //G4cout<<"scinNb= "<<scinNb<<" partId="<<partId<<G4endl;
178 }
179 else if( name=="physicalScinEcWest" ) {
180 partId = 2;
181 scinNb = number;
182 //G4cout<<"scinNb= "<<scinNb<<" partId="<<partId<<G4endl;
183 }
184 else if( name=="physicalScinEcEast" ) {
185 partId = 0;
186 scinNb = number;
187 //G4cout<<"scinNb= "<<scinNb<<" partId="<<partId<<G4endl;
188 }
189
190
191 //if construct by gdml files
192 else if( name=="logicalScinBr1" || name=="logicalScinBr2" ) {
193 partId = 1;
194 scinNb = (527-number)/3;
195 }
196 //copy number: east:1-95(scinNb:47-0), west:1-95(scinNb:47-0)
197 else if( name=="logicalScinEcEast" ) {
198 partId = 0;
199 scinNb = (95-number)/2;
200 }
201 else if( name=="logicalScinEcWest" ) {
202 partId = 2;
203 scinNb = (95-number)/2;
204 }
205 else { return false; }
206
207
208 if(name.contains("physical_gasLayer") || name.contains("logical_gasLayer"))
209 {
210 G4double zz = locPos.z()-0.5*mm+(24+3)*mm*6; //0.5mm is the offset between strip and gasLayer
211 if(zz<=0)
212 {
213 strip=0;
214 }
215 else if(zz>0 && zz<12*27*mm)
216 {
217 for(G4int i=0; i<12; i++)
218 {
219 if(zz>i*27*mm && zz<=(i+1)*27*mm)
220 {
221 strip = i;
222 //G4cout<<"Local z is "<<locPos.z()<<" strip number is: "<<strip<<" !!!"<<G4endl;
223 break;
224 }
225 }
226 }
227 else
228 {
229 strip=11;
230 }
231 if(strip<0) strip=0;
232 if(strip>11) strip=11;
233 }
234
235
236 newHit->SetPartId(partId);
237 newHit->SetScinNb(scinNb);
238 newHit->SetGapNb(gapNb);
239 newHit->SetLocPos(locPos);
240 newHit->SetModule(scinNb);
241 newHit->SetStrip(strip);
242 //newHit->SetIons(elIonPair->SampleNumberOfIonsAlongStep(aStep));
243 newHit->SetIons(SampleNumberOfIonsAlongStep(aStep,elIonPair));
244 //newHit->Draw();
245
246 G4int trackIndex, g4TrackId;
247 GetCurrentTrackIndex( trackIndex, g4TrackId );
248 newHit->SetTrackIndex( trackIndex );
249 //newHit->Print();
250 if( edep>0 ) {
251 m_besTofCollection->insert( newHit );
252 }
253
254 //for mc truth
255 if( m_besTofList ) {
256 G4int trackIndex, g4TrackId;
257 GetCurrentTrackIndex(trackIndex, g4TrackId);
258 newHit->SetTrackIndex(trackIndex);
259 if( m_trackIndex != trackIndex ) {
260 m_trackIndex = trackIndex;
261 //G4int size = m_trackIndexes.size();
262 //G4int flag=1;
263 //if(size>0)
264 //{
265 // for(G4int i=0;i<size;i++)
266 // if(m_trackIndexes[i] == trackIndex )
267 // {flag =0; break;}
268 //}
269 //if(flag && aStep->GetTrack()->GetTrackID()==g4TrackId )
270 G4int flag = 1;
271 G4int pdg = abs(aStep->GetTrack()->GetDefinition()->GetPDGEncoding());
272 if( pdg==12 || pdg==14 || pdg==16 ) { flag=0; }
273 if( flag && aStep->GetTrack()->GetTrackID()==g4TrackId ) {
274 m_trackIndexes.push_back(trackIndex);
275 BesTofHit* truHit = new BesTofHit();
276 *truHit = *newHit;
277 m_besTofList->insert(truHit);
278 }
279 }
280 }
281 if( edep<=0 ) { delete newHit; }
282
283 return true;
284
285}//Close Process Hits
286
287//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
288
289void BesTofSD::EndOfEvent( G4HCofThisEvent* HCE ) {
290 static G4int HCID=-1;
291 if( HCID<0 ) {
292 HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
293 }
294 HCE->AddHitsCollection(HCID,m_besTofCollection);
295}
296
297//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
298
299G4int BesTofSD::SampleNumberOfIonsAlongStep( const G4Step* step, G4ElectronIonPair* eipair ) {
300 G4double FanoFactor = m_fanoFactor; //0.2; //Is this value appropiate?
301 G4double meanion = eipair->MeanNumberOfIonsAlongStep(step);
302 G4double sig = std::sqrt(FanoFactor*meanion); //Be carefull : In originally Geant4 code this calculatin is wrong!
303 G4int nion = G4int(G4RandGauss::shoot(meanion,sig) + m_nionOff); //0.5);
304 return nion;
305}
306
307
G4THitsCollection< BesTofHit > BesTofHitsCollection
Definition: BesTofHit.hh:116
void GetCurrentTrackIndex(G4int &trackIndex, G4int &g4TrackId) const
void SetEvent(G4double event)
Definition: BesTofHit.hh:41
void SetModule(G4int module)
Definition: BesTofHit.hh:60
void SetPos(G4ThreeVector pos)
Definition: BesTofHit.hh:50
void SetDeltaT(G4double deltaT)
Definition: BesTofHit.hh:52
void SetCharge(G4double charge)
Definition: BesTofHit.hh:57
void SetPDGcode(G4int pdgcode)
Definition: BesTofHit.hh:58
void SetTrackIndex(G4int trackIndex)
Definition: BesTofHit.hh:43
void SetPDirection(G4ThreeVector pDirection)
Definition: BesTofHit.hh:53
void SetIons(G4int ions)
Definition: BesTofHit.hh:62
void SetPartId(G4int partId)
Definition: BesTofHit.hh:45
void SetLocPos(G4ThreeVector locPos)
Definition: BesTofHit.hh:63
void SetScinNb(G4int scinNb)
Definition: BesTofHit.hh:46
void SetGapNb(G4int gapNb)
Definition: BesTofHit.hh:64
void SetStrip(G4int strip)
Definition: BesTofHit.hh:61
void SetStepL(G4double stepL)
Definition: BesTofHit.hh:48
void SetTrackL(G4double length)
Definition: BesTofHit.hh:49
void SetTime(G4double time)
Definition: BesTofHit.hh:51
void SetEdep(G4double edep)
Definition: BesTofHit.hh:47
void SetMomentum(G4ThreeVector momentum)
Definition: BesTofHit.hh:54
void SetG4Index(G4int index)
Definition: BesTofHit.hh:44
void BeginOfTruthEvent(const G4Event *)
Definition: BesTofSD.cc:68
virtual void Initialize(G4HCofThisEvent *HCE)
Definition: BesTofSD.cc:62
virtual void EndOfEvent(G4HCofThisEvent *HCE)
Definition: BesTofSD.cc:289
BesTofSD(G4String name)
Definition: BesTofSD.cc:36
virtual ~BesTofSD()
Definition: BesTofSD.cc:57
void EndOfTruthEvent(const G4Event *)
Definition: BesTofSD.cc:75
G4int SampleNumberOfIonsAlongStep(const G4Step *, G4ElectronIonPair *)
Definition: BesTofSD.cc:299
virtual G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *)
Definition: BesTofSD.cc:86
float charge