CMS 3D CMS Logo

ConversionTools.cc
Go to the documentation of this file.
1 
2 #include <TMath.h>
12 
13 using namespace edm;
14 using namespace reco;
15 
16 
17 //--------------------------------------------------------------------------------------------------
18 bool ConversionTools::isGoodConversion(const Conversion &conv, const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
19 {
20 
21  //Check if a given conversion candidate passes the conversion selection cuts
22 
23  const reco::Vertex &vtx = conv.conversionVertex();
24 
25  //vertex validity
26  if (!vtx.isValid()) return false;
27 
28  //fit probability
29  if (TMath::Prob( vtx.chi2(), vtx.ndof() )<probMin) return false;
30 
31  //compute transverse decay length
33  double dbsx = vtx.x() - beamspot.x();
34  double dbsy = vtx.y() - beamspot.y();
35  double lxy = (mom.x()*dbsx + mom.y()*dbsy)/mom.rho();
36 
37  //transverse decay length
38  if ( lxy<lxyMin )
39  return false;
40 
41  //loop through daughters to check nhitsbeforevtx
42  for (std::vector<uint8_t>::const_iterator it = conv.nHitsBeforeVtx().begin(); it!=conv.nHitsBeforeVtx().end(); ++it) {
43  if ( (*it)>nHitsBeforeVtxMax ) return false;
44  }
45 
46  return true;
47 }
48 
49 //--------------------------------------------------------------------------------------------------
50 bool ConversionTools::matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch, bool allowAmbiguousGsfMatch)
51 {
52 
53  //check if a given GsfElectron matches a given conversion (no quality cuts applied)
54  //matching is always attempted through the gsf track ref, and optionally attempted through the
55  //closest ctf track ref
56 
57  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
58  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
59  if ( ele.reco::GsfElectron::gsfTrack().isNonnull() && ele.reco::GsfElectron::gsfTrack().id()==it->id() && ele.reco::GsfElectron::gsfTrack().key()==it->key()) return true;
60  else if ( allowCkfMatch && ele.reco::GsfElectron::closestCtfTrackRef().isNonnull() && ele.reco::GsfElectron::closestCtfTrackRef().id()==it->id() && ele.reco::GsfElectron::closestCtfTrackRef().key()==it->key() ) return true;
61  if (allowAmbiguousGsfMatch) {
63  if (tk->isNonnull() && tk->id()==it->id() && tk->key()==it->key()) return true;
64  }
65  }
66  }
67 
68  return false;
69 }
70 
71 //--------------------------------------------------------------------------------------------------
73 
74  //check if a given SuperCluster matches a given conversion (no quality cuts applied)
75  //matching is geometric between conversion momentum and vector joining conversion vertex
76  //to supercluster position
77 
78 
80 
81  const math::XYZPoint& scpos(sc.position());
83 
84 
85  math::XYZVector cscvector = scpos - cvtx;
86  float dR = reco::deltaR(mom,cscvector);
87  float dEta = mom.eta() - cscvector.eta();
88  float dPhi = reco::deltaPhi(mom.phi(),cscvector.phi());
89 
90  if (dR>dRMax) return false;
91  if (dEta>dEtaMax) return false;
92  if (dPhi>dPhiMax) return false;
93 
94  return true;
95 
96 }
97 
98 
99 //--------------------------------------------------------------------------------------------------
101 {
102 
103  //check if given track matches given conversion (matching by ref)
104 
105  if (trk.isNull()) return false;
106 
107  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
108  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
109  if (trk.id()==it->id() && trk.key()==it->key()) return true;
110  }
111 
112  return false;
113 }
114 
115 //--------------------------------------------------------------------------------------------------
117 {
118 
119  //check if given track matches given conversion (matching by ref)
120 
121  if (trk.isNull()) return false;
122 
123  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
124  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
125  if (trk.id()==it->id() && trk.key()==it->key()) return true;
126  }
127 
128  return false;
129 }
130 
131 //--------------------------------------------------------------------------------------------------
133 {
134 
135  //check if given track matches given conversion (matching by ref)
136 
137  if (trk.isNull()) return false;
138 
139  const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
140  for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
141  if (trk.id()==it->id() && trk.key()==it->key()) return true;
142  }
143 
144  return false;
145 }
146 
147 
148 //--------------------------------------------------------------------------------------------------
151  const math::XYZPoint &beamspot, bool allowCkfMatch, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
152 {
153  //check if a given electron candidate matches to at least one conversion candidate in the
154  //collection which also passes the selection cuts, optionally match with the closestckf track in
155  //in addition to just the gsf track (enabled in default arguments)
156 
157  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
158  if (!matchesConversion(ele, *it, allowCkfMatch)) continue;
159  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
160 
161  return true;
162  }
163 
164  return false;
165 
166 }
167 
168 //--------------------------------------------------------------------------------------------------
171  const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
172 {
173  //check if a given track matches to at least one conversion candidate in the
174  //collection which also passes the selection cuts
175 
176  if (trk.isNull()) return false;
177 
178  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
179  if (!matchesConversion(trk, *it)) continue;
180  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
181 
182  return true;
183  }
184 
185  return false;
186 
187 }
188 
189 //--------------------------------------------------------------------------------------------------
192  const math::XYZPoint &beamspot, float dRMax, float dEtaMax, float dPhiMax, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
193 {
194 
195  //check if a given SuperCluster matches to at least one conversion candidate in the
196  //collection which also passes the selection cuts
197 
198  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
199  if (!matchesConversion(sc, *it)) continue;
200  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
201 
202  return true;
203  }
204 
205  return false;
206 
207 }
208 
209 
210 //--------------------------------------------------------------------------------------------------
213  const math::XYZPoint &beamspot, bool allowCkfMatch, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
214 {
215  //check if a given electron candidate matches to at least one conversion candidate in the
216  //collection which also passes the selection cuts, optionally match with the closestckf track in
217  //in addition to just the gsf track (enabled in default arguments)
218  //If multiple conversions are found, returned reference corresponds to minimum
219  //conversion radius
220 
222 
223  double minRho = 999.;
224  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
225  float rho = it->conversionVertex().position().rho();
226  if (rho>minRho) continue;
227  if (!matchesConversion(ele, *it, allowCkfMatch)) continue;
228  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
229 
230  minRho = rho;
231  match = ConversionRef(convCol,it-convCol->begin());
232  }
233 
234  return match;
235 
236 }
237 
238 //--------------------------------------------------------------------------------------------------
241  const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
242 {
243  //check if a given track matches to at least one conversion candidate in the
244  //collection which also passes the selection cuts
245  //If multiple conversions are found, returned reference corresponds to minimum
246  //conversion radius
247 
249 
250  if (trk.isNull()) return match;
251 
252  double minRho = 999.;
253  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
254  float rho = it->conversionVertex().position().rho();
255  if (rho>minRho) continue;
256  if (!matchesConversion(trk, *it)) continue;
257  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
258 
259  minRho = rho;
260  match = ConversionRef(convCol,it-convCol->begin());
261  }
262 
263  return match;
264 
265 }
266 
267 //--------------------------------------------------------------------------------------------------
270  const math::XYZPoint &beamspot, float dRMax, float dEtaMax, float dPhiMax, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
271 {
272 
273  //check if a given SuperCluster matches to at least one conversion candidate in the
274  //collection which also passes the selection cuts
275  //If multiple conversions are found, returned reference corresponds to minimum
276  //conversion radius
277 
279 
280  double minRho = 999.;
281  for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
282  float rho = it->conversionVertex().position().rho();
283  if (rho>minRho) continue;
284  if (!matchesConversion(sc, *it, dRMax,dEtaMax,dPhiMax)) continue;
285  if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
286 
287  minRho = rho;
288  match = ConversionRef(convCol,it-convCol->begin());
289  }
290 
291  return match;
292 
293 }
294 
295 //--------------------------------------------------------------------------------------------------
297  const edm::Handle<reco::ConversionCollection> &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
298 {
299 
300  //check if a given SuperCluster matches to at least one GsfElectron having zero expected inner hits
301  //and not matching any conversion in the collection passing the quality cuts
302 
303  if (sc.isNull()) return false;
304 
305  for (GsfElectronCollection::const_iterator it = eleCol->begin(); it!=eleCol->end(); ++it) {
306  //match electron to supercluster
307  if (it->reco::GsfElectron::superCluster()!=sc) continue;
308 
309  //check expected inner hits
310  if (it->gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) > 0) continue;
311 
312  //check if electron is matching to a conversion
313  if (hasMatchedConversion(*it,convCol,beamspot,allowCkfMatch,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
314 
315 
316  return true;
317  }
318 
319  return false;
320 
321 
322 }
323 
324 
325 //--------------------------------------------------------------------------------------------------
327  const edm::Handle<reco::ConversionCollection> &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
328 {
329 
330  //check if a given SuperCluster matches to at least one GsfElectron having zero expected inner hits
331  //and not matching any conversion in the collection passing the quality cuts
332 
334 
335  if (sc.isNull()) return match;
336 
337  for (GsfElectronCollection::const_iterator it = eleCol->begin(); it!=eleCol->end(); ++it) {
338  //match electron to supercluster
339  if (it->reco::GsfElectron::superCluster()!=sc) continue;
340 
341  //check expected inner hits
342  if (it->gsfTrack()->hitPattern().numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS) > 0) continue;
343 
344  //check if electron is matching to a conversion
345  if (hasMatchedConversion(*it,convCol,beamspot,allowCkfMatch,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
346 
347 
348  match = GsfElectronRef(eleCol,it-eleCol->begin());
349  }
350 
351  return match;
352 
353 
354 }
const reco::Vertex & conversionVertex() const
returns the reco conversion vertex
Definition: Conversion.h:97
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:129
static bool matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch=true, bool allowAmbiguousGsfMatch=false)
GsfTrackRefVector::const_iterator ambiguousGsfTracksBegin() const
Definition: GsfElectron.h:709
static HepMC::IO_HEPEVT conv
double y() const
y coordinate
Definition: Vertex.h:113
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:68
static bool hasMatchedPromptElectron(const reco::SuperClusterRef &sc, const edm::Handle< reco::GsfElectronCollection > &eleCol, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
key_type key() const
Accessor for product key.
Definition: Ref.h:265
const Point & position() const
position
Definition: Vertex.h:109
ProductID id() const
Definition: RefToBase.h:242
math::XYZVectorF refittedPairMomentum() const
Conversion tracks momentum from the tracks refitted with vertex constraint.
Definition: Conversion.cc:248
ProductID id() const
Accessor for product ID.
Definition: Ref.h:259
GsfTrackRefVector::const_iterator ambiguousGsfTracksEnd() const
Definition: GsfElectron.h:710
size_t key() const
Definition: RefToBase.h:250
double chi2() const
chi-squares
Definition: Vertex.h:98
static bool hasMatchedConversion(const reco::GsfElectron &ele, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
bool isNull() const
Checks for null.
Definition: Ref.h:250
static bool isGoodConversion(const reco::Conversion &conv, const math::XYZPoint &beamspot, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=1)
edm::Ref< GsfElectronCollection > GsfElectronRef
reference to an object in a collection of GsfElectron objects
double ndof() const
Definition: Vertex.h:105
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
double x() const
x coordinate
Definition: Vertex.h:111
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
bool isNull() const
Checks for null.
Definition: RefToBase.h:328
fixed size matrix
HLT enums.
const std::vector< uint8_t > & nHitsBeforeVtx() const
Vector of the number of hits before the vertex along each track trajector.
Definition: Conversion.h:161
edm::Ref< ConversionCollection > ConversionRef
reference to an object in a collection of Conversion objects
Definition: ConversionFwd.h:15
ProductIndex id() const
Definition: ProductID.h:38
static reco::ConversionRef matchedConversion(const reco::GsfElectron &ele, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
static reco::GsfElectronRef matchedPromptElectron(const reco::SuperClusterRef &sc, const edm::Handle< reco::GsfElectronCollection > &eleCol, const edm::Handle< reco::ConversionCollection > &convCol, const math::XYZPoint &beamspot, bool allowCkfMatch=true, float lxyMin=2.0, float probMin=1e-6, unsigned int nHitsBeforeVtxMax=0)
std::vector< edm::RefToBase< reco::Track > > const & tracks() const
vector of track to base references
Definition: Conversion.cc:176
std::string match(BranchDescription const &a, BranchDescription const &b, std::string const &fileName)