CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackerOnlyConversionProducer.h
Go to the documentation of this file.
1 #ifndef RecoEgamma_EgammaPhotonProducers_TrackerOnlyConversionProducer_h
2 #define RecoEgamma_EgammaPhotonProducers_TrackerOnlyConversionProducer_h
3 
12 // user include files
15 
18 
21 
23 
24 //ECAL clusters
37 
46 
47 
48 //Tracker tracks
52 
53 
54 //photon data format
60 
63 
66 
68  public:
71 
72  void buildCollection( edm::Event& iEvent, const edm::EventSetup& iSetup,
73  //const reco::TrackRefVector& allTracks,
75  const std::multimap<double, reco::CaloClusterPtr>& basicClusterPtrs,
76  const reco::Vertex& the_pvtx,
77  reco::ConversionCollection & outputConvPhotonCollection);
78 
79  void buildCollection( edm::Event& iEvent, const edm::EventSetup& iSetup,
80  const reco::ConversionTrackCollection& allTracks,
81  const reco::CaloClusterPtr& basicClusterPtrs,
82  const reco::Vertex& the_pvtx,
83  reco::ConversionCollection & outputConvPhotonCollection);
84 
85  //track quality cut, returns pass or no
86  inline bool trackQualityFilter(const edm::RefToBase<reco::Track>& ref, bool isLeft);
87  inline bool trackD0Cut(const edm::RefToBase<reco::Track>& ref);
88  inline bool trackD0Cut(const edm::RefToBase<reco::Track>& ref, const reco::Vertex& the_pvtx);
89 
90  //track impact point at ECAL wall, returns validity to access position ew
91  bool getTrackImpactPosition(const reco::Track* tk_ref,
92  const TrackerGeometry* trackerGeom, const MagneticField* magField,
93  math::XYZPoint& ew);
94 
95  //distance at min approaching point, returns distance
96  // double getMinApproach(const edm::RefToBase<reco::Track>& ll, const edm::RefToBase<reco::Track>& rr,
97  // const MagneticField* magField);
98 
99  bool preselectTrackPair(const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r,
100  double& appDist);
101 
102  //cut-based selection, TODO remove global cut variables
104  const std::pair<edm::RefToBase<reco::Track>, reco::CaloClusterPtr>& rr);
105 
106  //kinematic vertex fitting, return true for valid vertex
107  bool checkVertex(const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r,
108  const MagneticField* magField,
109  reco::Vertex& the_vertex);
111  const TrackerGeometry* trackerGeom, const MagneticField* magField,
112  const reco::Vertex& the_vertex);
113 
114  //check the closest BC, returns true for found a BC
115  bool getMatchedBC(const std::multimap<double, reco::CaloClusterPtr>& bcMap,
116  const math::XYZPoint& trackImpactPosition,
117  reco::CaloClusterPtr& closestBC);
118 
119  bool getMatchedBC(const reco::CaloClusterPtrVector& bcMap,
120  const math::XYZPoint& trackImpactPosition,
121  reco::CaloClusterPtr& closestBC);
122 
123 
124 
125  private:
126 
127  virtual void produce(edm::Event&, const edm::EventSetup&);
128 
129  // ----------member data ---------------------------
130  std::string algoName_;
131 
133  typedef std::vector<Point> PointCollection;
134 
136 
140 
142 
143  bool usePvtx_;//if use primary vertices
144  std::string vertexProducer_;
146 
148 
149 
150  double halfWayEta_, halfWayPhi_;//halfway open angle to search in basic clusters
151  unsigned int maxNumOfTrackInPU_;
152 
153  double energyBC_;//1.5GeV for track BC selection
154  double energyTotalBC_;//5GeV for track pair BC selection
155  double d0Cut_;//0 for d0*charge cut
156  double dzCut_;//innerposition of z diff cut
157  double dEtaTkBC_, dPhiTkBC_;//0.06 0.6 for track and BC matching
158 
159  double maxChi2Left_, maxChi2Right_;//5. 5. for track chi2 quality
160  double minHitsLeft_, minHitsRight_;//5 2 for track hits quality
161 
162  double deltaCotTheta_, deltaPhi_, minApproachLow_, minApproachHigh_;//0.02 0.2 for track pair open angle and > -0.1 cm
163 
164 
165  double r_cut;//cross_r cut
166  double vtxChi2_;//vertex chi2 probablity cut
167 
168  bool allowSingleLeg_;//if single track conversion ?
169  bool rightBC_;//if right leg requires matching BC?
170 
171 };
172 
173 
174 inline const GeomDet * recHitDet( const TrackingRecHit & hit, const TrackingGeometry * geom ) {
175  return geom->idToDet( hit.geographicalId() );
176 }
177 
178 inline const BoundPlane & recHitSurface( const TrackingRecHit & hit, const TrackingGeometry * geom ) {
179  return recHitDet( hit, geom )->surface();
180 }
181 
182 inline LocalVector toLocal( const reco::Track::Vector & v, const Surface & s ) {
183  return s.toLocal( GlobalVector( v.x(), v.y(), v.z() ) );
184 }
185 
186 #endif
bool checkPhi(const edm::RefToBase< reco::Track > &tk_l, const edm::RefToBase< reco::Track > &tk_r, const TrackerGeometry *trackerGeom, const MagneticField *magField, const reco::Vertex &the_vertex)
std::vector< Conversion > ConversionCollection
collectin of Conversion objects
Definition: ConversionFwd.h:9
bool trackQualityFilter(const edm::RefToBase< reco::Track > &ref, bool isLeft)
std::vector< ConversionTrack > ConversionTrackCollection
collection of ConversionTracks
int iEvent
Definition: GenABIO.cc:243
LocalPoint toLocal(const GlobalPoint &gp) const
virtual void produce(edm::Event &, const edm::EventSetup &)
TrackerOnlyConversionProducer(const edm::ParameterSet &)
bool getTrackImpactPosition(const reco::Track *tk_ref, const TrackerGeometry *trackerGeom, const MagneticField *magField, math::XYZPoint &ew)
bool trackD0Cut(const edm::RefToBase< reco::Track > &ref)
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< double > > XYZPointD
point in space with cartesian internal representation
Definition: Point3D.h:9
const TransientTrackBuilder * thettbuilder_
virtual const GeomDet * idToDet(DetId) const =0
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
bool getMatchedBC(const std::multimap< double, reco::CaloClusterPtr > &bcMap, const math::XYZPoint &trackImpactPosition, reco::CaloClusterPtr &closestBC)
bool preselectTrackPair(const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r, double &appDist)
LocalVector toLocal(const reco::Track::Vector &v, const Surface &s)
bool checkTrackPair(const std::pair< edm::RefToBase< reco::Track >, reco::CaloClusterPtr > &ll, const std::pair< edm::RefToBase< reco::Track >, reco::CaloClusterPtr > &rr)
const GeomDet * recHitDet(const TrackingRecHit &hit, const TrackingGeometry *geom)
const BoundPlane & recHitSurface(const TrackingRecHit &hit, const TrackingGeometry *geom)
string s
Definition: asciidump.py:422
DetId geographicalId() const
bool checkVertex(const reco::TransientTrack &ttk_l, const reco::TransientTrack &ttk_r, const MagneticField *magField, reco::Vertex &the_vertex)
void buildCollection(edm::Event &iEvent, const edm::EventSetup &iSetup, const reco::ConversionTrackCollection &allTracks, const std::multimap< double, reco::CaloClusterPtr > &basicClusterPtrs, const reco::Vertex &the_pvtx, reco::ConversionCollection &outputConvPhotonCollection)
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:73
virtual const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
mathSSE::Vec4< T > v
Global3DVector GlobalVector
Definition: GlobalVector.h:10