00001 #include "FastSimulation/CaloRecHitsProducer/interface/EcalPreshowerRecHitsMaker.h"
00002 #include "SimDataFormats/CaloHit/interface/PCaloHitContainer.h"
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 #include "FWCore/Framework/interface/ESHandle.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "FWCore/Framework/interface/Event.h"
00007 #include "FWCore/Framework/interface/EventSetup.h"
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
00010 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
00011 #include "DataFormats/EcalDetId/interface/ESDetId.h"
00012 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00013 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00014 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00015 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00016
00017 #include "CLHEP/GenericFunctions/Erf.hh"
00018
00019 EcalPreshowerRecHitsMaker::EcalPreshowerRecHitsMaker(
00020 edm::ParameterSet const & p,
00021 const RandomEngine * myrandom)
00022 :
00023 random_(myrandom),
00024 myGaussianTailGenerator_(0)
00025 {
00026 edm::ParameterSet RecHitsParameters
00027 = p.getParameter<edm::ParameterSet>("ECALPreshower");
00028 noise_ = RecHitsParameters.getParameter<double>("Noise");
00029 threshold_ = RecHitsParameters.getParameter<double>("Threshold");
00030 inputCol_=RecHitsParameters.getParameter<edm::InputTag>("MixedSimHits");
00031
00032 initialized_=false;
00033
00034 Genfun::Erf myErf;
00035 if( noise_>0. ) {
00036 preshowerHotFraction_ = 0.5-0.5*myErf(threshold_/noise_/sqrt(2.));
00037 myGaussianTailGenerator_ = new GaussianTail(random_, noise_, threshold_);
00038 } else {
00039 preshowerHotFraction_ =0.;
00040 }
00041 }
00042
00043 EcalPreshowerRecHitsMaker::~EcalPreshowerRecHitsMaker()
00044 {
00045 initialized_=false;
00046 delete myGaussianTailGenerator_;
00047 }
00048
00049 void EcalPreshowerRecHitsMaker::clean()
00050 {
00051
00052 ecalsRecHits_.clear();
00053 }
00054
00055
00056 void EcalPreshowerRecHitsMaker::loadEcalPreshowerRecHits(edm::Event &iEvent,ESRecHitCollection & ecalHits)
00057 {
00058
00059
00060 if(ncells_==0) return;
00061 loadPCaloHits(iEvent);
00062 if( myGaussianTailGenerator_ ) noisify();
00063
00064 std::map<uint32_t,std::pair<float,bool> >::const_iterator
00065 it=ecalsRecHits_.begin();
00066 std::map<uint32_t,std::pair<float,bool> >::const_iterator
00067 itend=ecalsRecHits_.end();
00068
00069 for(;it!=itend;++it)
00070 {
00071
00072 if(it->second.second) continue;
00073
00074 if(it->second.first<threshold_) continue;
00075 ESDetId detid(it->first);
00076 ecalHits.push_back(EcalRecHit(detid,it->second.first,0.));
00077 }
00078 }
00079
00080 void EcalPreshowerRecHitsMaker::loadPCaloHits(const edm::Event & iEvent)
00081 {
00082
00083 clean();
00084
00085 edm::Handle<CrossingFrame<PCaloHit> > cf;
00086 iEvent.getByLabel(inputCol_,cf);
00087 std::auto_ptr<MixCollection<PCaloHit> > colcalo(new MixCollection<PCaloHit>(cf.product(),std::pair<int,int>(0,0) ));
00088
00089 MixCollection<PCaloHit>::iterator it=colcalo->begin();
00090 MixCollection<PCaloHit>::iterator itend=colcalo->end();
00091 for(;it!=itend;++it)
00092 {
00093 Fill(it->id(),it->energy(),ecalsRecHits_,it.getTrigger());
00094
00095 }
00096 }
00097
00098 void EcalPreshowerRecHitsMaker::init(const edm::EventSetup &es)
00099 {
00100 if(initialized_) return;
00101 ncells_=createVectorsOfCells(es);
00102 edm::LogInfo("CaloRecHitsProducer") << "Total number of cells in Preshower " << ncells_ << std::endl;
00103 initialized_=true;
00104 }
00105
00106 unsigned EcalPreshowerRecHitsMaker::createVectorsOfCells(const edm::EventSetup &es)
00107 {
00108 edm::ESHandle<CaloGeometry> pG;
00109 es.get<CaloGeometryRecord>().get(pG);
00110 unsigned total=0;
00111
00112 escells_.reserve(137728);
00113
00114 const CaloSubdetectorGeometry* geom=pG->getSubdetectorGeometry(DetId::Ecal,EcalPreshower);
00115 if(geom==0) return 0;
00116 const std::vector<DetId>& ids(geom->getValidDetIds(DetId::Ecal,EcalPreshower));
00117 for (std::vector<DetId>::const_iterator i=ids.begin(); i!=ids.end(); i++)
00118 {
00119 escells_.push_back(i->rawId());
00120 ++total;
00121 }
00122 return escells_.size();
00123 }
00124
00125 void EcalPreshowerRecHitsMaker::noisify()
00126 {
00127 if(ecalsRecHits_.size()<ncells_)
00128 {
00129
00130
00131 noisifySubdet(ecalsRecHits_,escells_,ncells_);
00132 }
00133 else
00134 edm::LogWarning("CaloRecHitsProducer") << "All HCAL(HB-HE) cells on ! " << std::endl;
00135 }
00136
00137
00138 void EcalPreshowerRecHitsMaker::noisifySubdet(std::map<uint32_t,std::pair<float,bool> >& theMap, const std::vector<uint32_t>& thecells, unsigned ncells)
00139 {
00140
00141 double mean = (double)(ncells-theMap.size())*preshowerHotFraction_;
00142 unsigned nps = random_->poissonShoot(mean);
00143
00144 unsigned ncell=0;
00145 unsigned cellindex=0;
00146 uint32_t cellnumber=0;
00147 std::map<uint32_t,std::pair<float,bool> >::const_iterator itcheck;
00148
00149 while(ncell < nps)
00150 {
00151 cellindex = (unsigned)(random_->flatShoot()*ncells);
00152 cellnumber = thecells[cellindex];
00153 itcheck=theMap.find(cellnumber);
00154 if(itcheck==theMap.end())
00155 {
00156 std::pair <float,bool> noisehit(myGaussianTailGenerator_->shoot(),
00157 false);
00158 theMap.insert(std::pair<uint32_t,std::pair<float,bool> >
00159 (cellnumber,noisehit));
00160 ++ncell;
00161 }
00162 }
00163
00164 }
00165
00166
00167
00168 void EcalPreshowerRecHitsMaker::Fill(uint32_t id,float energy, std::map<uint32_t,std::pair<float,bool> >& myHits,bool signal)
00169 {
00170
00171 if(signal)
00172 {
00173
00174
00175
00176 energy += random_->gaussShoot(0.,noise_);
00177 std::pair<float,bool> hit(energy,false);
00178
00179
00180 if(signal)
00181 myHits.insert(myHits.end(),std::pair<uint32_t,std::pair<float,bool> >(id,hit));
00182 }
00183 else
00184 {
00185 std::map<uint32_t,std::pair<float,bool> >::iterator itcheck=myHits.find(id);
00186 if(itcheck==myHits.end())
00187 {
00188 energy += random_->gaussShoot(0.,noise_);
00189 std::pair<float,bool> hit(energy,false);
00190 myHits.insert(std::pair<uint32_t,std::pair<float,bool> >(id,hit));
00191 }
00192 else
00193 {
00194 itcheck->second.first += energy;
00195 }
00196 }
00197 }
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213