CMS 3D CMS Logo

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