28 const double bFieldAtOrigin,
29 const double minFracSharedHits) {
31 std::vector<ConversionInfo>
temp =
getConversionInfos(*gsfElectron.
core(),ctftracks_h,gsftracks_h,bFieldAtOrigin,minFracSharedHits) ;
39 const double bFieldAtOrigin,
40 const double minFracSharedHits) {
42 std::vector<ConversionInfo>
temp =
getConversionInfos(gsfElectron,ctftracks_h,gsftracks_h,bFieldAtOrigin,minFracSharedHits) ;
52 const double bFieldAtOrigin,
53 const double minFracSharedHits) {
72 if(el_ctftrack.isNonnull() && el_ctftrack.id() != ctftracks_h.
id())
73 throw cms::Exception(
"ConversionFinderError") <<
"ProductID of ctf track collection does not match ProductID of electron's CTF track! \n";
74 if(el_gsftrack.isNonnull() && el_gsftrack.id() != gsftracks_h.
id())
75 throw cms::Exception(
"ConversionFinderError") <<
"ProductID of gsf track collection does not match ProductID of electron's GSF track! \n";
79 if(el_ctftrack.isNonnull() && gsfElectron.
ctfGsfOverlap() > minFracSharedHits)
80 el_ctftrack_p4 =
LorentzVector(el_ctftrack->px(), el_ctftrack->py(), el_ctftrack->pz(), el_ctftrack->p());
81 LorentzVector el_gsftrack_p4(el_gsftrack->px(), el_gsftrack->py(), el_gsftrack->pz(), el_gsftrack->p());
87 if(el_ctftrack.isNonnull() && gsfElectron.
ctfGsfOverlap() > minFracSharedHits)
88 ctfidx = static_cast<int>(el_ctftrack.key());
90 gsfidx =
static_cast<int>(el_gsftrack.key());
94 vector<ConversionInfo> v_candidatePartners;
101 for(TrackCollection::const_iterator ctftk = ctftracks->begin();
102 ctftk != ctftracks->end(); ctftk++, ctftk_i++) {
104 if((ctftk_i == ctfidx))
111 if(ctftk->ptError()/ctftk->pt() > 0.05)
113 if(ctftk->numberOfValidHits() < 5)
116 if(el_ctftrack.isNonnull() &&
118 fabs(ctftk_p4.Pt() - el_ctftrack->pt())/el_ctftrack->pt() < 0.2)
123 if(el_ctftrack.isNonnull() && gsfElectron.
ctfGsfOverlap() > minFracSharedHits &&
124 deltaR(el_ctftrack_p4, ctftk_p4) < 0.5 &&
125 (el_ctftrack->charge() + ctftk->charge() == 0) ) {
131 int deltaMissingHits = ctftk->trackerExpectedHitsInner().numberOfHits() - el_ctftrack->trackerExpectedHitsInner().numberOfHits();
141 v_candidatePartners.push_back(convInfo);
147 if(
deltaR(el_gsftrack_p4, ctftk_p4) < 0.5 &&
148 (el_gsftrack->charge() + ctftk->charge() == 0) &&
149 el_gsftrack->ptError()/el_gsftrack->pt() < 0.25) {
151 int deltaMissingHits = ctftk->trackerExpectedHitsInner().numberOfHits() - el_gsftrack->trackerExpectedHitsInner().numberOfHits();
155 convInfo.radiusOfConversion(),
156 convInfo.pointOfConversion(),
162 v_candidatePartners.push_back(convInfo);
169 for(GsfTrackCollection::const_iterator gsftk = gsftracks->begin();
170 gsftk != gsftracks->end(); gsftk++, gsftk_i++) {
173 if(gsfidx == gsftk_i)
179 if(gsftk->ptError()/gsftk->pt() > 0.5)
181 if(gsftk->numberOfValidHits() < 5)
184 if(fabs(gsftk->pt() - el_gsftrack->pt())/el_gsftrack->pt() < 0.25)
190 if(el_ctftrack.isNonnull() && gsfElectron.
ctfGsfOverlap() > minFracSharedHits &&
191 deltaR(el_ctftrack_p4, gsftk_p4) < 0.5 &&
192 (el_ctftrack->charge() + gsftk->charge() == 0)) {
194 int deltaMissingHits = gsftk->trackerExpectedHitsInner().numberOfHits() - el_ctftrack->trackerExpectedHitsInner().numberOfHits();
205 v_candidatePartners.push_back(convInfo);
210 if(
deltaR(el_gsftrack_p4, gsftk_p4) < 0.5 &&
211 (el_gsftrack->charge() + gsftk->charge() == 0) &&
212 (el_gsftrack->ptError()/el_gsftrack_p4.pt() < 0.5)) {
216 int deltaMissingHits = gsftk->trackerExpectedHitsInner().numberOfHits() - el_gsftrack->trackerExpectedHitsInner().numberOfHits();
226 v_candidatePartners.push_back(convInfo);
231 return v_candidatePartners;
239 const double bFieldAtOrigin) {
241 using namespace reco;
245 double elCurvature = -0.3*bFieldAtOrigin*(el_track->
charge()/el_tk_p4.pt())/100.;
246 double rEl = fabs(1./elCurvature);
247 double xEl = -1*(1./elCurvature - el_track->
d0())*
sin(el_tk_p4.phi());
248 double yEl = (1./elCurvature - el_track->
d0())*
cos(el_tk_p4.phi());
252 double candCurvature = -0.3*bFieldAtOrigin*(candPartnerTk->
charge()/cand_p4.pt())/100.;
253 double rCand = fabs(1./candCurvature);
254 double xCand = -1*(1./candCurvature - candPartnerTk->
d0())*
sin(cand_p4.phi());
255 double yCand = (1./candCurvature - candPartnerTk->
d0())*
cos(cand_p4.phi());
257 double d =
sqrt(
pow(xEl-xCand, 2) +
pow(yEl-yCand , 2));
258 double dist = d - (rEl + rCand);
259 double dcot = 1./
tan(el_tk_p4.theta()) - 1./
tan(cand_p4.theta());
262 double xa1 = xEl + (xCand-xEl) * rEl/d;
263 double xa2 = xCand + (xEl-xCand) * rCand/d;
264 double ya1 = yEl + (yCand-yEl) * rEl/d;
265 double ya2 = yCand + (yEl-yCand) * rCand/d;
267 double x=.5*(xa1+xa2);
268 double y=.5*(ya1+ya2);
270 double z = el_track->
dz() + rEl*el_track->
pz()*TMath::ACos(1-
pow(rconv,2)/(2.*
pow(rEl,2)))/el_track->
pt();
275 float tempsign = el_track->
px()*
x + el_track->
py()*
y;
276 tempsign = tempsign/fabs(tempsign);
277 rconv = tempsign*rconv;
301 if(v_convCandidates.size() == 1)
302 return v_convCandidates.at(0);
307 for(
unsigned int i = 1;
i < v_convCandidates.size();
i++) {
312 arbitratedConvInfo =
temp;
317 return arbitratedConvInfo;
326 if(v_convCandidates.size() == 0)
333 if(v_convCandidates.size() == 1)
334 return v_convCandidates.at(0);
336 vector<ConversionInfo> v_0;
337 vector<ConversionInfo> v_1;
338 vector<ConversionInfo> v_2;
339 vector<ConversionInfo> v_3;
341 for(
unsigned int i = 1;
i < v_convCandidates.size();
i++) {
344 if(temp.
flag() == 0) {
346 if(fabs(temp.
dist()) < 0.02 &&
347 fabs(temp.
dcot()) < 0.02 &&
360 if(temp.
flag() == 1) {
367 if(temp.
flag() == 2) {
375 if(temp.
flag() == 3) {
418 int trk1_q,
float trk1_d0,
420 int trk2_q,
float trk2_d0,
421 float bFieldAtOrigin) {
424 double tk1Curvature = -0.3*bFieldAtOrigin*(trk1_q/trk1_p4.pt())/100.;
425 double rTk1 = fabs(1./tk1Curvature);
426 double xTk1 = -1.*(1./tk1Curvature - trk1_d0)*
sin(trk1_p4.phi());
427 double yTk1 = (1./tk1Curvature - trk1_d0)*
cos(trk1_p4.phi());
430 double tk2Curvature = -0.3*bFieldAtOrigin*(trk2_q/trk2_p4.pt())/100.;
431 double rTk2 = fabs(1./tk2Curvature);
432 double xTk2 = -1.*(1./tk2Curvature - trk2_d0)*
sin(trk2_p4.phi());
433 double yTk2 = (1./tk2Curvature - trk2_d0)*
cos(trk2_p4.phi());
436 double dist =
sqrt(
pow(xTk1-xTk2, 2) +
pow(yTk1-yTk2 , 2));
437 dist = dist - (rTk1 + rTk2);
439 double dcot = 1./
tan(trk1_p4.theta()) - 1./
tan(trk2_p4.theta());
441 return std::make_pair(dist, dcot);
449 const double bFieldAtOrigin,
450 const double minFracSharedHits) {
453 using namespace reco;
462 if(el_ctftrack.isNonnull() && gsfElectron.
shFracInnerHits() > minFracSharedHits) {
463 ctfidx =
static_cast<int>(el_ctftrack.key());
475 LorentzVector el_tk_p4(el_track->px(), el_track->py(), el_track->pz(), el_track->p());
479 double mindcot = 9999.;
484 for(TrackCollection::const_iterator tk = ctftracks->begin();
485 tk != ctftracks->end(); tk++, tk_i++) {
493 double dR =
deltaR(el_tk_p4, tk_p4);
499 if(tk->charge() + el_track->charge() != 0)
502 double dcot = fabs(1./
tan(tk_p4.theta()) - 1./
tan(el_tk_p4.theta()));
510 if(!candCtfTrackRef.isNonnull())
519 double elCurvature = -0.3*bFieldAtOrigin*(el_track->charge()/el_tk_p4.pt())/100.;
520 double rEl = fabs(1./elCurvature);
521 double xEl = -1*(1./elCurvature - el_track->d0())*
sin(el_tk_p4.phi());
522 double yEl = (1./elCurvature - el_track->d0())*
cos(el_tk_p4.phi());
525 LorentzVector cand_p4 =
LorentzVector(candCtfTrackRef->px(), candCtfTrackRef->py(),candCtfTrackRef->pz(), candCtfTrackRef->p());
526 double candCurvature = -0.3*bFieldAtOrigin*(candCtfTrackRef->charge()/cand_p4.pt())/100.;
527 double rCand = fabs(1./candCurvature);
528 double xCand = -1*(1./candCurvature - candCtfTrackRef->d0())*
sin(cand_p4.phi());
529 double yCand = (1./candCurvature - candCtfTrackRef->d0())*
cos(cand_p4.phi());
531 double d =
sqrt(
pow(xEl-xCand, 2) +
pow(yEl-yCand , 2));
532 double dist = d - (rEl + rCand);
533 double dcot = 1./
tan(el_tk_p4.theta()) - 1./
tan(cand_p4.theta());
536 double xa1 = xEl + (xCand-xEl) * rEl/d;
537 double xa2 = xCand + (xEl-xCand) * rCand/d;
538 double ya1 = yEl + (yCand-yEl) * rEl/d;
539 double ya2 = yCand + (yEl-yCand) * rCand/d;
541 double x=.5*(xa1+xa2);
542 double y=.5*(ya1+ya2);
544 double z = el_track->dz() + rEl*el_track->pz()*TMath::ACos(1-
pow(rconv,2)/(2.*
pow(rEl,2)))/el_track->pt();
549 float tempsign = el_track->px()*
x + el_track->py()*
y;
550 tempsign = tempsign/fabs(tempsign);
551 rconv = tempsign*rconv;
553 int deltaMissingHits = -9999;
554 deltaMissingHits = candCtfTrackRef->trackerExpectedHitsInner().numberOfHits() - el_track->trackerExpectedHitsInner().numberOfHits();
edm::Ref< GsfTrackCollection > GsfTrackRef
persistent reference to a GsfTrack
double p() const
momentum vector magnitude
ConversionInfo findBestConversionMatch(const std::vector< ConversionInfo > &v_convCandidates)
double radiusOfConversion() const
double d0() const
dxy parameter in perigee convention (d0 = - dxy)
float ctfGsfOverlap() const
const reco::Track * getElectronTrack(const reco::GsfElectron &, const float minFracSharedHits=0.45)
ConversionInfo getConversionInfo(const reco::GsfElectronCore &, const edm::Handle< reco::TrackCollection > &ctftracks_h, const edm::Handle< reco::GsfTrackCollection > &gsftracks_h, const double bFieldAtOrigin, const double minFracSharedHits=0.45)
bool isFromConversion(const ConversionInfo &, double maxAbsDist=0.02, double maxAbsDcot=0.02)
Sin< T >::type sin(const T &t)
math::XYZTLorentzVector LorentzVector
const GsfTrackRef & gsfTrack() const
std::vector< Track > TrackCollection
collection of Tracks
double px() const
x coordinate of momentum vector
TrackRef closestCtfTrackRef() const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
bool isNonnull() const
Checks for non-null.
double pt() const
track transverse momentum
std::vector< GsfTrack > GsfTrackCollection
collection of GsfTracks
Cos< T >::type cos(const T &t)
Tan< T >::type tan(const T &t)
math::XYZPoint pointOfConversion() const
ConversionInfo arbitrateConversionPartnersbyR(const std::vector< ConversionInfo > &v_convCandidates)
int deltaMissingHits() const
double pz() const
z coordinate of momentum vector
std::vector< ConversionInfo > getConversionInfos(const reco::GsfElectronCore &, const edm::Handle< reco::TrackCollection > &ctftracks_h, const edm::Handle< reco::GsfTrackCollection > &gsftracks_h, const double bFieldAtOrigin, const double minFracSharedHits=0.45)
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
double deltaR(double eta1, double eta2, double phi1, double phi2)
float shFracInnerHits() const
GsfTrackRef gsfTrack() const
reference to a GsfTrack
XYZPointD XYZPoint
point in space with cartesian internal representation
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
T const * product() const
virtual GsfElectronCoreRef core() const
TrackRef ctfTrack() const
int charge() const
track electric charge
T const * get() const
Returns C++ pointer to the item.
Power< A, B >::type pow(const A &a, const B &b)
double py() const
y coordinate of momentum vector
math::PtEtaPhiELorentzVectorF LorentzVector