CMS 3D CMS Logo

HcalPulseShapes.cc
Go to the documentation of this file.
6 #include "CLHEP/Random/RandFlat.h"
8 
9 // #include "CalibCalorimetry/HcalAlgos/interface/HcalDbASCIIIO.h"
10 #include <cmath>
11 #include <iostream>
12 #include <fstream>
13 #include "TMath.h"
14 
16 : theDbService(nullptr),
17  theShapes()
18 {
19 /*
20 
21 Reco MC
22 -------------------------------------------------------------------------------------------
23 000 not used (reserved)
24 101 101 hpdShape_ HPD (original version)
25 102 102 =101 HPD BV 30 volts in HBP iphi54
26 103 123 hpdShape_v2,hpdShapeMC_v2 HPD (2011. oct version)
27 104 124 hpdBV30Shape_v2,hpdBV30ShapeMC_v2 HPD bv30 in HBP iph54
28 105 125 hpdShape_v2,hpdShapeMC_v2 HPD (2011.11.12 version)
29 201 201 siPMShapeHO_ SiPMs Zecotec shape (HO)
30 202 202 =201, SiPMs Hamamatsu shape (HO)
31 205 203 siPMShapeData2017_,siPMShapeMC2017_ SiPMs from Data, Hamamatsu shape (HE 2017)
32 207 206 siPMShapeData2018_,siPMShapeMC2018_ SiPMs from Data, Hamamatsu shape (HE 2018)
33 301 301 hfShape_ regular HF PMT shape
34 401 401 regular ZDC shape
35 -------------------------------------------------------------------------------------------
36 
37 */
38 
39 
40  float ts1, ts2, ts3, thpd, tpre, wd1, wd2, wd3;
41 
42  // HPD Shape Version 1 (used before CMSSW5, until Oct 2011)
43  ts1=8. ; ts2=10. ; ts3=29.3; thpd=4.0; tpre=9.0; wd1=2.0; wd2=0.7; wd3=1.0;
44  computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdShape_);
45  theShapes[101] = &hpdShape_;
46  theShapes[102] = theShapes[101];
47 
48  // HPD Shape Version 2 for CMSSW 5. Nov 2011 (RECO and MC separately)
49  ts1=8. ; ts2=10. ; ts3=25.0; thpd=4.0; tpre=9.0; wd1=2.0; wd2=0.7; wd3=1.0;
50  computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdShape_v2);
51  theShapes[103] = &hpdShape_v2;
52 
53  ts1=8. ; ts2=10. ; ts3=29.3; thpd=4.0; tpre=7.0; wd1=2.0; wd2=0.7; wd3=1.0;
54  computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdShapeMC_v2);
55  theShapes[123] = &hpdShapeMC_v2;
56 
57  // HPD Shape Version 3 for CMSSW 5. Nov 2011 (RECO and MC separately)
58  ts1=8. ; ts2=19. ; ts3=29.3; thpd=4.0; tpre=9.0; wd1=2.0; wd2=0.7; wd3=0.32;
59  computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdShape_v3);
60  theShapes[105] = &hpdShape_v3;
61 
62  ts1=8. ; ts2=10. ; ts3=22.3; thpd=4.0; tpre=7.0; wd1=2.0; wd2=0.7; wd3=1.0;
63  computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdShapeMC_v3);
64  theShapes[125] = &hpdShapeMC_v3;
65 
66  // HPD with Bias Voltage 30 volts, wider pulse. (HBPlus iphi54)
67 
68  ts1=8. ; ts2=12. ; ts3=31.7; thpd=9.0; tpre=9.0; wd1=2.0; wd2=0.7; wd3=1.0;
69  computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdBV30Shape_v2);
70  theShapes[104] = &hpdBV30Shape_v2;
71 
72  ts1=8. ; ts2=12. ; ts3=31.7; thpd=9.0; tpre=9.0; wd1=2.0; wd2=0.7; wd3=1.0;
73  computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdBV30ShapeMC_v2);
75 
76  // HF and SiPM
77 
84 
85  theShapes[201] = &siPMShapeHO_;
86  theShapes[202] = theShapes[201];
91  theShapes[301] = &hfShape_;
92  //theShapes[401] = new CaloCachedShapeIntegrator(&theZDCShape);
93 
94 }
95 
96 
98 }
99 
100 
102 {
103  edm::ESHandle<HcalDbService> conditions;
104  es.get<HcalDbRecord>().get(conditions);
105  theDbService = conditions.product();
106 }
107 
109 {
110  theDbService = conditions;
111 }
112 
113 //void HcalPulseShapes::computeHPDShape()
114 void HcalPulseShapes::computeHPDShape(float ts1, float ts2, float ts3, float thpd, float tpre,
115  float wd1, float wd2, float wd3, Shape &tmphpdShape_)
116 {
117 
118 // pulse shape time constants in ns
119 /*
120  const float ts1 = 8.; // scintillation time constants : 1,2,3
121  const float ts2 = 10.;
122  const float ts3 = 29.3;
123  const float thpd = 4.; // HPD current collection drift time
124  const float tpre = 9.; // preamp time constant (refit on TB04 data)
125 
126  const float wd1 = 2.; // relative weights of decay exponents
127  const float wd2 = 0.7;
128  const float wd3 = 1.;
129 */
130  // pulse shape components over a range of time 0 ns to 255 ns in 1 ns steps
131  unsigned int nbin = 256;
132  tmphpdShape_.setNBin(nbin);
133  std::vector<float> ntmp(nbin,0.0); // zeroing output pulse shape
134  std::vector<float> nth(nbin,0.0); // zeroing HPD drift shape
135  std::vector<float> ntp(nbin,0.0); // zeroing Binkley preamp shape
136  std::vector<float> ntd(nbin,0.0); // zeroing Scintillator decay shape
137 
138  unsigned int i,j,k;
139  float norm;
140 
141  // HPD starts at I and rises to 2I in thpd of time
142  norm=0.0;
143  for(j=0;j<thpd && j<nbin;j++){
144  nth[j] = 1.0 + ((float)j)/thpd;
145  norm += nth[j];
146  }
147  // normalize integrated current to 1.0
148  for(j=0;j<thpd && j<nbin;j++){
149  nth[j] /= norm;
150  }
151 
152  // Binkley shape over 6 time constants
153  norm=0.0;
154  for(j=0;j<6*tpre && j<nbin;j++){
155  ntp[j] = ((float)j)*exp(-((float)(j*j))/(tpre*tpre));
156  norm += ntp[j];
157  }
158  // normalize pulse area to 1.0
159  for(j=0;j<6*tpre && j<nbin;j++){
160  ntp[j] /= norm;
161  }
162 
163 // ignore stochastic variation of photoelectron emission
164 // <...>
165 
166 // effective tile plus wave-length shifter decay time over 4 time constants
167  unsigned int tmax = 6 * (int)ts3;
168 
169  norm=0.0;
170  for(j=0;j<tmax && j<nbin;j++){
171  ntd[j] = wd1 * exp(-((float)j)/ts1) +
172  wd2 * exp(-((float)j)/ts2) +
173  wd3 * exp(-((float)j)/ts3) ;
174  norm += ntd[j];
175  }
176  // normalize pulse area to 1.0
177  for(j=0;j<tmax && j<nbin;j++){
178  ntd[j] /= norm;
179  }
180 
181  unsigned int t1,t2,t3,t4;
182  for(i=0;i<tmax && i<nbin;i++){
183  t1 = i;
184  // t2 = t1 + top*rand;
185  // ignoring jitter from optical path length
186  t2 = t1;
187  for(j=0;j<thpd && j<nbin;j++){
188  t3 = t2 + j;
189  for(k=0;k<4*tpre && k<nbin;k++){ // here "4" is set deliberately,
190  t4 = t3 + k; // as in test fortran toy MC ...
191  if(t4<nbin){
192  unsigned int ntb=t4;
193  ntmp[ntb] += ntd[i]*nth[j]*ntp[k];
194  }
195  }
196  }
197  }
198 
199  // normalize for 1 GeV pulse height
200  norm = 0.;
201  for(i=0;i<nbin;i++){
202  norm += ntmp[i];
203  }
204 
205  for(i=0; i<nbin; i++){
206  ntmp[i] /= norm;
207  }
208 
209  for(i=0; i<nbin; i++){
210  tmphpdShape_.setShapeBin(i,ntmp[i]);
211  }
212 }
213 
215  // first create pulse shape over a range of time 0 ns to 255 ns in 1 ns steps
216  unsigned int nbin = 256;
217  hfShape_.setNBin(nbin);
218  std::vector<float> ntmp(nbin,0.0); //
219 
220  const float k0=0.7956; // shape parameters
221  const float p2=1.355;
222  const float p4=2.327;
223  const float p1=4.3; // position parameter
224 
225  float norm = 0.0;
226 
227  for(unsigned int j = 0; j < 25 && j < nbin; ++j){
228 
229  float r0 = j-p1;
230  float sigma0 = (r0<0) ? p2 : p2*p4;
231  r0 /= sigma0;
232  if(r0 < k0) ntmp[j] = exp(-0.5*r0*r0);
233  else ntmp[j] = exp(0.5*k0*k0-k0*r0);
234  norm += ntmp[j];
235  }
236  // normalize pulse area to 1.0
237  for(unsigned int j = 0; j < 25 && j < nbin; ++j){
238  ntmp[j] /= norm;
239  hfShape_.setShapeBin(j,ntmp[j]);
240  }
241 }
242 
243 
245 {
246 
247  //Combination of all phase scan data (May,Jul,Oct2017)
248  //runs: 294736-294740, 294929-294950, 298594-298598 and 305744-305758
249 
250  unsigned int nbin = 250;
251 
252  std::array<float, 250> nt {{
253  5.22174e-12, 7.04852e-10, 3.49584e-08, 7.78029e-07, 9.11847e-06, 6.39666e-05, 0.000297587, 0.000996661, 0.00256618, 0.00535396,
254  0.00944073, 0.0145521, 0.020145, 0.0255936, 0.0303632, 0.0341078, 0.0366849, 0.0381183, 0.0385392, 0.0381327,
255  0.0370956, 0.0356113, 0.0338366, 0.0318978, 0.029891, 0.0278866, 0.0259336, 0.0240643, 0.0222981, 0.0206453,
256  0.0191097, 0.0176902, 0.0163832, 0.0151829, 0.0140826, 0.0130752, 0.0121533, 0.01131, 0.0105382, 0.00983178,
257  0.00918467, 0.00859143, 0.00804709, 0.0075471, 0.00708733, 0.00666406, 0.00627393, 0.00591389, 0.00558122, 0.00527344,
258  0.00498834, 0.00472392, 0.00447837, 0.00425007, 0.00403754, 0.00383947, 0.00365465, 0.00348199, 0.00332052, 0.00316934,
259  0.00302764, 0.0028947, 0.00276983, 0.00265242, 0.00254193, 0.00243785, 0.00233971, 0.00224709, 0.0021596, 0.00207687,
260  0.0019986, 0.00192447, 0.00185421, 0.00178756, 0.0017243, 0.00166419, 0.00160705, 0.00155268, 0.00150093, 0.00145162,
261  0.00140461, 0.00135976, 0.00131696, 0.00127607, 0.00123699, 0.00119962, 0.00116386, 0.00112963, 0.00109683, 0.0010654,
262  0.00103526, 0.00100634, 0.000978578, 0.000951917, 0.000926299, 0.000901672, 0.000877987, 0.000855198, 0.00083326, 0.000812133,
263  0.000791778, 0.000772159, 0.000753242, 0.000734994, 0.000717384, 0.000700385, 0.000683967, 0.000668107, 0.000652779, 0.00063796,
264  0.000623629, 0.000609764, 0.000596346, 0.000583356, 0.000570777, 0.000558592, 0.000546785, 0.00053534, 0.000524243, 0.000513481,
265  0.00050304, 0.000492907, 0.000483072, 0.000473523, 0.000464248, 0.000455238, 0.000446483, 0.000437974, 0.0004297, 0.000421655,
266  0.00041383, 0.000406216, 0.000398807, 0.000391595, 0.000384574, 0.000377736, 0.000371076, 0.000364588, 0.000358266, 0.000352104,
267  0.000346097, 0.00034024, 0.000334528, 0.000328956, 0.00032352, 0.000318216, 0.000313039, 0.000307986, 0.000303052, 0.000298234,
268  0.000293528, 0.000288931, 0.000284439, 0.00028005, 0.000275761, 0.000271567, 0.000267468, 0.000263459, 0.000259538, 0.000255703,
269  0.000251951, 0.00024828, 0.000244688, 0.000241172, 0.00023773, 0.000234361, 0.000231061, 0.00022783, 0.000224666, 0.000221566,
270  0.000218528, 0.000215553, 0.000212636, 0.000209778, 0.000206977, 0.00020423, 0.000201537, 0.000198896, 0.000196307, 0.000193767,
271  0.000191275, 0.000188831, 0.000186432, 0.000184079, 0.000181769, 0.000179502, 0.000177277, 0.000175092, 0.000172947, 0.000170841,
272  0.000168772, 0.000166741, 0.000164745, 0.000162785, 0.000160859, 0.000158967, 0.000157108, 0.00015528, 0.000153484, 0.000151719,
273  0.000149984, 0.000148278, 0.000146601, 0.000144951, 0.000143329, 0.000141734, 0.000140165, 0.000138622, 0.000137104, 0.00013561,
274  0.000134141, 0.000132695, 0.000131272, 0.000129871, 0.000128493, 0.000127136, 0.000125801, 0.000124486, 0.000123191, 0.000121917,
275  0.000120662, 0.000119426, 0.000118209, 0.00011701, 0.000115829, 0.000114665, 0.000113519, 0.00011239, 0.000111278, 0.000110182,
276  0.000109102, 0.000108037, 0.000106988, 0.000105954, 0.000104935, 0.00010393, 0.000102939, 0.000101963, 0.000101, 0.000100051,
277  9.91146e-05, 9.81915e-05, 9.7281e-05, 9.63831e-05, 9.54975e-05, 9.46239e-05, 9.37621e-05, 9.2912e-05, 9.20733e-05, 9.12458e-05
278  }};
279 
281 
282  double norm = 0.;
283  for (unsigned int j = 0; j < nbin; ++j) {
284  norm += (nt[j]>0) ? nt[j] : 0.;
285  }
286 
287  for (unsigned int j = 0; j < nbin; ++j) {
288  nt[j] /= norm;
290  }
291 }
292 
293 
294 
296 {
297  //From Jay Lawhorn: derived from data Edward Laird phase scan may2017
298  //https://indico.cern.ch/event/641978/contributions/2604491/attachments/1468666/2271582/17-05-31-hcal-hep17-pulse-shape.pdf
299  //Run numbers are 294736-294740 and 294929-294950
300 
301  unsigned int nbin = 250;
302 
303  std::array<float, 250> nt {{
304  3.97958e-29, 1.11634e-22, 9.96106e-18, 6.25334e-14, 5.08863e-11, 8.59141e-09, 4.32285e-07, 8.56617e-06, 8.28549e-05, 0.000461447,
305  0.00168052, 0.00441395, 0.00901637, 0.0151806, 0.0220314, 0.028528, 0.0338471, 0.0375578, 0.0395985, 0.0401567,
306  0.0395398, 0.0380776, 0.0360669, 0.0337474, 0.0312984, 0.0288457, 0.0264721, 0.0242276, 0.0221393, 0.0202181,
307  0.0184647, 0.0168731, 0.0154335, 0.0141346, 0.0129639, 0.0119094, 0.0109594, 0.0101031, 0.0093305, 0.00863267,
308  0.0080015, 0.00742977, 0.00691107, 0.00643969, 0.00601059, 0.00561931, 0.00526188, 0.00493483, 0.00463505, 0.00435981,
309  0.00410667, 0.00387348, 0.00365832, 0.00345949, 0.00327547, 0.0031049, 0.00294656, 0.00279938, 0.00266237, 0.00253467,
310  0.00241548, 0.0023041, 0.00219989, 0.00210227, 0.00201072, 0.00192476, 0.00184397, 0.00176795, 0.00169634, 0.00162884,
311  0.00156512, 0.00150494, 0.00144803, 0.00139418, 0.00134317, 0.00129481, 0.00124894, 0.00120537, 0.00116398, 0.00112461,
312  0.00108715, 0.00105147, 0.00101747, 0.000985042, 0.000954096, 0.000924545, 0.000896308, 0.000869311, 0.000843482, 0.000818758,
313  0.000795077, 0.000772383, 0.000750623, 0.000729747, 0.00070971, 0.000690466, 0.000671977, 0.000654204, 0.00063711, 0.000620663,
314  0.000604831, 0.000589584, 0.000574894, 0.000560735, 0.000547081, 0.00053391, 0.0005212, 0.000508929, 0.000497078, 0.000485628,
315  0.000474561, 0.000463862, 0.000453514, 0.000443501, 0.000433811, 0.000424429, 0.000415343, 0.00040654, 0.00039801, 0.000389741,
316  0.000381722, 0.000373944, 0.000366398, 0.000359074, 0.000351964, 0.00034506, 0.000338353, 0.000331838, 0.000325505, 0.00031935,
317  0.000313365, 0.000307544, 0.000301881, 0.000296371, 0.000291009, 0.000285788, 0.000280705, 0.000275755, 0.000270932, 0.000266233,
318  0.000261653, 0.00025719, 0.000252837, 0.000248593, 0.000244454, 0.000240416, 0.000236475, 0.00023263, 0.000228876, 0.000225212,
319  0.000221633, 0.000218138, 0.000214724, 0.000211389, 0.00020813, 0.000204945, 0.000201831, 0.000198787, 0.000195811, 0.0001929,
320  0.000190053, 0.000187268, 0.000184543, 0.000181876, 0.000179266, 0.000176711, 0.00017421, 0.000171761, 0.000169363, 0.000167014,
321  0.000164713, 0.000162459, 0.00016025, 0.000158086, 0.000155964, 0.000153885, 0.000151847, 0.000149848, 0.000147888, 0.000145966,
322  0.000144081, 0.000142232, 0.000140418, 0.000138638, 0.000136891, 0.000135177, 0.000133494, 0.000131843, 0.000130221, 0.00012863,
323  0.000127066, 0.000125531, 0.000124023, 0.000122543, 0.000121088, 0.000119658, 0.000118254, 0.000116874, 0.000115518, 0.000114185,
324  0.000112875, 0.000111587, 0.000110321, 0.000109076, 0.000107851, 0.000106648, 0.000105464, 0.000104299, 0.000103154, 0.000102027,
325  0.000100918, 9.98271e-05, 9.87537e-05, 9.76974e-05, 9.66578e-05, 9.56346e-05, 9.46274e-05, 9.3636e-05, 9.26599e-05, 9.16989e-05,
326  9.07526e-05, 8.98208e-05, 8.89032e-05, 8.79995e-05, 8.71093e-05, 8.62325e-05, 8.53688e-05, 8.45179e-05, 8.36796e-05, 8.28536e-05,
327  8.20397e-05, 8.12376e-05, 8.04471e-05, 7.96681e-05, 7.89002e-05, 7.81433e-05, 7.73972e-05, 7.66616e-05, 7.59364e-05, 7.52213e-05,
328  7.45163e-05, 7.3821e-05, 7.31354e-05, 7.24592e-05, 7.17923e-05, 7.11345e-05, 7.04856e-05, 6.98455e-05, 6.9214e-05, 6.8591e-05
329  }};
330 
331 
333 
334  double norm = 0.;
335  for (unsigned int j = 0; j < nbin; ++j) {
336  norm += (nt[j]>0) ? nt[j] : 0.;
337  }
338 
339  for (unsigned int j = 0; j < nbin; ++j) {
340  nt[j] /= norm;
342  }
343 }
344 
346 {
347 
348  unsigned int nbin = 128;
349 
350 //From Jake Anderson: toy MC convolution of SiPM pulse + WLS fiber shape + SiPM nonlinear response
351  std::array<float, 128> nt {{
352  2.782980485851731e-6, 4.518134885954626e-5, 2.7689305197392056e-4, 9.18328418900969e-4, .002110072599166349, .003867856860331454, .006120046224897771, .008754774090536956,
353  0.0116469503358586, .01467007449455966, .01770489955229477, .02064621450689512, .02340678093764222, .02591874610854916, .02813325527435303, 0.0300189241965647,
354  .03155968107671164, .03275234052577155, .03360415306318798, .03413048377960748, .03435270899678218, .03429637464659661, .03398962975487166, .03346192884394954,
355  .03274298516247742, .03186195009136525, .03084679116113031, 0.0297238406141036, .02851748748929785, .02724998816332392, .02594137274487424, .02460942736731527,
356  .02326973510736116, .02193576080366117, 0.0206189674254987, .01932895378564653, 0.0180736052958666, .01685925112650875, 0.0156908225633535, .01457200857138456,
357  .01350540559602467, .01249265947824805, .01153459805300423, .01063135355597282, .009782474412011936, .008987026319784546, 0.00824368281357106, .007550805679909604,
358  .006906515742762193, .006308754629755056, .005755338185695127, .005244002229973356, .004772441359900532, .004338341490928299, .003939406800854143, 0.00357338171220501,
359  0.0032380685079891, .002931341133259233, .002651155690306086, .002395558090237333, .002162689279320922, .001950788415487319, .001758194329648101, .001583345567913682,
360  .001424779275191974, .001281129147671334, 0.00115112265163774, .001033577678808199, 9.273987838127585e-4, 8.315731274976846e-4, 7.451662302008696e-4, 6.673176219006913e-4,
361  5.972364609644049e-4, 5.341971801529036e-4, 4.775352065178378e-4, 4.266427928961177e-4, 3.8096498904225923e-4, 3.3999577417327287e-4, 3.032743659102713e-4, 2.703817158798329e-4,
362  2.4093719775272793e-4, 2.145954900503894e-4, 1.9104365317752797e-4, 1.6999839784346724e-4, 1.5120354022478893e-4, 1.3442763782650755e-4, 1.1946179895521507e-4, 1.0611765796993575e-4,
363  9.422550797617687e-5, 8.363258233342666e-5, 7.420147621931836e-5, 6.580869950304933e-5, 5.834335229919868e-5, 5.17059147771959e-5, 4.5807143072062634e-5, 4.0567063461299446e-5,
364  3.591405732740723e-5, 3.178402980354131e-5, 2.811965539165646e-5, 2.4869694240316126e-5, 2.1988373166730962e-5, 1.9434825899529382e-5, 1.717258740121378e-5, 1.5169137499243157e-5,
365  1.339548941011129e-5, 1.1825819079078403e-5, 1.0437131581057595e-5, 9.208961130078894e-6, 8.12310153137994e-6, 7.163364176588591e-6, 6.315360932244386e-6, 5.566309502463164e-6,
366  4.904859063429651e-6, 4.320934164082596e-6, 3.8055950719111903e-6, 3.350912911083174e-6, 2.9498580949517117e-6, 2.596200697612328e-6, 2.2844215378879293e-6, 2.0096328693141094e-6,
367  1.7675076766686654e-6, 1.5542166787225756e-6, 1.366372225473431e-6, 1.200978365778838e-6, 1.0553864128982371e-6, 9.272554464808518e-7, 8.145171945902259e-7, 7.153448381918271e-7
368  }};
369 
370  siPMShapeHO_.setNBin(nbin);
371 
372  double norm = 0.;
373  for (unsigned int j = 0; j < nbin; ++j) {
374  norm += (nt[j]>0) ? nt[j] : 0.;
375  }
376 
377  for (unsigned int j = 0; j < nbin; ++j) {
378  nt[j] /= norm;
380  }
381 }
382 
384 {
385  //numerical convolution of SiPM pulse + WLS fiber shape
386  std::vector<double> nt = convolve(nBinsSiPM_,analyticPulseShapeSiPMHE,Y11203);
387 
389 
390  //skip first bin, always 0
391  double norm = 0.;
392  for (unsigned int j = 1; j <= nBinsSiPM_; ++j) {
393  norm += (nt[j]>0) ? nt[j] : 0.;
394  }
395 
396  for (unsigned int j = 1; j <= nBinsSiPM_; ++j) {
397  nt[j] /= norm;
398  siPMShapeMC2017_.setShapeBin(j,nt[j]);
399  }
400 }
401 
403 {
404  //numerical convolution of SiPM pulse + WLS fiber shape
405  std::vector<double> nt = convolve(nBinsSiPM_,analyticPulseShapeSiPMHE,Y11206);
406 
408 
409  //Aligning 206 phase closer to 205 in order to have good reco agreement
410  int shift = -2;
411 
412  //skip first bin, always 0
413  double norm = 0.;
414  for (unsigned int j = std::max(1,-1*shift); j<=nBinsSiPM_; j++) {
415  norm += std::max(0., nt[j-shift]);
416  }
417  double normInv=1./norm;
418  for ( int j = 1; j<=nBinsSiPM_; j++) {
419  if ( j-shift>=0 ) {
420  nt[j-shift]*=normInv;
421  siPMShapeMC2018_.setShapeBin(j,nt[j-shift]);
422  }
423  else{
425  }
426  }
427 }
428 
430 HcalPulseShapes::getShape(int shapeType) const
431 {
432  ShapeMap::const_iterator shapeMapItr = theShapes.find(shapeType);
433  if(shapeMapItr == theShapes.end()) {
434  throw cms::Exception("HcalPulseShapes") << "unknown shapeType";
435  return hpdShape_; // should not return this, but...
436  } else {
437  return *(shapeMapItr->second);
438  }
439 }
440 
441 
443 HcalPulseShapes::shape(const HcalDetId & detId) const
444 {
445  if(!theDbService) {
446  return defaultShape(detId);
447  }
448  int shapeType = theDbService->getHcalMCParam(detId)->signalShape();
449 
450  ShapeMap::const_iterator shapeMapItr = theShapes.find(shapeType);
451  if(shapeMapItr == theShapes.end()) {
452  return defaultShape(detId);
453  } else {
454  return *(shapeMapItr->second);
455  }
456 }
457 
460 {
461  if(!theDbService) {
462  return defaultShape(detId);
463  }
464  int shapeType = theDbService->getHcalRecoParam(detId.rawId())->pulseShapeID();
465 
466  ShapeMap::const_iterator shapeMapItr = theShapes.find(shapeType);
467  if(shapeMapItr == theShapes.end()) {
468  return defaultShape(detId);
469  } else {
470  return *(shapeMapItr->second);
471  }
472 }
473 
474 
477 {
478  edm::LogWarning("HcalPulseShapes") << "Cannot find HCAL MC Params ";
479  HcalSubdetector subdet = detId.subdet();
480  switch(subdet) {
481  case HcalBarrel:
482  return hbShape();
483  case HcalEndcap:
484  return heShape();
485  case HcalForward:
486  return hfShape();
487  case HcalOuter:
488  //FIXME doesn't look for SiPMs
489  return hoShape(false);
490  default:
491  throw cms::Exception("HcalPulseShapes") << "unknown detId";
492  break;
493  }
494 }
495 
496 //SiPM helpers
497 
498 inline double gexp(double t, double A, double c, double t0, double s) {
499  static double const root2(sqrt(2));
500  return -A*0.5*exp(c*t+0.5*c*c*s*s-c*s)*(erf(-0.5*root2/s*(t-t0+c*s*s))-1);
501 }
502 
503 inline double onePulse(double t, double A, double sigma, double theta, double m) {
504  return (t<theta) ? 0 : A*TMath::LogNormal(t,sigma,theta,m);
505 }
506 
508  // HO SiPM pulse shape fit from Jake Anderson ca. 2013
509  double A1(0.08757), c1(-0.5257), t01(2.4013), s1(0.6721);
510  double A2(0.007598), c2(-0.1501), t02(6.9412), s2(0.8710);
511  return gexp(t,A1,c1,t01,s1) + gexp(t,A2,c2,t02,s2);
512 }
513 
515  // taken from fit to laser measurement taken by Iouri M. in Spring 2016.
516  double A1(5.204/6.94419), sigma1_shape(0.5387), theta1_loc(-0.3976), m1_scale(4.428);
517  double A2(1.855/6.94419), sigma2_shape(0.8132), theta2_loc(7.025), m2_scale(12.29);
518  return
519  onePulse(t,A1,sigma1_shape,theta1_loc,m1_scale) +
520  onePulse(t,A2,sigma2_shape,theta2_loc,m2_scale);
521 }
522 
523 double HcalPulseShapes::generatePhotonTime(CLHEP::HepRandomEngine* engine, unsigned int signalShape) {
524  if(signalShape==206) return generatePhotonTime206(engine);
525  else return generatePhotonTime203(engine);
526 }
527 
528 double HcalPulseShapes::generatePhotonTime203(CLHEP::HepRandomEngine* engine) {
529  double result(0.);
530  while (true) {
531  result = CLHEP::RandFlat::shoot(engine, HcalPulseShapes::Y11RANGE_);
532  if (CLHEP::RandFlat::shoot(engine, HcalPulseShapes::Y11MAX203_) < HcalPulseShapes::Y11203(result))
533  return result;
534  }
535 }
536 
537 double HcalPulseShapes::generatePhotonTime206(CLHEP::HepRandomEngine* engine) {
538  double result(0.);
539  while (true) {
540  result = CLHEP::RandFlat::shoot(engine, HcalPulseShapes::Y11RANGE_);
541  if (CLHEP::RandFlat::shoot(engine, HcalPulseShapes::Y11MAX206_) < HcalPulseShapes::Y11206(result))
542  return result;
543  }
544 }
545 
546 //Original scintillator+Y11 fit from Vasken's 2001 measurement
547 double HcalPulseShapes::Y11203(double t) {
548  return exp(-0.0635-0.1518*t + log(t)*2.528)/2485.9;
549 }
550 
551 //New scintillator+Y11 model from Vasken's 2017 measurement plus a Landau correction term
552 double HcalPulseShapes::Y11206(double t) {
553  //Shifting phase to have better comparison of digi shape with data
554  //If necessary, further digi phase adjustment can be done here:
555  //SimCalorimetry/HcalSimProducers/python/hcalSimParameters_cfi.py
556  //by changing "timePhase"
557  double shift = 7.2;
558 
559  //Fit From Deconvolved Data
560  double A,n,t0,fit;
561  A=0.104204; n=0.44064; t0=10.0186;
562  if(t>shift) fit = A*(1-exp(-(t-shift)/n))*exp(-(t-shift)/t0);
563  else fit = 0.0;
564 
565  //Correction Term
566  double norm,mpv,sigma,corTerm;
567  norm=0.0809882; mpv=0; sigma=20;
568  if(t>shift) corTerm = norm*TMath::Landau((t-shift),mpv,sigma);
569  else corTerm = 0.0;
570 
571  //Overall Y11
572  double frac = 0.11;
573  double val = (1-frac)*fit + frac*corTerm;
574 
575  if(val >= 0) return val;
576  else return 0.0;
577 }
const Shape & getShape(int shapeType) const
void setNBin(int n)
const Shape & heShape() const
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
void setShapeBin(int i, float f)
const Shape & hoShape(bool sipm=false) const
const HcalRecoParam * getHcalRecoParam(const HcalGenericDetId &fId) const
void beginRun(edm::EventSetup const &es)
static double generatePhotonTime206(CLHEP::HepRandomEngine *engine)
Geom::Theta< T > theta() const
static float Y11RANGE_
static const int nBinsSiPM_
const Shape & shapeForReco(const HcalDetId &detId) const
static double analyticPulseShapeSiPMHE(double t)
#define nullptr
double onePulse(double t, double A, double sigma, double theta, double m)
static float Y11MAX203_
const Shape & shape(const HcalDetId &detId) const
automatically figures out which shape to return
static float Y11MAX206_
uint32_t rawId() const
get the raw id
Definition: DetId.h:44
static std::vector< double > convolve(unsigned nbin, F1 f1, F2 f2)
auto const T2 &decltype(t1.eta()) t2
Definition: deltaR.h:16
T sqrt(T t)
Definition: SSEVec.h:18
double p4[4]
Definition: TauolaWrapper.h:92
static double Y11206(double t)
void computeSiPMShapeData2017()
HcalSubdetector
Definition: HcalAssistant.h:31
void computeSiPMShapeData2018()
const HcalMCParam * getHcalMCParam(const HcalGenericDetId &fId) const
double p2[4]
Definition: TauolaWrapper.h:90
const HcalDbService * theDbService
int nt
Definition: AMPTWrapper.h:32
int k[5][pyjets_maxn]
void computeSiPMShapeHE206()
static const double tmax[3]
const Shape & defaultShape(const HcalDetId &detId) const
in case of conditions problems
const Shape & hfShape() const
void computeHPDShape(float, float, float, float, float, float, float, float, Shape &)
static double generatePhotonTime(CLHEP::HepRandomEngine *engine, unsigned int signalShape)
const T & get() const
Definition: EventSetup.h:59
double gexp(double t, double A, double c, double t0, double s)
unsigned int signalShape() const
Definition: HcalMCParam.h:40
const Shape & hbShape() const
static double Y11203(double t)
double p1[4]
Definition: TauolaWrapper.h:89
static unsigned int const shift
static double generatePhotonTime203(CLHEP::HepRandomEngine *engine)
T const * product() const
Definition: ESHandle.h:86
void computeSiPMShapeHE203()
static double analyticPulseShapeSiPMHO(double t)