CMS 3D CMS Logo

hgcalDigitizer_cfi.py
Go to the documentation of this file.
1 import FWCore.ParameterSet.Config as cms
2 
3 from SimCalorimetry.HGCalSimProducers.hgcROCParameters_cfi import hgcROCSettings
4 
5 # Base configurations for HGCal digitizers
6 eV_per_eh_pair = 3.62
7 fC_per_ele = 1.6020506e-4
8 nonAgedCCEs = [1.0, 1.0, 1.0]
9 nonAgedNoises = [2100.0,2100.0,1600.0] #100,200,300 um (in electrons)
10 nonAgedNoises_v9 = [2000.0,2400.0,2000.0] # 120,200,300 um (in electrons)
11 thresholdTracksMIP = True
12 
13 ileakParam_600V = [0.993,-42.668]
14 ileakParam_800V = [0.996,-42.464]
15 HGCAL_ileakParam_toUse = cms.PSet(
16  ileakParam = cms.vdouble(ileakParam_600V)
17  )
18 
19 # line+log tdr 600V
20 cceParamFine_tdr600 = [1.5e+15, -3.00394e-17, 0.318083] #120
21 cceParamThin_tdr600 = [1.5e+15, -3.09878e-16, 0.211207] #200
22 cceParamThick_tdr600 = [6e+14, -7.96539e-16, 0.251751] #300
23 # line+log tdr 800V
24 cceParamFine_tdr800 = [4.2e+15, 2.35482e-18, 0.553187] #120
25 cceParamThin_tdr800 = [1.5e+15, -1.98109e-16, 0.280567] #200
26 cceParamThick_tdr800 = [6e+14, -5.24999e-16, 0.357616] #300
27 # line+log ttu 600V
28 cceParamFine_ttu600 = [1.5e+15, 9.98631e-18, 0.343774] #120
29 cceParamThin_ttu600 = [1.5e+15, -2.17083e-16, 0.304873] #200
30 cceParamThick_ttu600 = [6e+14, -8.01557e-16, 0.157375] #300
31 # line+log ttu 800V
32 cceParamFine_ttu800 = [1.5e+15, 3.35246e-17, 0.251679] #120
33 cceParamThin_ttu800 = [1.5e+15, -1.62096e-16, 0.293828] #200
34 cceParamThick_ttu800 = [6e+14, -5.95259e-16, 0.183929] #300
35 # line+log tdr 600V EPI
36 cceParamFine_epi600 = [3.5e+15, -9.73872e-19, 0.263812] #100
37 cceParamThin_epi600 = [1.5e+15, -3.09878e-16, 0.211207] #200
38 cceParamThick_epi600 = [6e+14, -7.96539e-16, 0.251751] #300
39 
40 HGCAL_cceParams_toUse = cms.PSet(
41  cceParamFine = cms.vdouble(cceParamFine_tdr600),
42  cceParamThin = cms.vdouble(cceParamThin_tdr600),
43  cceParamThick = cms.vdouble(cceParamThick_tdr600)
44  )
45 
46 HGCAL_noise_fC = cms.PSet(
47  scaleByDose = cms.bool(False),
48  scaleByDoseAlgo = cms.uint32(0),
49  doseMap = cms.string(""),
50  values = cms.vdouble( [x*fC_per_ele for x in nonAgedNoises] ), #100,200,300 um
51  )
52 
53 HGCAL_noise_heback = cms.PSet(
54  scaleByDose = cms.bool(False),
55  scaleByDoseAlgo = cms.uint32(0),
56  doseMap = cms.string(""), #empty dose map at begin-of-life
57  noise_MIP = cms.double(1./100.)
58  )
59 
60 HGCAL_chargeCollectionEfficiencies = cms.PSet(
61  values = cms.vdouble( nonAgedCCEs )
62  )
63 
64 HGCAL_noises = cms.PSet(
65  values = cms.vdouble([x for x in nonAgedNoises])
66  )
67 
68 # ECAL
69 hgceeDigitizer = cms.PSet(
70  accumulatorType = cms.string("HGCDigiProducer"),
71  hitCollection = cms.string("HGCHitsEE"),
72  digiCollection = cms.string("HGCDigisEE"),
73  NoiseGeneration_Method = cms.bool(True),
74  maxSimHitsAccTime = cms.uint32(100),
75  bxTime = cms.double(25),
76  eVPerEleHolePair = cms.double(eV_per_eh_pair),
77  tofDelay = cms.double(5),
78  geometryType = cms.uint32(0),
79  digitizationType = cms.uint32(0),
80  makeDigiSimLinks = cms.bool(False),
81  premixStage1 = cms.bool(False),
82  premixStage1MinCharge = cms.double(0),
83  premixStage1MaxCharge = cms.double(1e6),
84  useAllChannels = cms.bool(True),
85  verbosity = cms.untracked.uint32(0),
86  digiCfg = cms.PSet(
87  keV2fC = cms.double(0.044259), #1000 eV/3.62 (eV per e) / 6.24150934e3 (e per fC)
88  ileakParam = cms.PSet(refToPSet_ = cms.string("HGCAL_ileakParam_toUse")),
89  cceParams = cms.PSet(refToPSet_ = cms.string("HGCAL_cceParams_toUse")),
90  chargeCollectionEfficiencies = cms.PSet(refToPSet_ = cms.string("HGCAL_chargeCollectionEfficiencies")),
91  noise_fC = cms.PSet(refToPSet_ = cms.string("HGCAL_noise_fC")),
92  doTimeSamples = cms.bool(False),
93  thresholdFollowsMIP = cms.bool(thresholdTracksMIP),
94  feCfg = hgcROCSettings.clone()
95  )
96  )
97 
98 # HCAL front
99 hgchefrontDigitizer = cms.PSet(
100  accumulatorType = cms.string("HGCDigiProducer"),
101  hitCollection = cms.string("HGCHitsHEfront"),
102  digiCollection = cms.string("HGCDigisHEfront"),
103  NoiseGeneration_Method = cms.bool(True),
104  maxSimHitsAccTime = cms.uint32(100),
105  bxTime = cms.double(25),
106  tofDelay = cms.double(5),
107  geometryType = cms.uint32(0),
108  digitizationType = cms.uint32(0),
109  makeDigiSimLinks = cms.bool(False),
110  premixStage1 = cms.bool(False),
111  premixStage1MinCharge = cms.double(0),
112  premixStage1MaxCharge = cms.double(1e6),
113  useAllChannels = cms.bool(True),
114  verbosity = cms.untracked.uint32(0),
115  digiCfg = cms.PSet(
116  keV2fC = cms.double(0.044259), #1000 eV / 3.62 (eV per e) / 6.24150934e3 (e per fC)
117  ileakParam = cms.PSet(refToPSet_ = cms.string("HGCAL_ileakParam_toUse")),
118  cceParams = cms.PSet(refToPSet_ = cms.string("HGCAL_cceParams_toUse")),
119  chargeCollectionEfficiencies = cms.PSet(refToPSet_ = cms.string("HGCAL_chargeCollectionEfficiencies")),
120  noise_fC = cms.PSet(refToPSet_ = cms.string("HGCAL_noise_fC")),
121  doTimeSamples = cms.bool(False),
122  thresholdFollowsMIP = cms.bool(thresholdTracksMIP),
123  feCfg = hgcROCSettings.clone()
124  )
125 )
126 
127 # HCAL back
128 hgchebackDigitizer = cms.PSet(
129  accumulatorType = cms.string("HGCDigiProducer"),
130  hitCollection = cms.string("HcalHits"),
131  digiCollection = cms.string("HGCDigisHEback"),
132  NoiseGeneration_Method = cms.bool(True),
133  maxSimHitsAccTime = cms.uint32(100),
134  bxTime = cms.double(25),
135  tofDelay = cms.double(1),
136  geometryType = cms.uint32(0),
137  digitizationType = cms.uint32(1),
138  makeDigiSimLinks = cms.bool(False),
139  premixStage1 = cms.bool(False),
140  premixStage1MinCharge = cms.double(0),
141  premixStage1MaxCharge = cms.double(1e6),
142  useAllChannels = cms.bool(True),
143  verbosity = cms.untracked.uint32(0),
144  digiCfg = cms.PSet(
145  #0 empty digitizer, 1 calice digitizer, 2 realistic digitizer
146  algo = cms.uint32(2),
147  scaleByTileArea= cms.bool(False),
148  scaleBySipmArea= cms.bool(False),
149  sipmMap = cms.string("SimCalorimetry/HGCalSimProducers/data/sipmParams_geom-10.txt"),
150  noise = cms.PSet(refToPSet_ = cms.string("HGCAL_noise_heback")), #scales both for scint raddam and sipm dark current
151  keV2MIP = cms.double(1./675.0),
152  doTimeSamples = cms.bool(False),
153  nPEperMIP = cms.double(21.0),
154  nTotalPE = cms.double(7500),
155  xTalk = cms.double(0.01),
156  sdPixels = cms.double(1e-6), # this is additional photostatistics noise (as implemented), not sure why it's here...
157  thresholdFollowsMIP = cms.bool(thresholdTracksMIP),
158  feCfg = hgcROCSettings.clone(
159  adcNbits = 10, # standard ROC operations (was 2 bits more up to 11_0_0_pre12)
160  adcSaturation_fC = 68.75, # keep the adc LSB the same (i.e. set saturation one quarter value of pre12)
161  tdcSaturation_fC = 1000, # allow up to 1000 MIPs as a max range, including ToA mode
162  targetMIPvalue_ADC = 15, # to be used for HGCROC gain proposal
163  adcThreshold_fC = 0.5, # unchanged with respect to pre12
164  tdcOnset_fC = 55, # turn on TDC when 80% of the ADC range is reached (one quarter of pre12
165  # indicative at this point)
166  tdcForToAOnset_fC = cms.vdouble(12.,12.,12.), #turn ToA for 20% of the TDC threshold (indicative at this point)
167  )
168  )
169 )
170 
171 # HFNose
172 hfnoseDigitizer = cms.PSet(
173  accumulatorType = cms.string("HGCDigiProducer"),
174  hitCollection = cms.string("HFNoseHits"),
175  digiCollection = cms.string("HFNoseDigis"),
176  NoiseGeneration_Method = cms.bool(True),
177  maxSimHitsAccTime = cms.uint32(100),
178  bxTime = cms.double(25),
179  eVPerEleHolePair = cms.double(eV_per_eh_pair),
180  tofDelay = cms.double(5),
181  geometryType = cms.uint32(1),
182  digitizationType = cms.uint32(0),
183  makeDigiSimLinks = cms.bool(False),
184  premixStage1 = cms.bool(False),
185  premixStage1MinCharge = cms.double(0),
186  premixStage1MaxCharge = cms.double(1e6),
187  useAllChannels = cms.bool(True),
188  verbosity = cms.untracked.uint32(0),
189  digiCfg = cms.PSet(
190  keV2fC = cms.double(0.044259), #1000 eV/3.62 (eV per e) / 6.24150934e3 (e per fC)
191  ileakParam = cms.PSet(refToPSet_ = cms.string("HGCAL_ileakParam_toUse")),
192  cceParams = cms.PSet(refToPSet_ = cms.string("HGCAL_cceParams_toUse")),
193  chargeCollectionEfficiencies = cms.PSet(refToPSet_ = cms.string("HGCAL_chargeCollectionEfficiencies")),
194  noise_fC = cms.PSet(refToPSet_ = cms.string("HGCAL_noise_fC")),
195  doTimeSamples = cms.bool(False),
196  thresholdFollowsMIP = cms.bool(thresholdTracksMIP),
197  feCfg = hgcROCSettings.clone()
198  )
199  )
200 
201 # this bypasses the noise simulation
202 from Configuration.ProcessModifiers.premix_stage1_cff import premix_stage1
203 for _m in [hgceeDigitizer, hgchefrontDigitizer, hgchebackDigitizer, hfnoseDigitizer]:
204  premix_stage1.toModify(_m, premixStage1 = True)
205 
206 #function to set noise to aged HGCal
207 endOfLifeCCEs = [0.5, 0.5, 0.7]
208 endOfLifeNoises = [2400.0,2250.0,1750.0]
209 def HGCal_setEndOfLifeNoise(process,byDose=True,byDoseAlgo=0):
210  """includes all effects from radiation and gain choice"""
211  process=HGCal_setRealisticNoiseSi(process,byDose=byDose,byDoseAlgo=byDoseAlgo)
212  process=HGCal_setRealisticNoiseSci(process,byDose=byDose,byDoseAlgo=byDoseAlgo)
213  return process
214 
216  """include all effects except fluence impact on leakage current and CCE"""
217  #note: realistic electronics with Sci is not yet switched on
218  process=HGCal_setRealisticNoiseSi(process,byDose=True,byDoseAlgo=1)
219  return process
220 
221 def HGCal_setRealisticNoiseSi(process,byDose=True,byDoseAlgo=0):
222  process.HGCAL_noise_fC = cms.PSet(
223  scaleByDose = cms.bool(byDose),
224  scaleByDoseAlgo = cms.uint32(byDoseAlgo),
225  doseMap = cms.string("SimCalorimetry/HGCalSimProducers/data/doseParams_3000fb_fluka-3.5.15.9.txt"),
226  values = cms.vdouble( [x*fC_per_ele for x in endOfLifeNoises] ), #100,200,300 um
227  )
228  process.HGCAL_chargeCollectionEfficiencies = cms.PSet(
229  values = cms.vdouble(endOfLifeCCEs)
230  )
231  process.HGCAL_noises = cms.PSet(
232  values = cms.vdouble([x for x in endOfLifeNoises])
233  )
234  return process
235 
236 def HGCal_setRealisticNoiseSci(process,byDose=True,byDoseAlgo=0):
237  process.HGCAL_noise_heback = cms.PSet(
238  scaleByDose = cms.bool(byDose),
239  scaleByDoseAlgo = cms.uint32(byDoseAlgo),
240  doseMap = cms.string("SimCalorimetry/HGCalSimProducers/data/doseParams_3000fb_fluka-3.5.15.9.txt"),
241  noise_MIP = cms.double(1./5.) #uses noise map
242  )
243  return process
244 
245 
246 
247 def HGCal_disableNoise(process):
248  process.HGCAL_noise_fC = cms.PSet(
249  scaleByDose = cms.bool(False),
250  scaleByDoseAlgo = cms.uint32(0),
251  doseMap = cms.string(""),
252  values = cms.vdouble(0,0,0), #100,200,300 um
253  )
254  process.HGCAL_noise_heback = cms.PSet(
255  scaleByDose = cms.bool(False),
256  scaleByDoseAlgo = cms.uint32(0),
257  doseMap = cms.string(""),
258  noise_MIP = cms.double(0.) #zero noise
259  )
260  process.HGCAL_noises = cms.PSet(
261  values = cms.vdouble(0,0,0)
262  )
263  return process
264 
265 from Configuration.Eras.Modifier_phase2_hgcalV9_cff import phase2_hgcalV9
266 
267 phase2_hgcalV9.toModify( hgceeDigitizer,
268  geometryType = cms.uint32(1),
269 )
270 phase2_hgcalV9.toModify( hgchefrontDigitizer,
271  geometryType = cms.uint32(1),
272 )
273 phase2_hgcalV9.toModify( hgchebackDigitizer,
274  geometryType = cms.uint32(1),
275  hitCollection = cms.string("HGCHitsHEback"),
276 )
277 phase2_hgcalV9.toModify(HGCAL_noise_fC, values = [x*fC_per_ele for x in nonAgedNoises_v9])
278 phase2_hgcalV9.toModify(HGCAL_noises, values = [x for x in nonAgedNoises_v9])
hgcalDigitizer_cfi.HGCal_disableNoise
def HGCal_disableNoise(process)
Definition: hgcalDigitizer_cfi.py:247
hgcalDigitizer_cfi.HGCal_setRealisticNoiseSci
def HGCal_setRealisticNoiseSci(process, byDose=True, byDoseAlgo=0)
Definition: hgcalDigitizer_cfi.py:236
hgcalDigitizer_cfi.HGCal_setRealisticNoiseSi
def HGCal_setRealisticNoiseSi(process, byDose=True, byDoseAlgo=0)
Definition: hgcalDigitizer_cfi.py:221
hgcalDigitizer_cfi.HGCal_setEndOfLifeNoise
def HGCal_setEndOfLifeNoise(process, byDose=True, byDoseAlgo=0)
Definition: hgcalDigitizer_cfi.py:209
hgcalDigitizer_cfi.HGCal_setRealisticStartupNoise
def HGCal_setRealisticStartupNoise(process)
Definition: hgcalDigitizer_cfi.py:215