CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/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   threshold_ = RecHitsParameters.getParameter<std::vector<double> >("Threshold");
00061   doSaturation_ = RecHitsParameters.getParameter<bool>("EnableSaturation");
00062     
00063   refactor_ = RecHitsParameters.getParameter<double> ("Refactor");
00064   refactor_mean_ = RecHitsParameters.getParameter<double> ("Refactor_mean");
00065   hcalfileinpath_= RecHitsParameters.getParameter<std::string> ("fileNameHcal");  
00066   inputCol_=RecHitsParameters.getParameter<edm::InputTag>("MixedSimHits");
00067   nhbcells_=nhecells_=nhocells_=nhfcells_=0;
00068 
00069   if(det_==4)
00070     {
00071       hbhi_.reserve(2600);
00072       hehi_.reserve(2600);
00073     }
00074   else if (det_==5)
00075     hohi_.reserve(2200);
00076   else if (det_==6)
00077     hfhi_.reserve(1800);
00078 
00079   if(threshold_.size()!=noise_.size())
00080     {
00081       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;
00082       noise_.clear();
00083       noise_.push_back(0.);
00084     }
00085   else 
00086     nnoise_=noise_.size(); 
00087 
00088   //  edm::ParameterSet hcalparam = p2.getParameter<edm::ParameterSet>("hcalSimParam"); 
00089   //  myHcalSimParameterMap_ = new HcalSimParameterMap(hcalparam);
00090 
00091   // Computes the fraction of HCAL above the threshold
00092   Genfun::Erf myErf;
00093   hcalHotFraction_.resize(nnoise_,0.);
00094   myGaussianTailGenerators_.resize(nnoise_,0);
00095   if(noise_.size()>0) 
00096     {
00097       for(unsigned inoise=0;inoise<nnoise_;++inoise)
00098         {
00099           if(noise_[inoise]==0) 
00100             {
00101               hcalHotFraction_[inoise]=0.;
00102               continue;
00103             }
00104           else if(noise_[inoise]==-1) {
00105             noiseFromDb_=true;
00106             continue;
00107           }
00108           else
00109             {
00110               hcalHotFraction_.push_back(0.5-0.5*myErf(threshold_[inoise]/noise_[inoise]/sqrt(2.)));
00111               myGaussianTailGenerators_[inoise]=new GaussianTail(random_,noise_[inoise],threshold_[inoise]);
00112             }
00113         }   
00114     }  
00115 }
00116 
00117 HcalRecHitsMaker::~HcalRecHitsMaker()
00118 {
00119   clean();  
00120   if(myGaussianTailGenerators_.size()) 
00121     {
00122       for(unsigned igt=0; igt<myGaussianTailGenerators_.size();++igt)
00123         delete myGaussianTailGenerators_[igt];
00124     }
00125   myGaussianTailGenerators_.clear();
00126   theDetIds_.clear();
00127   hbhi_.clear();
00128   hehi_.clear();
00129   hohi_.clear();
00130   hfhi_.clear();
00131     
00132 }
00133 
00134 void HcalRecHitsMaker::init(const edm::EventSetup &es,bool doDigis,bool doMiscalib)
00135 {
00136   doDigis_=doDigis;
00137   doMiscalib_=doMiscalib;
00138 // needed for the noise simulation
00139   edm::ESHandle<HcalDbService> conditions;
00140   es.get<HcalDbRecord>().get(conditions);
00141   const HcalDbService * theDbService=conditions.product();
00142 
00143   // get the correction factors
00144   edm::ESHandle<HcalRespCorrs> rchandle;
00145   es.get<HcalRespCorrsRcd>().get(rchandle);
00146   myRespCorr= rchandle.product();
00147   
00148 
00149   if(!initialized_) 
00150     {     
00151       theDetIds_.resize(9201);
00152       unsigned ncells=createVectorsOfCells(es);
00153       edm::LogInfo("CaloRecHitsProducer") << "Total number of cells in HCAL " << ncells << std::endl;
00154       hcalRecHits_.resize(maxIndex_+1,0.);
00155       edm::LogInfo("CaloRecHitsProducer") << "Largest HCAL hashedindex" << maxIndex_ << std::endl;
00156 
00157       peds_.resize(9201);
00158       gains_.resize(9201);
00159       if(doSaturation_)
00160         sat_.resize(9201);
00161       if(noiseFromDb_)
00162         noisesigma_.resize(9201);
00163       
00164       
00165       
00166       miscalib_.resize(maxIndex_+1,1.);
00167       // Read from file ( a la HcalRecHitsRecalib.cc)
00168       // here read them from xml (particular to HCAL)
00169       CaloMiscalibMapHcal mapHcal;
00170       mapHcal.prefillMap();
00171       
00172 
00173       edm::FileInPath hcalfiletmp("CalibCalorimetry/CaloMiscalibTools/data/"+hcalfileinpath_);      
00174       std::string hcalfile=hcalfiletmp.fullPath();            
00175       MiscalibReaderFromXMLHcal hcalreader_(mapHcal);
00176       if(doMiscalib_ && !hcalfile.empty()) 
00177         {
00178           hcalreader_.parseXMLMiscalibFile(hcalfile);
00179           //      mapHcal.print();
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_[myDetId.hashed_index()]=refactor_mean_+(icalconst-refactor_mean_)*refactor_;
00187             }
00188         }
00189       
00190       
00191       // Open the histogram for the fC to ADC conversion
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           // Now fills the vector
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               //          std::cout << TPGFactor_[i] << std::endl;
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       //HB
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           //      std::cout << " NOISEHB " << theDetIds_[hbhi_[ic]].ieta() << " " << noisesigma_[hbhi_[ic]] << "  "<< std::endl;
00247             // 10000 (ADC scale) / 4. (to compute the mean) / 0.92  ADC/fC
00248   // *1.25 (only ~80% in 1ts Digi, while saturation applied to 4ts RecHit) 
00249           int ieta=theDetIds_[hbhi_[ic]].ieta();
00250           float tpgfactor=TPGFactor_[(ieta>0)?ieta+43:-ieta];
00251           mgain*=2500./0.92*tpgfactor ;// 10000 (ADC scale) / 4. (to compute the mean)
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       //HE
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           //      std::cout << " NOISEHE " << theDetIds_[hehi_[ic]].ieta() << " " << noisesigma_[hehi_[ic]] << "  "<< std::endl;
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       //HO
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           //      std::cout << " NOISEHO " << theDetIds_[hohi_[ic]].ieta() << " " << noisesigma_[hohi_[ic]] << "  "<< std::endl;
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       //HF
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           // additional 1/2 factor for the HF (Salavat)
00311           if(noiseFromDb_)
00312             {
00313               // computation from when the noise was taken 
00314               noisesigma_[hfhi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hfhi_[ic]])*mgain*0.25;
00315             }
00316           //      std::cout << " NOISEHF " << theDetIds_[hfhi_[ic]].ieta() << " " << noisesigma_[hfhi_[ic]] << "  "<< std::endl;
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   // clear the vector we don't need. It is a bit stupid 
00329   hcalRecHits_.resize(maxIndex_+1,0.);
00330 
00331 }
00332 
00333 
00334 // Get the PCaloHits from the event. They have to be stored in a map, because when
00335 // the pile-up is added thanks to the Mixing Module, the same cell can be present several times
00336 void HcalRecHitsMaker::loadPCaloHits(const edm::Event & iEvent)
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=detid.hashed_index();
00351       switch(detid.subdet())
00352         {
00353         case HcalBarrel: 
00354           {
00355             if(det_==4)
00356               Fill(hashedindex,it->energy(),firedCells_,noise_[0]);
00357           }
00358           break;
00359         case HcalEndcap: 
00360           {       
00361             if(det_==4)
00362               Fill(hashedindex,it->energy(),firedCells_,noise_[1]);
00363           }
00364           break;
00365         case HcalOuter: 
00366           {
00367             if(det_==5)
00368               Fill(hashedindex,it->energy(),firedCells_,noise_[0]);
00369           }
00370           break;                     
00371         case HcalForward: 
00372           {
00373             if(det_==6)
00374               Fill(hashedindex,it->energy(),firedCells_,noise_[0]);
00375           }
00376           break;
00377         default:
00378           edm::LogWarning("CaloRecHitsProducer") << "RecHit not registered\n";
00379           ;
00380         }
00381       ++counter;
00382     }
00383 }
00384 
00385 // Fills the collections. 
00386 void HcalRecHitsMaker::loadHcalRecHits(edm::Event &iEvent,HBHERecHitCollection& hbheHits, HBHEDigiCollection& hbheDigis)
00387 {
00388   loadPCaloHits(iEvent);
00389   noisify();
00390   hbheHits.reserve(firedCells_.size());
00391   if(doDigis_)
00392     {
00393       hbheDigis.reserve(firedCells_.size());
00394     }
00395   static HcalQIESample zeroSample(0,0,0,0);
00396   unsigned nhits=firedCells_.size();
00397   // HB and HE
00398 
00399   for(unsigned ihit=0;ihit<nhits;++ihit)
00400     {
00401       unsigned cellhashedindex=firedCells_[ihit];
00402       const HcalDetId& detid  = theDetIds_[cellhashedindex];
00403       unsigned subdet=(detid.subdet()==HcalBarrel) ? 0: 1;      
00404 
00405       // Check if it is above the threshold
00406       if(hcalRecHits_[cellhashedindex]<threshold_[subdet]) continue; 
00407       float energy=hcalRecHits_[cellhashedindex];
00408       // apply RespCorr only to the RecHit
00409       energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
00410       // poor man saturation
00411       if(energy>sat_[cellhashedindex]) 
00412         {
00413           //      std::cout << " Saturation " << energy << " " << sat_[cellhashedindex] << " HB " << std::endl;
00414           energy=sat_[cellhashedindex];
00415         }
00416       
00417 
00418       hbheHits.push_back(HBHERecHit(detid,energy,0.));      
00419       if(doDigis_)
00420         {
00421           HBHEDataFrame myDataFrame(detid);
00422           myDataFrame.setSize(2);
00423           double nfc=hcalRecHits_[cellhashedindex]/gains_[cellhashedindex]+peds_[cellhashedindex];
00424           int nadc=fCtoAdc(nfc);
00425           HcalQIESample qie(nadc, 0, 0, 0) ;
00426           myDataFrame.setSample(0,qie);
00427           myDataFrame.setSample(1,zeroSample);
00428           hbheDigis.push_back(myDataFrame);
00429         }
00430     }
00431 }
00432 
00433 
00434 // Fills the collections. 
00435 void HcalRecHitsMaker::loadHcalRecHits(edm::Event &iEvent, HORecHitCollection &hoHits, HODigiCollection & hoDigis)
00436 {
00437   loadPCaloHits(iEvent);
00438   noisify();
00439   hoHits.reserve(firedCells_.size());
00440   if(doDigis_)
00441     {
00442       hoDigis.reserve(firedCells_.size());
00443     }
00444   static HcalQIESample zeroSample(0,0,0,0);
00445 
00446   // HO
00447   unsigned nhits=firedCells_.size();
00448   for(unsigned ihit=0;ihit<nhits;++ihit)
00449     {
00450 
00451       unsigned cellhashedindex=firedCells_[ihit];
00452       // Check if it is above the threshold
00453       if(hcalRecHits_[cellhashedindex]<threshold_[0]) continue; 
00454 
00455       const HcalDetId&  detid=theDetIds_[cellhashedindex];
00456       int ieta = detid.ieta();
00457       
00458       // Force  only Ring#0 to remain
00459       if(ieta > 4 || ieta < -4 ) continue;
00460 
00461       float energy=hcalRecHits_[cellhashedindex];
00462       // apply RespCorr
00463       energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
00464 
00465       // poor man saturation
00466       if(energy>sat_[cellhashedindex]) energy=sat_[cellhashedindex];
00467 
00468       hoHits.push_back(HORecHit(detid,energy,0));
00469     }
00470 }
00471 
00472 // Fills the collections. 
00473 void HcalRecHitsMaker::loadHcalRecHits(edm::Event &iEvent,HFRecHitCollection &hfHits, HFDigiCollection& hfDigis)
00474 {
00475   loadPCaloHits(iEvent);
00476   noisify();
00477   hfHits.reserve(firedCells_.size());
00478   if(doDigis_)
00479     {
00480       hfDigis.reserve(firedCells_.size());
00481     }
00482   static HcalQIESample zeroSample(0,0,0,0);
00483 
00484   unsigned nhits=firedCells_.size();
00485   for(unsigned ihit=0;ihit<nhits;++ihit)
00486     {
00487       unsigned cellhashedindex=firedCells_[ihit];
00488       // Check if it is above the threshold
00489       if(hcalRecHits_[cellhashedindex]<threshold_[0]) continue; 
00490 
00491       float energy=hcalRecHits_[cellhashedindex];
00492 
00493       // apply RespCorr
00494       energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
00495 
00496       // poor man saturation
00497       if(energy>sat_[cellhashedindex]) energy=sat_[cellhashedindex];
00498 
00499       const HcalDetId & detid=theDetIds_[cellhashedindex];
00500       hfHits.push_back(HFRecHit(detid,energy,0.));      
00501       if(doDigis_)
00502         {
00503           HFDataFrame myDataFrame(detid);
00504           myDataFrame.setSize(1);
00505 
00506           double nfc= hcalRecHits_[cellhashedindex]/gains_[cellhashedindex]+peds_[cellhashedindex];
00507           int nadc=fCtoAdc(nfc/2.6);
00508           HcalQIESample qie(nadc, 0, 0, 0) ;
00509           myDataFrame.setSample(0,qie);
00510           hfDigis.push_back(myDataFrame);
00511         }
00512     }
00513 }
00514 
00515 
00516 // For a fast injection of the noise: the list of cell ids is stored
00517 unsigned HcalRecHitsMaker::createVectorsOfCells(const edm::EventSetup &es)
00518 {
00519     edm::ESHandle<CaloGeometry> pG;
00520     es.get<CaloGeometryRecord>().get(pG);     
00521     nhbcells_ = createVectorOfSubdetectorCells(*pG, HcalBarrel,  hbhi_);    
00522     nhecells_ = createVectorOfSubdetectorCells(*pG, HcalEndcap,  hehi_);
00523     nhocells_ = createVectorOfSubdetectorCells(*pG, HcalOuter,   hohi_);
00524     nhfcells_ = createVectorOfSubdetectorCells(*pG, HcalForward, hfhi_);    
00525 
00526     return nhbcells_+nhecells_+nhocells_+nhfcells_;
00527 }
00528 
00529 // list of the cellids for a given subdetector
00530 unsigned HcalRecHitsMaker::createVectorOfSubdetectorCells(const CaloGeometry& cg,int subdetn,std::vector<int>& cellsvec ) 
00531 {
00532 
00533   const CaloSubdetectorGeometry* geom=cg.getSubdetectorGeometry(DetId::Hcal,subdetn);  
00534   const std::vector<DetId>& ids=geom->getValidDetIds(DetId::Hcal,subdetn);  
00535 
00536   for (std::vector<DetId>::const_iterator i=ids.begin(); i!=ids.end(); i++) 
00537     {
00538       HcalDetId myDetId(*i);
00539       //      std::cout << myDetId << myHcalSimParameterMap_->simParameters(myDetId).simHitToPhotoelectrons() << std::endl;;
00540       //      std::cout << " hi " << hi << " " theDetIds_.size() << std::endl;
00541       unsigned hi=myDetId.hashed_index();
00542       theDetIds_[hi]=myDetId;
00543       //      std::cout << myDetId << " " << hi <<  std::endl;
00544       cellsvec.push_back(hi);      
00545 
00546       if(hi>maxIndex_)
00547         maxIndex_=hi;
00548     }
00549   return cellsvec.size();
00550 }
00551 
00552 // Takes a hit (from a PSimHit) and fills a map 
00553 void HcalRecHitsMaker::Fill(int id, float energy, std::vector<int>& theHits,float noise)
00554 {
00555   if(doMiscalib_) 
00556     energy*=miscalib_[id];
00557 
00558   if(noiseFromDb_)
00559     noise=noisesigma_[id];
00560 
00561   // Check if the RecHit exists
00562   if(hcalRecHits_[id]>0.)
00563     hcalRecHits_[id]+=energy;
00564   else
00565     {
00566       // the noise is injected only the first time
00567       hcalRecHits_[id]=energy + random_->gaussShoot(0.,noise);
00568       theHits.push_back(id);
00569     }
00570 }
00571 
00572 void HcalRecHitsMaker::noisify()
00573 {
00574   unsigned total=0;
00575   switch(det_)
00576     {
00577     case 4:
00578       {
00579         // do the HB
00580         if(noise_[0] != 0.) {
00581           total+=noisifySubdet(hcalRecHits_,firedCells_,hbhi_,nhbcells_,hcalHotFraction_[0],myGaussianTailGenerators_[0],noise_[0],threshold_[0]);
00582         }
00583         // do the HE
00584         if(noise_[1] != 0.) {    
00585           total+=noisifySubdet(hcalRecHits_,firedCells_,hehi_,nhecells_,hcalHotFraction_[1],myGaussianTailGenerators_[1],noise_[1],threshold_[1]);
00586         }
00587       }
00588       break;
00589     case 5:
00590       {
00591         // do the HO
00592         if(noise_[0] != 0.) {
00593           total+=noisifySubdet(hcalRecHits_,firedCells_,hohi_,nhocells_,hcalHotFraction_[0],myGaussianTailGenerators_[0],noise_[0],threshold_[0]);
00594         }
00595       }
00596       break;
00597     case 6:
00598       {
00599         // do the HF
00600         if(noise_[0] != 0.) {
00601           total+=noisifySubdet(hcalRecHits_,firedCells_,hfhi_,nhfcells_,hcalHotFraction_[0],myGaussianTailGenerators_[0],noise_[0],threshold_[0]);
00602         }
00603       }
00604       break;
00605     default:
00606       break;
00607     }
00608   edm::LogInfo("CaloRecHitsProducer") << "CaloRecHitsProducer : added noise in "<<  total << " HCAL cells "  << std::endl;
00609 }
00610 
00611 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)
00612 {
00613  // If the fraction of "hot " is small, use an optimized method to inject noise only in noisy cells. The 30% has not been tuned
00614   if(!noiseFromDb_ && hcalHotFraction==0.) return 0;
00615   if(hcalHotFraction<0.3 && !noiseFromDb_)
00616     {
00617       double mean = (double)(ncells-theHits.size())*hcalHotFraction;
00618       unsigned nhcal = random_->poissonShoot(mean);
00619   
00620       unsigned ncell=0;
00621       unsigned cellindex=0;
00622       uint32_t cellhashedindex=0;
00623       
00624       while(ncell < nhcal)
00625         {
00626           cellindex = (unsigned)(random_->flatShoot()*ncells);
00627           cellhashedindex = thecells[cellindex];
00628 
00629           if(hcalRecHits_[cellhashedindex]==0.) // new cell
00630             {
00631               hcalRecHits_[cellhashedindex]=myGT->shoot();
00632               theHits.push_back(cellhashedindex);
00633               ++ncell;
00634             }
00635         }
00636       return ncell;
00637     }
00638   else // otherwise, inject noise everywhere
00639     {
00640       uint32_t cellhashedindex=0;
00641       unsigned nhcal=thecells.size();
00642 
00643 
00644       for(unsigned ncell=0;ncell<nhcal;++ncell)
00645         {
00646           cellhashedindex = thecells[ncell];
00647           if(hcalRecHits_[cellhashedindex]==0.) // new cell
00648             {
00649               
00650               sigma=noisesigma_[cellhashedindex];
00651 
00652               double noise =random_->gaussShoot(0.,sigma);
00653               if(noise>threshold)
00654                 {
00655                   hcalRecHits_[cellhashedindex]=noise;              
00656                   theHits.push_back(cellhashedindex);
00657                 }
00658             }
00659         }
00660       return nhcal;
00661     }
00662   return 0;
00663 }
00664 
00665 void HcalRecHitsMaker::clean()
00666 {
00667   cleanSubDet(hcalRecHits_,firedCells_);
00668 }
00669 
00670 void HcalRecHitsMaker::cleanSubDet(std::vector<float>& hits,std::vector<int>& cells)
00671 {
00672   unsigned size=cells.size();
00673   // Reset the energies
00674   for(unsigned ic=0;ic<size;++ic)
00675     {
00676       hits[cells[ic]] = 0.;
00677     }
00678   // Clear the list of fired cells 
00679   cells.clear();
00680 }
00681 
00682 // fC to ADC conversion
00683 int HcalRecHitsMaker::fCtoAdc(double fc) const
00684 {
00685   if(fc<0.) return 0;
00686   if(fc>9985.) return 127;
00687   return fctoadc_[(unsigned)floor(fc)];
00688 }
00689 
00690 double HcalRecHitsMaker::noiseInfCfromDB(const HcalDbService * conditions,const HcalDetId & detId)
00691 {
00692   // method from Salavat
00693   // fetch pedestal widths (noise sigmas for all 4 CapID)
00694   const HcalPedestalWidth* pedWidth =
00695     conditions-> getPedestalWidth(detId); 
00696   double ssqq_1 = pedWidth->getSigma(0,0);
00697   double ssqq_2 = pedWidth->getSigma(1,1);
00698   double ssqq_3 = pedWidth->getSigma(2,2);
00699   double ssqq_4 = pedWidth->getSigma(3,3);
00700 
00701   // correction factors (hb,he,ho,hf)
00702   static float corrfac[4]={1.39,1.32,1.17,3.76};
00703 
00704   int sub   = detId.subdet();
00705 
00706   // HO: only Ring#0 matters 
00707   int ieta  = detId.ieta();
00708   if (sub == 3 && abs (ieta) > 4) return 0.;   
00709 
00710   // effective RecHits (for this particular detId) noise calculation :
00711   
00712   double sig_sq_mean =  0.25 * ( ssqq_1 + ssqq_2 + ssqq_3 + ssqq_4);
00713   
00714   // f - external parameter (currently set to 0.5 in the FullSim) !!!
00715   double f=0.5;
00716 
00717   double term  = sqrt (1. + sqrt(1. - 2.*f*f));
00718   double alpha = sqrt (0.5) * term;
00719   double beta  = sqrt (0.5) * f / term;
00720 
00721   //  double RMS1   = sqrt(sig_sq_mean) * sqrt (2.*beta*beta + alpha*alpha) ;
00722   double RMS4   = sqrt(sig_sq_mean) * sqrt (2.*beta*beta + 2.*(alpha-beta)*(alpha-beta) + 2.*(alpha-2.*beta)*(alpha-2.*beta)) ;
00723 
00724   double noise_rms_fC;
00725 
00726   //  if(sub == 4)  noise_rms_fC = RMS1;
00727   //  else          noise_rms_fC = RMS4;
00728   noise_rms_fC = RMS4;
00729 
00730   noise_rms_fC *= corrfac[sub-1];
00731 
00732   // to convert from above fC to GeV - multiply by gain (GeV/fC)        
00733   //  const HcalGain*  gain = conditions->getGain(detId); 
00734   //  double noise_rms_GeV = noise_rms_fC * gain->getValue(0); // Noise RMS (GeV) for detId channel
00735   return noise_rms_fC;
00736 }