CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/Calibration/HcalIsolatedTrackReco/src/HITRegionalPixelSeedGenerator.h

Go to the documentation of this file.
00001 #ifndef HITRegionalPixelSeedGenerator_h
00002 #define HITRegionalPixelSeedGenerator_h
00003 
00004 //
00005 // Class:           HITRegionalPixelSeedGenerator
00006 
00007 
00008 #include "FWCore/Framework/interface/EDProducer.h"
00009 #include "FWCore/Framework/interface/Event.h"
00010 #include "DataFormats/Common/interface/Handle.h"
00011 #include "FWCore/Framework/interface/EventSetup.h"
00012 #include "FWCore/Utilities/interface/InputTag.h"
00013 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
00014 #include "RecoTracker/TkTrackingRegions/interface/GlobalTrackingRegion.h"
00015 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
00016 #include "DataFormats/VertexReco/interface/Vertex.h"
00017 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00018 #include "DataFormats/TrackReco/interface/Track.h"
00019 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00020 #include "DataFormats/Math/interface/Vector3D.h"
00021 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
00022 #include "DataFormats/L1Trigger/interface/L1JetParticleFwd.h"
00023 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00024 
00025 // Math
00026 #include "Math/GenVector/VectorUtil.h"
00027 #include "Math/GenVector/PxPyPzE4D.h"
00028 
00029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00030 
00031 #include "DataFormats/Common/interface/Ref.h"
00032 #include "DataFormats/JetReco/interface/Jet.h"
00033 
00034 class HITRegionalPixelSeedGenerator : public TrackingRegionProducer {
00035  public:
00036   
00037   explicit HITRegionalPixelSeedGenerator(const edm::ParameterSet& conf_)
00038   {
00039     edm::LogInfo ("HITRegionalPixelSeedGenerator")<<"Enter the HITRegionalPixelSeedGenerator";
00040     
00041     edm::ParameterSet regionPSet = conf_.getParameter<edm::ParameterSet>("RegionPSet");
00042     
00043     ptmin=regionPSet.getParameter<double>("ptMin");
00044     originradius=regionPSet.getParameter<double>("originRadius");
00045     halflength=regionPSet.getParameter<double>("originHalfLength");
00046     vertexSrc=regionPSet.getParameter<std::string>("vertexSrc");
00047     etaCenter_=regionPSet.getParameter<double>("etaCenter");
00048     phiCenter_=regionPSet.getParameter<double>("phiCenter");
00049     deltaTrackEta = regionPSet.getParameter<double>("deltaEtaTrackRegion");
00050     deltaTrackPhi = regionPSet.getParameter<double>("deltaPhiTrackRegion");
00051     deltaL1JetEta = regionPSet.getParameter<double>("deltaEtaL1JetRegion");
00052     deltaL1JetPhi = regionPSet.getParameter<double>("deltaPhiL1JetRegion");
00053     tracksrc_ = regionPSet.getParameter<edm::InputTag>("trackSrc");
00054     isoTracksrc_ = regionPSet.getParameter<edm::InputTag>("isoTrackSrc");
00055     l1jetsrc_ = regionPSet.getParameter<edm::InputTag>("l1tjetSrc");
00056     usejets_ = regionPSet.getParameter<bool>("useL1Jets");
00057     usetracks_ = regionPSet.getParameter<bool>("useTracks");
00058     useIsoTracks_ = regionPSet.getParameter<bool>("useIsoTracks");
00059     fixedReg_ = regionPSet.getParameter<bool>("fixedReg");
00060   }
00061   
00062   virtual ~HITRegionalPixelSeedGenerator() {}
00063   
00064   
00065   virtual std::vector<TrackingRegion* > regions(const edm::Event& e, const edm::EventSetup& es) const 
00066     {
00067       std::vector<TrackingRegion* > result;
00068       float originz =0.;
00069       
00070       double deltaZVertex =  halflength;
00071       
00072       
00073       
00074       if (usetracks_)
00075         {
00076           edm::Handle<reco::TrackCollection> tracks;
00077           e.getByLabel(tracksrc_, tracks);
00078           
00079           edm::Handle<reco::VertexCollection> vertices;
00080           e.getByLabel(vertexSrc,vertices);
00081           const reco::VertexCollection vertCollection = *(vertices.product());
00082           reco::VertexCollection::const_iterator ci = vertCollection.begin();
00083           
00084           if(vertCollection.size() > 0) 
00085             {
00086               originz = ci->z();
00087             }
00088           else
00089             {
00090               deltaZVertex = 15.;
00091             }
00092       
00093           GlobalVector globalVector(0,0,1);
00094           if(tracks->size() == 0) return result;
00095           
00096           reco::TrackCollection::const_iterator itr = tracks->begin();
00097           for(;itr != tracks->end();itr++)
00098             {
00099               
00100               GlobalVector ptrVec((itr)->px(),(itr)->py(),(itr)->pz());
00101               globalVector = ptrVec;
00102               
00103               
00104               RectangularEtaPhiTrackingRegion* etaphiRegion = new  RectangularEtaPhiTrackingRegion(globalVector,
00105                                                                                                    GlobalPoint(0,0,originz), 
00106                                                                                                    ptmin,
00107                                                                                                    originradius,
00108                                                                                                    deltaZVertex,
00109                                                                                                    deltaTrackEta,
00110                                                                                                    deltaTrackPhi);
00111               result.push_back(etaphiRegion);
00112               
00113             }
00114         }
00115       
00116       if (useIsoTracks_)
00117         {
00118           edm::Handle<trigger::TriggerFilterObjectWithRefs> isotracks;
00119           e.getByLabel(isoTracksrc_, isotracks);
00120 
00121           std::vector< edm::Ref<reco::IsolatedPixelTrackCandidateCollection> > isoPixTrackRefs;
00122         
00123           isotracks->getObjects(trigger::TriggerTrack, isoPixTrackRefs);
00124           
00125           edm::Handle<reco::VertexCollection> vertices;
00126           e.getByLabel(vertexSrc,vertices);
00127           const reco::VertexCollection vertCollection = *(vertices.product());
00128           reco::VertexCollection::const_iterator ci = vertCollection.begin();
00129 
00130           if(vertCollection.size() > 0) 
00131             {
00132               originz = ci->z();
00133             }
00134           else
00135             {
00136               deltaZVertex = 15.;
00137             }
00138           
00139           GlobalVector globalVector(0,0,1);
00140           if(isoPixTrackRefs.size() == 0) return result;
00141           
00142           for(uint32_t p=0; p<isoPixTrackRefs.size(); p++)
00143             {
00144               GlobalVector ptrVec((isoPixTrackRefs[p]->track())->px(),(isoPixTrackRefs[p]->track())->py(),(isoPixTrackRefs[p]->track())->pz());
00145               globalVector = ptrVec;
00146               
00147               RectangularEtaPhiTrackingRegion* etaphiRegion = new  RectangularEtaPhiTrackingRegion(globalVector,
00148                                                                                                    GlobalPoint(0,0,originz),
00149                                                                                                    ptmin,
00150                                                                                                    originradius,
00151                                                                                                    deltaZVertex,
00152                                                                                                    deltaTrackEta,
00153                                                                                                    deltaTrackPhi);
00154               result.push_back(etaphiRegion);
00155             }
00156         }
00157       
00158       if (usejets_)
00159         {
00160           edm::Handle<l1extra::L1JetParticleCollection> jets;
00161           e.getByLabel(l1jetsrc_, jets);
00162           
00163           edm::Handle<reco::VertexCollection> vertices;
00164           e.getByLabel(vertexSrc,vertices);
00165           const reco::VertexCollection vertCollection = *(vertices.product());
00166           reco::VertexCollection::const_iterator ci = vertCollection.begin();
00167           if(vertCollection.size() > 0) 
00168             {
00169               originz = ci->z();
00170             }
00171           else
00172             {
00173               deltaZVertex = 15.;
00174             }
00175 
00176           GlobalVector globalVector(0,0,1);
00177           if(jets->size() == 0) return result;
00178           
00179           for (l1extra::L1JetParticleCollection::const_iterator iJet = jets->begin(); iJet != jets->end(); iJet++) 
00180             {
00181               GlobalVector jetVector(iJet->p4().x(), iJet->p4().y(), iJet->p4().z());
00182               GlobalPoint  vertex(0, 0, originz);
00183               
00184               RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion( jetVector,
00185                                                                                                    vertex,
00186                                                                                                    ptmin,
00187                                                                                                    originradius,
00188                                                                                                    deltaZVertex,
00189                                                                                                    deltaL1JetEta,
00190                                                                                                    deltaL1JetPhi );
00191               result.push_back(etaphiRegion);
00192             }
00193         }
00194       if (fixedReg_)
00195         {
00196           GlobalVector fixedVector(cos(phiCenter_)*sin(2*atan(exp(-etaCenter_))), sin(phiCenter_)*sin(2*atan(exp(-etaCenter_))), cos(2*atan(exp(-etaCenter_))));
00197           GlobalPoint  vertex(0, 0, originz);
00198           
00199           edm::Handle<reco::VertexCollection> vertices;
00200           e.getByLabel(vertexSrc,vertices);
00201           const reco::VertexCollection vertCollection = *(vertices.product());
00202           reco::VertexCollection::const_iterator ci = vertCollection.begin();
00203           if(vertCollection.size() > 0) 
00204             {
00205               originz = ci->z();
00206             }
00207           else
00208             {
00209               deltaZVertex = 15.;
00210             }
00211           
00212           RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion( fixedVector,
00213                                                                                                vertex,
00214                                                                                                ptmin,
00215                                                                                                originradius,
00216                                                                                                deltaZVertex,
00217                                                                                                deltaL1JetEta,
00218                                                                                                deltaL1JetPhi );
00219           result.push_back(etaphiRegion);
00220         }
00221       
00222       return result;
00223     }
00224   
00225   
00226  private:
00227   edm::ParameterSet conf_;
00228   
00229   float ptmin;
00230   float originradius;
00231   float halflength;
00232   double etaCenter_;
00233   double phiCenter_;
00234   float deltaTrackEta;
00235   float deltaTrackPhi;
00236   float deltaL1JetEta;
00237   float deltaL1JetPhi;
00238   edm::InputTag tracksrc_;
00239   edm::InputTag isoTracksrc_;
00240   std::string vertexSrc;
00241   edm::InputTag l1jetsrc_;
00242   bool usejets_;
00243   bool usetracks_;
00244   bool fixedReg_;
00245   bool useIsoTracks_;
00246 
00247 };
00248 
00249 #endif
00250 
00251 
00252