CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/FastSimulation/CaloRecHitsProducer/src/HcalRecHitsMaker.cc

Go to the documentation of this file.
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 
00027 #include "TFile.h"
00028 #include "TGraph.h"
00029 #include "TROOT.h"
00030 #include <fstream>
00031 
00032 class RandomEngine;
00033 
00034 bool HcalRecHitsMaker::initialized_ = false; 
00035 std::vector<float> HcalRecHitsMaker::peds_;
00036 std::vector<int> HcalRecHitsMaker::fctoadc_;
00037 std::vector<float> HcalRecHitsMaker::sat_;
00038 std::vector<float> HcalRecHitsMaker::gains_;
00039 std::vector<float> HcalRecHitsMaker::noisesigma_;
00040 std::vector<float> HcalRecHitsMaker::TPGFactor_;
00041 std::vector<float> HcalRecHitsMaker::miscalib_;
00042 std::vector<HcalDetId> HcalRecHitsMaker::theDetIds_;
00043 std::vector<int> HcalRecHitsMaker::hbhi_;
00044 std::vector<int> HcalRecHitsMaker::hehi_;
00045 std::vector<int> HcalRecHitsMaker::hohi_;
00046 std::vector<int> HcalRecHitsMaker::hfhi_;
00047 unsigned HcalRecHitsMaker::maxIndex_  = 0 ; 
00048 
00049 HcalRecHitsMaker::HcalRecHitsMaker(edm::ParameterSet const & p, int det,
00050                                    const RandomEngine * myrandom)
00051   :
00052   det_(det),
00053   doDigis_(false),
00054   noiseFromDb_(false),
00055   random_(myrandom)
00056   //,myHcalSimParameterMap_(0)
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   //  edm::ParameterSet hcalparam = p2.getParameter<edm::ParameterSet>("hcalSimParam"); 
00090   //  myHcalSimParameterMap_ = new HcalSimParameterMap(hcalparam);
00091 
00092   // Computes the fraction of HCAL above the threshold
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 // needed for the noise simulation
00136   edm::ESHandle<HcalDbService> conditions;
00137   es.get<HcalDbRecord>().get(conditions);
00138   const HcalDbService * theDbService=conditions.product();
00139 
00140   // get the correction factors
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       // Read from file ( a la HcalRecHitsRecalib.cc)
00164       // here read them from xml (particular to HCAL)
00165       CaloMiscalibMapHcal mapHcal;
00166       mapHcal.prefillMap();
00167       
00168 
00169       edm::FileInPath hcalfiletmp("CalibCalorimetry/CaloMiscalibTools/data/"+hcalfileinpath_);      
00170       std::string hcalfile=hcalfiletmp.fullPath();            
00171       if(doMiscalib_ && !hcalfile.empty()) 
00172         {
00173           MiscalibReaderFromXMLHcal hcalreader_(mapHcal);
00174 
00175           hcalreader_.parseXMLMiscalibFile(hcalfile);
00176           //      mapHcal.print();
00177           std::map<uint32_t,float>::const_iterator it=mapHcal.get().begin();
00178           std::map<uint32_t,float>::const_iterator itend=mapHcal.get().end();
00179           for(;it!=itend;++it)
00180             {
00181               HcalDetId myDetId(it->first);
00182               float icalconst=it->second;
00183               miscalib_[myDetId.hashed_index()]=refactor_mean_+(icalconst-refactor_mean_)*refactor_;
00184             }
00185         }
00186       
00187       
00188       // Open the histogram for the fC to ADC conversion
00189       gROOT->cd();
00190       edm::FileInPath myDataFile("FastSimulation/CaloRecHitsProducer/data/adcvsfc.root");
00191       TFile * myFile = new TFile(myDataFile.fullPath().c_str(),"READ");
00192       TGraph * myGraf = (TGraph*)myFile->Get("adcvsfc");
00193       unsigned size=myGraf->GetN();
00194       fctoadc_.resize(10000);
00195       unsigned p_index=0;
00196       fctoadc_[0]=0;
00197       int prev_nadc=0;
00198       int nadc=0;
00199       for(unsigned ibin=0;ibin<size;++ibin)
00200         {
00201           double x,y;
00202           myGraf->GetPoint(ibin,x,y);
00203           int index=(int)floor(x);
00204           if(index<0||index>=10000) continue;
00205           prev_nadc=nadc;
00206           nadc=(int)y;
00207           // Now fills the vector
00208           for(unsigned ivec=p_index;ivec<(unsigned)index;++ivec)
00209             {
00210               fctoadc_[ivec] = prev_nadc;
00211             }
00212           p_index = index;
00213         }
00214       myFile->Close();
00215       gROOT->cd();
00216       edm::FileInPath myTPGFilePath("CalibCalorimetry/HcalTPGAlgos/data/RecHit-TPG-calib.dat");
00217       TPGFactor_.resize(87,1.2);
00218       std::ifstream  myTPGFile(myTPGFilePath.fullPath().c_str(),ifstream::in);
00219       if(myTPGFile)
00220         {
00221           float gain;
00222           myTPGFile >> gain;
00223           for(unsigned i=0;i<86;++i)
00224             {
00225               myTPGFile >> TPGFactor_[i] ;
00226               //          std::cout << TPGFactor_[i] << std::endl;
00227             }
00228         }
00229       else
00230         {
00231           std::cout << " Unable to open CalibCalorimetry/HcalTPGAlgos/data/RecHit-TPG-calib.dat" << std::endl;
00232           std::cout <<  " Using a constant 1.2 factor " << std::endl;
00233         }
00234       //HB
00235       for(unsigned ic=0;ic<nhbcells_;++ic)
00236         {
00237           float gain = theDbService->getGain(theDetIds_[hbhi_[ic]])->getValue(0);
00238           float mgain=0.;
00239           for(unsigned ig=0;ig<4;++ig)
00240             mgain+=theDbService->getGain(theDetIds_[hbhi_[ic]])->getValue(ig);
00241           if(noiseFromDb_)
00242             noisesigma_[hbhi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hbhi_[ic]])*mgain*0.25;
00243           //      std::cout << " NOISEHB " << theDetIds_[hbhi_[ic]].ieta() << " " << noisesigma_[hbhi_[ic]] << "  "<< std::endl;
00244             // 10000 (ADC scale) / 4. (to compute the mean) / 0.92  ADC/fC
00245   // *1.25 (only ~80% in 1ts Digi, while saturation applied to 4ts RecHit) 
00246           int ieta=theDetIds_[hbhi_[ic]].ieta();
00247           float tpgfactor=TPGFactor_[(ieta>0)?ieta+43:-ieta];
00248           mgain*=2500./0.92*tpgfactor ;// 10000 (ADC scale) / 4. (to compute the mean)
00249           sat_[hbhi_[ic]]=(doSaturation_)?mgain:99999.;
00250       
00251           peds_[hbhi_[ic]]=theDbService->getPedestal(theDetIds_[hbhi_[ic]])->getValue(0);
00252 
00253           gain*=tpgfactor;
00254           gains_[hbhi_[ic]]=gain;
00255         }
00256       //HE
00257 
00258       for(unsigned ic=0;ic<nhecells_;++ic)
00259         {
00260           float gain= theDbService->getGain(theDetIds_[hehi_[ic]])->getValue(0);
00261           int ieta=theDetIds_[hehi_[ic]].ieta();
00262           float mgain=0.;
00263           for(unsigned ig=0;ig<4;++ig)
00264             mgain+=theDbService->getGain(theDetIds_[hehi_[ic]])->getValue(ig);
00265           if(noiseFromDb_)
00266             noisesigma_[hehi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hehi_[ic]])*mgain*0.25;
00267       
00268           //      std::cout << " NOISEHE " << theDetIds_[hehi_[ic]].ieta() << " " << noisesigma_[hehi_[ic]] << "  "<< std::endl;
00269           float tpgfactor=TPGFactor_[(ieta>0)?ieta+44:-ieta+1];
00270           mgain*=2500./0.92*tpgfactor;
00271           sat_[hehi_[ic]]=(doSaturation_)?mgain:99999.;
00272       
00273           gain*=tpgfactor;
00274           peds_[hehi_[ic]]=theDbService->getPedestal(theDetIds_[hehi_[ic]])->getValue(0);
00275           gains_[hehi_[ic]]=gain;
00276         }
00277       //HO
00278 
00279       for(unsigned ic=0;ic<nhocells_;++ic)
00280         {
00281           float ped=theDbService->getPedestal(theDetIds_[hohi_[ic]])->getValue(0);
00282           float gain=theDbService->getGain(theDetIds_[hohi_[ic]])->getValue(0);
00283           float mgain=0.;
00284           for(unsigned ig=0;ig<4;++ig)
00285             mgain+=theDbService->getGain(theDetIds_[hohi_[ic]])->getValue(ig);
00286           if(noiseFromDb_)
00287             noisesigma_[hohi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hohi_[ic]])*mgain*0.25;
00288           //      std::cout << " NOISEHO " << theDetIds_[hohi_[ic]].ieta() << " " << noisesigma_[hohi_[ic]] << "  "<< std::endl;
00289           int ieta=HcalDetId(hohi_[ic]).ieta();
00290           float tpgfactor=TPGFactor_[(ieta>0)?ieta+43:-ieta];
00291           mgain*=2500./0.92*tpgfactor;
00292           sat_[hohi_[ic]]=(doSaturation_)?mgain:99999.;
00293 
00294           gain*=tpgfactor;
00295           peds_[hohi_[ic]]=ped;
00296           gains_[hohi_[ic]]=gain;
00297         }
00298       //HF
00299 
00300       for(unsigned ic=0;ic<nhfcells_;++ic)
00301         {
00302           float ped=theDbService->getPedestal(theDetIds_[hfhi_[ic]])->getValue(0);
00303           float gain=theDbService->getGain(theDetIds_[hfhi_[ic]])->getValue(0);
00304           float mgain=0.;
00305           for(unsigned ig=0;ig<4;++ig)
00306             mgain+=theDbService->getGain(theDetIds_[hfhi_[ic]])->getValue(ig);
00307           // additional 1/2 factor for the HF (Salavat)
00308           if(noiseFromDb_)
00309             {
00310               // computation from when the noise was taken 
00311               noisesigma_[hfhi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hfhi_[ic]])*mgain*0.25;
00312             }
00313           //      std::cout << " NOISEHF " << theDetIds_[hfhi_[ic]].ieta() << " " << noisesigma_[hfhi_[ic]] << "  "<< std::endl;
00314       
00315           mgain*=2500./0.36;
00316           sat_[hfhi_[ic]]=(doSaturation_)?mgain:99999.;
00317           int ieta=theDetIds_[hfhi_[ic]].ieta();
00318           gain*=TPGFactor_[(ieta>0)?ieta+45:-ieta+2];
00319           peds_[hfhi_[ic]]=ped;
00320           gains_[hfhi_[ic]]=gain;
00321         }
00322       initialized_=true; 
00323     }
00324   
00325   // clear the vector we don't need. It is a bit stupid 
00326   hcalRecHits_.resize(maxIndex_+1,0.);
00327 
00328 }
00329 
00330 
00331 // Get the PCaloHits from the event. They have to be stored in a map, because when
00332 // the pile-up is added thanks to the Mixing Module, the same cell can be present several times
00333 void HcalRecHitsMaker::loadPCaloHits(const edm::Event & iEvent)
00334 {
00335   clean();
00336 
00337   edm::Handle<CrossingFrame<PCaloHit> > cf;
00338   iEvent.getByLabel(inputCol_,cf);
00339   std::auto_ptr<MixCollection<PCaloHit> > colcalo(new MixCollection<PCaloHit>(cf.product(),std::pair<int,int>(0,0) ));
00340 
00341   MixCollection<PCaloHit>::iterator it=colcalo->begin();;
00342   MixCollection<PCaloHit>::iterator itend=colcalo->end();
00343   unsigned counter=0;
00344   for(;it!=itend;++it)
00345     {
00346       HcalDetId detid(it->id());
00347       int hashedindex=detid.hashed_index();
00348 
00349       // apply ToF correction
00350       int time_slice=0; // temporary
00351       double fTOF=1.;
00352       if (detid.subdet()==HcalForward) fTOF = (time_slice==0) ? 1. : 0.;        
00353       else fTOF = fractionOOT(time_slice);
00354 
00355       switch(detid.subdet())
00356         {
00357         case HcalBarrel: 
00358           {
00359             if(det_==4)
00360               Fill(hashedindex,fTOF*(it->energy()),firedCells_,noise_[0],corrfac_[0]);
00361           }
00362           break;
00363         case HcalEndcap: 
00364           {       
00365             if(det_==4)
00366               Fill(hashedindex,fTOF*(it->energy()),firedCells_,noise_[1],corrfac_[1]);
00367           }
00368           break;
00369         case HcalOuter: 
00370           {
00371             if(det_==5)
00372               Fill(hashedindex,fTOF*(it->energy()),firedCells_,noise_[0],corrfac_[0]);
00373           }
00374           break;                     
00375         case HcalForward: 
00376           {
00377             if(det_==6 && time_slice==0) // skip the HF hit if out-of-time
00378               Fill(hashedindex,it->energy(),firedCells_,noise_[0],corrfac_[0]);
00379           }
00380           break;
00381         default:
00382           edm::LogWarning("CaloRecHitsProducer") << "RecHit not registered\n";
00383           ;
00384         }
00385       ++counter;
00386     }
00387 }
00388 
00389 // Fills the collections. 
00390 void HcalRecHitsMaker::loadHcalRecHits(edm::Event &iEvent,HBHERecHitCollection& hbheHits, HBHEDigiCollection& hbheDigis)
00391 {
00392   loadPCaloHits(iEvent);
00393   noisify();
00394   hbheHits.reserve(firedCells_.size());
00395   if(doDigis_)
00396     {
00397       hbheDigis.reserve(firedCells_.size());
00398     }
00399   static HcalQIESample zeroSample(0,0,0,0);
00400   unsigned nhits=firedCells_.size();
00401   // HB and HE
00402 
00403   for(unsigned ihit=0;ihit<nhits;++ihit)
00404     {
00405       unsigned cellhashedindex=firedCells_[ihit];
00406       const HcalDetId& detid  = theDetIds_[cellhashedindex];
00407       unsigned subdet=(detid.subdet()==HcalBarrel) ? 0: 1;      
00408 
00409       float energy=hcalRecHits_[cellhashedindex];
00410       // Check if it is above the threshold
00411       if(energy<threshold_[subdet]) continue; 
00412       // apply RespCorr only to the RecHit
00413       energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
00414       // poor man saturation
00415       if(energy>sat_[cellhashedindex]) 
00416         {
00417           //      std::cout << " Saturation " << energy << " " << sat_[cellhashedindex] << " HB " << std::endl;
00418           energy=sat_[cellhashedindex];
00419         }
00420       
00421 
00422       hbheHits.push_back(HBHERecHit(detid,energy,0.));      
00423       if(doDigis_)
00424         {
00425           HBHEDataFrame myDataFrame(detid);
00426           myDataFrame.setSize(2);
00427           double nfc=hcalRecHits_[cellhashedindex]/gains_[cellhashedindex]+peds_[cellhashedindex];
00428           int nadc=fCtoAdc(nfc);
00429           HcalQIESample qie(nadc, 0, 0, 0) ;
00430           myDataFrame.setSample(0,qie);
00431           myDataFrame.setSample(1,zeroSample);
00432           hbheDigis.push_back(myDataFrame);
00433         }
00434     }
00435 }
00436 
00437 
00438 // Fills the collections. 
00439 void HcalRecHitsMaker::loadHcalRecHits(edm::Event &iEvent, HORecHitCollection &hoHits, HODigiCollection & hoDigis)
00440 {
00441   loadPCaloHits(iEvent);
00442   noisify();
00443   hoHits.reserve(firedCells_.size());
00444   if(doDigis_)
00445     {
00446       hoDigis.reserve(firedCells_.size());
00447     }
00448   static HcalQIESample zeroSample(0,0,0,0);
00449 
00450   // HO
00451   unsigned nhits=firedCells_.size();
00452   for(unsigned ihit=0;ihit<nhits;++ihit)
00453     {
00454 
00455       unsigned cellhashedindex=firedCells_[ihit];
00456       // Check if it is above the threshold
00457       if(hcalRecHits_[cellhashedindex]<threshold_[0]) continue; 
00458 
00459       const HcalDetId&  detid=theDetIds_[cellhashedindex];
00460       int ieta = detid.ieta();
00461       
00462       // Force  only Ring#0 to remain
00463       if(ieta > 4 || ieta < -4 ) continue;
00464 
00465       float energy=hcalRecHits_[cellhashedindex];
00466       // apply RespCorr
00467       energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
00468 
00469       // poor man saturation
00470       if(energy>sat_[cellhashedindex]) energy=sat_[cellhashedindex];
00471 
00472       hoHits.push_back(HORecHit(detid,energy,0));
00473     }
00474 }
00475 
00476 // Fills the collections. 
00477 void HcalRecHitsMaker::loadHcalRecHits(edm::Event &iEvent,HFRecHitCollection &hfHits, HFDigiCollection& hfDigis)
00478 {
00479   loadPCaloHits(iEvent);
00480   noisify();
00481   hfHits.reserve(firedCells_.size());
00482   if(doDigis_)
00483     {
00484       hfDigis.reserve(firedCells_.size());
00485     }
00486   static HcalQIESample zeroSample(0,0,0,0);
00487 
00488   unsigned nhits=firedCells_.size();
00489   for(unsigned ihit=0;ihit<nhits;++ihit)
00490     {
00491       unsigned cellhashedindex=firedCells_[ihit];
00492       // Check if it is above the threshold
00493       if(hcalRecHits_[cellhashedindex]<threshold_[0]) continue; 
00494 
00495       float energy=hcalRecHits_[cellhashedindex];
00496 
00497       // apply RespCorr
00498       energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
00499 
00500       // poor man saturation
00501       if(energy>sat_[cellhashedindex]) energy=sat_[cellhashedindex];
00502 
00503       const HcalDetId & detid=theDetIds_[cellhashedindex];
00504       hfHits.push_back(HFRecHit(detid,energy,0.));      
00505       if(doDigis_)
00506         {
00507           HFDataFrame myDataFrame(detid);
00508           myDataFrame.setSize(1);
00509 
00510           double nfc= hcalRecHits_[cellhashedindex]/gains_[cellhashedindex]+peds_[cellhashedindex];
00511           int nadc=fCtoAdc(nfc/2.6);
00512           HcalQIESample qie(nadc, 0, 0, 0) ;
00513           myDataFrame.setSample(0,qie);
00514           hfDigis.push_back(myDataFrame);
00515         }
00516     }
00517 }
00518 
00519 
00520 // For a fast injection of the noise: the list of cell ids is stored
00521 unsigned HcalRecHitsMaker::createVectorsOfCells(const edm::EventSetup &es)
00522 {
00523     edm::ESHandle<CaloGeometry> pG;
00524     es.get<CaloGeometryRecord>().get(pG);     
00525     nhbcells_ = createVectorOfSubdetectorCells(*pG, HcalBarrel,  hbhi_);    
00526     nhecells_ = createVectorOfSubdetectorCells(*pG, HcalEndcap,  hehi_);
00527     nhocells_ = createVectorOfSubdetectorCells(*pG, HcalOuter,   hohi_);
00528     nhfcells_ = createVectorOfSubdetectorCells(*pG, HcalForward, hfhi_);    
00529 
00530     return nhbcells_+nhecells_+nhocells_+nhfcells_;
00531 }
00532 
00533 // list of the cellids for a given subdetector
00534 unsigned HcalRecHitsMaker::createVectorOfSubdetectorCells(const CaloGeometry& cg,int subdetn,std::vector<int>& cellsvec ) 
00535 {
00536 
00537   const CaloSubdetectorGeometry* geom=cg.getSubdetectorGeometry(DetId::Hcal,subdetn);  
00538   const std::vector<DetId>& ids=geom->getValidDetIds(DetId::Hcal,subdetn);  
00539 
00540   for (std::vector<DetId>::const_iterator i=ids.begin(); i!=ids.end(); i++) 
00541     {
00542       HcalDetId myDetId(*i);
00543       //      std::cout << myDetId << myHcalSimParameterMap_->simParameters(myDetId).simHitToPhotoelectrons() << std::endl;;
00544       //      std::cout << " hi " << hi << " " theDetIds_.size() << std::endl;
00545       unsigned hi=myDetId.hashed_index();
00546       theDetIds_[hi]=myDetId;
00547       //      std::cout << myDetId << " " << hi <<  std::endl;
00548       cellsvec.push_back(hi);      
00549 
00550       if(hi>maxIndex_)
00551         maxIndex_=hi;
00552     }
00553   return cellsvec.size();
00554 }
00555 
00556 // Takes a hit (from a PSimHit) and fills a map 
00557 void HcalRecHitsMaker::Fill(int id, float energy, std::vector<int>& theHits,float noise,float correctionfactor)
00558 {
00559   if(doMiscalib_) 
00560     energy*=miscalib_[id];
00561 
00562   if(noiseFromDb_)
00563     noise=noisesigma_[id]*correctionfactor;
00564 
00565   // Check if the RecHit exists
00566   if(hcalRecHits_[id]>0.)
00567     hcalRecHits_[id]+=energy;
00568   else
00569     {
00570       // the noise is injected only the first time
00571       hcalRecHits_[id]=energy + random_->gaussShoot(0.,noise);
00572       theHits.push_back(id);
00573     }
00574 }
00575 
00576 void HcalRecHitsMaker::noisify()
00577 {
00578   unsigned total=0;
00579   switch(det_)
00580     {
00581     case 4:
00582       {
00583         // do the HB
00584         if(noise_[0] != 0.) {
00585           total+=noisifySubdet(hcalRecHits_,firedCells_,hbhi_,nhbcells_,hcalHotFraction_[0],myGaussianTailGenerators_[0],noise_[0],threshold_[0],corrfac_[0]);
00586         }
00587         // do the HE
00588         if(noise_[1] != 0.) {    
00589           total+=noisifySubdet(hcalRecHits_,firedCells_,hehi_,nhecells_,hcalHotFraction_[1],myGaussianTailGenerators_[1],noise_[1],threshold_[1],corrfac_[1]);
00590         }
00591       }
00592       break;
00593     case 5:
00594       {
00595         // do the HO
00596         if(noise_[0] != 0.) {
00597           total+=noisifySubdet(hcalRecHits_,firedCells_,hohi_,nhocells_,hcalHotFraction_[0],myGaussianTailGenerators_[0],noise_[0],threshold_[0],corrfac_[0]);
00598         }
00599       }
00600       break;
00601     case 6:
00602       {
00603         // do the HF
00604         if(noise_[0] != 0.) {
00605           total+=noisifySubdet(hcalRecHits_,firedCells_,hfhi_,nhfcells_,hcalHotFraction_[0],myGaussianTailGenerators_[0],noise_[0],threshold_[0],corrfac_[0]);
00606         }
00607       }
00608       break;
00609     default:
00610       break;
00611     }
00612   edm::LogInfo("CaloRecHitsProducer") << "CaloRecHitsProducer : added noise in "<<  total << " HCAL cells "  << std::endl;
00613 }
00614 
00615 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)
00616 {
00617  // If the fraction of "hot " is small, use an optimized method to inject noise only in noisy cells. The 30% has not been tuned
00618   if(!noiseFromDb_ && hcalHotFraction==0.) return 0;
00619   if(hcalHotFraction<0.3 && !noiseFromDb_)
00620     {
00621       double mean = (double)(ncells-theHits.size())*hcalHotFraction;
00622       unsigned nhcal = random_->poissonShoot(mean);
00623   
00624       unsigned ncell=0;
00625       unsigned cellindex=0;
00626       uint32_t cellhashedindex=0;
00627       
00628       while(ncell < nhcal)
00629         {
00630           cellindex = (unsigned)(random_->flatShoot()*ncells);
00631           cellhashedindex = thecells[cellindex];
00632 
00633           if(hcalRecHits_[cellhashedindex]==0.) // new cell
00634             {
00635               hcalRecHits_[cellhashedindex]=myGT->shoot();
00636               theHits.push_back(cellhashedindex);
00637               ++ncell;
00638             }
00639         }
00640       return ncell;
00641     }
00642   else // otherwise, inject noise everywhere
00643     {
00644       uint32_t cellhashedindex=0;
00645       unsigned nhcal=thecells.size();
00646 
00647 
00648       for(unsigned ncell=0;ncell<nhcal;++ncell)
00649         {
00650           cellhashedindex = thecells[ncell];
00651           if(hcalRecHits_[cellhashedindex]==0.) // new cell
00652             {
00653               
00654               sigma=noisesigma_[cellhashedindex]*correctionfactor;
00655 
00656               double noise =random_->gaussShoot(0.,sigma);
00657               if(noise>threshold)
00658                 {
00659                   hcalRecHits_[cellhashedindex]=noise;              
00660                   theHits.push_back(cellhashedindex);
00661                 }
00662             }
00663         }
00664       return nhcal;
00665     }
00666   return 0;
00667 }
00668 
00669 void HcalRecHitsMaker::clean()
00670 {
00671   cleanSubDet(hcalRecHits_,firedCells_);
00672 }
00673 
00674 void HcalRecHitsMaker::cleanSubDet(std::vector<float>& hits,std::vector<int>& cells)
00675 {
00676   unsigned size=cells.size();
00677   // Reset the energies
00678   for(unsigned ic=0;ic<size;++ic)
00679     {
00680       hits[cells[ic]] = 0.;
00681     }
00682   // Clear the list of fired cells 
00683   cells.clear();
00684 }
00685 
00686 // fC to ADC conversion
00687 int HcalRecHitsMaker::fCtoAdc(double fc) const
00688 {
00689   if(fc<0.) return 0;
00690   if(fc>9985.) return 127;
00691   return fctoadc_[(unsigned)floor(fc)];
00692 }
00693 
00694 double HcalRecHitsMaker::noiseInfCfromDB(const HcalDbService * conditions,const HcalDetId & detId)
00695 {
00696   // method from Salavat
00697   // fetch pedestal widths (noise sigmas for all 4 CapID)
00698   const HcalPedestalWidth* pedWidth =
00699     conditions-> getPedestalWidth(detId); 
00700   double ssqq_1 = pedWidth->getSigma(0,0);
00701   double ssqq_2 = pedWidth->getSigma(1,1);
00702   double ssqq_3 = pedWidth->getSigma(2,2);
00703   double ssqq_4 = pedWidth->getSigma(3,3);
00704 
00705   // correction factors (hb,he,ho,hf)
00706   //  static float corrfac[4]={1.39,1.32,1.17,3.76};
00707   //  static float corrfac[4]={0.93,0.88,1.17,2.51}; // divided HB, HE, HF by 1.5 (to take into account the halving of the number of time slices), HO did not change
00708 
00709   int sub   = detId.subdet();
00710 
00711   // HO: only Ring#0 matters 
00712   int ieta  = detId.ieta();
00713   if (sub == 3 && abs (ieta) > 4) return 0.;   
00714 
00715   // effective RecHits (for this particular detId) noise calculation :
00716   
00717   double sig_sq_mean =  0.25 * ( ssqq_1 + ssqq_2 + ssqq_3 + ssqq_4);
00718   
00719   // f - external parameter (currently set to 0.5 in the FullSim) !!!
00720   double f=0.5;
00721 
00722   double term  = sqrt (1. + sqrt(1. - 2.*f*f));
00723   double alpha = sqrt (0.5) * term;
00724   double beta  = sqrt (0.5) * f / term;
00725 
00726   //  double RMS1   = sqrt(sig_sq_mean) * sqrt (2.*beta*beta + alpha*alpha) ;
00727   double RMS4   = sqrt(sig_sq_mean) * sqrt (2.*beta*beta + 2.*(alpha-beta)*(alpha-beta) + 2.*(alpha-2.*beta)*(alpha-2.*beta)) ;
00728 
00729   double noise_rms_fC;
00730 
00731   //  if(sub == 4)  noise_rms_fC = RMS1;
00732   //  else          noise_rms_fC = RMS4;
00733   noise_rms_fC = RMS4;
00734 
00735 
00736   // to convert from above fC to GeV - multiply by gain (GeV/fC)        
00737   //  const HcalGain*  gain = conditions->getGain(detId); 
00738   //  double noise_rms_GeV = noise_rms_fC * gain->getValue(0); // Noise RMS (GeV) for detId channel
00739   return noise_rms_fC;
00740 }
00741 
00742 // fraction of energy collected as a function of ToF (for out-of-time particles; use case is out-of-time pileup)
00743 double HcalRecHitsMaker::fractionOOT(int time_slice)// in units of 25 ns; 0 means in-time
00744 {
00745   double SF = 100./88.; // to normalize to in-time signal (88% is in the reco window)
00746   if (time_slice==-4) {
00747     return 0.02*SF;
00748   } else if (time_slice==-3) {
00749     return 0.06*SF;
00750   } else if (time_slice==-2) {
00751     return 0.19*SF;
00752   } else if (time_slice==-1) {
00753     return 0.24*SF;
00754   } else if (time_slice==0) {
00755     return 1.;
00756   } else if (time_slice==1) {
00757     return 0.70*SF;
00758   } else {
00759     return 0.;
00760   }
00761 
00762 }