CMS 3D CMS Logo

Conversion.h
Go to the documentation of this file.
1 #ifndef EgammaCandidates_Conversion_h
2 #define EgammaCandidates_Conversion_h
3 
12 #include <bitset>
21 
22 
23 namespace reco {
24  class Conversion {
25  public:
26 
30  mixed=3,
31  pflow=4,
32  algoSize=5};
33 
43 
44  static const std::string algoNames[];
45 
46  // Default constructor
47  Conversion();
48 
50  const std::vector<edm::RefToBase<reco::Track> >& tr,
51  const std::vector<math::XYZPointF>& trackPositionAtEcal ,
52  const reco::Vertex & convVtx,
53  const std::vector<reco::CaloClusterPtr> & matchingBC,
54  const float DCA,
55  const std::vector<math::XYZPointF> & innPoint,
56  const std::vector<math::XYZVectorF> & trackPin ,
57  const std::vector<math::XYZVectorF> & trackPout,
58  const std::vector<uint8_t>& nHitsBeforeVtx,
59  const std::vector<Measurement1DFloat> & dlClosestHitToVtx,
60  uint8_t nSharedHits,
61  const float mva,
63 
64 
66  const std::vector<reco::TrackRef>& tr,
67  const std::vector<math::XYZPointF>& trackPositionAtEcal ,
68  const reco::Vertex & convVtx,
69  const std::vector<reco::CaloClusterPtr> & matchingBC,
70  const float DCA,
71  const std::vector<math::XYZPointF> & innPoint,
72  const std::vector<math::XYZVectorF> & trackPin ,
73  const std::vector<math::XYZVectorF> & trackPout,
74  const float mva,
76 
77 
78 
79 
81  const std::vector<reco::TrackRef>& tr,
82  const reco::Vertex & convVtx,
84 
86  const std::vector<edm::RefToBase<reco::Track> >& tr,
87  const reco::Vertex & convVtx,
89 
90 
92  Conversion * clone() const;
96  std::vector<edm::RefToBase<reco::Track> > const& tracks() const ;
98  const reco::Vertex & conversionVertex() const { return theConversionVertex_ ; }
100  bool isConverted() const;
102  unsigned int nTracks() const {return tracks().size(); }
104  double MVAout() const { return theMVAout_;}
106  std::vector<float> const oneLegMVA() {return theOneLegMVA_;}
108  double pairInvariantMass() const;
110  double pairCotThetaSeparation() const;
120  double EoverP() const;
124  double EoverPrefittedTracks() const;
125  // Dist of minimum approach between tracks
127  // deltaPhi tracks at innermost point
128  double dPhiTracksAtVtx() const;
129  // deltaPhi tracks at ECAl
130  double dPhiTracksAtEcal() const;
131  // deltaEta tracks at ECAl
132  double dEtaTracksAtEcal() const;
133 
134  //impact parameter and decay length computed with respect to given beamspot or vertex
135  //computed from refittedPairMomentum
136 
137  //transverse impact parameter
138  double dxy(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
139  //longitudinal impact parameter
140  double dz(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
141  //transverse decay length
142  double lxy(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
143  //longitudinal decay length
144  double lz(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const;
145  //z position of intersection with beamspot in rz plane (possible tilt of beamspot is neglected)
146  double zOfPrimaryVertexFromTracks(const math::XYZPoint& myBeamSpot = math::XYZPoint()) const { return dz(myBeamSpot) + myBeamSpot.z(); }
147 
150  const std::vector<math::XYZPointF> & ecalImpactPosition() const {return thePositionAtEcal_;}
151  // pair of BC matching a posteriori the tracks
152  const std::vector<reco::CaloClusterPtr>& bcMatchingWithTracks() const { return theMatchingBCs_;}
154  std::vector<double> tracksSigned_d0() const ;
156  const std::vector<math::XYZPointF>& tracksInnerPosition() const {return theTrackInnerPosition_;}
158  const std::vector<math::XYZVectorF>& tracksPout() const {return theTrackPout_;}
160  const std::vector<math::XYZVectorF>& tracksPin() const {return theTrackPin_;}
162  const std::vector<uint8_t> &nHitsBeforeVtx() const { return nHitsBeforeVtx_; }
164  const std::vector<Measurement1DFloat> &dlClosestHitToVtx() const { return dlClosestHitToVtx_; }
166  uint8_t nSharedHits() const { return nSharedHits_; }
167 
168 
170  void setMVAout(const float& mva) { theMVAout_=mva;}
172  void setOneLegMVA(const std::vector<float>& mva) { theOneLegMVA_=mva;}
173  // Set the ptr to the Super cluster if not set in the constructor
176  void setConversionAlgorithm(const ConversionAlgorithm a, bool set=true) { if (set) algorithm_=a; else algorithm_=undefined;}
177  ConversionAlgorithm algo() const ;
178  std::string algoName() const;
181 
182  bool quality(ConversionQuality q) const { return (qualityMask_ & (1<<q))>>q; }
183  void setQuality(ConversionQuality q, bool b);
184 
185 
186 
187  private:
188 
192  std::vector<edm::RefToBase<reco::Track> > trackToBaseRefs_;
194  std::vector<math::XYZPointF> thePositionAtEcal_;
198  std::vector<reco::CaloClusterPtr> theMatchingBCs_;
200  std::vector<math::XYZPointF> theTrackInnerPosition_;
202  std::vector<math::XYZVectorF> theTrackPin_;
204  std::vector<math::XYZVectorF> theTrackPout_;
206  std::vector<uint8_t> nHitsBeforeVtx_;
208  std::vector<Measurement1DFloat> dlClosestHitToVtx_;
210  std::vector<float> theOneLegMVA_;
214  float theMVAout_;
215  uint16_t qualityMask_;
217  uint8_t nSharedHits_;
219  uint8_t algorithm_;
220 
221 
222  };
223 
226  }
227 
228 
230 
231  switch(algorithm_)
232  {
233  case undefined: return "undefined";
234  case ecalSeeded: return "ecalSeeded";
235  case trackerOnly: return "trackerOnly";
236  case mixed: return "mixed";
237  case pflow: return "pflow";
238 
239  }
240  return "undefined";
241  }
242 
244  if(int(a) < int(algoSize) && int(a)>0) return algoNames[int(a)];
245  return "undefined";
246  }
247 
249  if (b)//regular OR if setting value to true
250  qualityMask_ |= (1<<q) ;
251  else // doing "half-XOR" if unsetting value
252  qualityMask_ &= (~(1<<q));
253 
254  }
255 
256 }
257 
258 #endif
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:98
static const std::string algoNames[]
Definition: Conversion.h:44
std::vector< math::XYZVectorF > theTrackPout_
P_out of tracks.
Definition: Conversion.h:204
ConversionAlgorithm algo() const
Definition: Conversion.h:224
bool quality(ConversionQuality q) const
Definition: Conversion.h:182
void setQuality(ConversionQuality q, bool b)
Definition: Conversion.h:248
double EoverP() const
Definition: Conversion.cc:259
const std::vector< Measurement1DFloat > & dlClosestHitToVtx() const
Vector of signed decay length with uncertainty from nearest hit on track to the conversion vtx positi...
Definition: Conversion.h:164
reco::CaloClusterPtrVector caloCluster_
vector pointer to a/multiple seed CaloCluster(s)
Definition: Conversion.h:190
double zOfPrimaryVertexFromTracks(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.h:146
double distOfMinimumApproach() const
Definition: Conversion.h:126
std::vector< double > tracksSigned_d0() const
signed transverse impact parameter for each track
Definition: Conversion.cc:303
double lxy(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.cc:392
static ConversionAlgorithm algoByName(const std::string &name)
Definition: Conversion.cc:163
double pairCotThetaSeparation() const
Delta cot(Theta) where Theta is the angle in the (y,z) plane between the two tracks. Original tracks are used.
Definition: Conversion.cc:209
void setConversionAlgorithm(const ConversionAlgorithm a, bool set=true)
Conversion Track algorithm/provenance.
Definition: Conversion.h:176
double pairInvariantMass() const
if nTracks=2 returns the pair invariant mass. Original tracks are used here
Definition: Conversion.cc:190
math::XYZVectorF refittedPairMomentum() const
Conversion tracks momentum from the tracks refitted with vertex constraint.
Definition: Conversion.cc:248
uint8_t nSharedHits_
number of shared hits between tracks
Definition: Conversion.h:217
double EoverPrefittedTracks() const
Definition: Conversion.cc:281
const std::vector< math::XYZVectorF > & tracksPout() const
Vector of track momentum measured at the outermost hit.
Definition: Conversion.h:158
double MVAout() const
get the value of the TMVA output
Definition: Conversion.h:104
uint16_t qualityMask_
Definition: Conversion.h:215
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float > > XYZVectorF
spatial vector with cartesian internal representation
Definition: Vector3D.h:17
const std::vector< reco::CaloClusterPtr > & bcMatchingWithTracks() const
Definition: Conversion.h:152
std::vector< float > theOneLegMVA_
vectors of TMVA outputs from pflow for one leg conversions
Definition: Conversion.h:210
std::string algoName() const
Definition: Conversion.h:229
std::vector< math::XYZPointF > theTrackInnerPosition_
P_in of tracks.
Definition: Conversion.h:200
std::vector< float > const oneLegMVA()
get the MVS output from PF for one leg conversions
Definition: Conversion.h:106
unsigned int nTracks() const
Number of tracks= 0,1,2.
Definition: Conversion.h:102
Conversion * clone() const
returns a clone of the candidate
Definition: Conversion.cc:171
math::XYZTLorentzVectorF refittedPair4Momentum() const
Conversion track pair 4-momentum from the tracks refitted with vertex constraint. ...
Definition: Conversion.cc:235
const std::vector< math::XYZVectorF > & tracksPin() const
Vector of track momentum measured at the innermost hit.
Definition: Conversion.h:160
uint8_t algorithm_
conversion algorithm/provenance
Definition: Conversion.h:219
const std::vector< math::XYZPointF > & tracksInnerPosition() const
Vector containing the position of the innermost hit of each track.
Definition: Conversion.h:156
math::XYZVectorF pairMomentum() const
Conversion tracks momentum from the tracks inner momentum.
Definition: Conversion.cc:223
const std::vector< math::XYZPointF > & ecalImpactPosition() const
Definition: Conversion.h:150
uint8_t nSharedHits() const
number of shared hits btw the two track
Definition: Conversion.h:166
double dEtaTracksAtEcal() const
Definition: Conversion.cc:351
double dz(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.cc:380
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< math::XYZVectorF > theTrackPin_
P_in of tracks.
Definition: Conversion.h:202
float theMVAout_
TMVA output.
Definition: Conversion.h:214
double b
Definition: hdecay.h:120
std::vector< reco::CaloClusterPtr > theMatchingBCs_
Clusters mathing the tracks (these are not the seeds)
Definition: Conversion.h:198
double lz(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.cc:406
void setMVAout(const float &mva)
set the value of the TMVA output
Definition: Conversion.h:170
float theMinDistOfApproach_
Distance of min approach of the two tracks.
Definition: Conversion.h:212
std::vector< math::XYZPointF > thePositionAtEcal_
position at the ECAl surface of the track extrapolation
Definition: Conversion.h:194
fixed size matrix
void setOneLegMVA(const std::vector< float > &mva)
set the MVS output from PF for one leg conversions
Definition: Conversion.h:172
double a
Definition: hdecay.h:121
std::vector< edm::RefToBase< reco::Track > > trackToBaseRefs_
vector Track RefToBase
Definition: Conversion.h:192
const std::vector< uint8_t > & nHitsBeforeVtx() const
Vector of the number of hits before the vertex along each track trajector.
Definition: Conversion.h:162
std::vector< Measurement1DFloat > dlClosestHitToVtx_
signed decay length and uncertainty from nearest hit on track to conversion vertex ...
Definition: Conversion.h:208
double dPhiTracksAtEcal() const
Definition: Conversion.cc:326
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:22
reco::Vertex theConversionVertex_
Fitted Kalman conversion vertex.
Definition: Conversion.h:196
bool isConverted() const
Bool flagging objects having track size >0.
Definition: Conversion.cc:182
std::vector< uint8_t > nHitsBeforeVtx_
number of hits before the vertex on each trackerOnly
Definition: Conversion.h:206
reco::CaloClusterPtrVector caloCluster() const
Pointer to CaloCluster (foe Egamma Conversions it points to a SuperCluster)
Definition: Conversion.h:94
void setMatchingSuperCluster(const reco::CaloClusterPtrVector &sc)
Definition: Conversion.h:174
double dPhiTracksAtVtx() const
Definition: Conversion.cc:313
std::vector< edm::RefToBase< reco::Track > > const & tracks() const
vector of track to base references
Definition: Conversion.cc:176
double dxy(const math::XYZPoint &myBeamSpot=math::XYZPoint()) const
Definition: Conversion.cc:368