Belle II Software  release-05-02-19
NoKickRTSel Class Reference

This class implement some methods useful for the application of cuts evaluated in NoKickCutsEval module. More...

#include <NoKickRTSel.h>

Inheritance diagram for NoKickRTSel:
Collaboration diagram for NoKickRTSel:

Public Member Functions

 NoKickRTSel (const std::string &fileName, bool outputHisto)
 Constructor with input file for use specific cuts file and allows validation.
 
 NoKickRTSel ()
 Empty Constructor that uses the defaults cuts file.
 
void initNoKickRTSel ()
 Inizialize the class cleaning the member vectors.
 
void hitXPBuilder (const RecoTrack &track)
 this method build a vector of hitXP from a track. More...
 
void hit8TrackBuilder (const RecoTrack &track)
 this metod build a vector of hitXP from a track selecting the first hit on each layer of VXD (8 hit for SVD only, counting overlaps). More...
 
bool trackSelector (const RecoTrack &track)
 This method return true if every segment (see segmentSelector) of the input track respects the cuts contraints. More...
 
bool segmentSelector (hitXP hit1, hitXP hit2, std::vector< double > selCut, NoKickCuts::EParameters par, bool is0=false)
 This method return true if a couple of hits resects the cuts constraints. More...
 
bool globalCut (const std::vector< hitXP > &track8)
 This method make some global cuts on the tracks (layer 3 and 6 required, d0 and z0 inside beam pipe). More...
 
void initHistoNoKick (bool outHisto)
 This metod initialize some validation histograms of the Training Sample Selection. More...
 
void produceHistoNoKick ()
 This method produce the validation histograms (to be used the endrun combined with the filling in trackSelector method)
 
 ClassDef (NoKickRTSel, 1)
 Making this class a ROOT class.
 

Public Attributes

std::vector< hitXPm_hitXP
 vector of hit, to convert the track
 
std::set< hitXP, hitXP::timeComparem_setHitXP
 set of hit to order the hit in time
 
std::vector< hitXPm_8hitTrack
 vector of selected hit
 
NoKickCuts m_trackCuts
 auxiliary member to apply the cuts
 
double m_pmax = 10.
 range analyzed with cuts
 
int m_numberOfCuts
 number of catastrophic interaction for each track
 
bool m_outputFlag
 true=produce validation output
 
TFile * m_noKickOutputTFile
 validartion output TFile
 
TH1F * m_momSel
 histogram of selected tracks
 
TH1F * m_momCut
 histrogram of cutted tracks
 
TH1F * m_momEff
 histogram for efficiency
 
TH1F * m_PDGIDCut
 histogram for PDGID of cutted track
 
TH1F * m_PDGIDSel
 histogram for PDGID of selected track
 
TH1F * m_PDGIDEff
 histogram for efficiency for each PDGID
 
TH1F * m_nCutHit
 histogram for number of cutted hist per track
 
bool m_isCutted
 Indicator if cut is applied.
 
double m_pMag
 momentum magnitut
 
double m_pt
 transverse momentum
 
double m_pdgID
 pdg Code
 
int m_Ncuts
 number of times the cut is applied on a particle
 
TTree * m_noKickTree
 TTree to which the information is written.
 

Detailed Description

This class implement some methods useful for the application of cuts evaluated in NoKickCutsEval module.

Use and auxiliary class (NoKickCuts) that cointains the cuts used in selection.

Definition at line 43 of file NoKickRTSel.h.

Member Function Documentation

◆ globalCut()

bool globalCut ( const std::vector< hitXP > &  track8)

This method make some global cuts on the tracks (layer 3 and 6 required, d0 and z0 inside beam pipe).

Return false if this filter fails. input (the selected hit of the track)

Definition at line 113 of file NoKickRTSel.cc.

114 {
115  int flagd0 = 1;
116  int flagz0 = 1;
117  int lay3 = 0;
118  int lay4 = 0;
119  int lay5 = 0;
120  int lay6 = 0;
121  for (hitXP XP : track8) {
122  if (XP.getSensorLayer() == 3) lay3 = 1;
123  if (XP.getSensorLayer() == 4) lay4 = 1;
124  if (XP.getSensorLayer() == 5) lay5 = 1;
125  if (XP.getSensorLayer() == 6) lay6 = 1;
126  // if (fabs(XP.getD0Entry()) > 1.) flagd0 = 0;
127  // if (fabs(XP.getZ0Entry()) > 1.) flagz0 = 0;
128  }
129  int N_lay = lay3 + lay4 + lay5 + lay6;
130  if (N_lay >= 3) N_lay = 1;
131  else N_lay = 0;
132  int flagTot = flagd0 * flagz0 * N_lay;
133  if (flagTot == 1) return true;
134  else return false;
135 }

◆ hit8TrackBuilder()

void hit8TrackBuilder ( const RecoTrack track)

this metod build a vector of hitXP from a track selecting the first hit on each layer of VXD (8 hit for SVD only, counting overlaps).

The ouput is the member of the class. input (one reconstructed track)

Definition at line 98 of file NoKickRTSel.cc.

◆ hitXPBuilder()

void hitXPBuilder ( const RecoTrack track)

this method build a vector of hitXP from a track.

The ouput is the member of the class. input (one reconstructed track)

Definition at line 23 of file NoKickRTSel.cc.

◆ initHistoNoKick()

void initHistoNoKick ( bool  outHisto)
inline

This metod initialize some validation histograms of the Training Sample Selection.

The input boolean allows the initialization, otherwise the method is empty (no validation)

Definition at line 128 of file NoKickRTSel.h.

◆ segmentSelector()

bool segmentSelector ( hitXP  hit1,
hitXP  hit2,
std::vector< double >  selCut,
NoKickCuts::EParameters  par,
bool  is0 = false 
)

This method return true if a couple of hits resects the cuts constraints.

input (first hit, second hit, selected cut to apply, track parameter, it is first hit the IP?)

Definition at line 137 of file NoKickRTSel.cc.

◆ trackSelector()

bool trackSelector ( const RecoTrack track)

This method return true if every segment (see segmentSelector) of the input track respects the cuts contraints.

input (one reconstructed track)

Definition at line 196 of file NoKickRTSel.cc.


The documentation for this class was generated from the following files:
Belle2::hitXP
This class collects some information of a TrueHit, using SVDCLuster and MCParticle information too.
Definition: hitXP.h:42