CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/RecoParticleFlow/PFTracking/plugins/PFTrackProducer.cc

Go to the documentation of this file.
00001 #include <memory>
00002 #include "RecoParticleFlow/PFTracking/interface/PFTrackProducer.h"
00003 #include "RecoParticleFlow/PFTracking/interface/PFTrackTransformer.h"
00004 #include "DataFormats/ParticleFlowReco/interface/PFRecTrack.h"
00005 #include "DataFormats/ParticleFlowReco/interface/PFRecTrackFwd.h"
00006 #include "DataFormats/TrackReco/interface/Track.h"
00007 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00008 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00009 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00010 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
00011 #include "DataFormats/EgammaReco/interface/ElectronSeedFwd.h"
00012 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00013 #include "FWCore/Framework/interface/ESHandle.h"
00014 #include "MagneticField/Engine/interface/MagneticField.h"
00015 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00016 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00017 #include "DataFormats/MuonReco/interface/Muon.h"
00018 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00019 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00020 #include "TrackingTools/IPTools/interface/IPTools.h"
00021 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00022 #include "DataFormats/VertexReco/interface/Vertex.h"
00023 
00024 using namespace std;
00025 using namespace edm;
00026 using namespace reco;
00027 PFTrackProducer::PFTrackProducer(const ParameterSet& iConfig):
00028   pfTransformer_(0)
00029 {
00030   produces<reco::PFRecTrackCollection>();
00031   
00032   tracksContainers_ = 
00033     iConfig.getParameter< vector < InputTag > >("TkColList");
00034   
00035   useQuality_   = iConfig.getParameter<bool>("UseQuality");
00036   
00037   gsfTrackLabel_ = iConfig.getParameter<InputTag>
00038     ("GsfTrackModuleLabel");  
00039 
00040   trackQuality_=reco::TrackBase::qualityByName(iConfig.getParameter<std::string>("TrackQuality"));
00041   
00042   muonColl_ = iConfig.getParameter< InputTag >("MuColl");
00043   
00044   trajinev_ = iConfig.getParameter<bool>("TrajInEvents");
00045   
00046   gsfinev_ = iConfig.getParameter<bool>("GsfTracksInEvents");
00047   vtx_h=iConfig.getParameter<edm::InputTag>("PrimaryVertexLabel");
00048 }
00049 
00050 PFTrackProducer::~PFTrackProducer()
00051 {
00052   delete pfTransformer_;
00053 }
00054 
00055 void
00056 PFTrackProducer::produce(Event& iEvent, const EventSetup& iSetup)
00057 {
00058   
00059   //create the empty collections 
00060   auto_ptr< reco::PFRecTrackCollection > 
00061     PfTrColl (new reco::PFRecTrackCollection);
00062   
00063   //read track collection
00064   Handle<GsfTrackCollection> gsftrackcoll;
00065   bool foundgsf = iEvent.getByLabel(gsfTrackLabel_,gsftrackcoll);
00066   GsfTrackCollection gsftracks;
00067   //Get PV for STIP calculation, if there is none then take the dummy  
00068   Handle<reco::VertexCollection> vertex;
00069   iEvent.getByLabel(vtx_h, vertex);
00070   reco::Vertex dummy;
00071   const reco::Vertex* pv=&dummy;  
00072   if (vertex.isValid()) 
00073     {
00074       pv = &*vertex->begin();    
00075     } 
00076   else 
00077     { // create a dummy PV
00078       reco::Vertex::Error e;
00079       e(0, 0) = 0.0015 * 0.0015;
00080       e(1, 1) = 0.0015 * 0.0015;
00081       e(2, 2) = 15. * 15.;
00082       reco::Vertex::Point p(0, 0, 0);
00083       dummy = reco::Vertex(p, e, 0, 0, 0);   
00084     } 
00085   
00086   edm::ESHandle<TransientTrackBuilder> builder;
00087   iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", builder);
00088   TransientTrackBuilder thebuilder = *(builder.product());
00089   
00090   
00091   if(gsfinev_) {
00092     if(!foundgsf )
00093       LogError("PFTrackProducer")
00094         <<" cannot get GsfTracks (probably in HI events): "
00095         << " please set GsfTracksInEvents = False in RecoParticleFlow/PFTracking/python/pfTrack_cfi.py" << endl;
00096     else
00097       gsftracks  = *(gsftrackcoll.product());
00098   }  
00099   
00100   // read muon collection
00101   Handle< reco::MuonCollection > recMuons;
00102   iEvent.getByLabel(muonColl_, recMuons);
00103   
00104   
00105   for (unsigned int istr=0; istr<tracksContainers_.size();istr++){
00106     
00107     //Track collection
00108     Handle<reco::TrackCollection> tkRefCollection;
00109     iEvent.getByLabel(tracksContainers_[istr], tkRefCollection);
00110     reco::TrackCollection  Tk=*(tkRefCollection.product());
00111     
00112     vector<Trajectory> Tj(0);
00113     if(trajinev_) {
00114       //Trajectory collection
00115       Handle<vector<Trajectory> > tjCollection;
00116       bool found = iEvent.getByLabel(tracksContainers_[istr], tjCollection);
00117         if(!found )
00118           LogError("PFTrackProducer")
00119             <<" cannot get Trajectories of: "
00120             <<  tracksContainers_[istr]
00121             << " please set TrajInEvents = False in RecoParticleFlow/PFTracking/python/pfTrack_cfi.py" << endl;
00122         
00123         Tj =*(tjCollection.product());
00124     }
00125     
00126     
00127     for(unsigned int i=0;i<Tk.size();i++){
00128       
00129       reco::TrackRef trackRef(tkRefCollection, i);
00130         
00131       if (useQuality_ &&
00132           (!(Tk[i].quality(trackQuality_)))){
00133         
00134         bool isMuCandidate = false;
00135         
00136         //TrackRef trackRef(tkRefCollection, i);
00137         
00138         if(recMuons.isValid() ) {
00139           for(unsigned j=0;j<recMuons->size(); j++) {
00140             reco::MuonRef muonref( recMuons, j );
00141             if (muonref->track().isNonnull()) 
00142               if( muonref->track() == trackRef && muonref->isGlobalMuon()){
00143                 isMuCandidate=true;
00144                 //cout<<" SAVING TRACK "<<endl;
00145                 break;
00146               }
00147           }
00148         }
00149 
00150         if(!isMuCandidate)
00151           {
00152             continue;     
00153           }
00154         
00155       }
00156            
00157       // find the pre-id kf track
00158       bool preId = false;
00159       if(foundgsf) {
00160         for (unsigned int igsf=0; igsf<gsftracks.size();igsf++) {
00161           GsfTrackRef gsfTrackRef(gsftrackcoll, igsf);
00162           if (gsfTrackRef->seedRef().isNull()) continue;
00163           ElectronSeedRef ElSeedRef= gsfTrackRef->extra()->seedRef().castTo<ElectronSeedRef>();
00164           if (ElSeedRef->ctfTrack().isNonnull()) {
00165             if(ElSeedRef->ctfTrack() == trackRef) preId = true;
00166           }
00167         }
00168       }
00169       if(preId) {
00170         // Set PFRecTrack of type KF_ElCAND
00171         reco::PFRecTrack pftrack( trackRef->charge(), 
00172                                   reco::PFRecTrack::KF_ELCAND, 
00173                                   i, trackRef );
00174         
00175         bool valid = false;
00176         if(trajinev_) {
00177           valid = pfTransformer_->addPoints( pftrack, *trackRef, Tj[i]);
00178         }
00179         else {
00180           Trajectory FakeTraj;
00181           valid = pfTransformer_->addPoints( pftrack, *trackRef, FakeTraj);
00182         }
00183         if(valid) {
00184           //calculate STIP
00185           double stip=-999;
00186           const reco::PFTrajectoryPoint& atECAL=pftrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance);
00187           if(atECAL.isValid()) //if track extrapolates to ECAL
00188             {
00189               GlobalVector direction(pftrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().x(),
00190                                      pftrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().y(), 
00191                                      pftrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().z());
00192               stip = IPTools::signedTransverseImpactParameter(thebuilder.build(*trackRef), direction, *pv).second.significance();
00193             }
00194           pftrack.setSTIP(stip);
00195           PfTrColl->push_back(pftrack);
00196         }               
00197       }
00198       else {
00199         reco::PFRecTrack pftrack( trackRef->charge(), 
00200                                   reco::PFRecTrack::KF, 
00201                                   i, trackRef );
00202         bool valid = false;
00203         if(trajinev_) {
00204           valid = pfTransformer_->addPoints( pftrack, *trackRef, Tj[i]);
00205         }
00206         else {
00207           Trajectory FakeTraj;
00208           valid = pfTransformer_->addPoints( pftrack, *trackRef, FakeTraj);
00209         }
00210         
00211         if(valid) {
00212           double stip=-999;
00213           const reco::PFTrajectoryPoint& atECAL=pftrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance);
00214           if(atECAL.isValid())
00215             {
00216               GlobalVector direction(pftrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().x(),
00217                                      pftrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().y(), 
00218                                      pftrack.extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance).position().z());
00219               stip = IPTools::signedTransverseImpactParameter(thebuilder.build(*trackRef), direction, *pv).second.significance();
00220             }
00221           pftrack.setSTIP(stip);
00222           PfTrColl->push_back(pftrack);
00223         }
00224       }
00225     }
00226   }
00227   iEvent.put(PfTrColl);
00228 }
00229 
00230 // ------------ method called once each job just before starting event loop  ------------
00231 void 
00232 PFTrackProducer::beginRun(edm::Run& run,
00233                           const EventSetup& iSetup)
00234 {
00235   ESHandle<MagneticField> magneticField;
00236   iSetup.get<IdealMagneticFieldRecord>().get(magneticField);
00237   pfTransformer_= new PFTrackTransformer(math::XYZVector(magneticField->inTesla(GlobalPoint(0,0,0))));
00238   if(!trajinev_)
00239     pfTransformer_->OnlyProp();
00240 }
00241 
00242 // ------------ method called once each job just after ending the event loop  ------------
00243 void 
00244 PFTrackProducer::endRun() {
00245   delete pfTransformer_;
00246 }