ATLAS Offline Software
Loading...
Searching...
No Matches
LUCID_DigitizationToolBox.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3*/
4
5
8#include <algorithm>
9#include <utility>
10
11#include "CLHEP/Random/RandFlat.h"
12#include "CLHEP/Random/RandGaussQ.h"
13#include "CLHEP/Random/RandPoissonQ.h"
14#include "GaudiKernel/MsgStream.h"
15#include "GaudiKernel/Bootstrap.h"
16#include "GaudiKernel/ITHistSvc.h"
17
18//--------------------------------------------------------------------------
20 m_tubeID(0),
21 m_npe (0),
22 m_npeGas(0),
23 m_npePmt(0),
24 m_qdc (0),
25 m_tdcTot(0),
26 m_tdc (0)
27{
28 m_digitContainer = nullptr;
29 m_tubeInfo = nullptr;
30
31 m_numTubes = 0;
33 m_qdcPedestal = 0;
37
39 m_AmpFactor = 0;
40 m_Q1bin = 0;
41 m_NoiseCharge = 0;
42 m_numDyinodes = 0;
46
48}
49
50//--------------------------------------------------------------------------
52 double qdcChannelsPerPE,
53 double qdcPedestal,
54 double qdcFedNoiseFactor,
55 double tdcPmtNoiseFactor,
56 double tdcFedNoiseFactor,
57 double TotalPMTgain,
58 double AmpFactor,
59 double Q1bin,
60 int NoiseCharge,
61 int numDyinodes,
62 double dynodeGammaFactor,
63 std::vector<double> pmtSmearing,
64 std::vector<double> pmtScaling,
65 std::vector<double> gasScaling,
66 std::vector<double> npeThreshold,
67 bool fillRootTree):
68 m_numTubes (numTubes),
69 m_qdcChannelsPerPE (qdcChannelsPerPE),
70 m_qdcPedestal (qdcPedestal),
71 m_qdcFedNoiseFactor (qdcFedNoiseFactor),
72 m_tdcPmtNoiseFactor (tdcPmtNoiseFactor),
73 m_tdcFedNoiseFactor (tdcFedNoiseFactor),
74 m_TotalPMTgain (TotalPMTgain),
75 m_AmpFactor (AmpFactor),
76 m_Q1bin (Q1bin),
77 m_NoiseCharge (NoiseCharge),
78 m_numDyinodes (numDyinodes),
79 m_dynodeGammaFactor (dynodeGammaFactor),
80 m_pmtSmearing (std::move(pmtSmearing)),
81 m_pmtScaling (std::move(pmtScaling)),
82 m_gasScaling (std::move(gasScaling)),
83 m_npeThreshold (std::move(npeThreshold)),
84 m_fillRootTree (fillRootTree),
85 m_tubeID (0),
86 m_npe (0),
87 m_npeGas (0),
88 m_npePmt (0),
89 m_qdc (0),
90 m_tdcTot (0),
91 m_tdc (0)
92{
93 m_gainPerDynode = std::pow(m_TotalPMTgain/m_dynodeGammaFactor, 1.0/static_cast<double>(m_numDyinodes));
94 m_ChargeToQdcFactor = (double)1.6e-19 * m_AmpFactor/m_Q1bin;
95
96 m_digitContainer = nullptr;
97 m_tubeInfo = nullptr;
98}
99
100//--------------------------------------------------------------------------
102 SmartIF<ITHistSvc> digitHistSvc{Gaudi::svcLocator()->service("THistSvc")};
103 if (!digitHistSvc) { return StatusCode::FAILURE; }
104
105 m_tubeInfo = new TTree("t", "LUCID_LUMI_SUMMARY");
106
107 StatusCode sc = digitHistSvc->regTree("/AANT/LUCID_LUMI_SUMMARY", m_tubeInfo);
108
109 if (sc.isFailure()) return StatusCode::FAILURE;
110
111 m_tubeInfo->Branch("tubeID", &m_tubeID, "tubeID/s");
112 m_tubeInfo->Branch("npe" , &m_npe , "npe/D");
113 m_tubeInfo->Branch("npeGas", &m_npeGas, "npeGas/s");
114 m_tubeInfo->Branch("npePmt", &m_npePmt, "npePmt/s");
115 m_tubeInfo->Branch("qdc" , &m_qdc , "qdc/s");
116 m_tubeInfo->Branch("tdc" , &m_tdc , "tdc/s");
117
118 return StatusCode::SUCCESS;
119}
120
121//--------------------------------------------------------------------------
123
124 assert(x >= INT_MIN-0.5);
125 assert(x <= INT_MAX+0.5);
126
127 if (x >= 0) return static_cast<int>(x+0.5);
128
129 return static_cast<int>(x-0.5);
130}
131
132//--------------------------------------------------------------------------
133StatusCode LUCID_DigitizationToolBox::createAndStoreDigit(unsigned short tubeID, CLHEP::HepRandomEngine* rndEngine) {
134
135 double npePmt = m_npePmt * m_pmtScaling[tubeID] * CLHEP::RandGaussQ::shoot(rndEngine, 1, m_pmtSmearing[tubeID]);
136 double npeGas = m_npeGas * m_gasScaling[tubeID];
137 double npeTot = npeGas + npePmt;
138 double qdcSmeared = DynodeChainSimulation(npeTot, rndEngine);
139
140 m_qdc = roundoff(qdcSmeared);
141
143
144 double tdcTotSmeared = 0;
145
146 tdcTotSmeared = CLHEP::RandGaussQ::shoot(rndEngine, m_tdcTot , m_tdcPmtNoiseFactor);
147 tdcTotSmeared = CLHEP::RandGaussQ::shoot(rndEngine, tdcTotSmeared, m_tdcFedNoiseFactor);
148 tdcTotSmeared = (npeTot) ? tdcTotSmeared/npeTot : tdcTotSmeared;
149 tdcTotSmeared = (tdcTotSmeared < 0) ? 0 : tdcTotSmeared;
150
151 m_tdc = static_cast<unsigned int>(tdcTotSmeared);
152
153 if (m_fillRootTree) m_tubeInfo->Fill();
154
155 m_digitContainer->push_back(new LUCID_Digit(m_tubeID,
156 m_npe,
157 m_npeGas,
158 m_npePmt,
159 m_qdc,
160 m_tdc,
161 m_npe > m_npeThreshold[tubeID]));
162 return StatusCode::SUCCESS;
163}
164
166StatusCode LUCID_DigitizationToolBox::fillDigitContainer(TimedHitCollection<LUCID_SimHit>& thpclucid, CLHEP::HepRandomEngine* rndEngine) {
167
168 for (m_tubeID = 0; m_tubeID < m_numTubes; m_tubeID++) {
169
170 if (m_tubeID >= 16 && m_tubeID <= 19) continue; // skip fiber tubes on sideA
171 if (m_tubeID >= 36 && m_tubeID <= 39) continue; // skip fiber tubes on sideC
172
174
175 TimedHitCollection<LUCID_SimHit> thpc = thpclucid;
177
178 while (thpc.nextDetectorElement(i, e)) for (hitIt = i; hitIt != e; ++hitIt) {
179
180 if (m_tubeID != (*hitIt)->GetTubeID()) continue;
181
182 if (!(*hitIt)->isDetected(rndEngine)) continue; // applying quantum efficiency
183
184 if ((*hitIt)->GetGenVolume() == 1) m_npeGas++;
185 if ((*hitIt)->GetGenVolume() == 2) m_npeGas++;
186 if ((*hitIt)->GetGenVolume() == 3) m_npePmt++;
187
188 m_tdcTot += (*hitIt)->GetPreStepTime();
189 }
190
191 if (createAndStoreDigit(m_tubeID, rndEngine).isFailure()) return StatusCode::FAILURE;
192 }
193
194 return StatusCode::SUCCESS;
195
196}
197
199StatusCode LUCID_DigitizationToolBox::fillDigitContainer(LUCID_SimHitCollection* thpclucid, CLHEP::HepRandomEngine* rndEngine) {
200
201 for (m_tubeID = 0; m_tubeID < m_numTubes; m_tubeID++) {
202
203 if (m_tubeID >= 16 && m_tubeID <= 19) continue; // skip fiber tubes on sideA
204 if (m_tubeID >= 36 && m_tubeID <= 39) continue; // skip fiber tubes on sideC
205
207
209 LUCID_SimHitCollection::const_iterator hitItE = thpclucid->end();
210
211 for (; hitIt != hitItE; ++hitIt) {
212
213 if (m_tubeID != (*hitIt).GetTubeID()) continue;
214
215 if (!(*hitIt).isDetected(rndEngine)) continue; // applying quantum efficiency
216
217 if ((*hitIt).GetGenVolume() == 1) m_npeGas++;
218 if ((*hitIt).GetGenVolume() == 2) m_npeGas++;
219 if ((*hitIt).GetGenVolume() == 3) m_npePmt++;
220
221 m_tdcTot += (*hitIt).GetPreStepTime();
222 }
223
224 if (createAndStoreDigit(m_tubeID, rndEngine).isFailure()) return StatusCode::FAILURE;
225 }
226
227 return StatusCode::SUCCESS;
228}
229
230
231//--------------------------------------------------------------------------
232StatusCode LUCID_DigitizationToolBox::recordContainers(const ServiceHandle<StoreGateSvc>& digitsStore, const std::string& key_digitCnt) {
233
235
236 StatusCode sc = digitsStore->record(m_digitContainer, key_digitCnt);
237
238 return sc;
239}
240
241//--------------------------------------------------------------------------
242double LUCID_DigitizationToolBox::DynodeChainSimulation(double npe, CLHEP::HepRandomEngine* rndEngine) const{
243
244 double smearedQDC = 0;
245
246 if (npe == 0) return CLHEP::RandGaussQ::shoot(rndEngine, npe, m_ChargeToQdcFactor * m_NoiseCharge) + m_qdcPedestal;
247
248 for (int i_pe=1; i_pe <= m_numDyinodes; i_pe++ ) {
249
250 if (i_pe == 1) smearedQDC = DynodeGainSmearing(npe * m_gainPerDynode * m_dynodeGammaFactor, rndEngine);
251 else smearedQDC = DynodeGainSmearing(smearedQDC * m_gainPerDynode , rndEngine);
252 }
253
254 smearedQDC = smearedQDC*m_ChargeToQdcFactor;
255
256 return static_cast<double>(CLHEP::RandGaussQ::shoot(rndEngine, smearedQDC, m_NoiseCharge * m_ChargeToQdcFactor)) + m_qdcPedestal;
257}
258
259//--------------------------------------------------------------------------
260double LUCID_DigitizationToolBox::DynodeGainSmearing(double npe, CLHEP::HepRandomEngine* rndEngine){
261
262 if (npe < 20) return CLHEP::RandPoissonQ::shoot(rndEngine, npe);
263 else return CLHEP::RandGaussQ::shoot (rndEngine, npe, std::sqrt(npe));
264}
static Double_t sc
AtlasHitsVector< LUCID_SimHit > LUCID_SimHitCollection
the preferred mechanism to access information from the different event stores in a pileup job.
#define x
CONT::const_iterator const_iterator
const_iterator begin() const
const_iterator end() const
LUCID_DigitContainer * m_digitContainer
static double DynodeGainSmearing(double npe, CLHEP::HepRandomEngine *rndEngine)
StatusCode fillDigitContainer(TimedHitCollection< LUCID_SimHit > &, CLHEP::HepRandomEngine *)
double DynodeChainSimulation(double npe, CLHEP::HepRandomEngine *rndEngine) const
StatusCode recordContainers(const ServiceHandle< StoreGateSvc > &, const std::string &)
static unsigned int roundoff(double x)
StatusCode createAndStoreDigit(unsigned short tubeID, CLHEP::HepRandomEngine *rndEngine)
bool nextDetectorElement(const_iterator &b, const_iterator &e)
sets an iterator range with the hits of current detector element returns a bool when done
TimedVector::const_iterator const_iterator
STL namespace.