ATLAS Offline Software
Loading...
Searching...
No Matches
MinBiasScintillatorSD.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5//************************************************************
6//
7// Class MinBiasScintillatorSD
8// Sensitive detector for the Minimum Bias Scintillator
9//
10// Author: A. Solodkov
11// April 01, 2006.
12//
13//************************************************************
14
15#ifndef MINBIASSCINTILLATOR_MINBIASSCINTILLATORSD_H
16#define MINBIASSCINTILLATOR_MINBIASSCINTILLATORSD_H
17
18// Base class
19#include "G4VSensitiveDetector.hh"
20
21// For the hits
24
26
27// STL headers
28#include <vector>
29
30// G4 needed classes
31class G4Step;
32class G4TouchableHistory;
33
34class TileTBID;
35
36class MinBiasScintillatorSD: public G4VSensitiveDetector {
37 public:
38 MinBiasScintillatorSD(const G4String& name, const std::string& hitCollectionName,
39 const MinBiasScintSDOptions& opts);
41
42 virtual void Initialize(G4HCofThisEvent*) override final;
43
44 virtual G4bool ProcessHits(G4Step*, G4TouchableHistory*) override final;
45
47
48 private:
49 // Options for the SD configuration
51
52 // The hits collection
54
56
57 static const int N_SIDE = 2;
58 static const int N_PHI = 8;
59 static const int N_ETA = 2;
60 static const int N_CELLS = N_SIDE * N_PHI * N_ETA;
61 static const int N_DIST = N_SIDE * N_ETA;
62
63 inline int cell_index(int side, int phi, int eta) const {
64 return (side * N_PHI + phi) * N_ETA + eta;
65 }
66 inline int dist_index(int side, int eta) const {
67 return side * N_ETA + eta;
68 }
69
73 std::vector<int> m_numberOfHitsInCell;
74
77 std::vector<TileSimHit*> m_tempSimHit;
78
83
86 double m_deltaT;
87 std::vector<double> m_deltaTvec;
88
91 double m_timeCut;
92
97
99
100 G4double BirkLaw(const G4Step* aStep) const;
101
105
106 inline double deltaT(double time) const {
107 unsigned int i = 0;
108 double delta = m_deltaTvec[i++];
109 while (i < m_deltaTvec.size()) {
110 if (time > m_deltaTvec[i++] && time < m_deltaTvec[i++]) break;
111 delta = m_deltaTvec[i++];
112 }
113 return delta;
114 }
115};
116
117#endif //MINBIASSCINTILLATOR_MINBIASSCINTILLATORSD_H
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
Handle class for recording to StoreGate.
double m_deltaT
granularity in time for hits
~MinBiasScintillatorSD()=default
double deltaT(double time) const
function to provide correct deltaT bin width for given time
std::vector< TileSimHit * > m_tempSimHit
local temporary TileSimHit for each cell (size=nCell).
const MinBiasScintSDOptions m_options
G4double BirkLaw(const G4Step *aStep) const
MinBiasScintillatorSD(const G4String &name, const std::string &hitCollectionName, const MinBiasScintSDOptions &opts)
int cell_index(int side, int phi, int eta) const
Identifier m_channelID[N_CELLS]
array to cache look-ups of Identifiers for each cell (size=nCell).
std::vector< int > m_numberOfHitsInCell
number of contributions to the energy in each cell (size=nCell).
int dist_index(int side, int eta) const
double m_lateHitTime
time for hits which are above m_timeCut threshold it is equal to m_tilesize_deltaT - m_deltaT
virtual void Initialize(G4HCofThisEvent *) override final
double m_timeCut
max allowed time for hits
std::vector< double > m_deltaTvec
SG::WriteHandle< TileHitVector > m_HitColl
virtual G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
Helper class for TileCal offline identifiers of ancillary testbeam detectors and MBTS.