CMS 3D CMS Logo

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