Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:46

0001 #ifndef ElectronSeedGenerator_H
0002 #define ElectronSeedGenerator_H
0003 
0004 /** \class ElectronSeedGenerator
0005 
0006  * Class to generate the trajectory seed from two hits in
0007  *  the pixel detector which have been found compatible with
0008  *  an ECAL cluster.
0009  *
0010  * \author U.Berthon, C.Charlot, LLR Palaiseau
0011  *
0012  * \version   1st Version May 30, 2006
0013  *
0014  ************************************************************/
0015 
0016 #include "RecoEgamma/EgammaElectronAlgos/interface/PixelHitMatcher.h"
0017 
0018 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
0019 #include "TrackingTools/PatternTools/interface/TrajectoryStateUpdator.h"
0020 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
0021 
0022 #include "RecoTracker/TransientTrackingRecHit/interface/TSiPixelRecHit.h"
0023 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0024 
0025 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
0026 #include "DataFormats/EgammaReco/interface/ElectronSeedFwd.h"
0027 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
0028 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
0029 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0030 #include "DataFormats/VertexReco/interface/Vertex.h"
0031 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0032 
0033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0034 #include "FWCore/Framework/interface/ConsumesCollector.h"
0035 #include "FWCore/Framework/interface/EventSetup.h"
0036 #include "FWCore/Framework/interface/Event.h"
0037 #include "FWCore/Framework/interface/ESHandle.h"
0038 #include "FWCore/Framework/interface/ESWatcher.h"
0039 #include "FWCore/Framework/interface/EDConsumerBase.h"
0040 
0041 #include "TrackingTools/KalmanUpdators/interface/KFUpdator.h"
0042 
0043 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0044 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0045 
0046 class ElectronSeedGenerator {
0047 public:
0048   struct Tokens {
0049     edm::EDGetTokenT<std::vector<reco::Vertex> > token_vtx;
0050     edm::EDGetTokenT<reco::BeamSpot> token_bs;
0051   };
0052 
0053   typedef edm::OwnVector<TrackingRecHit> PRecHitContainer;
0054   typedef TransientTrackingRecHit::ConstRecHitPointer ConstRecHitPointer;
0055   typedef TransientTrackingRecHit::RecHitPointer RecHitPointer;
0056   typedef TransientTrackingRecHit::RecHitContainer RecHitContainer;
0057 
0058   ElectronSeedGenerator(const edm::ParameterSet&, const Tokens&, edm::ConsumesCollector&&);
0059 
0060   void setupES(const edm::EventSetup& setup);
0061   void run(edm::Event&,
0062            const reco::SuperClusterRefVector&,
0063            const std::vector<const TrajectorySeedCollection*>& seedsV,
0064            reco::ElectronSeedCollection&);
0065 
0066 private:
0067   void seedsFromThisCluster(edm::Ref<reco::SuperClusterCollection> seedCluster,
0068                             reco::BeamSpot const& beamSpot,
0069                             std::vector<reco::Vertex> const* vertices,
0070                             reco::ElectronSeedCollection& out);
0071 
0072   const bool dynamicPhiRoad_;
0073 
0074   const edm::EDGetTokenT<std::vector<reco::Vertex> > verticesTag_;
0075   const edm::EDGetTokenT<reco::BeamSpot> beamSpotTag_;
0076   const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> magFieldToken_;
0077   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> trackerGeometryToken_;
0078   edm::ESWatcher<IdealMagneticFieldRecord> magneticFieldWatcher_;
0079   edm::ESWatcher<TrackerDigiGeometryRecord> trackerGeometryWatcher_;
0080 
0081   const float lowPtThresh_;
0082   const float highPtThresh_;
0083   const float nSigmasDeltaZ1_;     // first z window size if not using the reco vertex
0084   const float deltaZ1WithVertex_;  // first z window size when using the reco vertex
0085   const float sizeWindowENeg_;
0086 
0087   const float deltaPhi1Low_;
0088   const float deltaPhi1High_;
0089 
0090   // so that deltaPhi1 = dPhi1Coef1_ + dPhi1Coef2_/clusterEnergyT
0091   const double dPhi1Coef2_;
0092   const double dPhi1Coef1_;
0093 
0094   const std::vector<const TrajectorySeedCollection*>* initialSeedCollectionVector_ = nullptr;
0095 
0096   const bool useRecoVertex_;
0097 
0098   const float deltaPhi2B_;
0099   const float deltaPhi2F_;
0100 
0101   const float phiMin2B_;
0102   const float phiMin2F_;
0103   const float phiMax2B_;
0104   const float phiMax2F_;
0105 
0106   PixelHitMatcher matcher_;
0107 };
0108 
0109 #endif  // ElectronSeedGenerator_H