00001
00002 #include "RecoJets/JetProducers/interface/PileUpSubtractor.h"
00003
00004 #include "DataFormats/Common/interface/View.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008
00009 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00010 #include "DataFormats/Candidate/interface/Candidate.h"
00011 #include "DataFormats/Math/interface/deltaR.h"
00012 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00013 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00014
00015 #include <map>
00016 using namespace std;
00017
00018 PileUpSubtractor::PileUpSubtractor(const edm::ParameterSet& iConfig) :
00019 reRunAlgo_ (iConfig.getUntrackedParameter<bool>("reRunAlgo",false)),
00020 doAreaFastjet_ (iConfig.getParameter<bool> ("doAreaFastjet")),
00021 doRhoFastjet_ (iConfig.getParameter<bool> ("doRhoFastjet")),
00022 jetPtMin_(iConfig.getParameter<double> ("jetPtMin")),
00023 nSigmaPU_(iConfig.getParameter<double>("nSigmaPU")),
00024 radiusPU_(iConfig.getParameter<double>("radiusPU")),
00025 geo_(0)
00026 {
00027 if (iConfig.exists("puPtMin"))
00028 puPtMin_=iConfig.getParameter<double> ("puPtMin");
00029 else{
00030 puPtMin_=10;
00031 edm::LogWarning("MisConfiguration")<<"the parameter puPtMin is now necessary for PU substraction. setting it to "<<puPtMin_;
00032 }
00033 if ( doAreaFastjet_ || doRhoFastjet_ ) {
00034
00035 double ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
00036
00037 int activeAreaRepeats = iConfig.getParameter<int> ("Active_Area_Repeats");
00038
00039 double ghostArea = iConfig.getParameter<double> ("GhostArea");
00040 fjActiveArea_ = ActiveAreaSpecPtr(new fastjet::ActiveAreaSpec(ghostEtaMax,
00041 activeAreaRepeats,
00042 ghostArea));
00043 fjRangeDef_ = RangeDefPtr( new fastjet::RangeDefinition(ghostEtaMax) );
00044 }
00045 }
00046
00047 void PileUpSubtractor::reset(std::vector<edm::Ptr<reco::Candidate> >& input,
00048 std::vector<fastjet::PseudoJet>& towers,
00049 std::vector<fastjet::PseudoJet>& output){
00050
00051 inputs_ = &input;
00052 fjInputs_ = &towers;
00053 fjJets_ = &output;
00054 fjOriginalInputs_ = (*fjInputs_);
00055 for (unsigned int i = 0; i < fjInputs_->size(); ++i){
00056 fjOriginalInputs_[i].set_user_index((*fjInputs_)[i].user_index());
00057 }
00058
00059 }
00060
00061 void PileUpSubtractor::setAlgorithm(ClusterSequencePtr& algorithm){
00062 fjClusterSeq_ = algorithm;
00063 }
00064
00065 void PileUpSubtractor::setupGeometryMap(edm::Event& iEvent,const edm::EventSetup& iSetup)
00066 {
00067
00068 LogDebug("PileUpSubtractor")<<"The subtractor setting up geometry...\n";
00069
00070 if(geo_ == 0) {
00071 edm::ESHandle<CaloGeometry> pG;
00072 iSetup.get<CaloGeometryRecord>().get(pG);
00073 geo_ = pG.product();
00074 std::vector<DetId> alldid = geo_->getValidDetIds();
00075
00076 int ietaold = -10000;
00077 ietamax_ = -10000;
00078 ietamin_ = 10000;
00079 for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
00080 if( (*did).det() == DetId::Hcal ){
00081 HcalDetId hid = HcalDetId(*did);
00082 if( (hid).depth() == 1 ) {
00083 allgeomid_.push_back(*did);
00084
00085 if((hid).ieta() != ietaold){
00086 ietaold = (hid).ieta();
00087 geomtowers_[(hid).ieta()] = 1;
00088 if((hid).ieta() > ietamax_) ietamax_ = (hid).ieta();
00089 if((hid).ieta() < ietamin_) ietamin_ = (hid).ieta();
00090 }
00091 else{
00092 geomtowers_[(hid).ieta()]++;
00093 }
00094 }
00095 }
00096 }
00097 }
00098
00099 for (int i = ietamin_; i<ietamax_+1; i++) {
00100 ntowersWithJets_[i] = 0;
00101 }
00102 }
00103
00104
00105
00106
00107 void PileUpSubtractor::calculatePedestal( vector<fastjet::PseudoJet> const & coll )
00108 {
00109 LogDebug("PileUpSubtractor")<<"The subtractor calculating pedestals...\n";
00110
00111 map<int,double> emean2;
00112 map<int,int> ntowers;
00113
00114 int ietaold = -10000;
00115 int ieta0 = -100;
00116
00117
00118
00119 for(int i = ietamin_; i < ietamax_+1; i++)
00120 {
00121 emean_[i] = 0.;
00122 emean2[i] = 0.;
00123 esigma_[i] = 0.;
00124 ntowers[i] = 0;
00125 }
00126
00127 for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin (),
00128 fjInputsEnd = coll.end();
00129 input_object != fjInputsEnd; ++input_object) {
00130
00131 const reco::CandidatePtr & originalTower=(*inputs_)[ input_object->user_index()];
00132 ieta0 = ieta( originalTower );
00133 double Original_Et = originalTower->et();
00134
00135 if( ieta0-ietaold != 0 )
00136 {
00137 emean_[ieta0] = emean_[ieta0]+Original_Et;
00138 emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
00139 ntowers[ieta0] = 1;
00140 ietaold = ieta0;
00141 }
00142 else
00143 {
00144 emean_[ieta0] = emean_[ieta0]+Original_Et;
00145 emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
00146 ntowers[ieta0]++;
00147 }
00148
00149 }
00150
00151 for(map<int,int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++)
00152 {
00153
00154 int it = (*gt).first;
00155
00156 double e1 = (*(emean_.find(it))).second;
00157 double e2 = (*emean2.find(it)).second;
00158 int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
00159
00160 LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" number of towers : "<<nt<<" e1 : "<<e1<<" e2 : "<<e2<<"\n";
00161
00162 if(nt > 0) {
00163 emean_[it] = e1/nt;
00164 double eee = e2/nt - e1*e1/(nt*nt);
00165 if(eee<0.) eee = 0.;
00166 esigma_[it] = nSigmaPU_*sqrt(eee);
00167 }
00168 else
00169 {
00170 emean_[it] = 0.;
00171 esigma_[it] = 0.;
00172 }
00173 LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" Pedestals : "<<emean_[it]<<" "<<esigma_[it]<<"\n";
00174 }
00175 }
00176
00177
00178
00179
00180
00181 void PileUpSubtractor::subtractPedestal(vector<fastjet::PseudoJet> & coll)
00182 {
00183
00184 LogDebug("PileUpSubtractor")<<"The subtractor subtracting pedestals...\n";
00185
00186 int it = -100;
00187 int ip = -100;
00188
00189 for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
00190 fjInputsEnd = coll.end();
00191 input_object != fjInputsEnd; ++input_object) {
00192
00193 reco::CandidatePtr const & itow = (*inputs_)[ input_object->user_index() ];
00194
00195 it = ieta( itow );
00196 ip = iphi( itow );
00197
00198 double etnew = itow->et() - (*(emean_.find(it))).second - (*(esigma_.find(it))).second;
00199 float mScale = etnew/input_object->Et();
00200 if(etnew < 0.) mScale = 0.;
00201
00202 math::XYZTLorentzVectorD towP4(input_object->px()*mScale, input_object->py()*mScale,
00203 input_object->pz()*mScale, input_object->e()*mScale);
00204
00205 int index = input_object->user_index();
00206 input_object->reset ( towP4.px(),
00207 towP4.py(),
00208 towP4.pz(),
00209 towP4.energy() );
00210 input_object->set_user_index(index);
00211 }
00212 }
00213
00214 void PileUpSubtractor::calculateOrphanInput(vector<fastjet::PseudoJet> & orphanInput)
00215 {
00216
00217 LogDebug("PileUpSubtractor")<<"The subtractor calculating orphan input...\n";
00218
00219 (*fjInputs_) = fjOriginalInputs_;
00220
00221 vector<int> jettowers;
00222 vector<pair<int,int> > excludedTowers;
00223
00224 vector <fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
00225 fjJetsEnd = fjJets_->end();
00226
00227 for (; pseudojetTMP != fjJetsEnd ; ++pseudojetTMP) {
00228 if(pseudojetTMP->perp() < puPtMin_) continue;
00229
00230
00231 for(vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++)
00232 {
00233 double dr = reco::deltaR(geo_->getPosition((DetId)(*im)),(*pseudojetTMP));
00234 vector<pair<int,int> >::const_iterator exclude = find(excludedTowers.begin(),excludedTowers.end(),pair<int,int>(im->ieta(),im->iphi()));
00235 if( dr < radiusPU_ && exclude == excludedTowers.end()) {
00236 ntowersWithJets_[(*im).ieta()]++;
00237 excludedTowers.push_back(pair<int,int>(im->ieta(),im->iphi()));
00238 }
00239 }
00240
00241 vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(),
00242 fjInputsEnd = fjInputs_->end();
00243
00244 for (; it != fjInputsEnd; ++it ) {
00245 int index = it->user_index();
00246 int ie = ieta((*inputs_)[index]);
00247 int ip = iphi((*inputs_)[index]);
00248
00249 vector<pair<int,int> >::const_iterator exclude = find(excludedTowers.begin(),excludedTowers.end(),pair<int,int>(ie,ip));
00250 if(exclude != excludedTowers.end()) {
00251 jettowers.push_back(index);
00252 }
00253 }
00254 }
00255
00256
00257
00258
00259 for(vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(),
00260 fjInputsEnd = fjInputs_->end(); it != fjInputsEnd; ++it ) {
00261
00262 int index = it->user_index();
00263 vector<int>::const_iterator itjet = find(jettowers.begin(),jettowers.end(),index);
00264 if( itjet == jettowers.end() ){
00265
00266 const reco::CandidatePtr& originalTower = (*inputs_)[index];
00267 fastjet::PseudoJet orphan(originalTower->px(),originalTower->py(),originalTower->pz(),originalTower->energy());
00268 orphan.set_user_index(index);
00269
00270 orphanInput.push_back(orphan);
00271 }
00272 }
00273 }
00274
00275
00276 void PileUpSubtractor::offsetCorrectJets()
00277 {
00278
00279 LogDebug("PileUpSubtractor")<<"The subtractor correcting jets...\n";
00280 jetOffset_.clear();
00281 using namespace reco;
00282
00283
00284
00285
00286
00287 jetOffset_.reserve(fjJets_->size());
00288
00289 vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
00290 jetsEnd = fjJets_->end();
00291 for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
00292
00293 int ijet = pseudojetTMP - fjJets_->begin();
00294 jetOffset_[ijet] = 0;
00295
00296 std::vector<fastjet::PseudoJet> towers =
00297 sorted_by_pt(fjClusterSeq_->constituents(*pseudojetTMP));
00298
00299 double newjetet = 0.;
00300 for(vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(),
00301 towEnd = towers.end();
00302 ito != towEnd;
00303 ++ito)
00304 {
00305 const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
00306 int it = ieta( originalTower );
00307 double Original_Et = originalTower->et();
00308 double etnew = Original_Et - (*emean_.find(it)).second - (*esigma_.find(it)).second;
00309 if(etnew < 0.) etnew = 0;
00310 newjetet = newjetet + etnew;
00311 jetOffset_[ijet] += Original_Et - etnew;
00312 }
00313
00314 double mScale = newjetet/pseudojetTMP->Et();
00315 LogDebug("PileUpSubtractor")<<"pseudojetTMP->Et() : "<<pseudojetTMP->Et()<<"\n";
00316 LogDebug("PileUpSubtractor")<<"newjetet : "<<newjetet<<"\n";
00317 LogDebug("PileUpSubtractor")<<"jetOffset_[ijet] : "<<jetOffset_[ijet]<<"\n";
00318 LogDebug("PileUpSubtractor")<<"pseudojetTMP->Et() - jetOffset_[ijet] : "<<pseudojetTMP->Et() - jetOffset_[ijet]<<"\n";
00319 LogDebug("PileUpSubtractor")<<"Scale is : "<<mScale<<"\n";
00320
00321 int cshist = pseudojetTMP->cluster_hist_index();
00322 pseudojetTMP->reset(pseudojetTMP->px()*mScale, pseudojetTMP->py()*mScale,
00323 pseudojetTMP->pz()*mScale, pseudojetTMP->e()*mScale);
00324 pseudojetTMP->set_cluster_hist_index(cshist);
00325
00326 }
00327 }
00328
00329 double PileUpSubtractor::getCone(double cone, double eta, double phi, double& et, double& pu){
00330 pu = 0;
00331
00332 for(vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++){
00333 if( im->depth() != 1 ) continue;
00334 const GlobalPoint& point = geo_->getPosition((DetId)(*im));
00335 double dr = reco::deltaR(point.eta(),point.phi(),eta,phi);
00336 if( dr < cone){
00337 pu += (*emean_.find(im->ieta())).second+(*esigma_.find(im->ieta())).second;
00338 }
00339 }
00340
00341 return pu;
00342 }
00343
00344 double PileUpSubtractor::getMeanAtTower(const reco::CandidatePtr & in) const{
00345 int it = ieta(in);
00346 return (*emean_.find(it)).second;
00347 }
00348
00349 double PileUpSubtractor::getSigmaAtTower(const reco::CandidatePtr & in) const {
00350 int it = ieta(in);
00351 return (*esigma_.find(it)).second;
00352 }
00353
00354 double PileUpSubtractor::getPileUpAtTower(const reco::CandidatePtr & in) const {
00355 int it = ieta(in);
00356 return (*emean_.find(it)).second + (*esigma_.find(it)).second;
00357 }
00358
00359 int PileUpSubtractor::getN(const reco::CandidatePtr & in) const {
00360 int it = ieta(in);
00361
00362 int n = (*(geomtowers_.find(it))).second - (*(ntowersWithJets_.find(it))).second;
00363 return n;
00364
00365 }
00366
00367 int PileUpSubtractor::getNwithJets(const reco::CandidatePtr & in) const {
00368 int it = ieta(in);
00369 int n = (*(ntowersWithJets_.find(it))).second;
00370 return n;
00371
00372 }
00373
00374
00375 int PileUpSubtractor::ieta(const reco::CandidatePtr & in) const {
00376 int it = 0;
00377 const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00378 if(ctc){
00379 it = ctc->id().ieta();
00380 } else
00381 {
00382 throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
00383 }
00384 return it;
00385 }
00386
00387 int PileUpSubtractor::iphi(const reco::CandidatePtr & in) const {
00388 int it = 0;
00389 const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00390 if(ctc){
00391 it = ctc->id().iphi();
00392 } else
00393 {
00394 throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
00395 }
00396 return it;
00397 }
00398
00399
00400 #include "FWCore/PluginManager/interface/PluginFactory.h"
00401 EDM_REGISTER_PLUGINFACTORY(PileUpSubtractorFactory,"PileUpSubtractorFactory");
00402
00403
00404