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