CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/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 ( doAreaFastjet_ || doRhoFastjet_ ) {
00028       // default Ghost_EtaMax should be 5
00029       double ghostEtaMax = iConfig.getParameter<double>("Ghost_EtaMax");
00030       // default Active_Area_Repeats 1
00031       int    activeAreaRepeats = iConfig.getParameter<int> ("Active_Area_Repeats");
00032       // default GhostArea 0.01
00033       double ghostArea = iConfig.getParameter<double> ("GhostArea");
00034       fjActiveArea_ =  ActiveAreaSpecPtr(new fastjet::ActiveAreaSpec(ghostEtaMax,
00035                                                                      activeAreaRepeats,
00036                                                                      ghostArea));
00037       fjRangeDef_ = RangeDefPtr( new fastjet::RangeDefinition(ghostEtaMax) );
00038    } 
00039 }
00040 
00041 void PileUpSubtractor::reset(std::vector<edm::Ptr<reco::Candidate> >& input,
00042                              std::vector<fastjet::PseudoJet>& towers,
00043                              std::vector<fastjet::PseudoJet>& output){
00044   
00045   inputs_ = &input;
00046   fjInputs_ = &towers;
00047   fjJets_ = &output;
00048   fjOriginalInputs_ = (*fjInputs_);
00049   for (unsigned int i = 0; i < fjInputs_->size(); ++i){
00050     fjOriginalInputs_[i].set_user_index((*fjInputs_)[i].user_index());
00051   }
00052   
00053 }
00054 
00055 void PileUpSubtractor::setAlgorithm(ClusterSequencePtr& algorithm){
00056   fjClusterSeq_ = algorithm;
00057 }
00058 
00059 void PileUpSubtractor::setupGeometryMap(edm::Event& iEvent,const edm::EventSetup& iSetup)
00060 {
00061 
00062   LogDebug("PileUpSubtractor")<<"The subtractor setting up geometry...\n";
00063 
00064   if(geo_ == 0) {
00065     edm::ESHandle<CaloGeometry> pG;
00066     iSetup.get<CaloGeometryRecord>().get(pG);
00067     geo_ = pG.product();
00068     std::vector<DetId> alldid =  geo_->getValidDetIds();
00069 
00070     int ietaold = -10000;
00071     ietamax_ = -10000;
00072     ietamin_ = 10000;
00073     for(std::vector<DetId>::const_iterator did=alldid.begin(); did != alldid.end(); did++){
00074       if( (*did).det() == DetId::Hcal ){
00075         HcalDetId hid = HcalDetId(*did);
00076         if( (hid).depth() == 1 ) {
00077           allgeomid_.push_back(*did);
00078 
00079           if((hid).ieta() != ietaold){
00080             ietaold = (hid).ieta();
00081             geomtowers_[(hid).ieta()] = 1;
00082             if((hid).ieta() > ietamax_) ietamax_ = (hid).ieta();
00083             if((hid).ieta() < ietamin_) ietamin_ = (hid).ieta();
00084           }
00085           else{
00086             geomtowers_[(hid).ieta()]++;
00087           }
00088         }
00089       }
00090     }
00091   }
00092 
00093   for (int i = ietamin_; i<ietamax_+1; i++) {
00094     ntowersWithJets_[i] = 0;
00095   }
00096 }
00097 
00098 //
00099 // Calculate mean E and sigma from jet collection "coll".  
00100 //
00101 void PileUpSubtractor::calculatePedestal( vector<fastjet::PseudoJet> const & coll )
00102 {
00103   LogDebug("PileUpSubtractor")<<"The subtractor calculating pedestals...\n";
00104 
00105   map<int,double> emean2;
00106   map<int,int> ntowers;
00107     
00108   int ietaold = -10000;
00109   int ieta0 = -100;
00110    
00111   // Initial values for emean_, emean2, esigma_, ntowers
00112 
00113   for(int i = ietamin_; i < ietamax_+1; i++)
00114     {
00115       emean_[i] = 0.;
00116       emean2[i] = 0.;
00117       esigma_[i] = 0.;
00118       ntowers[i] = 0;
00119     }
00120     
00121   for (vector<fastjet::PseudoJet>::const_iterator input_object = coll.begin (),
00122          fjInputsEnd = coll.end();  
00123        input_object != fjInputsEnd; ++input_object) {
00124 
00125      const reco::CandidatePtr & originalTower=(*inputs_)[ input_object->user_index()];
00126     ieta0 = ieta( originalTower );
00127     double Original_Et = originalTower->et();
00128 
00129   if( ieta0-ietaold != 0 )
00130       {
00131         emean_[ieta0] = emean_[ieta0]+Original_Et;
00132         emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
00133         ntowers[ieta0] = 1;
00134         ietaold = ieta0;
00135       }
00136   else
00137       {
00138         emean_[ieta0] = emean_[ieta0]+Original_Et;
00139         emean2[ieta0] = emean2[ieta0]+Original_Et*Original_Et;
00140         ntowers[ieta0]++;
00141       }
00142 
00143   }
00144 
00145   for(map<int,int>::const_iterator gt = geomtowers_.begin(); gt != geomtowers_.end(); gt++)    
00146     {
00147 
00148       int it = (*gt).first;
00149        
00150       double e1 = (*(emean_.find(it))).second;
00151       double e2 = (*emean2.find(it)).second;
00152       int nt = (*gt).second - (*(ntowersWithJets_.find(it))).second;
00153 
00154       LogDebug("PileUpSubtractor")<<" ieta : "<<it<<" number of towers : "<<nt<<" e1 : "<<e1<<" e2 : "<<e2<<"\n";
00155         
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   int ip = -100;
00182         
00183   for (vector<fastjet::PseudoJet>::iterator input_object = coll.begin (),
00184          fjInputsEnd = coll.end(); 
00185        input_object != fjInputsEnd; ++input_object) {
00186     
00187      reco::CandidatePtr const & itow =  (*inputs_)[ input_object->user_index() ];
00188     
00189     it = ieta( itow );
00190     ip = iphi( itow );
00191     
00192     double etnew = itow->et() - (*(emean_.find(it))).second - (*(esigma_.find(it))).second;
00193     float mScale = etnew/input_object->Et(); 
00194     if(etnew < 0.) mScale = 0.;
00195     
00196     math::XYZTLorentzVectorD towP4(input_object->px()*mScale, input_object->py()*mScale,
00197                                    input_object->pz()*mScale, input_object->e()*mScale);
00198     
00199     int index = input_object->user_index();
00200     input_object->reset ( towP4.px(),
00201                           towP4.py(),
00202                           towP4.pz(),
00203                           towP4.energy() );
00204     input_object->set_user_index(index);
00205   }
00206 }
00207 
00208 void PileUpSubtractor::calculateOrphanInput(vector<fastjet::PseudoJet> & orphanInput) 
00209 {
00210 
00211   LogDebug("PileUpSubtractor")<<"The subtractor calculating orphan input...\n";
00212 
00213   (*fjInputs_) = fjOriginalInputs_;
00214 
00215   vector<int> jettowers; // vector of towers indexed by "user_index"
00216   vector<pair<int,int> >  excludedTowers; // vector of excluded ieta, iphi values
00217 
00218   vector <fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
00219     fjJetsEnd = fjJets_->end();
00220     
00221   for (; pseudojetTMP != fjJetsEnd ; ++pseudojetTMP) {
00222 
00223     vector<fastjet::PseudoJet> newtowers;
00224     // find towers within radiusPU_ of this jet
00225     for(vector<HcalDetId>::const_iterator im = allgeomid_.begin(); im != allgeomid_.end(); im++)
00226       {
00227           double dr = reco::deltaR(geo_->getPosition((DetId)(*im)),(*pseudojetTMP));
00228           if( dr < radiusPU_) {
00229             ntowersWithJets_[(*im).ieta()]++;     
00230             excludedTowers.push_back(pair<int,int>(im->ieta(),im->iphi()));
00231           }
00232       }
00233     
00234     vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(),
00235       fjInputsEnd = fjInputs_->end();
00236       
00237     for (; it != fjInputsEnd; ++it ) {
00238       
00239       int index = it->user_index();
00240       int ie = ieta((*inputs_)[index]);
00241       int ip = iphi((*inputs_)[index]);
00242       
00243       vector<pair<int,int> >::const_iterator exclude = find(excludedTowers.begin(),excludedTowers.end(),pair<int,int>(ie,ip));
00244       if(exclude != excludedTowers.end()) {
00245         newtowers.push_back(*it);
00246         jettowers.push_back(index);
00247       } //dr < 0.5
00248     } // initial input collection  
00249   } // pseudojets
00250   
00251   //
00252   // Create a new collections from the towers not included in jets 
00253   //
00254   for(vector<fastjet::PseudoJet>::const_iterator it = fjInputs_->begin(),
00255         fjInputsEnd = fjInputs_->end(); it != fjInputsEnd; ++it ) {
00256 
00257     int index = it->user_index();
00258     vector<int>::const_iterator itjet = find(jettowers.begin(),jettowers.end(),index);
00259     if( itjet == jettowers.end() ){
00260 
00261       const reco::CandidatePtr& originalTower = (*inputs_)[index];
00262       fastjet::PseudoJet orphan(originalTower->px(),originalTower->py(),originalTower->pz(),originalTower->energy());
00263       orphan.set_user_index(index);
00264 
00265       orphanInput.push_back(orphan); 
00266     }
00267   }
00268 }
00269 
00270 
00271 void PileUpSubtractor::offsetCorrectJets() 
00272 {
00273 
00274   LogDebug("PileUpSubtractor")<<"The subtractor correcting jets...\n";
00275   jetOffset_.clear();
00276   using namespace reco;
00277   
00278   //    
00279   // Reestimate energy of jet (energy of jet with initial map)
00280   //
00281 
00282   jetOffset_.reserve(fjJets_->size());
00283   
00284   vector<fastjet::PseudoJet>::iterator pseudojetTMP = fjJets_->begin (),
00285     jetsEnd = fjJets_->end();
00286   for (; pseudojetTMP != jetsEnd; ++pseudojetTMP) {
00287     
00288     int ijet = pseudojetTMP - fjJets_->begin();
00289     jetOffset_[ijet] = 0;
00290     
00291     std::vector<fastjet::PseudoJet> towers =
00292       sorted_by_pt(fjClusterSeq_->constituents(*pseudojetTMP));
00293     
00294     double newjetet = 0.;       
00295     for(vector<fastjet::PseudoJet>::const_iterator ito = towers.begin(),
00296           towEnd = towers.end(); 
00297         ito != towEnd; 
00298         ++ito)
00299       {
00300         const reco::CandidatePtr& originalTower = (*inputs_)[ito->user_index()];
00301         int it = ieta( originalTower );
00302         double Original_Et = originalTower->et();
00303         double etnew = Original_Et - (*emean_.find(it)).second - (*esigma_.find(it)).second; 
00304         if(etnew < 0.) etnew = 0;
00305         newjetet = newjetet + etnew;
00306         jetOffset_[ijet] += Original_Et - etnew;
00307       }
00308     
00309     double mScale = newjetet/pseudojetTMP->Et();
00310     LogDebug("PileUpSubtractor")<<"pseudojetTMP->Et() : "<<pseudojetTMP->Et()<<"\n";
00311     LogDebug("PileUpSubtractor")<<"newjetet : "<<newjetet<<"\n";
00312     LogDebug("PileUpSubtractor")<<"jetOffset_[ijet] : "<<jetOffset_[ijet]<<"\n";
00313     LogDebug("PileUpSubtractor")<<"pseudojetTMP->Et() - jetOffset_[ijet] : "<<pseudojetTMP->Et() - jetOffset_[ijet]<<"\n";
00314     LogDebug("PileUpSubtractor")<<"Scale is : "<<mScale<<"\n";
00315     
00316     int cshist = pseudojetTMP->cluster_hist_index();
00317     pseudojetTMP->reset(pseudojetTMP->px()*mScale, pseudojetTMP->py()*mScale,
00318                         pseudojetTMP->pz()*mScale, pseudojetTMP->e()*mScale);
00319     pseudojetTMP->set_cluster_hist_index(cshist);
00320     
00321   }
00322 }
00323 
00324 double PileUpSubtractor::getMeanAtTower(const reco::CandidatePtr & in) const{
00325   int it = ieta(in);
00326   return (*emean_.find(it)).second;
00327 }
00328 
00329 double PileUpSubtractor::getSigmaAtTower(const reco::CandidatePtr & in) const {
00330    int it = ieta(in);
00331    return (*esigma_.find(it)).second;
00332 }
00333 
00334 double PileUpSubtractor::getPileUpAtTower(const reco::CandidatePtr & in) const {
00335   int it = ieta(in);
00336   return (*emean_.find(it)).second + (*esigma_.find(it)).second;
00337 }
00338 
00339 int PileUpSubtractor::ieta(const reco::CandidatePtr & in) const {
00340   int it = 0;
00341   const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00342   if(ctc){
00343     it = ctc->id().ieta();
00344   } else
00345     {
00346       throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
00347     }
00348   return it;
00349 }
00350 
00351 int PileUpSubtractor::iphi(const reco::CandidatePtr & in) const {
00352   int it = 0;
00353   const CaloTower* ctc = dynamic_cast<const CaloTower*>(in.get());
00354   if(ctc){
00355     it = ctc->id().iphi();
00356   } else
00357     {
00358       throw cms::Exception("Invalid Constituent") << "CaloJet constituent is not of CaloTower type";
00359     }
00360   return it;
00361 }
00362 
00363 
00364 #include "FWCore/PluginManager/interface/PluginFactory.h"
00365 EDM_REGISTER_PLUGINFACTORY(PileUpSubtractorFactory,"PileUpSubtractorFactory");
00366 
00367 
00368