Fork me on GitHub

Changeset 9458496b in git for modules


Ignore:
Timestamp:
Oct 30, 2017, 7:58:06 PM (7 years ago)
Author:
Kevin Pedro <kpedro88@…>
Branches:
ImprovedOutputFile, Timing, dual_readout, llp, master
Children:
ca9e119
Parents:
e40b9cf
Message:

optimize track counting algorithm

File:
1 edited

Legend:

Unmodified
Added
Removed
  • modules/TrackCountingBTagging.cc

    re40b9cf r9458496b  
    117117    fItTrackInputArray->Reset();
    118118    count = 0;
    119     while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
     119    // stop once we have enough tracks
     120    while((track = static_cast<Candidate*>(fItTrackInputArray->Next())) and count < fNtracks)
    120121    {
    121122      const TLorentzVector &trkMomentum = track->Momentum;
     123      tpt = trkMomentum.Pt();
     124      if(tpt < fPtMin) continue;
    122125     
     126      d0 = TMath::Abs(track->D0);
     127      if(d0 > fIPmax) continue;
     128
    123129      dr = jetMomentum.DeltaR(trkMomentum);
    124       tpt = trkMomentum.Pt();
     130      if(dr > fDeltaR) continue;
     131
    125132      xd = track->Xd;
    126133      yd = track->Yd;
    127134      zd = track->Zd;
    128       d0 = TMath::Abs(track->D0);
    129135      dd0 = TMath::Abs(track->ErrorD0);
    130136      dz = TMath::Abs(track->DZ);
    131137      ddz = TMath::Abs(track->ErrorDZ);
    132138
    133       if(tpt < fPtMin) continue;
    134       if(dr > fDeltaR) continue;
    135       if(d0 > fIPmax) continue;
    136 
    137139      if(fUse3D){
    138140        sign = (jpx*xd + jpy*yd + jpz*zd > 0.0) ? 1 : -1;
    139         //add transvers and longitudinal significances in quadrature
     141        //add transverse and longitudinal significances in quadrature
    140142        sip = sign * TMath::Sqrt( TMath::Power(d0 / dd0, 2) + TMath::Power(dz / ddz, 2) );
    141143      }
Note: See TracChangeset for help on using the changeset viewer.