13#include "ReadBoostRoot.hh"
14#include "G4HCofThisEvent.hh"
16#include "G4ThreeVector.hh"
17#include "G4SDManager.hh"
19#include "G4UnitsTable.hh"
21#include "GaudiKernel/ISvcLocator.h"
22#include "GaudiKernel/Bootstrap.h"
23#include "GaudiKernel/MsgStream.h"
24#include "G4LossTableManager.hh"
25#include "G4ElectronIonPair.hh"
30 collectionName.insert(
"BesTofHitsCollection");
31 collectionName.insert(
"BesTofHitsList");
32 elIonPair =
new G4ElectronIonPair();
50 m_trackIndexes.clear();
57 HLID = G4SDManager::GetSDMpointer()->GetCollectionID( collectionName[1] );
59 G4HCofThisEvent* HCE = evt->GetHCofThisEvent();
60 HCE->AddHitsCollection( HLID, m_besTofList );
67 G4double chg=aStep->GetTrack()->GetDefinition()->GetPDGCharge();
68 G4double edep = aStep->GetTotalEnergyDeposit();
69 G4double stepL=aStep->GetStepLength();
70 G4double deltaT=aStep->GetDeltaTime();
71 G4StepPoint* preStep = aStep->GetPreStepPoint();
72 G4ThreeVector pDirection=preStep->GetMomentumDirection();
73 G4String particleName = aStep->GetTrack()->GetDefinition()->GetParticleName();
74 G4Material* scinMaterial = aStep->GetTrack()->GetMaterial();
75 G4double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
76 G4int pdgcode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
78 if( chg==0 && edep==0 && stepL==0 ) {
return false; }
81 G4int trackId = aStep->GetTrack()->GetTrackID();
84 newHit->
SetG4Index(aStep->GetTrack()->GetTrackID());
88 newHit->
SetTrackL(aStep->GetTrack()->GetTrackLength());
89 G4ThreeVector pos=preStep->GetPosition();
91 G4double globalTime=preStep->GetGlobalTime();
101 G4ThreeVector locPos(0,0,0);
102 G4TouchableHistory* theTouchable = (G4TouchableHistory*)( preStep->GetTouchable() );
113 name = theTouchable->GetVolume(0)->GetName();
115 G4int partId=-1, scinNb=-1, gapNb=-1, number=-1;
117 gapNb = theTouchable->GetReplicaNumber(0);
118 number = theTouchable->GetReplicaNumber(2);
121 if(name.contains(
"physical_gasLayer"))
123 locPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(pos);
124 number = theTouchable->GetReplicaNumber(3);
127 G4String name1 = theTouchable->GetVolume(4)->GetName();
128 if(name1 ==
"physicalEcTofEast") partId=3;
129 else if(name1 ==
"physicalEcTofWest") partId=4;
133 else if(name==
"logical_gasLayer")
135 locPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(pos);
136 number = theTouchable->GetReplicaNumber(3);
139 G4String name1 = theTouchable->GetVolume(4)->GetName();
140 if(name1 ==
"logicalEcTofEast") partId=3;
141 else if(name1 ==
"logicalEcTofWest") partId=4;
147 else if( name==
"physicalScinBr1" ) {
152 else if( name==
"physicalScinBr2" ) {
157 else if( name==
"physicalScinEcWest" ) {
162 else if( name==
"physicalScinEcEast" ) {
170 else if( name==
"logicalScinBr1" || name==
"logicalScinBr2" ) {
172 scinNb = (527-number)/3;
175 else if( name==
"logicalScinEcEast" ) {
177 scinNb = (95-number)/2;
179 else if( name==
"logicalScinEcWest" ) {
181 scinNb = (95-number)/2;
183 else {
return false; }
186 if(name.contains(
"physical_gasLayer") || name.contains(
"logical_gasLayer"))
188 G4double zz = locPos.z()-0.5*mm+(24+3)*mm*6;
193 else if(zz>0 && zz<12*27*mm)
195 for(G4int i=0; i<12; i++)
197 if(zz>i*27*mm && zz<=(i+1)*27*mm)
210 if(strip>11) strip=11;
224 G4int trackIndex, g4TrackId;
229 m_besTofCollection->insert( newHit );
234 G4int trackIndex, g4TrackId;
237 if( m_trackIndex != trackIndex ) {
238 m_trackIndex = trackIndex;
249 G4int pdg =
abs(aStep->GetTrack()->GetDefinition()->GetPDGEncoding());
250 if( pdg==12 || pdg==14 || pdg==16 ) {
flag=0; }
251 if(
flag && aStep->GetTrack()->GetTrackID()==g4TrackId ) {
252 m_trackIndexes.push_back(trackIndex);
255 m_besTofList->insert(truHit);
259 if( edep<=0 ) {
delete newHit; }
268 static G4int HCID=-1;
270 HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
272 HCE->AddHitsCollection(HCID,m_besTofCollection);
278 G4double FanoFactor = 0.2;
279 G4double meanion = eipair->MeanNumberOfIonsAlongStep(step);
280 G4double sig = std::sqrt(FanoFactor*meanion);
281 G4int nion = G4int(G4RandGauss::shoot(meanion,sig) + 0.5);
G4THitsCollection< BesTofHit > BesTofHitsCollection
void GetCurrentTrackIndex(G4int &trackIndex, G4int &g4TrackId) const
void SetModule(G4int module)
void SetPos(G4ThreeVector pos)
void SetDeltaT(G4double deltaT)
void SetCharge(G4double charge)
void SetPDGcode(G4int pdgcode)
void SetTrackIndex(G4int trackIndex)
void SetPDirection(G4ThreeVector pDirection)
void SetPartId(G4int partId)
void SetLocPos(G4ThreeVector locPos)
void SetScinNb(G4int scinNb)
void SetGapNb(G4int gapNb)
void SetStrip(G4int strip)
void SetStepL(G4double stepL)
void SetTrackL(G4double length)
void SetTime(G4double time)
void SetEdep(G4double edep)
void SetMomentum(G4ThreeVector momentum)
void SetG4Index(G4int index)
void BeginOfTruthEvent(const G4Event *)
virtual void Initialize(G4HCofThisEvent *HCE)
virtual void EndOfEvent(G4HCofThisEvent *HCE)
void EndOfTruthEvent(const G4Event *)
G4int SampleNumberOfIonsAlongStep(const G4Step *, G4ElectronIonPair *)
virtual G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *)