00001 #include "CalibCalorimetry/HcalAlgos/interface/HcalPulseShapes.h"
00002 #include "CondFormats/HcalObjects/interface/HcalMCParam.h"
00003 #include "CondFormats/HcalObjects/interface/HcalMCParams.h"
00004 #include "CondFormats/DataRecord/interface/HcalMCParamsRcd.h"
00005 #include "CondFormats/HcalObjects/interface/HcalRecoParam.h"
00006 #include "CondFormats/HcalObjects/interface/HcalRecoParams.h"
00007 #include "CondFormats/DataRecord/interface/HcalRecoParamsRcd.h"
00008 #include "FWCore/Utilities/interface/Exception.h"
00009 #include "FWCore/Framework/interface/ESHandle.h"
00010 #include "FWCore/Framework/interface/EventSetup.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 #include "Geometry/CaloTopology/interface/HcalTopology.h"
00013 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00014
00015
00016 #include <cmath>
00017
00018 #include <iostream>
00019 #include <fstream>
00020
00021 HcalPulseShapes::HcalPulseShapes()
00022 : theMCParams(0),
00023 theTopology(0),
00024 theRecoParams(0),
00025 theShapes()
00026 {
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046 float ts1, ts2, ts3, thpd, tpre, wd1, wd2, wd3;
00047
00048
00049 ts1=8. ; ts2=10. ; ts3=29.3; thpd=4.0; tpre=9.0; wd1=2.0; wd2=0.7; wd3=1.0;
00050 computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdShape_);
00051 theShapes[101] = &hpdShape_;
00052 theShapes[102] = theShapes[101];
00053
00054
00055 ts1=8. ; ts2=10. ; ts3=25.0; thpd=4.0; tpre=9.0; wd1=2.0; wd2=0.7; wd3=1.0;
00056 computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdShape_v2);
00057 theShapes[103] = &hpdShape_v2;
00058
00059 ts1=8. ; ts2=10. ; ts3=29.3; thpd=4.0; tpre=7.0; wd1=2.0; wd2=0.7; wd3=1.0;
00060 computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdShapeMC_v2);
00061 theShapes[123] = &hpdShapeMC_v2;
00062
00063
00064 ts1=8. ; ts2=19. ; ts3=29.3; thpd=4.0; tpre=9.0; wd1=2.0; wd2=0.7; wd3=0.32;
00065 computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdShape_v3);
00066 theShapes[105] = &hpdShape_v3;
00067
00068 ts1=8. ; ts2=10. ; ts3=22.3; thpd=4.0; tpre=7.0; wd1=2.0; wd2=0.7; wd3=1.0;
00069 computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdShapeMC_v3);
00070 theShapes[125] = &hpdShapeMC_v3;
00071
00072
00073
00074 ts1=8. ; ts2=12. ; ts3=31.7; thpd=9.0; tpre=9.0; wd1=2.0; wd2=0.7; wd3=1.0;
00075 computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdBV30Shape_v2);
00076 theShapes[104] = &hpdBV30Shape_v2;
00077
00078 ts1=8. ; ts2=12. ; ts3=31.7; thpd=9.0; tpre=9.0; wd1=2.0; wd2=0.7; wd3=1.0;
00079 computeHPDShape(ts1,ts2,ts3,thpd,tpre,wd1,wd2,wd3, hpdBV30ShapeMC_v2);
00080 theShapes[124] = &hpdBV30ShapeMC_v2;
00081
00082
00083
00084 computeHFShape();
00085 computeSiPMShape();
00086
00087 theShapes[201] = &siPMShape_;
00088 theShapes[202] = theShapes[201];
00089 theShapes[301] = &hfShape_;
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101 }
00102
00103
00104 HcalPulseShapes::~HcalPulseShapes() {
00105 if (theMCParams) delete theMCParams;
00106 if (theRecoParams) delete theRecoParams;
00107 if (theTopology) delete theTopology;
00108 }
00109
00110
00111 void HcalPulseShapes::beginRun(edm::EventSetup const & es)
00112 {
00113 edm::ESHandle<HcalMCParams> p;
00114 es.get<HcalMCParamsRcd>().get(p);
00115 theMCParams = new HcalMCParams(*p.product());
00116
00117 edm::ESHandle<HcalTopology> htopo;
00118 es.get<IdealGeometryRecord>().get(htopo);
00119 theTopology=new HcalTopology(*htopo);
00120 theMCParams->setTopo(theTopology);
00121
00122 edm::ESHandle<HcalRecoParams> q;
00123 es.get<HcalRecoParamsRcd>().get(q);
00124 theRecoParams = new HcalRecoParams(*q.product());
00125 theRecoParams->setTopo(theTopology);
00126
00127
00128
00129
00130 }
00131
00132
00133 void HcalPulseShapes::endRun()
00134 {
00135 if (theMCParams) delete theMCParams;
00136 if (theRecoParams) delete theRecoParams;
00137 if (theTopology) delete theTopology;
00138
00139
00140 theMCParams = 0;
00141 theRecoParams = 0;
00142 theTopology = 0;
00143 }
00144
00145
00146
00147 void HcalPulseShapes::computeHPDShape(float ts1, float ts2, float ts3, float thpd, float tpre,
00148 float wd1, float wd2, float wd3, Shape &tmphpdShape_)
00149 {
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172 unsigned int nbin = 256;
00173 tmphpdShape_.setNBin(nbin);
00174 std::vector<float> ntmp(nbin,0.0);
00175 std::vector<float> nth(nbin,0.0);
00176 std::vector<float> ntp(nbin,0.0);
00177 std::vector<float> ntd(nbin,0.0);
00178
00179 unsigned int i,j,k;
00180 float norm;
00181
00182
00183 norm=0.0;
00184 for(j=0;j<thpd && j<nbin;j++){
00185 nth[j] = 1.0 + ((float)j)/thpd;
00186 norm += nth[j];
00187 }
00188
00189 for(j=0;j<thpd && j<nbin;j++){
00190 nth[j] /= norm;
00191 }
00192
00193
00194 norm=0.0;
00195 for(j=0;j<6*tpre && j<nbin;j++){
00196 ntp[j] = ((float)j)*exp(-((float)(j*j))/(tpre*tpre));
00197 norm += ntp[j];
00198 }
00199
00200 for(j=0;j<6*tpre && j<nbin;j++){
00201 ntp[j] /= norm;
00202 }
00203
00204
00205
00206
00207
00208 unsigned int tmax = 6 * (int)ts3;
00209
00210 norm=0.0;
00211 for(j=0;j<tmax && j<nbin;j++){
00212 ntd[j] = wd1 * exp(-((float)j)/ts1) +
00213 wd2 * exp(-((float)j)/ts2) +
00214 wd3 * exp(-((float)j)/ts3) ;
00215 norm += ntd[j];
00216 }
00217
00218 for(j=0;j<tmax && j<nbin;j++){
00219 ntd[j] /= norm;
00220 }
00221
00222 unsigned int t1,t2,t3,t4;
00223 for(i=0;i<tmax && i<nbin;i++){
00224 t1 = i;
00225
00226
00227 t2 = t1;
00228 for(j=0;j<thpd && j<nbin;j++){
00229 t3 = t2 + j;
00230 for(k=0;k<4*tpre && k<nbin;k++){
00231 t4 = t3 + k;
00232 if(t4<nbin){
00233 unsigned int ntb=t4;
00234 ntmp[ntb] += ntd[i]*nth[j]*ntp[k];
00235 }
00236 }
00237 }
00238 }
00239
00240
00241 norm = 0.;
00242 for(i=0;i<nbin;i++){
00243 norm += ntmp[i];
00244 }
00245
00246
00247 for(i=0; i<nbin; i++){
00248 ntmp[i] /= norm;
00249
00250 }
00251
00252 for(i=0; i<nbin; i++){
00253 tmphpdShape_.setShapeBin(i,ntmp[i]);
00254 }
00255 }
00256
00257 void HcalPulseShapes::computeHFShape() {
00258
00259 unsigned int nbin = 256;
00260 hfShape_.setNBin(nbin);
00261 std::vector<float> ntmp(nbin,0.0);
00262
00263 const float k0=0.7956;
00264 const float p2=1.355;
00265 const float p4=2.327;
00266 const float p1=4.3;
00267
00268 float norm = 0.0;
00269
00270 for(unsigned int j = 0; j < 25 && j < nbin; ++j){
00271
00272 float r0 = j-p1;
00273 float sigma0 = (r0<0) ? p2 : p2*p4;
00274 r0 /= sigma0;
00275 if(r0 < k0) ntmp[j] = exp(-0.5*r0*r0);
00276 else ntmp[j] = exp(0.5*k0*k0-k0*r0);
00277 norm += ntmp[j];
00278 }
00279
00280 for(unsigned int j = 0; j < 25 && j < nbin; ++j){
00281 ntmp[j] /= norm;
00282 hfShape_.setShapeBin(j,ntmp[j]);
00283 }
00284 }
00285
00286
00287 void HcalPulseShapes::computeSiPMShape()
00288 {
00289 unsigned int nbin = 512;
00290 siPMShape_.setNBin(nbin);
00291 std::vector<float> nt(nbin,0.0);
00292
00293 double norm = 0.;
00294 for (unsigned int j = 0; j < nbin; ++j) {
00295 if (j <= 31.)
00296 nt[j] = 2.15*j;
00297 else if ((j > 31) && (j <= 96))
00298 nt[j] = 102.1 - 1.12*j;
00299 else
00300 nt[j] = 0.0076*j - 6.4;
00301 norm += (nt[j]>0) ? nt[j] : 0.;
00302 }
00303
00304 for (unsigned int j = 0; j < nbin; ++j) {
00305 nt[j] /= norm;
00306 siPMShape_.setShapeBin(j,nt[j]);
00307 }
00308 }
00309
00310
00311 const HcalPulseShapes::Shape &
00312 HcalPulseShapes::getShape(int shapeType) const
00313 {
00314
00315
00316
00317
00318 ShapeMap::const_iterator shapeMapItr = theShapes.find(shapeType);
00319 if(shapeMapItr == theShapes.end()) {
00320 throw cms::Exception("HcalPulseShapes") << "unknown shapeType";
00321 return hpdShape_;
00322 } else {
00323 return *(shapeMapItr->second);
00324 }
00325 }
00326
00327
00328 const HcalPulseShapes::Shape &
00329 HcalPulseShapes::shape(const HcalDetId & detId) const
00330 {
00331 if(!theMCParams) {
00332 return defaultShape(detId);
00333 }
00334 int shapeType = theMCParams->getValues(detId)->signalShape();
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346
00347 ShapeMap::const_iterator shapeMapItr = theShapes.find(shapeType);
00348 if(shapeMapItr == theShapes.end()) {
00349 return defaultShape(detId);
00350 } else {
00351 return *(shapeMapItr->second);
00352 }
00353 }
00354
00355 const HcalPulseShapes::Shape &
00356 HcalPulseShapes::shapeForReco(const HcalDetId & detId) const
00357 {
00358 if(!theRecoParams) {
00359 return defaultShape(detId);
00360 }
00361 int shapeType = theRecoParams->getValues(detId.rawId())->pulseShapeID();
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375 ShapeMap::const_iterator shapeMapItr = theShapes.find(shapeType);
00376 if(shapeMapItr == theShapes.end()) {
00377 return defaultShape(detId);
00378 } else {
00379 return *(shapeMapItr->second);
00380 }
00381 }
00382
00383
00384 const HcalPulseShapes::Shape &
00385 HcalPulseShapes::defaultShape(const HcalDetId & detId) const
00386 {
00387 edm::LogWarning("HcalPulseShapes") << "Cannot find HCAL MC Params ";
00388 HcalSubdetector subdet = detId.subdet();
00389 switch(subdet) {
00390 case HcalBarrel:
00391 return hbShape();
00392 case HcalEndcap:
00393 return heShape();
00394 case HcalForward:
00395 return hfShape();
00396 case HcalOuter:
00397
00398 return hoShape(false);
00399 default:
00400 throw cms::Exception("HcalPulseShapes") << "unknown detId";
00401 break;
00402 }
00403 }
00404