CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/RecoTauTag/HLTProducers/src/TauRegionalPixelSeedGenerator.h

Go to the documentation of this file.
00001 #ifndef TauRegionalPixelSeedGenerator_h
00002 #define TauRegionalPixelSeedGenerator_h
00003 
00004 //
00005 // Class:           TauRegionalPixelSeedGenerator
00006 
00007 
00008 #include "FWCore/Framework/interface/EDProducer.h"
00009 #include "FWCore/Framework/interface/Event.h"
00010 #include "FWCore/Framework/interface/EventSetup.h"
00011 #include "FWCore/Utilities/interface/InputTag.h"
00012 #include "DataFormats/Common/interface/Handle.h"
00013 #include "DataFormats/VertexReco/interface/Vertex.h"
00014 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00015 #include "DataFormats/TrackReco/interface/Track.h"
00016 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00017 #include "DataFormats/Math/interface/Vector3D.h"
00018 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
00019 #include "RecoTracker/TkTrackingRegions/interface/GlobalTrackingRegion.h"
00020 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
00021 // Math
00022 #include "Math/GenVector/VectorUtil.h"
00023 #include "Math/GenVector/PxPyPzE4D.h"
00024 
00025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00026 
00027 #include "DataFormats/Common/interface/Ref.h"
00028 #include "DataFormats/JetReco/interface/Jet.h"
00029 
00030 
00031 class TauRegionalPixelSeedGenerator : public TrackingRegionProducer {
00032   public:
00033     
00034     explicit TauRegionalPixelSeedGenerator(const edm::ParameterSet& conf_){
00035       edm::LogInfo ("TauRegionalPixelSeedGenerator") << "Enter the TauRegionalPixelSeedGenerator";
00036 
00037       edm::ParameterSet regionPSet = conf_.getParameter<edm::ParameterSet>("RegionPSet");
00038 
00039       m_ptMin        = regionPSet.getParameter<double>("ptMin");
00040       m_originRadius = regionPSet.getParameter<double>("originRadius");
00041       m_halfLength   = regionPSet.getParameter<double>("originHalfLength");
00042       m_deltaEta     = regionPSet.getParameter<double>("deltaEtaRegion");
00043       m_deltaPhi     = regionPSet.getParameter<double>("deltaPhiRegion");
00044       m_jetSrc       = regionPSet.getParameter<edm::InputTag>("JetSrc");
00045       m_vertexSrc    = regionPSet.getParameter<edm::InputTag>("vertexSrc");
00046       if (regionPSet.exists("searchOpt")){
00047         m_searchOpt    = regionPSet.getParameter<bool>("searchOpt");
00048       }
00049       else{
00050         m_searchOpt = false;
00051       }
00052       m_measurementTracker ="";
00053       m_howToUseMeasurementTracker=0;
00054       if (regionPSet.exists("measurementTrackerName")){
00055         m_measurementTracker = regionPSet.getParameter<std::string>("measurementTrackerName");
00056         if (regionPSet.exists("howToUseMeasurementTracker")){
00057           m_howToUseMeasurementTracker = regionPSet.getParameter<double>("howToUseMeasurementTracker");
00058         }
00059       }
00060     }
00061   
00062     virtual ~TauRegionalPixelSeedGenerator() {}
00063     
00064 
00065     virtual std::vector<TrackingRegion* > regions(const edm::Event& e, const edm::EventSetup& es) const {
00066       std::vector<TrackingRegion* > result;
00067 
00068       //      double originZ;
00069       double deltaZVertex, deltaRho;
00070         GlobalPoint vertex;
00071       // get the primary vertex
00072       edm::Handle<reco::VertexCollection> h_vertices;
00073       e.getByLabel(m_vertexSrc, h_vertices);
00074       const reco::VertexCollection & vertices = * h_vertices;
00075       if (not vertices.empty()) {
00076 //        originZ      = vertices.front().z();
00077         GlobalPoint myTmp(vertices.at(0).position().x(),vertices.at(0).position().y(), vertices.at(0).position().z());
00078           vertex = myTmp;
00079           deltaZVertex = m_halfLength;
00080           deltaRho = m_originRadius;
00081       } else {
00082   //      originZ      =  0.;
00083           GlobalPoint myTmp(0.,0.,0.);
00084           vertex = myTmp;
00085           deltaRho = 1.;
00086          deltaZVertex = 15.;
00087       }
00088       
00089       // get the jet direction
00090       edm::Handle<edm::View<reco::Candidate> > h_jets;
00091       e.getByLabel(m_jetSrc, h_jets);
00092       
00093       for (unsigned int iJet =0; iJet < h_jets->size(); ++iJet)
00094         {
00095           const reco::Candidate & myJet = (*h_jets)[iJet];
00096           GlobalVector jetVector(myJet.momentum().x(),myJet.momentum().y(),myJet.momentum().z());
00097 //          GlobalPoint  vertex(0, 0, originZ);
00098           RectangularEtaPhiTrackingRegion* etaphiRegion = new RectangularEtaPhiTrackingRegion( jetVector,
00099                                                                                                vertex,
00100                                                                                                m_ptMin,
00101                                                                                                deltaRho,
00102                                                                                                deltaZVertex,
00103                                                                                                m_deltaEta,
00104                                                                                                m_deltaPhi,
00105                                                                                                m_howToUseMeasurementTracker,
00106                                                                                                true,
00107                                                                                                m_measurementTracker,
00108                                                                                                m_searchOpt);
00109           result.push_back(etaphiRegion);
00110       }
00111 
00112       return result;
00113     }
00114   
00115  private:
00116   edm::ParameterSet conf_;
00117 
00118   float m_ptMin;
00119   float m_originRadius;
00120   float m_halfLength;
00121   float m_deltaEta;
00122   float m_deltaPhi;
00123   edm::InputTag m_jetSrc;
00124   edm::InputTag m_vertexSrc;
00125   std::string m_measurementTracker;
00126   double m_howToUseMeasurementTracker;
00127   bool m_searchOpt;
00128 };
00129 
00130 #endif