Loading...
Searching...
No Matches
CalorimeterSD.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//
26//
27/// \file B4/B4c/src/CalorimeterSD.cc
28/// \brief Implementation of the B4c::CalorimeterSD class
29
30#include "CalorimeterSD.hh"
31#include "G4HCofThisEvent.hh"
32#include "G4Step.hh"
33#include "G4ThreeVector.hh"
34#include "G4SDManager.hh"
35#include "G4ios.hh"
36
37namespace B4c
38{
39
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
41
43 const G4String& hitsCollectionName,
44 G4int nofCells)
46 fNofCells(nofCells)
47{
48 collectionName.insert(hitsCollectionName);
49}
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
54{
55 // Create hits collection
57 = new CalorHitsCollection(SensitiveDetectorName, collectionName[0]);
58
59 // Add this collection in hce
60 auto hcID
61 = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
62 hce->AddHitsCollection( hcID, fHitsCollection );
63
64 // Create hits
65 // fNofCells for cells + one more for total sums
66 for (G4int i=0; i<fNofCells+1; i++ ) {
67 fHitsCollection->insert(new CalorHit());
68 }
69}
70
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
72
75{
76 // energy deposit
77 auto edep = step->GetTotalEnergyDeposit();
78
79 // step length
80 G4double stepLength = 0.;
81 if ( step->GetTrack()->GetDefinition()->GetPDGCharge() != 0. ) {
82 stepLength = step->GetStepLength();
83 }
84
85 if ( edep==0. && stepLength == 0. ) return false;
86
87 auto touchable = (step->GetPreStepPoint()->GetTouchable());
88
89 // Get calorimeter cell id
90 auto layerNumber = touchable->GetReplicaNumber(1);
91
92 // Get hit accounting data for this cell
93 auto hit = (*fHitsCollection)[layerNumber];
94 if ( ! hit ) {
95 G4ExceptionDescription msg;
96 msg << "Cannot access hit " << layerNumber;
97 G4Exception("CalorimeterSD::ProcessHits()",
98 "MyCode0004", FatalException, msg);
99 }
100
101 // Get hit for total accounting
102 auto hitTotal
103 = (*fHitsCollection)[fHitsCollection->entries()-1];
104
105 // Add values
106 hit->Add(edep, stepLength);
107 hitTotal->Add(edep, stepLength);
108
109 return true;
110}
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113
115{
116 if ( verboseLevel>1 ) {
117 auto nofHits = fHitsCollection->entries();
118 G4cout
119 << G4endl
120 << "-------->Hits Collection: in this event they are " << nofHits
121 << " hits in the tracker chambers: " << G4endl;
122 for ( std::size_t i=0; i<nofHits; ++i ) (*fHitsCollection)[i]->Print();
123 }
124}
125
126//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
127
128}
Definition of the B4c::CalorimeterSD class.
Calorimeter hit class.
void Initialize(G4HCofThisEvent *hitCollection) override
CalorHitsCollection * fHitsCollection
void EndOfEvent(G4HCofThisEvent *hitCollection) override
CalorimeterSD(const G4String &name, const G4String &hitsCollectionName, G4int nofCells)
G4bool ProcessHits(G4Step *step, G4TouchableHistory *history) override
G4THitsCollection< CalorHit > CalorHitsCollection
Definition CalorHit.hh:80

Applications | User Support | Publications | Collaboration