Go to the documentation of this file.00001 #ifndef TrackingRegionsFromBeamSpotAndL2Tau_h
00002 #define TrackingRegionsFromBeamSpotAndL2Tau_h
00003
00004
00005
00006
00007
00008
00009
00010 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
00011 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
00012
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/Utilities/interface/InputTag.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "DataFormats/Common/interface/Handle.h"
00018 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00019 #include "DataFormats/Candidate/interface/Candidate.h"
00020
00024 class TrackingRegionsFromBeamSpotAndL2Tau : public TrackingRegionProducer
00025 {
00026 public:
00027
00028 explicit TrackingRegionsFromBeamSpotAndL2Tau(const edm::ParameterSet& conf)
00029 {
00030 edm::LogInfo ("TrackingRegionsFromBeamSpotAndL2Tau") << "Enter the TrackingRegionsFromBeamSpotAndL2Tau";
00031
00032 edm::ParameterSet regionPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
00033
00034 m_ptMin = regionPSet.getParameter<double>("ptMin");
00035 m_originRadius = regionPSet.getParameter<double>("originRadius");
00036 m_originHalfLength = regionPSet.getParameter<double>("originHalfLength");
00037 m_deltaEta = regionPSet.getParameter<double>("deltaEta");
00038 m_deltaPhi = regionPSet.getParameter<double>("deltaPhi");
00039 m_jetSrc = regionPSet.getParameter<edm::InputTag>("JetSrc");
00040 m_jetMinPt = regionPSet.getParameter<double>("JetMinPt");
00041 m_jetMaxEta = regionPSet.getParameter<double>("JetMaxEta");
00042 m_jetMaxN = regionPSet.getParameter<int>("JetMaxN");
00043 m_beamSpotTag = regionPSet.getParameter<edm::InputTag>("beamSpot");
00044 m_precise = regionPSet.getParameter<bool>("precise");
00045
00046 if (regionPSet.exists("searchOpt")) m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
00047 else m_searchOpt = false;
00048
00049 m_measurementTrackerName ="";
00050 m_whereToUseMeasurementTracker=0;
00051 if (regionPSet.exists("measurementTrackerName"))
00052 {
00053 m_measurementTrackerName = regionPSet.getParameter<std::string>("measurementTrackerName");
00054 if (regionPSet.exists("whereToUseMeasurementTracker"))
00055 m_whereToUseMeasurementTracker = regionPSet.getParameter<double>("whereToUseMeasurementTracker");
00056 }
00057 }
00058
00059 virtual ~TrackingRegionsFromBeamSpotAndL2Tau() {}
00060
00061
00062 virtual std::vector<TrackingRegion* > regions(const edm::Event& e, const edm::EventSetup& es) const
00063 {
00064 std::vector<TrackingRegion* > result;
00065
00066
00067 edm::Handle<reco::BeamSpot> bsHandle;
00068 e.getByLabel( m_beamSpotTag, bsHandle);
00069 if(!bsHandle.isValid()) return result;
00070 const reco::BeamSpot & bs = *bsHandle;
00071 GlobalPoint origin(bs.x0(), bs.y0(), bs.z0());
00072
00073
00074 edm::Handle< reco::CandidateView > objects;
00075 e.getByLabel( m_jetSrc, objects );
00076 size_t n_objects = objects->size();
00077 if (n_objects == 0) return result;
00078
00079
00080
00081
00082 int n_regions = 0;
00083 for (size_t i =0; i < n_objects && n_regions < m_jetMaxN; ++i)
00084 {
00085 const reco::Candidate & jet = (*objects)[i];
00086 if ( jet.pt() < m_jetMinPt || std::abs(jet.eta()) > m_jetMaxEta ) continue;
00087
00088 GlobalVector direction(jet.momentum().x(), jet.momentum().y(), jet.momentum().z());
00089
00090 RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion(
00091 direction,
00092 origin,
00093 m_ptMin,
00094 m_originRadius,
00095 m_originHalfLength,
00096 m_deltaEta,
00097 m_deltaPhi,
00098 m_whereToUseMeasurementTracker,
00099 m_precise,
00100 m_measurementTrackerName,
00101 m_searchOpt
00102 );
00103 result.push_back(etaphiRegion);
00104 ++n_regions;
00105 }
00106
00107 return result;
00108 }
00109
00110 private:
00111
00112 float m_ptMin;
00113 float m_originRadius;
00114 float m_originHalfLength;
00115 float m_deltaEta;
00116 float m_deltaPhi;
00117 edm::InputTag m_jetSrc;
00118 float m_jetMinPt;
00119 float m_jetMaxEta;
00120 int m_jetMaxN;
00121 std::string m_measurementTrackerName;
00122 float m_whereToUseMeasurementTracker;
00123 bool m_searchOpt;
00124 edm::InputTag m_beamSpotTag;
00125 bool m_precise;
00126 };
00127
00128 #endif