CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/HLTrigger/btau/src/L3MumuTrackingRegion.h

Go to the documentation of this file.
00001 #ifndef HLTrigger_btau_L3MumuTrackingRegion_H 
00002 #define HLTrigger_btau_L3MumuTrackingRegion_H 
00003 
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
00006 #include "RecoTracker/TkTrackingRegions/interface/GlobalTrackingRegion.h"
00007 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
00008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00009 #include "FWCore/Utilities/interface/InputTag.h"
00010 #include "DataFormats/VertexReco/interface/Vertex.h"
00011 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00012 #include "DataFormats/TrackReco/interface/Track.h"
00013 
00014 
00015 class L3MumuTrackingRegion : public TrackingRegionProducer {
00016 
00017 public:
00018 
00019   L3MumuTrackingRegion(const edm::ParameterSet& cfg) { 
00020 
00021     edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
00022 
00023     theVertexSrc   = regionPSet.getParameter<std::string>("vertexSrc");
00024     theInputTrkSrc = regionPSet.getParameter<edm::InputTag>("TrkSrc");
00025 
00026     useVtxTks = regionPSet.getParameter<bool>("UseVtxTks");
00027 
00028     thePtMin              = regionPSet.getParameter<double>("ptMin");
00029     theOriginRadius       = regionPSet.getParameter<double>("originRadius");
00030     theOriginHalfLength   = regionPSet.getParameter<double>("originHalfLength");
00031     theOriginZPos  = regionPSet.getParameter<double>("vertexZDefault");
00032 
00033     theDeltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
00034     theDeltaPhi =  regionPSet.getParameter<double>("deltaPhiRegion");
00035     if (regionPSet.exists("searchOpt")){
00036       m_searchOpt    = regionPSet.getParameter<bool>("searchOpt");
00037     }
00038     else{
00039       m_searchOpt = false;
00040     }
00041     m_measurementTracker ="";
00042     m_howToUseMeasurementTracker=0;
00043     if (regionPSet.exists("measurementTrackerName")){
00044       m_measurementTracker = regionPSet.getParameter<std::string>("measurementTrackerName");
00045       if (regionPSet.exists("howToUseMeasurementTracker")){
00046         m_howToUseMeasurementTracker = regionPSet.getParameter<double>("howToUseMeasurementTracker");
00047       }
00048     }
00049   }   
00050 
00051   virtual ~L3MumuTrackingRegion(){}
00052 
00053   virtual std::vector<TrackingRegion* > regions(const edm::Event& ev, 
00054       const edm::EventSetup& es) const {
00055 
00056     std::vector<TrackingRegion* > result;
00057 
00058     // optional constraint for vertex
00059     // get highest Pt pixel vertex (if existing)
00060     double deltaZVertex =  theOriginHalfLength;
00061     double originz = theOriginZPos;
00062     if (theVertexSrc.length()>1) {
00063       edm::Handle<reco::VertexCollection> vertices;
00064       ev.getByLabel(theVertexSrc,vertices);
00065       const reco::VertexCollection vertCollection = *(vertices.product());
00066       reco::VertexCollection::const_iterator ci = vertCollection.begin();
00067       if (vertCollection.size()>0) {
00068             originz = ci->z();
00069       } else {
00070             originz = theOriginZPos;
00071             deltaZVertex = 15.;
00072       }
00073       if (useVtxTks) {
00074 
00075         for(ci=vertCollection.begin();ci!=vertCollection.end();ci++)
00076           for (reco::Vertex::trackRef_iterator trackIt =  ci->tracks_begin();trackIt !=  ci->tracks_end();trackIt++){
00077             reco::TrackRef iTrk =  (*trackIt).castTo<reco::TrackRef>() ;
00078             GlobalVector dirVector((iTrk)->px(),(iTrk)->py(),(iTrk)->pz());
00079             result.push_back(
00080                              new RectangularEtaPhiTrackingRegion( dirVector, GlobalPoint(0,0,float(ci->z())),
00081                                                                   thePtMin, theOriginRadius, deltaZVertex, theDeltaEta, theDeltaPhi,
00082                                                                   m_howToUseMeasurementTracker,
00083                                                                   true,
00084                                                                   m_measurementTracker,
00085                                                                   m_searchOpt) );
00086           }
00087         return result;
00088       }
00089     }
00090 
00091     edm::Handle<reco::TrackCollection> trks;
00092     ev.getByLabel(theInputTrkSrc, trks);
00093 
00094     for(reco::TrackCollection::const_iterator iTrk = trks->begin();iTrk != trks->end();iTrk++) {
00095       GlobalVector dirVector((iTrk)->px(),(iTrk)->py(),(iTrk)->pz());
00096       result.push_back( 
00097           new RectangularEtaPhiTrackingRegion( dirVector, GlobalPoint(0,0,float(originz)), 
00098                                                thePtMin, theOriginRadius, deltaZVertex, theDeltaEta, theDeltaPhi,
00099                                                m_howToUseMeasurementTracker,
00100                                                true,
00101                                                m_measurementTracker,
00102                                                m_searchOpt) );
00103     }
00104     return result;
00105   }
00106 
00107 private:
00108 
00109   std::string theVertexSrc;
00110   edm::InputTag theInputTrkSrc;
00111 
00112   bool useVtxTks;
00113 
00114   double thePtMin; 
00115   double theOriginRadius; 
00116   double theOriginHalfLength; 
00117   double theOriginZPos;
00118 
00119   double theDeltaEta; 
00120   double theDeltaPhi;
00121   std::string m_measurementTracker;
00122   double m_howToUseMeasurementTracker;
00123   bool m_searchOpt;
00124 };
00125 
00126 #endif 
00127