Go to the documentation of this file.00001 #ifndef CandidateSeededTrackingRegionsProducer_h
00002 #define CandidateSeededTrackingRegionsProducer_h
00003
00004
00005
00006 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
00007 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
00008
00009 #include "FWCore/Framework/interface/Event.h"
00010 #include "FWCore/Framework/interface/EventSetup.h"
00011 #include "FWCore/Utilities/interface/InputTag.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 #include "DataFormats/Common/interface/Handle.h"
00014 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00015 #include "DataFormats/VertexReco/interface/Vertex.h"
00016 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00017 #include "DataFormats/Candidate/interface/Candidate.h"
00018
00044 class CandidateSeededTrackingRegionsProducer : public TrackingRegionProducer
00045 {
00046 public:
00047
00048 typedef enum {BEAM_SPOT_FIXED, BEAM_SPOT_SIGMA, VERTICES_FIXED, VERTICES_SIGMA } Mode;
00049
00050 explicit CandidateSeededTrackingRegionsProducer(const edm::ParameterSet& conf)
00051 {
00052 edm::ParameterSet regPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
00053
00054
00055 std::string modeString = regPSet.getParameter<std::string>("mode");
00056 if (modeString == "BeamSpotFixed") m_mode = BEAM_SPOT_FIXED;
00057 else if (modeString == "BeamSpotSigma") m_mode = BEAM_SPOT_SIGMA;
00058 else if (modeString == "VerticesFixed") m_mode = VERTICES_FIXED;
00059 else if (modeString == "VerticesSigma") m_mode = VERTICES_SIGMA;
00060 else edm::LogError ("CandidateSeededTrackingRegionsProducer")<<"Unknown mode string: "<<modeString;
00061
00062
00063 m_input = regPSet.getParameter<edm::InputTag>("input");
00064 m_maxNRegions = regPSet.getParameter<int>("maxNRegions");
00065 m_beamSpot = regPSet.getParameter<edm::InputTag>("beamSpot");
00066 m_vertexCollection = edm::InputTag();
00067 m_maxNVertices = 1;
00068 if (m_mode == VERTICES_FIXED || m_mode == VERTICES_SIGMA)
00069 {
00070 m_vertexCollection = regPSet.getParameter<edm::InputTag>("vertexCollection");
00071 m_maxNVertices = regPSet.getParameter<int>("maxNVertices");
00072 }
00073
00074
00075 m_ptMin = regPSet.getParameter<double>("ptMin");
00076 m_originRadius = regPSet.getParameter<double>("originRadius");
00077 m_zErrorBeamSpot = regPSet.getParameter<double>("zErrorBeamSpot");
00078 m_deltaEta = regPSet.getParameter<double>("deltaEta");
00079 m_deltaPhi = regPSet.getParameter<double>("deltaPhi");
00080 m_precise = regPSet.getParameter<bool>("precise");
00081 m_measurementTrackerName = "";
00082 m_whereToUseMeasurementTracker = 0;
00083 if (regPSet.exists("measurementTrackerName"))
00084 {
00085 m_measurementTrackerName = regPSet.getParameter<std::string>("measurementTrackerName");
00086 if (regPSet.exists("whereToUseMeasurementTracker"))
00087 m_whereToUseMeasurementTracker = regPSet.getParameter<double>("whereToUseMeasurementTracker");
00088 }
00089 m_searchOpt = false;
00090 if (regPSet.exists("searchOpt")) m_searchOpt = regPSet.getParameter<bool>("searchOpt");
00091
00092
00093 if (m_mode == VERTICES_SIGMA) m_nSigmaZVertex = regPSet.getParameter<double>("nSigmaZVertex");
00094 if (m_mode == VERTICES_FIXED) m_zErrorVetex = regPSet.getParameter<double>("zErrorVetex");
00095 m_nSigmaZBeamSpot = -1.;
00096 if (m_mode == BEAM_SPOT_SIGMA)
00097 {
00098 m_nSigmaZBeamSpot = regPSet.getParameter<double>("nSigmaZBeamSpot");
00099 if (m_nSigmaZBeamSpot < 0.)
00100 edm::LogError ("CandidateSeededTrackingRegionsProducer")<<"nSigmaZBeamSpot must be positive for BeamSpotSigma mode!";
00101 }
00102 }
00103
00104 virtual ~CandidateSeededTrackingRegionsProducer() {}
00105
00106
00107 virtual std::vector<TrackingRegion* > regions(const edm::Event& e, const edm::EventSetup& es) const
00108 {
00109 std::vector<TrackingRegion* > result;
00110
00111
00112 edm::Handle< reco::CandidateView > objects;
00113 e.getByLabel( m_input, objects );
00114 size_t n_objects = objects->size();
00115 if (n_objects == 0) return result;
00116
00117
00118 edm::Handle< reco::BeamSpot > bs;
00119 e.getByLabel( m_beamSpot, bs );
00120 if( !bs.isValid() ) return result;
00121
00122
00123 GlobalPoint default_origin( bs->x0(), bs->y0(), bs->z0() );
00124
00125
00126 std::vector< std::pair< GlobalPoint, float > > origins;
00127
00128
00129 if (m_mode == BEAM_SPOT_FIXED || m_mode == BEAM_SPOT_SIGMA)
00130 {
00131 origins.push_back( std::make_pair(
00132 default_origin,
00133 (m_mode == BEAM_SPOT_FIXED) ? m_zErrorBeamSpot : m_nSigmaZBeamSpot*bs->sigmaZ()
00134 ));
00135 }
00136 else if (m_mode == VERTICES_FIXED || m_mode == VERTICES_SIGMA)
00137 {
00138 edm::Handle< reco::VertexCollection > vertices;
00139 e.getByLabel( m_vertexCollection, vertices );
00140 int n_vert = 0;
00141 for (reco::VertexCollection::const_iterator v = vertices->begin(); v != vertices->end() && n_vert < m_maxNVertices; ++v)
00142 {
00143 if ( v->isFake() || !v->isValid() ) continue;
00144
00145 origins.push_back( std::make_pair(
00146 GlobalPoint( v->x(), v->y(), v->z() ),
00147 (m_mode == VERTICES_FIXED) ? m_zErrorVetex : m_nSigmaZVertex*v->zError()
00148 ));
00149 ++n_vert;
00150 }
00151
00152 if (origins.empty())
00153 {
00154 origins.push_back( std::make_pair(
00155 default_origin,
00156 (m_nSigmaZBeamSpot > 0.) ? m_nSigmaZBeamSpot*bs->z0Error() : m_zErrorBeamSpot
00157 ));
00158 }
00159 }
00160
00161
00162
00163 int n_regions = 0;
00164 for(size_t i = 0; i < n_objects && n_regions < m_maxNRegions; ++i )
00165 {
00166 const reco::Candidate & object = (*objects)[i];
00167 GlobalVector direction( object.momentum().x(), object.momentum().y(), object.momentum().z() );
00168
00169 for (size_t j=0; j<origins.size() && n_regions < m_maxNRegions; ++j)
00170 {
00171 result.push_back( new RectangularEtaPhiTrackingRegion(
00172 direction,
00173 origins[j].first,
00174 m_ptMin,
00175 m_originRadius,
00176 origins[j].second,
00177 m_deltaEta,
00178 m_deltaPhi,
00179 m_whereToUseMeasurementTracker,
00180 m_precise,
00181 m_measurementTrackerName,
00182 m_searchOpt
00183 ));
00184 ++n_regions;
00185 }
00186 }
00187
00188 edm::LogInfo ("CandidateSeededTrackingRegionsProducer") << "produced "<<n_regions<<" regions";
00189
00190 return result;
00191 }
00192
00193 private:
00194
00195 Mode m_mode;
00196
00197 edm::InputTag m_input;
00198 int m_maxNRegions;
00199 edm::InputTag m_beamSpot;
00200 edm::InputTag m_vertexCollection;
00201 int m_maxNVertices;
00202
00203 float m_ptMin;
00204 float m_originRadius;
00205 float m_zErrorBeamSpot;
00206 float m_deltaEta;
00207 float m_deltaPhi;
00208 bool m_precise;
00209 std::string m_measurementTrackerName;
00210 float m_whereToUseMeasurementTracker;
00211 bool m_searchOpt;
00212
00213 float m_nSigmaZVertex;
00214 float m_zErrorVetex;
00215 float m_nSigmaZBeamSpot;
00216 };
00217
00218 #endif