CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoJets/JetProducers/src/PileUpSubtractor.cc

Go to the documentation of this file.
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       // default Ghost_EtaMax should be 5
00035       double ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
00036       // default Active_Area_Repeats 1
00037       int    activeAreaRepeats = iConfig.getParameter<int> ("Active_Area_Repeats");
00038       // default GhostArea 0.01
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::setDefinition(JetDefPtr const & jetDef){
00062   fjJetDefinition_ = JetDefPtr( new fastjet::JetDefinition( *jetDef ) );
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 // Calculate mean E and sigma from jet collection "coll".  
00106 //
00107 void PileUpSubtractor::calculatePedestal( vector<fastjet::PseudoJet> const & coll )
00108 {
00109   LogDebug("PileUpSubtractor")<<"The subtractor calculating pedestals...\n";
00110   map<int,double> emean2;
00111   map<int,int> ntowers;
00112     
00113   int ietaold = -10000;
00114   int ieta0 = -100;
00115    
00116   // Initial values for emean_, emean2, esigma_, ntowers
00117 
00118   for(int i = ietamin_; i < ietamax_+1; i++)
00119     {
00120       emean_[i] = 0.;
00121       emean2[i] = 0.;
00122       esigma_[i] = 0.;
00123       ntowers[i] = 0;
00124     }
00125     
00126   for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin (),
00127          fjInputsEnd = coll.end();  
00128        input_object != fjInputsEnd; ++input_object) {
00129      const reco::CandidatePtr & originalTower=(*inputs_)[ input_object->user_index()];
00130     ieta0 = ieta( originalTower );
00131     double Original_Et = originalTower->et();
00132   if( ieta0-ietaold != 0 )
00133       {
00134         emean_[ieta0] = emean_[ieta0]+Original_Et;
00135         emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
00136         ntowers[ieta0] = 1;
00137         ietaold = ieta0;
00138       }
00139   else
00140       {
00141         emean_[ieta0] = emean_[ieta0]+Original_Et;
00142         emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
00143         ntowers[ieta0]++;
00144       }
00145   }
00146 
00147   for(map<int,int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++)    
00148     {
00149       int it = (*gt).first;
00150        
00151       double e1 = (*(emean_.find(it))).second;
00152       double e2 = (*emean2.find(it)).second;
00153       int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
00154 
00155       LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" number of towers : "<<nt<<" e1 : "<<e1<<" e2 : "<<e2<<"\n";
00156       if(nt > 0) {
00157         emean_[it] = e1/nt;
00158         double eee = e2/nt - e1*e1/(nt*nt);    
00159         if(eee<0.) eee = 0.;
00160         esigma_[it] = nSigmaPU_*sqrt(eee);
00161       }
00162       else
00163         {
00164           emean_[it] = 0.;
00165           esigma_[it] = 0.;
00166         }
00167       LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" Pedestals : "<<emean_[it]<<"  "<<esigma_[it]<<"\n";
00168     }
00169 }
00170 
00171 
00172 //
00173 // Subtract mean and sigma from fjInputs_
00174 //    
00175 void PileUpSubtractor::subtractPedestal(vector<fastjet::PseudoJet> & coll)
00176 {
00177 
00178   LogDebug("PileUpSubtractor")<<"The subtractor subtracting pedestals...\n";
00179 
00180   int it = -100;
00181   for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
00182          fjInputsEnd = coll.end(); 
00183        input_object != fjInputsEnd; ++input_object) {
00184     
00185      reco::CandidatePtr const & itow =  (*inputs_)[ input_object->user_index() ];
00186     
00187     it = ieta( itow );
00188     
00189     double etnew = itow->et() - (*(emean_.find(it))).second - (*(esigma_.find(it))).second;
00190     float mScale = etnew/input_object->Et(); 
00191     if(etnew < 0.) mScale = 0.;
00192     
00193     math::XYZTLorentzVectorD towP4(input_object->px()*mScale, input_object->py()*mScale,
00194                                    input_object->pz()*mScale, input_object->e()*mScale);
00195     
00196     int index = input_object->user_index();
00197     input_object->reset_momentum ( towP4.px(),
00198                                    towP4.py(),
00199                                    towP4.pz(),
00200                                    towP4.energy() );
00201     input_object->set_user_index(index);
00202   }
00203 }
00204 
00205 void PileUpSubtractor::calculateOrphanInput(vector<fastjet::PseudoJet> & orphanInput) 
00206 {
00207 
00208   LogDebug("PileUpSubtractor")<<"The subtractor calculating orphan input...\n";
00209 
00210   (*fjInputs_) = fjOriginalInputs_;
00211 
00212   vector<int> jettowers; // vector of towers indexed by "user_index"
00213   vector<pair<int,int> >  excludedTowers; // vector of excluded ieta, iphi values
00214 
00215   vector <fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
00216     fjJetsEnd = fjJets_->end();
00217   for (; pseudojetTMP != fjJetsEnd ; ++pseudojetTMP) {
00218     if(pseudojetTMP->perp() < puPtMin_) continue;
00219 
00220     // find towers within radiusPU_ of this jet
00221     for(vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++)
00222       {
00223         double dr = reco::deltaR(geo_->getPosition((DetId)(*im)),(*pseudojetTMP));
00224         vector<pair<int,int> >::const_iterator exclude = find(excludedTowers.begin(),excludedTowers.end(),pair<int,int>(im->ieta(),im->iphi()));
00225         if( dr < radiusPU_ && exclude == excludedTowers.end()) {
00226           ntowersWithJets_[(*im).ieta()]++;     
00227           excludedTowers.push_back(pair<int,int>(im->ieta(),im->iphi()));
00228         }
00229       }
00230     vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(),
00231       fjInputsEnd = fjInputs_->end();
00232       
00233     for (; it != fjInputsEnd; ++it ) {
00234       int index = it->user_index();
00235       int ie = ieta((*inputs_)[index]);
00236       int ip = iphi((*inputs_)[index]);
00237       vector<pair<int,int> >::const_iterator exclude = find(excludedTowers.begin(),excludedTowers.end(),pair<int,int>(ie,ip));
00238       if(exclude != excludedTowers.end()) {
00239         jettowers.push_back(index);
00240       } //dr < radiusPU_
00241     } // initial input collection  
00242   } // pseudojets
00243   
00244   //
00245   // Create a new collections from the towers not included in jets 
00246   //
00247   for(vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(),
00248         fjInputsEnd = fjInputs_->end(); it != fjInputsEnd; ++it ) {
00249     int index = it->user_index();
00250     vector<int>::const_iterator itjet = find(jettowers.begin(),jettowers.end(),index);
00251     if( itjet == jettowers.end() ){
00252       const reco::CandidatePtr& originalTower = (*inputs_)[index];
00253       fastjet::PseudoJet orphan(originalTower->px(),originalTower->py(),originalTower->pz(),originalTower->energy());
00254       orphan.set_user_index(index);
00255 
00256       orphanInput.push_back(orphan); 
00257     }
00258   }
00259 }
00260 
00261 
00262 void PileUpSubtractor::offsetCorrectJets() 
00263 {
00264   LogDebug("PileUpSubtractor")<<"The subtractor correcting jets...\n";
00265   jetOffset_.clear();
00266   using namespace reco;
00267   
00268   //    
00269   // Reestimate energy of jet (energy of jet with initial map)
00270   //
00271   jetOffset_.reserve(fjJets_->size());
00272   vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
00273     jetsEnd = fjJets_->end();
00274   for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
00275     int ijet = pseudojetTMP - fjJets_->begin();
00276     jetOffset_[ijet] = 0;
00277     
00278     std::vector<fastjet::PseudoJet> towers =
00279       fastjet::sorted_by_pt( pseudojetTMP->constituents() );
00280     double newjetet = 0.;       
00281     for(vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(),
00282           towEnd = towers.end(); 
00283         ito != towEnd; 
00284         ++ito)
00285       {
00286         const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
00287         int it = ieta( originalTower );
00288         double Original_Et = originalTower->et();
00289         double etnew = Original_Et - (*emean_.find(it)).second - (*esigma_.find(it)).second; 
00290         if(etnew < 0.) etnew = 0;
00291         newjetet = newjetet + etnew;
00292         jetOffset_[ijet] += Original_Et - etnew;
00293       }
00294     double mScale = newjetet/pseudojetTMP->Et();
00295     LogDebug("PileUpSubtractor")<<"pseudojetTMP->Et() : "<<pseudojetTMP->Et()<<"\n";
00296     LogDebug("PileUpSubtractor")<<"newjetet : "<<newjetet<<"\n";
00297     LogDebug("PileUpSubtractor")<<"jetOffset_[ijet] : "<<jetOffset_[ijet]<<"\n";
00298     LogDebug("PileUpSubtractor")<<"pseudojetTMP->Et() - jetOffset_[ijet] : "<<pseudojetTMP->Et() - jetOffset_[ijet]<<"\n";
00299     LogDebug("PileUpSubtractor")<<"Scale is : "<<mScale<<"\n";
00300     int cshist = pseudojetTMP->cluster_hist_index();
00301     pseudojetTMP->reset_momentum(pseudojetTMP->px()*mScale, pseudojetTMP->py()*mScale,
00302                                  pseudojetTMP->pz()*mScale, pseudojetTMP->e()*mScale);
00303     pseudojetTMP->set_cluster_hist_index(cshist);
00304     
00305   }
00306 }
00307 
00308 double PileUpSubtractor::getCone(double cone, double eta, double phi, double& et, double& pu){
00309   pu = 0;
00310   
00311   for(vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++){
00312      if( im->depth() != 1 ) continue;    
00313     const GlobalPoint& point = geo_->getPosition((DetId)(*im));
00314     double dr = reco::deltaR(point.eta(),point.phi(),eta,phi);      
00315     if( dr < cone){
00316       pu += (*emean_.find(im->ieta())).second+(*esigma_.find(im->ieta())).second;
00317     }
00318   }
00319   
00320   return pu;
00321 }
00322 
00323 double PileUpSubtractor::getMeanAtTower(const reco::CandidatePtr & in) const{
00324   int it = ieta(in);
00325   return (*emean_.find(it)).second;
00326 }
00327 
00328 double PileUpSubtractor::getSigmaAtTower(const reco::CandidatePtr & in) const {
00329    int it = ieta(in);
00330    return (*esigma_.find(it)).second;
00331 }
00332 
00333 double PileUpSubtractor::getPileUpAtTower(const reco::CandidatePtr & in) const {
00334   int it = ieta(in);
00335   return (*emean_.find(it)).second + (*esigma_.find(it)).second;
00336 }
00337 
00338 int PileUpSubtractor::getN(const reco::CandidatePtr & in) const {
00339    int it = ieta(in);
00340    
00341    int n = (*(geomtowers_.find(it))).second - (*(ntowersWithJets_.find(it))).second;
00342    return n;
00343 
00344 }
00345 
00346 int PileUpSubtractor::getNwithJets(const reco::CandidatePtr & in) const {
00347    int it = ieta(in);
00348    int n = (*(ntowersWithJets_.find(it))).second;
00349    return n;
00350 
00351 }
00352 
00353 
00354 int PileUpSubtractor::ieta(const reco::CandidatePtr & in) const {
00355   int it = 0;
00356   const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00357   if(ctc){
00358     it = ctc->id().ieta();
00359   } else
00360     {
00361       throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
00362     }
00363   return it;
00364 }
00365 
00366 int PileUpSubtractor::iphi(const reco::CandidatePtr & in) const {
00367   int it = 0;
00368   const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00369   if(ctc){
00370     it = ctc->id().iphi();
00371   } else
00372     {
00373       throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
00374     }
00375   return it;
00376 }
00377 
00378 
00379 #include "FWCore/PluginManager/interface/PluginFactory.h"
00380 EDM_REGISTER_PLUGINFACTORY(PileUpSubtractorFactory,"PileUpSubtractorFactory");
00381 
00382 
00383