CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/RecoJets/JetPlusTracks/plugins/JetPlusTrackProducerAA.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    JetPlusTrack
00004 // Class:      JetPlusTrack
00005 // 
00013 //
00014 // Original Author:  Olga Kodolova,40 R-A12,+41227671273,
00015 //         Created:  Fri Feb 19 10:14:02 CET 2010
00016 // $Id: JetPlusTrackProducerAA.cc,v 1.8 2011/07/01 08:16:15 kodolova Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
00022 #include <memory>
00023 
00024 // user include files
00025 #include "FWCore/Framework/interface/Frameworkfwd.h"
00026 #include "FWCore/Framework/interface/EDProducer.h"
00027 
00028 #include "FWCore/Framework/interface/Event.h"
00029 #include "FWCore/Framework/interface/MakerMacros.h"
00030 
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 
00033 #include "RecoJets/JetPlusTracks/plugins/JetPlusTrackProducerAA.h"
00034 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00035 #include "DataFormats/JetReco/interface/CaloJet.h"
00036 #include "DataFormats/JetReco/interface/JPTJetCollection.h"
00037 #include "DataFormats/JetReco/interface/JPTJet.h"
00038 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00039 #include "DataFormats/TrackReco/interface/Track.h"
00040 #include "DataFormats/JetReco/interface/Jet.h"
00041 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00042 #include "DataFormats/VertexReco/interface/Vertex.h"
00043 #include "DataFormats/Math/interface/deltaPhi.h"
00044 #include "DataFormats/Math/interface/deltaR.h"
00045 
00046 //=>
00047 #include "RecoJets/JetAssociationAlgorithms/interface/JetTracksAssociationXtrpCalo.h"
00048 #include "DataFormats/GeometrySurface/interface/Cylinder.h"
00049 #include "DataFormats/GeometrySurface/interface/Plane.h"
00050 #include "DataFormats/Math/interface/deltaR.h"
00051 #include "DataFormats/Math/interface/Vector3D.h"
00052 #include "MagneticField/Engine/interface/MagneticField.h"
00053 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
00054 #include "TrackingTools/TrajectoryParametrization/interface/GlobalTrajectoryParameters.h"
00055 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00056 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00057 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00058 #include "DataFormats/DetId/interface/DetId.h"
00059 #include "DataFormats/JetReco/interface/CaloJet.h"
00060 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00061 //=>
00062 
00063 #include <string>
00064 
00065 using namespace std;
00066 using namespace jpt;
00067 
00068 //
00069 // constants, enums and typedefs
00070 //
00071 
00072 
00073 //
00074 // static data member definitions
00075 //
00076 
00077 //
00078 // constructors and destructor
00079 //
00080 JetPlusTrackProducerAA::JetPlusTrackProducerAA(const edm::ParameterSet& iConfig)
00081 {
00082    //register your products
00083    src = iConfig.getParameter<edm::InputTag>("src");
00084    alias = iConfig.getUntrackedParameter<string>("alias");
00085    mTracks = iConfig.getParameter<edm::InputTag> ("tracks");
00086    srcPVs_ = iConfig.getParameter<edm::InputTag>("srcPVs"); 
00087    vectorial_ = iConfig.getParameter<bool>("VectorialCorrection");
00088    useZSP = iConfig.getParameter<bool>("UseZSP");
00089    std::string tq = iConfig.getParameter<std::string>("TrackQuality");
00090    trackQuality_ = reco::TrackBase::qualityByName(tq);
00091    mConeSize = iConfig.getParameter<double> ("coneSize");
00092 //=>
00093    mExtrapolations = iConfig.getParameter<edm::InputTag> ("extrapolations");
00094 //=>
00095    mJPTalgo  = new JetPlusTrackCorrector(iConfig);
00096    if(useZSP) mZSPalgo  = new ZSPJPTJetCorrector(iConfig);
00097 
00098    produces<reco::JPTJetCollection>().setBranchAlias(alias); 
00099      
00100 }
00101 
00102 
00103 JetPlusTrackProducerAA::~JetPlusTrackProducerAA()
00104 {
00105  
00106    // do anything here that needs to be done at desctruction time
00107    // (e.g. close files, deallocate resources etc.)
00108 
00109 }
00110 
00111 
00112 //
00113 // member functions
00114 //
00115 
00116 // ------------ method called to produce the data  ------------
00117 void
00118 JetPlusTrackProducerAA::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00119 {
00120    using namespace edm; 
00121 
00122 
00123 //  std::cout<<" RecoJets::JetPlusTrackProducerAA::produce "<<std::endl;
00124 
00125 // get stuff from Event
00126   edm::Handle <edm::View <reco::CaloJet> > jets_h;
00127   iEvent.getByLabel (src, jets_h);
00128   
00129   edm::Handle <reco::TrackCollection> tracks_h;
00130   iEvent.getByLabel (mTracks, tracks_h);
00131   
00132   std::vector <reco::TrackRef> fTracks;
00133   fTracks.reserve (tracks_h->size());
00134   for (unsigned i = 0; i < tracks_h->size(); ++i) {
00135              fTracks.push_back (reco::TrackRef (tracks_h, i));
00136   } 
00137 
00138 //=>
00139   edm::Handle <std::vector<reco::TrackExtrapolation> > extrapolations_h;
00140   iEvent.getByLabel (mExtrapolations, extrapolations_h);
00141 
00142 //  std::cout<<"JetPlusTrackProducerAA::produce, extrapolations_h="<<extrapolations_h->size()<<std::endl;  
00143 //=>
00144 
00145   std::auto_ptr<reco::JPTJetCollection> pOut(new reco::JPTJetCollection());
00146   
00147   reco::JPTJetCollection tmpColl;
00148 
00149   for (unsigned i = 0; i < jets_h->size(); ++i) {
00150 
00151    const reco::CaloJet* oldjet = &(*(jets_h->refAt(i)));
00152    
00153    reco::CaloJet corrected = *oldjet; 
00154    
00155 // ZSP corrections    
00156 
00157    double factorZSP = 1.;
00158    if(useZSP) factorZSP = mZSPalgo->correction(corrected, iEvent, iSetup);
00159 
00160 //   std::cout << " UseZSP = "<<useZSP<<std::endl;
00161 
00162 
00163    corrected.scaleEnergy (factorZSP);
00164 
00165 // JPT corrections 
00166 
00167    double scaleJPT = 1.; 
00168 
00169    math::XYZTLorentzVector p4;
00170    
00171    if ( !vectorial_ ) {
00172             
00173    scaleJPT = mJPTalgo->correction ( corrected, *oldjet, iEvent, iSetup );
00174    p4 = math::XYZTLorentzVector( corrected.px()*scaleJPT, 
00175                                  corrected.py()*scaleJPT,
00176                                  corrected.pz()*scaleJPT, 
00177                                  corrected.energy()*scaleJPT );
00178    } else {
00179    scaleJPT = mJPTalgo->correction( corrected, *oldjet, iEvent, iSetup, p4 );
00180   }         
00181 
00182 // Construct JPTJet constituent
00183   jpt::MatchedTracks pions;
00184   jpt::MatchedTracks muons;
00185   jpt::MatchedTracks elecs;
00186   
00187   bool ok = mJPTalgo->matchTracks( *oldjet, 
00188                                     iEvent, 
00189                                     iSetup,
00190                                      pions, 
00191                                      muons, 
00192                                      elecs );
00193 
00194    
00195   reco::JPTJet::Specific specific;
00196 
00197   if(ok) {
00198 //    std::cout<<" Size of Pion in-in "<<pions.inVertexInCalo_.size()<<" in-out "<<pions.inVertexOutOfCalo_.size()
00199 //             <<" out-in "<<pions.outOfVertexInCalo_.size()<<" Oldjet "<<oldjet->et()<<" factorZSP "<<factorZSP
00200 //             <<"  "<<corrected.et()<<" scaleJPT "<<scaleJPT<<" after JPT "<<p4.pt()<<std::endl;
00201     
00202 
00203     specific.pionsInVertexInCalo = pions.inVertexInCalo_;
00204     specific.pionsInVertexOutCalo = pions.inVertexOutOfCalo_;
00205     specific.pionsOutVertexInCalo = pions.outOfVertexInCalo_;
00206     specific.muonsInVertexInCalo = muons.inVertexInCalo_;
00207     specific.muonsInVertexOutCalo = muons.inVertexOutOfCalo_;
00208     specific.muonsOutVertexInCalo = muons.outOfVertexInCalo_;
00209     specific.elecsInVertexInCalo = elecs.inVertexInCalo_;
00210     specific.elecsInVertexOutCalo = elecs.inVertexOutOfCalo_;
00211     specific.elecsOutVertexInCalo = elecs.outOfVertexInCalo_;
00212   }
00213 
00214 // Fill JPT Specific
00215     edm::RefToBase<reco::Jet> myjet = (edm::RefToBase<reco::Jet>)jets_h->refAt(i);
00216     specific.theCaloJetRef = myjet;
00217     specific.mZSPCor = factorZSP;
00218     specific.mResponseOfChargedWithEff = (float)mJPTalgo->getResponseOfChargedWithEff();
00219     specific.mResponseOfChargedWithoutEff = (float)mJPTalgo->getResponseOfChargedWithoutEff();
00220     specific.mSumPtOfChargedWithEff = (float)mJPTalgo->getSumPtWithEff();
00221     specific.mSumPtOfChargedWithoutEff = (float)mJPTalgo->getSumPtWithoutEff();
00222     specific.mSumEnergyOfChargedWithEff = (float)mJPTalgo->getSumEnergyWithEff();
00223     specific.mSumEnergyOfChargedWithoutEff = (float)mJPTalgo->getSumEnergyWithoutEff();
00224     specific.mChargedHadronEnergy = (float)mJPTalgo->getSumEnergyWithoutEff();
00225 // Fill Charged Jet shape parameters
00226    double deR2Tr = 0.;
00227    double deEta2Tr = 0.;
00228    double dePhi2Tr = 0.;
00229    double Zch = 0.;
00230    double Pout2 = 0.;
00231    double Pout = 0.;
00232    double denominator_tracks = 0.;
00233    int ntracks = 0;
00234 
00235    for( reco::TrackRefVector::const_iterator it = pions.inVertexInCalo_.begin(); it != pions.inVertexInCalo_.end(); it++) {
00236     double deR = deltaR((*it)->eta(), (*it)->phi(), p4.eta(), p4.phi());
00237     double deEta = (*it)->eta() - p4.eta();
00238     double dePhi = deltaPhi((*it)->phi(), p4.phi());
00239      if((**it).ptError()/(**it).pt() < 0.1) {
00240        deR2Tr   =  deR2Tr + deR*deR*(*it)->pt();
00241        deEta2Tr = deEta2Tr + deEta*deEta*(*it)->pt();
00242        dePhi2Tr = dePhi2Tr + dePhi*dePhi*(*it)->pt();
00243        denominator_tracks = denominator_tracks + (*it)->pt();
00244        Zch    = Zch + (*it)->pt();
00245        
00246        Pout2 = Pout2 + (**it).p()*(**it).p() - (Zch*p4.P())*(Zch*p4.P());
00247        ntracks++;
00248      }
00249    }
00250    for( reco::TrackRefVector::const_iterator it = muons.inVertexInCalo_.begin(); it != muons.inVertexInCalo_.end(); it++) {
00251     double deR = deltaR((*it)->eta(), (*it)->phi(), p4.eta(), p4.phi());
00252     double deEta = (*it)->eta() - p4.eta();
00253     double dePhi = deltaPhi((*it)->phi(), p4.phi());
00254      if((**it).ptError()/(**it).pt() < 0.1) {
00255        deR2Tr   =  deR2Tr + deR*deR*(*it)->pt();
00256        deEta2Tr = deEta2Tr + deEta*deEta*(*it)->pt();
00257        dePhi2Tr = dePhi2Tr + dePhi*dePhi*(*it)->pt();
00258        denominator_tracks = denominator_tracks + (*it)->pt();
00259        Zch    = Zch + (*it)->pt();
00260        
00261        Pout2 = Pout2 + (**it).p()*(**it).p() - (Zch*p4.P())*(Zch*p4.P());
00262        ntracks++;
00263      }
00264    }
00265    for( reco::TrackRefVector::const_iterator it = elecs.inVertexInCalo_.begin(); it != elecs.inVertexInCalo_.end(); it++) {
00266     double deR = deltaR((*it)->eta(), (*it)->phi(), p4.eta(), p4.phi());
00267     double deEta = (*it)->eta() - p4.eta();
00268     double dePhi = deltaPhi((*it)->phi(), p4.phi());
00269      if((**it).ptError()/(**it).pt() < 0.1) {
00270        deR2Tr   =  deR2Tr + deR*deR*(*it)->pt();
00271        deEta2Tr = deEta2Tr + deEta*deEta*(*it)->pt();
00272        dePhi2Tr = dePhi2Tr + dePhi*dePhi*(*it)->pt();
00273        denominator_tracks = denominator_tracks + (*it)->pt();
00274        Zch    = Zch + (*it)->pt();
00275        
00276        Pout2 = Pout2 + (**it).p()*(**it).p() - (Zch*p4.P())*(Zch*p4.P());
00277        ntracks++;
00278      }
00279    }
00280    for( reco::TrackRefVector::const_iterator it = pions.inVertexOutOfCalo_.begin(); it != pions.inVertexOutOfCalo_.end(); it++) { 
00281      Zch    =  Zch + (*it)->pt();
00282    }
00283    for( reco::TrackRefVector::const_iterator it = muons.inVertexOutOfCalo_.begin(); it != muons.inVertexOutOfCalo_.end(); it++) { 
00284      Zch    =  Zch + (*it)->pt();
00285    }
00286    for( reco::TrackRefVector::const_iterator it = elecs.inVertexOutOfCalo_.begin(); it != elecs.inVertexOutOfCalo_.end(); it++) { 
00287      Zch    =  Zch + (*it)->pt();
00288    }
00289 
00290      if(mJPTalgo->getSumPtForBeta()> 0.) Zch = Zch/mJPTalgo->getSumPtForBeta();
00291 
00292  //    std::cout<<" Zch "<< Zch<<" "<<mJPTalgo->getSumPtForBeta()<<std::endl;
00293 
00294         if(ntracks > 0) {
00295           Pout   = sqrt(fabs(Pout2))/ntracks;
00296         }
00297 
00298 
00299           if (denominator_tracks!=0){
00300             deR2Tr  = deR2Tr/denominator_tracks;
00301             deEta2Tr= deEta2Tr/denominator_tracks;
00302             dePhi2Tr= dePhi2Tr/denominator_tracks;
00303           }
00304 
00305       specific.R2momtr = deR2Tr;
00306       specific.Eta2momtr = deEta2Tr;
00307       specific.Phi2momtr = dePhi2Tr;
00308       specific.Pout = Pout;
00309       specific.Zch = Zch;
00310 
00311 // Create JPT jet
00312 
00313    reco::Particle::Point vertex_=reco::Jet::Point(0,0,0);
00314    
00315 // If we add primary vertex
00316    edm::Handle<reco::VertexCollection> pvCollection;
00317    iEvent.getByLabel(srcPVs_,pvCollection);
00318    if ( pvCollection.isValid() && pvCollection->size()>0 ) vertex_=pvCollection->begin()->position();
00319  
00320    reco::JPTJet fJet(p4, vertex_, specific, corrected.getJetConstituents()); 
00321   // fJet.printJet();
00322 
00323 // Temporarily collection before correction for background
00324  
00325    tmpColl.push_back(fJet);     
00326 
00327   }
00328 
00329 //=======================================================================================================>
00330 // Correction for background
00331 
00332   reco::TrackRefVector trBgOutOfCalo;
00333   reco::TrackRefVector trBgOutOfVertex = calculateBGtracksJet(tmpColl,fTracks,extrapolations_h,trBgOutOfCalo);
00334 
00335 //===> Area without Jets 
00336     std::map<reco::JPTJetCollection::iterator, double> AreaNonJet;
00337     
00338     for(reco::JPTJetCollection::iterator ij1=tmpColl.begin(); ij1!=tmpColl.end(); ij1++) 
00339     {
00340       int nj1 = 1;
00341       for(reco::JPTJetCollection::iterator ij2=tmpColl.begin(); ij2!=tmpColl.end(); ij2++) 
00342       {
00343        if(ij2 == ij1) continue;
00344        if(fabs((*ij1).eta() - (*ij2).eta()) > 0.5 ) continue;
00345        nj1++;
00346 
00347       }
00348 
00349       AreaNonJet[ij1] = 4*M_PI*mConeSize - nj1*4*mConeSize*mConeSize;
00350 
00351 //      std::cout<<"+++AreaNonJet[ij1]="<<AreaNonJet[ij1]<<" nj1="<<nj1<<std::endl;
00352     }
00353 
00354 //===>
00355 
00356 //  std::cout<<" The size of BG tracks: trBgOutOfVertex= "<<trBgOutOfVertex.size()
00357 //           <<" trBgOutOfCalo= "<<trBgOutOfCalo.size()<<std::endl;
00358 //
00359 //  std::cout<<" The size of JPT jet collection "<<tmpColl.size()<<std::endl;
00360   
00361   for(reco::JPTJetCollection::iterator ij=tmpColl.begin(); ij!=tmpColl.end(); ij++)
00362   {    
00363 // Correct JPTjet for background tracks
00364 
00365    const reco::TrackRefVector pioninin  = (*ij).getPionsInVertexInCalo();
00366    const reco::TrackRefVector pioninout = (*ij).getPionsInVertexOutCalo();
00367    
00368    double ja = (AreaNonJet.find(ij))->second;
00369 
00370 //    std::cout<<"+++ ja="<<ja<<" pioninout="<<pioninout.size()<<std::endl;
00371 
00372    double factorPU = mJPTalgo->correctAA(*ij,trBgOutOfVertex,mConeSize,pioninin,pioninout,ja,trBgOutOfCalo);
00373 
00374    (*ij).scaleEnergy (factorPU);
00375    
00376 //   std::cout<<" FactorPU "<<factorPU<<std::endl;
00377    
00378 // Output module
00379     pOut->push_back(*ij);
00380     
00381 //    std::cout<<" New JPT energy "<<(*ij).et()<<" "<<(*ij).pt()<<" "<<(*ij).eta()<<" "<<(*ij).phi()<<std::endl;
00382     
00383   }
00384   
00385    iEvent.put(pOut);
00386    
00387 }
00388 // -----------------------------------------------
00389 // ------------ calculateBGtracksJet  ------------
00390 // ------------ Tracks not included in jets ------
00391 // -----------------------------------------------
00392 reco::TrackRefVector  JetPlusTrackProducerAA::calculateBGtracksJet(reco::JPTJetCollection& fJets, std::vector <reco::TrackRef>& fTracks,
00393                                                                 edm::Handle <std::vector<reco::TrackExtrapolation> > & extrapolations_h, 
00394                                                                                                    reco::TrackRefVector& trBgOutOfCalo){ 
00395 
00396   
00397   reco::TrackRefVector trBgOutOfVertex;
00398   
00399   for (unsigned t = 0; t < fTracks.size(); ++t) {
00400 
00401      int track_bg = 0;
00402      
00403     // if(!(*fTracks[t]).quality(trackQuality_))
00404     // {
00405       // cout<<"BG, BAD trackQuality, ptBgV="<<fTracks[t]->pt()<<" etaBgV = "<<fTracks[t]->eta()<<" phiBgV = "<<fTracks[t]->phi()<<endl;
00406       // continue;
00407     // }
00408 
00409     const reco::Track* track = &*(fTracks[t]);
00410     double trackEta = track->eta();
00411     double trackPhi = track->phi();
00412 
00413   //  std::cout<<"++++++++++++++++>  track="<<t<<" trackEta="<<trackEta<<" trackPhi="<<trackPhi
00414   //           <<" coneSize="<<mConeSize<<std::endl;
00415    
00416    //loop on jets
00417     for (unsigned j = 0; j < fJets.size(); ++j) {
00418 
00419      const reco::Jet* jet = &(fJets[j]);
00420      double jetEta = jet->eta();
00421      double jetPhi = jet->phi();
00422 
00423     //  std::cout<<"-jet="<<j<<" jetEt ="<<jet->pt()
00424     //  <<" jetE="<<jet->energy()<<" jetEta="<<jetEta<<" jetPhi="<<jetPhi<<std::endl;
00425 
00426       if(fabs(jetEta - trackEta) < mConeSize) {
00427        double dphiTrackJet = fabs(trackPhi - jetPhi);
00428        if(dphiTrackJet > M_PI) dphiTrackJet = 2*M_PI - dphiTrackJet;
00429 
00430        if(dphiTrackJet < mConeSize) 
00431         {
00432          track_bg = 1;
00433      //    std::cout<<"===>>>> Track inside jet at vertex, track_bg="<< track_bg <<" track="<<t<<" jet="<<j
00434       //            <<" trackEta="<<trackEta<<" trackPhi="<<trackPhi
00435       //            <<" jetEta="<<jetEta<<" jetPhi="<<jetPhi<<std::endl;
00436         }
00437       }      
00438     } //jets
00439 
00440     if( track_bg == 0 ) 
00441      {
00442        trBgOutOfVertex.push_back (fTracks[t]);
00443     
00444 //       std::cout<<"------Track outside jet at vertex, track_bg="<< track_bg<<" track="<<t
00445 //               <<" trackEta="<<trackEta<<" trackPhi="<<trackPhi <<std::endl;    
00446      }
00447 
00448   } //tracks    
00449 
00450 //=====> Propagate BG tracks to calo 
00451     int nValid = 0;
00452     for ( std::vector<reco::TrackExtrapolation>::const_iterator xtrpBegin = extrapolations_h->begin(),
00453           xtrpEnd = extrapolations_h->end(), ixtrp = xtrpBegin;
00454           ixtrp != xtrpEnd; ++ixtrp ) {
00455 
00456 //    std::cout<<"JetPlusTrackProducerAA::calculateBGtracksJet: initial track pt= "<<ixtrp->track()->pt()
00457 //             <<" eta= "<<ixtrp->track()->eta()<<" phi="<<ixtrp->track()->phi()
00458 //             <<" Valid? "<<ixtrp->isValid().at(0)<<std::endl;
00459 
00460           //if( ixtrp->isValid().at(0) == 0 ) continue;
00461           //in DF change in 4.2, all entries are valid.
00462           nValid++;
00463 
00464           reco::TrackRefVector::iterator it = find(trBgOutOfVertex.begin(),trBgOutOfVertex.end(),(*ixtrp).track() );
00465 
00466           if ( it != trBgOutOfVertex.end() ){
00467              trBgOutOfCalo.push_back (*it);
00468 
00469 //          std::cout<<"+++trBgOutOfCalo, pt= "<<ixtrp->track()->pt()<<" eta= "<<ixtrp->track()->eta()<<" phi="<<ixtrp->track()->phi()
00470 //                   <<" Valid? "<<ixtrp->isValid().at(0)<<std::endl;
00471           }
00472 
00473     }
00474 
00475 //     std::cout<<"calculateBGtracksJet, trBgOutOfCalo="<<trBgOutOfCalo.size()
00476 //              <<" trBgOutOfVertex="<<trBgOutOfVertex.size()<<" nValid="<<nValid<<endl;
00477 //=====>
00478 
00479   return trBgOutOfVertex;
00480 }
00481 // ------------ method called once each job just before starting event loop  ------------
00482 void 
00483 JetPlusTrackProducerAA::beginJob()
00484 {
00485 }
00486 
00487 // ------------ method called once each job just after ending the event loop  ------------
00488 void 
00489 JetPlusTrackProducerAA::endJob() {
00490 }
00491 
00492 //define this as a plug-in
00493 //DEFINE_FWK_MODULE(JetPlusTrackProducerAA);