Fork me on GitHub

Ignore:
Timestamp:
Apr 17, 2014, 10:54:17 AM (11 years ago)
Author:
Pavel Demin
Message:

several minor corrections

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/modules/TrackCountingBTagging.cc

    r1367 r1372  
    22/** \class TrackCountingBTagging
    33 *
    4  *  btagging algorithm based on counting tracks with large impact parameterDetermines origin of jet,
     4 *  b-tagging algorithm based on counting tracks with large impact parameter
    55 *
    66 *  $Date: 2014-03-27 12:39:14 +0200 (Fri, 27 March 2014) $
     
    2525#include "TLorentzVector.h"
    2626
    27 #include <algorithm> 
     27#include <algorithm>
    2828#include <stdexcept>
    2929#include <iostream>
     
    3131
    3232using namespace std;
    33 
    3433
    3534//------------------------------------------------------------------------------
     
    5049void TrackCountingBTagging::Init()
    5150{
    52  
    5351  fBitNumber = GetInt("BitNumber", 0);
    54  
     52
    5553  fPtMin     = GetDouble("TrackPtMin", 1.0);
    5654  fDeltaR    = GetDouble("DeltaR", 0.3);
    5755  fIPmax     = GetDouble("TrackIPMax", 2.0);
    58  
     56
    5957  fSigMin    = GetDouble("SigMin", 6.5);
    60   fNtracks   = GetInt("Ntracks", 3); 
    61    
     58  fNtracks   = GetInt("Ntracks", 3);
     59
    6260  // import input array(s)
    6361
    6462  fTrackInputArray = ImportArray(GetString("TrackInputArray", "Calorimeter/eflowTracks"));
    6563  fItTrackInputArray = fTrackInputArray->MakeIterator();
    66  
     64
    6765  fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
    6866  fItJetInputArray = fJetInputArray->MakeIterator();
     
    8280{
    8381  Candidate *jet, *track;
    84  
     82
    8583  Double_t jpx, jpy;
    8684  Double_t dr, tpx, tpy, tpt;
    8785  Double_t xd, yd, dxy, ddxy, ip, sip;
    88  
     86
    8987  Int_t sign;
    90  
     88
    9189  Int_t count;
    92  
     90
    9391  // loop over all input jets
    9492  fItJetInputArray->Reset();
     
    9896    jpx = jetMomentum.Px();
    9997    jpy = jetMomentum.Py();
    100    
     98
    10199    // loop over all input tracks
    102100    fItTrackInputArray->Reset();
     
    104102    while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
    105103    {
    106      const TLorentzVector &trkMomentum = track->Momentum;
    107      
    108      dr = jetMomentum.DeltaR(trkMomentum);
    109      
    110      tpt = trkMomentum.Pt();
    111      tpx = trkMomentum.Px();
    112      tpy = trkMomentum.Py();
    113    
    114      xd   = track->Xd;
    115      yd   = track->Yd;
    116      dxy  = TMath::Abs(track->Dxy);
    117      ddxy = track->SDxy;
    118      
    119      if( tpt < fPtMin ) continue;
    120      if( dr  > fDeltaR ) continue;
    121      if( dxy > fIPmax ) continue;
    122        
    123      sign  = ( jpx*xd + jpy*yd > 0.0 ) ? 1 : -1;
    124        
    125      ip  = sign*dxy;
    126      sip = ip / TMath::Abs(ddxy);
    127      
    128      if( sip > fSigMin ) count++;
    129      
     104      const TLorentzVector &trkMomentum = track->Momentum;
     105
     106      dr = jetMomentum.DeltaR(trkMomentum);
     107
     108      tpt = trkMomentum.Pt();
     109      tpx = trkMomentum.Px();
     110      tpy = trkMomentum.Py();
     111
     112      xd   = track->Xd;
     113      yd   = track->Yd;
     114      dxy  = TMath::Abs(track->Dxy);
     115      ddxy = track->SDxy;
     116
     117      if(tpt < fPtMin) continue;
     118      if(dr  > fDeltaR) continue;
     119      if(dxy > fIPmax) continue;
     120
     121      sign  = (jpx*xd + jpy*yd > 0.0) ? 1 : -1;
     122
     123      ip  = sign*dxy;
     124      sip = ip / TMath::Abs(ddxy);
     125
     126      if(sip > fSigMin) count++;
    130127    }
    131      
     128
    132129    // set BTag flag to true if count >= Ntracks
    133     jet->BTag |= ( count >= fNtracks ) << fBitNumber;
     130    jet->BTag |= (count >= fNtracks) << fBitNumber;
    134131  }
    135132}
Note: See TracChangeset for help on using the changeset viewer.