Loading...
Searching...
No Matches
Par03SensitiveDetector.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
27#include "Par03Hit.hh"
28
29#include "G4HCofThisEvent.hh"
30#include "G4TouchableHistory.hh"
31#include "G4Track.hh"
32#include "G4Step.hh"
33#include "G4SDManager.hh"
34
37{
38 collectionName.insert("hits");
39}
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
41
43 G4int aNumRho, G4int aNumPhi)
45 , fCellNoZ(aNumLayers)
46 , fCellNoRho(aNumRho)
47 , fCellNoPhi(aNumPhi)
48{
49 collectionName.insert("hits");
50}
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57
59{
61 new Par03HitsCollection(SensitiveDetectorName, collectionName[0]);
62 if(fHitCollectionID < 0)
63 {
65 G4SDManager::GetSDMpointer()->GetCollectionID(fHitsCollection);
66 }
67 aHCE->AddHitsCollection(fHitCollectionID, fHitsCollection);
68
69 // fill calorimeter hits with zero energy deposition
70 for(G4int iphi = 0; iphi < fCellNoPhi; iphi++)
71 for(G4int irho = 0; irho < fCellNoRho; irho++)
72 for(G4int iz = 0; iz < fCellNoZ; iz++)
73 {
74 auto hit = new Par03Hit();
75 fHitsCollection->insert(hit);
76 }
77}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80
82{
83 G4double edep = aStep->GetTotalEnergyDeposit();
84 if(edep == 0.)
85 return true;
86
87 auto aTouchable =
88 (G4TouchableHistory*) (aStep->GetPreStepPoint()->GetTouchable());
89
90 auto hit = RetrieveAndSetupHit(aTouchable);
91
92 // Add energy deposit from G4Step
93 hit->AddEdep(edep);
94
95 // Fill time information from G4Step
96 // If it's already filled, choose hit with earliest global time
97 if(hit->GetTime() == -1 ||
98 hit->GetTime() > aStep->GetTrack()->GetGlobalTime())
99 hit->SetTime(aStep->GetTrack()->GetGlobalTime());
100
101 // Set hit type to full simulation (only if hit is not already marked as fast
102 // sim)
103 if(hit->GetType() != 1)
104 hit->SetType(0);
105
106 return true;
107}
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110
112 const G4FastTrack* aTrack,
113 G4TouchableHistory* aTouchable)
114{
115 G4double edep = aHit->GetEnergy();
116 if(edep == 0.)
117 return true;
118
119 auto hit = RetrieveAndSetupHit(aTouchable);
120
121 // Add energy deposit from G4FastHit
122 hit->AddEdep(edep);
123
124 // Fill time information from G4FastTrack
125 // If it's already filled, choose hit with earliest global time
126 if(hit->GetTime() == -1 ||
127 hit->GetTime() > aTrack->GetPrimaryTrack()->GetGlobalTime())
128 {
129 hit->SetTime(aTrack->GetPrimaryTrack()->GetGlobalTime());
130 }
131
132 // Set hit type to fast simulation (even if hit was already marked as full
133 // sim, overwrite it)
134 hit->SetType(1);
135
136 return true;
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140
142 G4TouchableHistory* aTouchable)
143{
144 G4int rhoNo = aTouchable->GetCopyNumber(0); // cell
145 G4int phiNo = aTouchable->GetCopyNumber(1); // segment
146 G4int zNo = aTouchable->GetCopyNumber(2); // layer
147
148 std::size_t hitID = fCellNoRho * fCellNoZ * phiNo + fCellNoZ * rhoNo + zNo;
149
150 if(hitID >= fHitsCollection->entries())
151 {
152 G4Exception(
153 "Par03SensitiveDetector::RetrieveAndSetupHit()", "InvalidSetup",
154 FatalException,
155 "Size of hit collection in Par03SensitiveDetector is smaller than the "
156 "number of cells created in Par03DetectorConstruction!");
157 }
158 Par03Hit* hit = (*fHitsCollection)[hitID];
159
160 if(hit->GetRhoId() < 0)
161 {
162 hit->SetRhoId(rhoNo);
163 hit->SetPhiId(phiNo);
164 hit->SetZid(zNo);
165 hit->SetLogV(aTouchable->GetVolume(0)->GetLogicalVolume());
166 G4AffineTransform transform = aTouchable->GetHistory()->GetTopTransform();
167 hit->SetRot(transform.NetRotation());
168 transform.Invert();
169 hit->SetPos(transform.NetTranslation());
170 }
171 return hit;
172}
G4THitsCollection< Par03Hit > Par03HitsCollection
Definition Par03Hit.hh:132
G4double GetEnergy() const
Hit class to store energy deposited in the sensitive detector.
G4int GetRhoId() const
Get rho id of the cell in the readout segmentation.
Definition Par03Hit.hh:93
void SetPhiId(G4int aPhi)
Set phi id of the cell in the readout segmentation.
Definition Par03Hit.hh:95
void SetZid(G4int aZ)
Set Z id of the cell in the readout segmentation.
Definition Par03Hit.hh:87
void SetRhoId(G4int aRho)
Set Rho id of the cell in the readout segmentation.
Definition Par03Hit.hh:91
void SetRot(G4RotationMatrix aXYZ)
Set rotation.
Definition Par03Hit.hh:77
void SetLogV(G4LogicalVolume *aLogVol)
Definition Par03Hit.hh:107
void SetPos(G4ThreeVector aXYZ)
Set position.
Definition Par03Hit.hh:73
G4int fCellNoPhi
Number of readout cells along azimuthal angle.
G4int fHitCollectionID
ID of collection of hits.
virtual ~Par03SensitiveDetector()
Par03Hit * RetrieveAndSetupHit(G4TouchableHistory *aTouchable)
Process energy deposit - common part for full and fast simulation It is invoked from ProcessHits() me...
G4int fCellNoZ
Number of readout cells along z axis.
virtual void Initialize(G4HCofThisEvent *HCE) final
Create hit collection.
Par03HitsCollection * fHitsCollection
Collection of hits.
G4int fCellNoRho
Number of readout cells along radius of cylinder.
virtual G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *aROhist) final
Process energy deposit from the full simulation.

Applications | User Support | Publications | Collaboration