00001 #include "FastSimulation/CaloRecHitsProducer/interface/HcalRecHitsMaker.h"
00002 #include "SimDataFormats/CaloHit/interface/PCaloHitContainer.h"
00003
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/Framework/interface/EventSetup.h"
00008 #include "FWCore/Framework/interface/Event.h"
00009 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
00010 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
00011 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00012 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00013 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00014 #include "DataFormats/HcalDigi/interface/HBHEDataFrame.h"
00015 #include "CLHEP/GenericFunctions/Erf.hh"
00016 #include "CalibFormats/HcalObjects/interface/HcalTPGRecord.h"
00017 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00018 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00019 #include "CondFormats/HcalObjects/interface/HcalGains.h"
00020 #include "CondFormats/HcalObjects/interface/HcalPedestal.h"
00021 #include "SimCalorimetry/HcalSimAlgos/interface/HcalSimParameterMap.h"
00022 #include "CalibCalorimetry/CaloMiscalibTools/interface/MiscalibReaderFromXMLHcal.h"
00023 #include "CalibCalorimetry/CaloMiscalibTools/interface/CaloMiscalibMapHcal.h"
00024 #include "CondFormats/HcalObjects/interface/HcalRespCorrs.h"
00025 #include "CondFormats/DataRecord/interface/HcalRespCorrsRcd.h"
00026 #include "Geometry/CaloTopology/interface/HcalTopology.h"
00027
00028 #include "TFile.h"
00029 #include "TGraph.h"
00030 #include "TROOT.h"
00031 #include <fstream>
00032 #include <iostream>
00033
00034 class RandomEngine;
00035
00036 bool HcalRecHitsMaker::initialized_ = false;
00037 bool HcalRecHitsMaker::initializedHB_ = false;
00038 bool HcalRecHitsMaker::initializedHE_ = false;
00039 bool HcalRecHitsMaker::initializedHO_ = false;
00040 bool HcalRecHitsMaker::initializedHF_ = false;
00041 std::vector<float> HcalRecHitsMaker::peds_;
00042 std::vector<int> HcalRecHitsMaker::fctoadc_;
00043 std::vector<float> HcalRecHitsMaker::sat_;
00044 std::vector<float> HcalRecHitsMaker::gains_;
00045 std::vector<float> HcalRecHitsMaker::noisesigma_;
00046 std::vector<float> HcalRecHitsMaker::TPGFactor_;
00047 std::vector<float> HcalRecHitsMaker::miscalib_;
00048 std::vector<HcalDetId> HcalRecHitsMaker::theDetIds_;
00049 std::vector<int> HcalRecHitsMaker::hbhi_;
00050 std::vector<int> HcalRecHitsMaker::hehi_;
00051 std::vector<int> HcalRecHitsMaker::hohi_;
00052 std::vector<int> HcalRecHitsMaker::hfhi_;
00053 unsigned HcalRecHitsMaker::maxIndex_ = 0 ;
00054
00055 HcalRecHitsMaker::HcalRecHitsMaker(edm::ParameterSet const & p, int det,
00056 const RandomEngine * myrandom)
00057 :
00058 det_(det),
00059 doDigis_(false),
00060 random_(myrandom)
00061 {
00062 edm::ParameterSet RecHitsParameters=p.getParameter<edm::ParameterSet>("HCAL");
00063 noise_ = RecHitsParameters.getParameter<std::vector<double> >("Noise");
00064 corrfac_ = RecHitsParameters.getParameter<std::vector<double> >("NoiseCorrectionFactor");
00065 threshold_ = RecHitsParameters.getParameter<std::vector<double> >("Threshold");
00066 doSaturation_ = RecHitsParameters.getParameter<bool>("EnableSaturation");
00067
00068 refactor_ = RecHitsParameters.getParameter<double> ("Refactor");
00069 refactor_mean_ = RecHitsParameters.getParameter<double> ("Refactor_mean");
00070 hcalfileinpath_= RecHitsParameters.getParameter<std::string> ("fileNameHcal");
00071 inputCol_=RecHitsParameters.getParameter<edm::InputTag>("MixedSimHits");
00072 nhbcells_=nhecells_=nhocells_=nhfcells_=0;
00073
00074 if(det_==4)
00075 {
00076 hbhi_.reserve(2600);
00077 hehi_.reserve(2600);
00078 }
00079 else if (det_==5)
00080 hohi_.reserve(2200);
00081 else if (det_==6)
00082 hfhi_.reserve(1800);
00083
00084 nnoise_=noise_.size();
00085 if(threshold_.size()!=noise_.size())
00086 {
00087 edm::LogWarning("CaloRecHitsProducer") << " WARNING : HCAL Noise simulation, the number of parameters should be the same for the noise and the thresholds. Disabling the noise simulation " << std::endl;
00088 noise_.resize(nnoise_,0);
00089 }
00090
00091
00092
00093
00094
00095 Genfun::Erf myErf;
00096 noiseFromDb_.resize(nnoise_,false);
00097 hcalHotFraction_.resize(nnoise_,0.);
00098 myGaussianTailGenerators_.resize(nnoise_,0);
00099 if(noise_.size()>0)
00100 {
00101 for(unsigned inoise=0;inoise<nnoise_;++inoise)
00102 {
00103 if(noise_[inoise]==0) {
00104 hcalHotFraction_[inoise]=0.;
00105 continue;
00106 } else if(noise_[inoise]==-1) {
00107 if(det_==4) noiseFromDb_[inoise]=true;
00108 else if(det_==5) noiseFromDb_[inoise]=true;
00109 else if(det_==6) noiseFromDb_[inoise]=true;
00110 continue;
00111 } else {
00112 hcalHotFraction_[inoise] = 0.5-0.5*myErf(threshold_[inoise]/noise_[inoise]/sqrt(2.));
00113 myGaussianTailGenerators_[inoise]=new GaussianTail(random_,noise_[inoise],threshold_[inoise]);
00114 }
00115 }
00116 }
00117 }
00118
00119 HcalRecHitsMaker::~HcalRecHitsMaker()
00120 {
00121 clean();
00122 if(myGaussianTailGenerators_.size())
00123 {
00124 for(unsigned igt=0; igt<myGaussianTailGenerators_.size();++igt)
00125 delete myGaussianTailGenerators_[igt];
00126 }
00127 myGaussianTailGenerators_.clear();
00128 theDetIds_.clear();
00129 hbhi_.clear();
00130 hehi_.clear();
00131 hohi_.clear();
00132 hfhi_.clear();
00133
00134 }
00135
00136 void HcalRecHitsMaker::init(const edm::EventSetup &es,bool doDigis,bool doMiscalib)
00137 {
00138 doDigis_=doDigis;
00139 doMiscalib_=doMiscalib;
00140
00141 edm::ESHandle<HcalDbService> conditions;
00142 es.get<HcalDbRecord>().get(conditions);
00143 const HcalDbService * theDbService=conditions.product();
00144
00145
00146 edm::ESHandle<HcalRespCorrs> rchandle;
00147 es.get<HcalRespCorrsRcd>().get(rchandle);
00148 myRespCorr= rchandle.product();
00149
00150
00151 if(!initialized_)
00152 {
00153 theDetIds_.resize(9201);
00154 unsigned ncells=createVectorsOfCells(es);
00155 edm::LogInfo("CaloRecHitsProducer") << "Total number of cells in HCAL " << ncells << std::endl;
00156 hcalRecHits_.resize(maxIndex_+1,0.);
00157 edm::LogInfo("CaloRecHitsProducer") << "Largest HCAL hashedindex" << maxIndex_ << std::endl;
00158
00159 peds_.resize(9201);
00160 gains_.resize(9201);
00161 sat_.resize(9201);
00162 noisesigma_.resize(9201);
00163
00164 miscalib_.resize(maxIndex_+1,1.);
00165
00166
00167 CaloMiscalibMapHcal mapHcal;
00168
00169 edm::ESHandle<HcalTopology> topo;
00170 es.get<IdealGeometryRecord>().get( topo );
00171 mapHcal.prefillMap(*topo);
00172
00173 edm::FileInPath hcalfiletmp("CalibCalorimetry/CaloMiscalibTools/data/"+hcalfileinpath_);
00174 std::string hcalfile=hcalfiletmp.fullPath();
00175 if(doMiscalib_ && !hcalfile.empty())
00176 {
00177 MiscalibReaderFromXMLHcal hcalreader_(mapHcal);
00178
00179 hcalreader_.parseXMLMiscalibFile(hcalfile);
00180
00181 std::map<uint32_t,float>::const_iterator it=mapHcal.get().begin();
00182 std::map<uint32_t,float>::const_iterator itend=mapHcal.get().end();
00183 for(;it!=itend;++it)
00184 {
00185 HcalDetId myDetId(it->first);
00186 float icalconst=it->second;
00187 miscalib_[topo->detId2denseId(myDetId)]=refactor_mean_+(icalconst-refactor_mean_)*refactor_;
00188 }
00189 }
00190
00191
00192 gROOT->cd();
00193 edm::FileInPath myDataFile("FastSimulation/CaloRecHitsProducer/data/adcvsfc.root");
00194 TFile * myFile = new TFile(myDataFile.fullPath().c_str(),"READ");
00195 TGraph * myGraf = (TGraph*)myFile->Get("adcvsfc");
00196 unsigned size=myGraf->GetN();
00197 fctoadc_.resize(10000);
00198 unsigned p_index=0;
00199 fctoadc_[0]=0;
00200 int prev_nadc=0;
00201 int nadc=0;
00202 for(unsigned ibin=0;ibin<size;++ibin)
00203 {
00204 double x,y;
00205 myGraf->GetPoint(ibin,x,y);
00206 int index=(int)floor(x);
00207 if(index<0||index>=10000) continue;
00208 prev_nadc=nadc;
00209 nadc=(int)y;
00210
00211 for(unsigned ivec=p_index;ivec<(unsigned)index;++ivec)
00212 {
00213 fctoadc_[ivec] = prev_nadc;
00214 }
00215 p_index = index;
00216 }
00217 myFile->Close();
00218 gROOT->cd();
00219 edm::FileInPath myTPGFilePath("CalibCalorimetry/HcalTPGAlgos/data/RecHit-TPG-calib.dat");
00220 TPGFactor_.resize(87,1.2);
00221 std::ifstream myTPGFile(myTPGFilePath.fullPath().c_str(),ifstream::in);
00222 if(myTPGFile)
00223 {
00224 float gain;
00225 myTPGFile >> gain;
00226 for(unsigned i=0;i<86;++i)
00227 {
00228 myTPGFile >> TPGFactor_[i] ;
00229
00230 }
00231 }
00232 else
00233 {
00234 std::cout << " Unable to open CalibCalorimetry/HcalTPGAlgos/data/RecHit-TPG-calib.dat" << std::endl;
00235 std::cout << " Using a constant 1.2 factor " << std::endl;
00236 }
00237
00238 initialized_=true;
00239 }
00240
00241 if(!initializedHB_ && det_==4){
00242 nhbcells_ = hbhi_.size();
00243
00244 for(unsigned ic=0;ic<nhbcells_;++ic)
00245 {
00246 float gain = theDbService->getGain(theDetIds_[hbhi_[ic]])->getValue(0);
00247 float mgain=0.;
00248 for(unsigned ig=0;ig<4;++ig)
00249 mgain+=theDbService->getGain(theDetIds_[hbhi_[ic]])->getValue(ig);
00250 if(noiseFromDb_[0])
00251 noisesigma_[hbhi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hbhi_[ic]])*mgain*0.25;
00252
00253
00254
00255 int ieta=theDetIds_[hbhi_[ic]].ieta();
00256 float tpgfactor=TPGFactor_[(ieta>0)?ieta+43:-ieta];
00257 mgain*=2500./0.92*tpgfactor ;
00258 sat_[hbhi_[ic]]=(doSaturation_)?mgain:99999.;
00259
00260 peds_[hbhi_[ic]]=theDbService->getPedestal(theDetIds_[hbhi_[ic]])->getValue(0);
00261
00262 gain*=tpgfactor;
00263 gains_[hbhi_[ic]]=gain;
00264 }
00265
00266 initializedHB_=true;
00267 }
00268
00269 if(!initializedHE_ && det_==4){
00270 nhecells_ = hehi_.size();
00271
00272 for(unsigned ic=0;ic<nhecells_;++ic)
00273 {
00274 float gain= theDbService->getGain(theDetIds_[hehi_[ic]])->getValue(0);
00275 int ieta=theDetIds_[hehi_[ic]].ieta();
00276 float mgain=0.;
00277 for(unsigned ig=0;ig<4;++ig)
00278 mgain+=theDbService->getGain(theDetIds_[hehi_[ic]])->getValue(ig);
00279 if(noiseFromDb_[1])
00280 noisesigma_[hehi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hehi_[ic]])*mgain*0.25;
00281
00282 float tpgfactor=TPGFactor_[(ieta>0)?ieta+44:-ieta+1];
00283 mgain*=2500./0.92*tpgfactor;
00284 sat_[hehi_[ic]]=(doSaturation_)?mgain:99999.;
00285
00286 gain*=tpgfactor;
00287 peds_[hehi_[ic]]=theDbService->getPedestal(theDetIds_[hehi_[ic]])->getValue(0);
00288 gains_[hehi_[ic]]=gain;
00289 }
00290
00291 initializedHE_=true;
00292 }
00293
00294 if(!initializedHO_ && det_==5){
00295 nhocells_ = hohi_.size();
00296
00297 for(unsigned ic=0;ic<nhocells_;++ic)
00298 {
00299 float ped=theDbService->getPedestal(theDetIds_[hohi_[ic]])->getValue(0);
00300 float gain=theDbService->getGain(theDetIds_[hohi_[ic]])->getValue(0);
00301 float mgain=0.;
00302 for(unsigned ig=0;ig<4;++ig)
00303 mgain+=theDbService->getGain(theDetIds_[hohi_[ic]])->getValue(ig);
00304 if(noiseFromDb_[0])
00305 noisesigma_[hohi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hohi_[ic]])*mgain*0.25;
00306
00307 int ieta=HcalDetId(hohi_[ic]).ieta();
00308 float tpgfactor=TPGFactor_[(ieta>0)?ieta+43:-ieta];
00309 mgain*=2500./0.92*tpgfactor;
00310 sat_[hohi_[ic]]=(doSaturation_)?mgain:99999.;
00311
00312 gain*=tpgfactor;
00313 peds_[hohi_[ic]]=ped;
00314 gains_[hohi_[ic]]=gain;
00315 }
00316
00317 initializedHO_=true;
00318 }
00319
00320 if(!initializedHF_ && det_==6){
00321 nhfcells_ = hfhi_.size();
00322
00323 for(unsigned ic=0;ic<nhfcells_;++ic)
00324 {
00325 float ped=theDbService->getPedestal(theDetIds_[hfhi_[ic]])->getValue(0);
00326 float gain=theDbService->getGain(theDetIds_[hfhi_[ic]])->getValue(0);
00327 float mgain=0.;
00328 for(unsigned ig=0;ig<4;++ig)
00329 mgain+=theDbService->getGain(theDetIds_[hfhi_[ic]])->getValue(ig);
00330
00331 if(noiseFromDb_[0])
00332 {
00333
00334 noisesigma_[hfhi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hfhi_[ic]])*mgain*0.25;
00335 }
00336
00337
00338 mgain*=2500./0.36;
00339 sat_[hfhi_[ic]]=(doSaturation_)?mgain:99999.;
00340 int ieta=theDetIds_[hfhi_[ic]].ieta();
00341 gain*=TPGFactor_[(ieta>0)?ieta+45:-ieta+2];
00342 peds_[hfhi_[ic]]=ped;
00343 gains_[hfhi_[ic]]=gain;
00344 }
00345
00346 initializedHF_=true;
00347 }
00348
00349
00350 hcalRecHits_.resize(maxIndex_+1,0.);
00351
00352 }
00353
00354
00355
00356
00357 void HcalRecHitsMaker::loadPCaloHits(const edm::Event & iEvent, const HcalTopology& topo)
00358 {
00359 clean();
00360
00361 edm::Handle<CrossingFrame<PCaloHit> > cf;
00362 iEvent.getByLabel(inputCol_,cf);
00363 std::auto_ptr<MixCollection<PCaloHit> > colcalo(new MixCollection<PCaloHit>(cf.product(),std::pair<int,int>(0,0) ));
00364
00365 MixCollection<PCaloHit>::iterator it=colcalo->begin();;
00366 MixCollection<PCaloHit>::iterator itend=colcalo->end();
00367 unsigned counter=0;
00368 for(;it!=itend;++it)
00369 {
00370 HcalDetId detid(it->id());
00371 int hashedindex=topo.detId2denseId(detid);
00372
00373
00374 int time_slice=0;
00375 double fTOF=1.;
00376 if (detid.subdet()==HcalForward) fTOF = (time_slice==0) ? 1. : 0.;
00377 else fTOF = fractionOOT(time_slice);
00378
00379 switch(detid.subdet())
00380 {
00381 case HcalBarrel:
00382 {
00383 if(det_==4)
00384 Fill(hashedindex,fTOF*(it->energy()),firedCells_,noise_[0],corrfac_[0]);
00385 }
00386 break;
00387 case HcalEndcap:
00388 {
00389 if(det_==4)
00390 Fill(hashedindex,fTOF*(it->energy()),firedCells_,noise_[1],corrfac_[1]);
00391 }
00392 break;
00393 case HcalOuter:
00394 {
00395 if(det_==5)
00396 Fill(hashedindex,fTOF*(it->energy()),firedCells_,noise_[0],corrfac_[0]);
00397 }
00398 break;
00399 case HcalForward:
00400 {
00401 if(det_==6 && time_slice==0)
00402 Fill(hashedindex,it->energy(),firedCells_,noise_[0],corrfac_[0]);
00403 }
00404 break;
00405 default:
00406 edm::LogWarning("CaloRecHitsProducer") << "RecHit not registered\n";
00407 ;
00408 }
00409 ++counter;
00410 }
00411 }
00412
00413
00414 void HcalRecHitsMaker::loadHcalRecHits(edm::Event &iEvent,const HcalTopology& topo, HBHERecHitCollection& hbheHits, HBHEDigiCollection& hbheDigis)
00415 {
00416 loadPCaloHits(iEvent,topo);
00417 noisify();
00418 hbheHits.reserve(firedCells_.size());
00419 if(doDigis_)
00420 {
00421 hbheDigis.reserve(firedCells_.size());
00422 }
00423 static HcalQIESample zeroSample(0,0,0,0);
00424 unsigned nhits=firedCells_.size();
00425
00426
00427 for(unsigned ihit=0;ihit<nhits;++ihit)
00428 {
00429 unsigned cellhashedindex=firedCells_[ihit];
00430 const HcalDetId& detid = theDetIds_[cellhashedindex];
00431 unsigned subdet=(detid.subdet()==HcalBarrel) ? 0: 1;
00432
00433 float energy=hcalRecHits_[cellhashedindex];
00434
00435 if(energy<threshold_[subdet]) continue;
00436
00437 energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
00438
00439 if(energy>sat_[cellhashedindex])
00440 {
00441
00442 energy=sat_[cellhashedindex];
00443 }
00444
00445
00446 hbheHits.push_back(HBHERecHit(detid,energy,0.));
00447 if(doDigis_)
00448 {
00449 HBHEDataFrame myDataFrame(detid);
00450 myDataFrame.setSize(2);
00451 double nfc=hcalRecHits_[cellhashedindex]/gains_[cellhashedindex]+peds_[cellhashedindex];
00452 int nadc=fCtoAdc(nfc);
00453 HcalQIESample qie(nadc, 0, 0, 0) ;
00454 myDataFrame.setSample(0,qie);
00455 myDataFrame.setSample(1,zeroSample);
00456 hbheDigis.push_back(myDataFrame);
00457 }
00458 }
00459 }
00460
00461
00462
00463 void HcalRecHitsMaker::loadHcalRecHits(edm::Event &iEvent,const HcalTopology& topo, HORecHitCollection &hoHits, HODigiCollection & hoDigis)
00464 {
00465 loadPCaloHits(iEvent,topo);
00466 noisify();
00467 hoHits.reserve(firedCells_.size());
00468 if(doDigis_)
00469 {
00470 hoDigis.reserve(firedCells_.size());
00471 }
00472 static HcalQIESample zeroSample(0,0,0,0);
00473
00474
00475 unsigned nhits=firedCells_.size();
00476 for(unsigned ihit=0;ihit<nhits;++ihit)
00477 {
00478
00479 unsigned cellhashedindex=firedCells_[ihit];
00480
00481 if(hcalRecHits_[cellhashedindex]<threshold_[0]) continue;
00482
00483 const HcalDetId& detid=theDetIds_[cellhashedindex];
00484 int ieta = detid.ieta();
00485
00486
00487 if(ieta > 4 || ieta < -4 ) continue;
00488
00489 float energy=hcalRecHits_[cellhashedindex];
00490
00491 energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
00492
00493
00494 if(energy>sat_[cellhashedindex]) energy=sat_[cellhashedindex];
00495
00496 hoHits.push_back(HORecHit(detid,energy,0));
00497 }
00498 }
00499
00500
00501 void HcalRecHitsMaker::loadHcalRecHits(edm::Event &iEvent,const HcalTopology& topo,HFRecHitCollection &hfHits, HFDigiCollection& hfDigis)
00502 {
00503 loadPCaloHits(iEvent,topo);
00504 noisify();
00505 hfHits.reserve(firedCells_.size());
00506 if(doDigis_)
00507 {
00508 hfDigis.reserve(firedCells_.size());
00509 }
00510 static HcalQIESample zeroSample(0,0,0,0);
00511
00512 unsigned nhits=firedCells_.size();
00513 for(unsigned ihit=0;ihit<nhits;++ihit)
00514 {
00515 unsigned cellhashedindex=firedCells_[ihit];
00516
00517 if(hcalRecHits_[cellhashedindex]<threshold_[0]) continue;
00518
00519 float energy=hcalRecHits_[cellhashedindex];
00520
00521
00522 energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
00523
00524
00525 if(energy>sat_[cellhashedindex]) energy=sat_[cellhashedindex];
00526
00527 const HcalDetId & detid=theDetIds_[cellhashedindex];
00528 hfHits.push_back(HFRecHit(detid,energy,0.));
00529 if(doDigis_)
00530 {
00531 HFDataFrame myDataFrame(detid);
00532 myDataFrame.setSize(1);
00533
00534 double nfc= hcalRecHits_[cellhashedindex]/gains_[cellhashedindex]+peds_[cellhashedindex];
00535 int nadc=fCtoAdc(nfc/2.6);
00536 HcalQIESample qie(nadc, 0, 0, 0) ;
00537 myDataFrame.setSample(0,qie);
00538 hfDigis.push_back(myDataFrame);
00539 }
00540 }
00541 }
00542
00543
00544
00545 unsigned HcalRecHitsMaker::createVectorsOfCells(const edm::EventSetup &es)
00546 {
00547 edm::ESHandle<CaloGeometry> pG;
00548 es.get<CaloGeometryRecord>().get(pG);
00549 edm::ESHandle<HcalTopology> topo;
00550 es.get<IdealGeometryRecord>().get( topo );
00551
00552 nhbcells_ = createVectorOfSubdetectorCells(*pG, *topo, HcalBarrel, hbhi_);
00553 nhecells_ = createVectorOfSubdetectorCells(*pG, *topo, HcalEndcap, hehi_);
00554 nhocells_ = createVectorOfSubdetectorCells(*pG, *topo, HcalOuter, hohi_);
00555 nhfcells_ = createVectorOfSubdetectorCells(*pG, *topo, HcalForward, hfhi_);
00556
00557 return nhbcells_+nhecells_+nhocells_+nhfcells_;
00558 }
00559
00560
00561 unsigned HcalRecHitsMaker::createVectorOfSubdetectorCells(const CaloGeometry& cg, const HcalTopology& topo, int subdetn,std::vector<int>& cellsvec )
00562 {
00563
00564 const CaloSubdetectorGeometry* geom=cg.getSubdetectorGeometry(DetId::Hcal,subdetn);
00565 const std::vector<DetId>& ids=geom->getValidDetIds(DetId::Hcal,subdetn);
00566
00567 for (std::vector<DetId>::const_iterator i=ids.begin(); i!=ids.end(); i++)
00568 {
00569 HcalDetId myDetId(*i);
00570
00571
00572 unsigned int hi=topo.detId2denseId(myDetId);
00573 theDetIds_[hi]=myDetId;
00574
00575 cellsvec.push_back(hi);
00576
00577 if(hi>maxIndex_)
00578 maxIndex_=hi;
00579 }
00580 return cellsvec.size();
00581 }
00582
00583
00584 void HcalRecHitsMaker::Fill(int id, float energy, std::vector<int>& theHits,float noise,float correctionfactor)
00585 {
00586 if(doMiscalib_)
00587 energy*=miscalib_[id];
00588
00589 if(noise==-1)
00590 noise=noisesigma_[id]*correctionfactor;
00591
00592
00593 if(hcalRecHits_[id]>0.)
00594 hcalRecHits_[id]+=energy;
00595 else
00596 {
00597
00598 hcalRecHits_[id]=energy + random_->gaussShoot(0.,noise);
00599 theHits.push_back(id);
00600 }
00601 }
00602
00603 void HcalRecHitsMaker::noisify()
00604 {
00605 unsigned total=0;
00606 switch(det_)
00607 {
00608 case 4:
00609 {
00610
00611 if(noise_[0] != 0.) {
00612 total+=noisifySubdet(hcalRecHits_,firedCells_,hbhi_,nhbcells_,hcalHotFraction_[0],myGaussianTailGenerators_[0],noise_[0],threshold_[0],corrfac_[0]);
00613 }
00614
00615 if(noise_[1] != 0.) {
00616 total+=noisifySubdet(hcalRecHits_,firedCells_,hehi_,nhecells_,hcalHotFraction_[1],myGaussianTailGenerators_[1],noise_[1],threshold_[1],corrfac_[1]);
00617 }
00618 }
00619 break;
00620 case 5:
00621 {
00622
00623 if(noise_[0] != 0.) {
00624 total+=noisifySubdet(hcalRecHits_,firedCells_,hohi_,nhocells_,hcalHotFraction_[0],myGaussianTailGenerators_[0],noise_[0],threshold_[0],corrfac_[0]);
00625 }
00626 }
00627 break;
00628 case 6:
00629 {
00630
00631 if(noise_[0] != 0.) {
00632 total+=noisifySubdet(hcalRecHits_,firedCells_,hfhi_,nhfcells_,hcalHotFraction_[0],myGaussianTailGenerators_[0],noise_[0],threshold_[0],corrfac_[0]);
00633 }
00634 }
00635 break;
00636 default:
00637 break;
00638 }
00639 edm::LogInfo("CaloRecHitsProducer") << "CaloRecHitsProducer : added noise in "<< total << " HCAL cells " << std::endl;
00640 }
00641
00642 unsigned HcalRecHitsMaker::noisifySubdet(std::vector<float>& theMap, std::vector<int>& theHits, const std::vector<int>& thecells, unsigned ncells, double hcalHotFraction,const GaussianTail *myGT,double noise_,double threshold,double correctionfactor)
00643 {
00644
00645 if(noise_!=-1 && hcalHotFraction==0.) return 0;
00646 if(hcalHotFraction<0.3 && noise_!=-1)
00647 {
00648 double mean = (double)(ncells-theHits.size())*hcalHotFraction;
00649 unsigned nhcal = random_->poissonShoot(mean);
00650
00651 unsigned ncell=0;
00652 unsigned cellindex=0;
00653 uint32_t cellhashedindex=0;
00654
00655 while(ncell < nhcal)
00656 {
00657 cellindex = (unsigned)(random_->flatShoot()*ncells);
00658 cellhashedindex = thecells[cellindex];
00659
00660 if(hcalRecHits_[cellhashedindex]==0.)
00661 {
00662 hcalRecHits_[cellhashedindex]=myGT->shoot();
00663 theHits.push_back(cellhashedindex);
00664 ++ncell;
00665 }
00666 }
00667 return ncell;
00668 }
00669 else
00670 {
00671 uint32_t cellhashedindex=0;
00672 unsigned nhcal=thecells.size();
00673
00674 double sigma = noise_;
00675 for(unsigned ncell=0;ncell<nhcal;++ncell)
00676 {
00677 cellhashedindex = thecells[ncell];
00678 if(hcalRecHits_[cellhashedindex]==0.)
00679 {
00680
00681
00682 if(noise_==-1) sigma=noisesigma_[cellhashedindex]*correctionfactor;
00683
00684 double noise =random_->gaussShoot(0.,sigma);
00685 if(noise>threshold)
00686 {
00687 hcalRecHits_[cellhashedindex]=noise;
00688 theHits.push_back(cellhashedindex);
00689 }
00690 }
00691 }
00692 return nhcal;
00693 }
00694 return 0;
00695 }
00696
00697 void HcalRecHitsMaker::clean()
00698 {
00699 cleanSubDet(hcalRecHits_,firedCells_);
00700 }
00701
00702 void HcalRecHitsMaker::cleanSubDet(std::vector<float>& hits,std::vector<int>& cells)
00703 {
00704 unsigned size=cells.size();
00705
00706 for(unsigned ic=0;ic<size;++ic)
00707 {
00708 hits[cells[ic]] = 0.;
00709 }
00710
00711 cells.clear();
00712 }
00713
00714
00715 int HcalRecHitsMaker::fCtoAdc(double fc) const
00716 {
00717 if(fc<0.) return 0;
00718 if(fc>9985.) return 127;
00719 return fctoadc_[(unsigned)floor(fc)];
00720 }
00721
00722 double HcalRecHitsMaker::noiseInfCfromDB(const HcalDbService * conditions,const HcalDetId & detId)
00723 {
00724
00725
00726 const HcalPedestalWidth* pedWidth =
00727 conditions-> getPedestalWidth(detId);
00728 double ssqq_1 = pedWidth->getSigma(0,0);
00729 double ssqq_2 = pedWidth->getSigma(1,1);
00730 double ssqq_3 = pedWidth->getSigma(2,2);
00731 double ssqq_4 = pedWidth->getSigma(3,3);
00732
00733
00734
00735
00736
00737 int sub = detId.subdet();
00738
00739
00740 int ieta = detId.ieta();
00741 if (sub == 3 && abs (ieta) > 4) return 0.;
00742
00743
00744
00745 double sig_sq_mean = 0.25 * ( ssqq_1 + ssqq_2 + ssqq_3 + ssqq_4);
00746
00747
00748 double f=0.5;
00749
00750 double term = sqrt (1. + sqrt(1. - 2.*f*f));
00751 double alpha = sqrt (0.5) * term;
00752 double beta = sqrt (0.5) * f / term;
00753
00754
00755 double RMS4 = sqrt(sig_sq_mean) * sqrt (2.*beta*beta + 2.*(alpha-beta)*(alpha-beta) + 2.*(alpha-2.*beta)*(alpha-2.*beta)) ;
00756
00757 double noise_rms_fC;
00758
00759
00760
00761 noise_rms_fC = RMS4;
00762
00763
00764
00765
00766
00767 return noise_rms_fC;
00768 }
00769
00770
00771 double HcalRecHitsMaker::fractionOOT(int time_slice)
00772 {
00773 double SF = 100./88.;
00774 if (time_slice==-4) {
00775 return 0.02*SF;
00776 } else if (time_slice==-3) {
00777 return 0.06*SF;
00778 } else if (time_slice==-2) {
00779 return 0.19*SF;
00780 } else if (time_slice==-1) {
00781 return 0.24*SF;
00782 } else if (time_slice==0) {
00783 return 1.;
00784 } else if (time_slice==1) {
00785 return 0.70*SF;
00786 } else {
00787 return 0.;
00788 }
00789
00790 }