00001
00002
00003 #include <TMath.h>
00004 #include "RecoEgamma/EgammaTools/interface/ConversionTools.h"
00005 #include "DataFormats/EgammaCandidates/interface/Conversion.h"
00006 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00007 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00008 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00009 #include "DataFormats/TrackReco/interface/Track.h"
00010 #include "DataFormats/Common/interface/RefToPtr.h"
00011 #include "DataFormats/Math/interface/deltaPhi.h"
00012 #include "DataFormats/Math/interface/deltaR.h"
00013
00014 using namespace edm;
00015 using namespace reco;
00016
00017
00018
00019 bool ConversionTools::isGoodConversion(const Conversion &conv, const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
00020 {
00021
00022
00023
00024 const reco::Vertex &vtx = conv.conversionVertex();
00025
00026
00027 if (!vtx.isValid()) return false;
00028
00029
00030 if (TMath::Prob( vtx.chi2(), vtx.ndof() )<probMin) return false;
00031
00032
00033 math::XYZVector mom(conv.refittedPairMomentum());
00034 double dbsx = vtx.x() - beamspot.x();
00035 double dbsy = vtx.y() - beamspot.y();
00036 double lxy = (mom.x()*dbsx + mom.y()*dbsy)/mom.rho();
00037
00038
00039 if ( lxy<lxyMin )
00040 return false;
00041
00042
00043 for (std::vector<uint8_t>::const_iterator it = conv.nHitsBeforeVtx().begin(); it!=conv.nHitsBeforeVtx().end(); ++it) {
00044 if ( (*it)>nHitsBeforeVtxMax ) return false;
00045 }
00046
00047 return true;
00048 }
00049
00050
00051 bool ConversionTools::matchesConversion(const reco::GsfElectron &ele, const reco::Conversion &conv, bool allowCkfMatch)
00052 {
00053
00054
00055
00056
00057
00058 const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
00059 for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
00060 if ( ele.gsfTrack().isNonnull() && ele.gsfTrack().id()==it->id() && ele.gsfTrack().key()==it->key()) return true;
00061 else if ( allowCkfMatch && ele.closestCtfTrackRef().isNonnull() && ele.closestCtfTrackRef().id()==it->id() && ele.closestCtfTrackRef().key()==it->key() ) return true;
00062 }
00063
00064 return false;
00065 }
00066
00067
00068 bool ConversionTools::matchesConversion(const reco::SuperCluster &sc, const reco::Conversion &conv, float dRMax, float dEtaMax, float dPhiMax) {
00069
00070
00071
00072
00073
00074
00075 math::XYZVector mom(conv.refittedPairMomentum());
00076
00077 math::XYZPoint scpos(sc.position());
00078 math::XYZPoint cvtx(conv.conversionVertex().position());
00079
00080
00081 math::XYZVector cscvector = scpos - cvtx;
00082 float dR = reco::deltaR(mom,cscvector);
00083 float dEta = mom.eta() - cscvector.eta();
00084 float dPhi = reco::deltaPhi(mom.phi(),cscvector.phi());
00085
00086 if (dR>dRMax) return false;
00087 if (dEta>dEtaMax) return false;
00088 if (dPhi>dPhiMax) return false;
00089
00090 return true;
00091
00092 }
00093
00094
00095
00096 bool ConversionTools::matchesConversion(const edm::RefToBase<reco::Track> &trk, const reco::Conversion &conv)
00097 {
00098
00099
00100
00101 if (trk.isNull()) return false;
00102
00103 const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
00104 for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
00105 if (trk.id()==it->id() && trk.key()==it->key()) return true;
00106 }
00107
00108 return false;
00109 }
00110
00111
00112 bool ConversionTools::matchesConversion(const reco::TrackRef &trk, const reco::Conversion &conv)
00113 {
00114
00115
00116
00117 if (trk.isNull()) return false;
00118
00119 const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
00120 for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
00121 if (trk.id()==it->id() && trk.key()==it->key()) return true;
00122 }
00123
00124 return false;
00125 }
00126
00127
00128 bool ConversionTools::matchesConversion(const reco::GsfTrackRef &trk, const reco::Conversion &conv)
00129 {
00130
00131
00132
00133 if (trk.isNull()) return false;
00134
00135 const std::vector<edm::RefToBase<reco::Track> > &convTracks = conv.tracks();
00136 for (std::vector<edm::RefToBase<reco::Track> >::const_iterator it=convTracks.begin(); it!=convTracks.end(); ++it) {
00137 if (trk.id()==it->id() && trk.key()==it->key()) return true;
00138 }
00139
00140 return false;
00141 }
00142
00143
00144
00145 bool ConversionTools::hasMatchedConversion(const reco::GsfElectron &ele,
00146 const edm::Handle<reco::ConversionCollection> &convCol,
00147 const math::XYZPoint &beamspot, bool allowCkfMatch, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
00148 {
00149
00150
00151
00152
00153 for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
00154 if (!matchesConversion(ele, *it, allowCkfMatch)) continue;
00155 if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
00156
00157 return true;
00158 }
00159
00160 return false;
00161
00162 }
00163
00164
00165 bool ConversionTools::hasMatchedConversion(const reco::TrackRef &trk,
00166 const edm::Handle<reco::ConversionCollection> &convCol,
00167 const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
00168 {
00169
00170
00171
00172 if (trk.isNull()) return false;
00173
00174 for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
00175 if (!matchesConversion(trk, *it)) continue;
00176 if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
00177
00178 return true;
00179 }
00180
00181 return false;
00182
00183 }
00184
00185
00186 bool ConversionTools::hasMatchedConversion(const reco::SuperCluster &sc,
00187 const edm::Handle<reco::ConversionCollection> &convCol,
00188 const math::XYZPoint &beamspot, float dRMax, float dEtaMax, float dPhiMax, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
00189 {
00190
00191
00192
00193
00194 for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
00195 if (!matchesConversion(sc, *it)) continue;
00196 if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
00197
00198 return true;
00199 }
00200
00201 return false;
00202
00203 }
00204
00205
00206
00207 reco::ConversionRef ConversionTools::matchedConversion(const reco::GsfElectron &ele,
00208 const edm::Handle<reco::ConversionCollection> &convCol,
00209 const math::XYZPoint &beamspot, bool allowCkfMatch, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
00210 {
00211
00212
00213
00214
00215
00216
00217 ConversionRef match;
00218
00219 double minRho = 999.;
00220 for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
00221 float rho = it->conversionVertex().position().rho();
00222 if (rho>minRho) continue;
00223 if (!matchesConversion(ele, *it, allowCkfMatch)) continue;
00224 if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
00225
00226 minRho = rho;
00227 match = ConversionRef(convCol,it-convCol->begin());
00228 }
00229
00230 return match;
00231
00232 }
00233
00234
00235 reco::ConversionRef ConversionTools::matchedConversion(const reco::TrackRef &trk,
00236 const edm::Handle<reco::ConversionCollection> &convCol,
00237 const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
00238 {
00239
00240
00241
00242
00243
00244 ConversionRef match;
00245
00246 if (trk.isNull()) return match;
00247
00248 double minRho = 999.;
00249 for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
00250 float rho = it->conversionVertex().position().rho();
00251 if (rho>minRho) continue;
00252 if (!matchesConversion(trk, *it)) continue;
00253 if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
00254
00255 minRho = rho;
00256 match = ConversionRef(convCol,it-convCol->begin());
00257 }
00258
00259 return match;
00260
00261 }
00262
00263
00264 reco::ConversionRef ConversionTools::matchedConversion(const reco::SuperCluster &sc,
00265 const edm::Handle<reco::ConversionCollection> &convCol,
00266 const math::XYZPoint &beamspot, float dRMax, float dEtaMax, float dPhiMax, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
00267 {
00268
00269
00270
00271
00272
00273
00274 ConversionRef match;
00275
00276 double minRho = 999.;
00277 for (ConversionCollection::const_iterator it = convCol->begin(); it!=convCol->end(); ++it) {
00278 float rho = it->conversionVertex().position().rho();
00279 if (rho>minRho) continue;
00280 if (!matchesConversion(sc, *it, dRMax,dEtaMax,dPhiMax)) continue;
00281 if (!isGoodConversion(*it,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
00282
00283 minRho = rho;
00284 match = ConversionRef(convCol,it-convCol->begin());
00285 }
00286
00287 return match;
00288
00289 }
00290
00291
00292 bool ConversionTools::hasMatchedPromptElectron(const reco::SuperClusterRef &sc, const edm::Handle<reco::GsfElectronCollection> &eleCol,
00293 const edm::Handle<reco::ConversionCollection> &convCol, const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
00294 {
00295
00296
00297
00298
00299 if (sc.isNull()) return false;
00300
00301 for (GsfElectronCollection::const_iterator it = eleCol->begin(); it!=eleCol->end(); ++it) {
00302
00303 if (it->superCluster()!=sc) continue;
00304
00305
00306 if (it->gsfTrack()->trackerExpectedHitsInner().numberOfHits()>0) continue;
00307
00308
00309 if (hasMatchedConversion(*it,convCol,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
00310
00311
00312 return true;
00313 }
00314
00315 return false;
00316
00317
00318 }
00319
00320
00321
00322 reco::GsfElectronRef ConversionTools::matchedPromptElectron(const reco::SuperClusterRef &sc, const edm::Handle<reco::GsfElectronCollection> &eleCol,
00323 const edm::Handle<reco::ConversionCollection> &convCol, const math::XYZPoint &beamspot, float lxyMin, float probMin, unsigned int nHitsBeforeVtxMax)
00324 {
00325
00326
00327
00328
00329 GsfElectronRef match;
00330
00331 if (sc.isNull()) return match;
00332
00333 for (GsfElectronCollection::const_iterator it = eleCol->begin(); it!=eleCol->end(); ++it) {
00334
00335 if (it->superCluster()!=sc) continue;
00336
00337
00338 if (it->gsfTrack()->trackerExpectedHitsInner().numberOfHits()>0) continue;
00339
00340
00341 if (hasMatchedConversion(*it,convCol,beamspot,lxyMin,probMin,nHitsBeforeVtxMax)) continue;
00342
00343
00344 match = GsfElectronRef(eleCol,it-eleCol->begin());
00345 }
00346
00347 return match;
00348
00349
00350 }