CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoTracker/TkTrackingRegions/plugins/GlobalTrackingRegionWithVerticesProducer.h

Go to the documentation of this file.
00001 #ifndef RecoTracker_TkTrackingRegions_GlobalTrackingRegionWithVerticesProducer_H 
00002 #define RecoTracker_TkTrackingRegions_GlobalTrackingRegionWithVerticesProducer_H
00003 
00004 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
00005 #include "RecoTracker/TkTrackingRegions/interface/GlobalTrackingRegion.h"
00006 #include "FWCore/Framework/interface/Event.h"
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008 
00009 #include "DataFormats/VertexReco/interface/Vertex.h"
00010 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00011 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00012 #include "FWCore/Utilities/interface/InputTag.h"
00013 
00014 class GlobalTrackingRegionWithVerticesProducer : public TrackingRegionProducer
00015 {
00016 public:
00017 
00018   GlobalTrackingRegionWithVerticesProducer(const edm::ParameterSet& cfg)
00019   { 
00020     edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
00021 
00022     thePtMin            = regionPSet.getParameter<double>("ptMin");
00023     theOriginRadius     = regionPSet.getParameter<double>("originRadius");
00024     theNSigmaZ          = regionPSet.getParameter<double>("nSigmaZ");
00025     theBeamSpotTag      = regionPSet.getParameter<edm::InputTag>("beamSpot");
00026     thePrecise          = regionPSet.getParameter<bool>("precise"); 
00027 
00028     theSigmaZVertex     = regionPSet.getParameter<double>("sigmaZVertex");
00029     theFixedError       = regionPSet.getParameter<double>("fixedError");
00030 
00031     theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices");
00032     theUseFixedError    = regionPSet.getParameter<bool>("useFixedError");
00033     vertexCollName      = regionPSet.getParameter<edm::InputTag>("VertexCollection");
00034   }   
00035 
00036   virtual ~GlobalTrackingRegionWithVerticesProducer(){}
00037 
00038   virtual std::vector<TrackingRegion* > regions
00039     (const edm::Event& ev, const edm::EventSetup&) const
00040   {
00041     std::vector<TrackingRegion* > result;
00042 
00043     GlobalPoint theOrigin;
00044     edm::Handle<reco::BeamSpot> bsHandle;
00045     ev.getByLabel( theBeamSpotTag, bsHandle);
00046     double bsSigmaZ;
00047     if(bsHandle.isValid()) {
00048       const reco::BeamSpot & bs = *bsHandle; 
00049       bsSigmaZ = theNSigmaZ*bs.sigmaZ();
00050       theOrigin = GlobalPoint(bs.x0(), bs.y0(), bs.z0());
00051     }else{
00052       throw cms::Exception("Seeding") << "ERROR: input beamSpot is not valid in GlobalTrackingRegionWithVertices";
00053     }
00054 
00055     if(theUseFoundVertices)
00056     {
00057       edm::Handle<reco::VertexCollection> vertexCollection;
00058       ev.getByLabel(vertexCollName,vertexCollection);
00059 
00060       for(reco::VertexCollection::const_iterator iV=vertexCollection->begin(); iV != vertexCollection->end() ; iV++) {
00061           if (iV->isFake() || !iV->isValid()) continue;
00062           GlobalPoint theOrigin_       = GlobalPoint(iV->x(),iV->y(),iV->z());
00063           double theOriginHalfLength_ = (theUseFixedError ? theFixedError : (iV->zError())*theSigmaZVertex); 
00064           result.push_back( new GlobalTrackingRegion(thePtMin, theOrigin_, theOriginRadius, theOriginHalfLength_, thePrecise) );
00065       }
00066       
00067       if (result.empty()) {
00068         result.push_back( new GlobalTrackingRegion(thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise) );
00069       }
00070     }
00071     else
00072     {
00073       result.push_back(
00074         new GlobalTrackingRegion(thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise) );
00075     }
00076 
00077     return result;
00078   }
00079 
00080 private:
00081   double thePtMin; 
00082   double theOriginRadius; 
00083   double theNSigmaZ;
00084   edm::InputTag theBeamSpotTag;
00085 
00086   double theSigmaZVertex;
00087   double theFixedError;
00088   bool thePrecise;
00089   
00090   bool theUseFoundVertices;
00091   bool theUseFixedError;
00092   edm::InputTag vertexCollName;
00093 
00094 
00095 };
00096 
00097 #endif