CMS 3D CMS Logo

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/ParameterSet/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     thePtMin              = regionPSet.getParameter<double>("ptMin");
00027     theOriginRadius       = regionPSet.getParameter<double>("originRadius");
00028     theOriginHalfLength   = regionPSet.getParameter<double>("originHalfLength");
00029     theOriginZPos  = regionPSet.getParameter<double>("vertexZDefault");
00030 
00031     theDeltaEta = regionPSet.getParameter<double>("deltaEtaRegion");
00032     theDeltaPhi =  regionPSet.getParameter<double>("deltaPhiRegion");
00033   }   
00034 
00035   virtual ~L3MumuTrackingRegion(){}
00036 
00037   virtual std::vector<TrackingRegion* > regions(const edm::Event& ev, 
00038       const edm::EventSetup& es) const {
00039 
00040     std::vector<TrackingRegion* > result;
00041 
00042     // optional constraint for vertex
00043     // get highest Pt pixel vertex (if existing)
00044     double deltaZVertex =  theOriginHalfLength;
00045     double originz = theOriginZPos;
00046     if (theVertexSrc.length()>1) {
00047       edm::Handle<reco::VertexCollection> vertices;
00048       ev.getByLabel(theVertexSrc,vertices);
00049       const reco::VertexCollection vertCollection = *(vertices.product());
00050       reco::VertexCollection::const_iterator ci = vertCollection.begin();
00051       if (vertCollection.size()>0) {
00052             originz = ci->z();
00053       } else {
00054             originz = theOriginZPos;
00055             deltaZVertex = 15.;
00056       }
00057     }
00058 
00059     edm::Handle<reco::TrackCollection> trks;
00060     ev.getByLabel(theInputTrkSrc, trks);
00061 
00062     for(reco::TrackCollection::const_iterator iTrk = trks->begin();iTrk != trks->end();iTrk++) {
00063       GlobalVector dirVector((iTrk)->px(),(iTrk)->py(),(iTrk)->pz());
00064       result.push_back( 
00065           new RectangularEtaPhiTrackingRegion( dirVector, GlobalPoint(0,0,float(originz)), 
00066           thePtMin, theOriginRadius, deltaZVertex, theDeltaEta, theDeltaPhi) );
00067     }
00068     return result;
00069   }
00070 
00071 private:
00072 
00073   std::string theVertexSrc;
00074   edm::InputTag theInputTrkSrc;
00075 
00076   double thePtMin; 
00077   double theOriginRadius; 
00078   double theOriginHalfLength; 
00079   bool   theVertexZconstrained;
00080   double theOriginZPos;
00081 
00082   double theDeltaEta; 
00083   double theDeltaPhi;
00084 };
00085 
00086 #endif 
00087 

Generated on Tue Jun 9 17:37:39 2009 for CMSSW by  doxygen 1.5.4