Belle II Software  release-08-01-10
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 33 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 111 of file NoKickRTSel.cc.

112 {
113  int flagd0 = 1;
114  int flagz0 = 1;
115  int lay3 = 0;
116  int lay4 = 0;
117  int lay5 = 0;
118  int lay6 = 0;
119  for (hitXP XP : track8) {
120  if (XP.getSensorLayer() == 3) lay3 = 1;
121  if (XP.getSensorLayer() == 4) lay4 = 1;
122  if (XP.getSensorLayer() == 5) lay5 = 1;
123  if (XP.getSensorLayer() == 6) lay6 = 1;
124  // if (fabs(XP.getD0Entry()) > 1.) flagd0 = 0;
125  // if (fabs(XP.getZ0Entry()) > 1.) flagz0 = 0;
126  }
127  int N_lay = lay3 + lay4 + lay5 + lay6;
128  if (N_lay >= 3) N_lay = 1;
129  else N_lay = 0;
130  int flagTot = flagd0 * flagz0 * N_lay;
131  if (flagTot == 1) return true;
132  else return false;
133 }
This class collects some information of a TrueHit, using SVDCLuster and MCParticle information too.
Definition: hitXP.h:32

◆ 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 96 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 21 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 118 of file NoKickRTSel.h.

119  {
120  if (outHisto) {
121  m_noKickOutputTFile = new TFile("TrackSelection_NoKick.root", "RECREATE");
122  m_momSel = new TH1F("m_momSel", "m_momSel", 100, 0, 4);
123  m_momCut = new TH1F("m_momCut", "m_momCut", 100, 0, 4);
124  m_momEff = new TH1F("m_momEff", "m_momEff", 100, 0, 4);
125 
126  m_PDGIDSel = new TH1F("m_PDGIDSel", "m_PDGIDSel", 6000, -3000, 3000);
127  m_PDGIDCut = new TH1F("m_PDGIDCut", "m_PDGIDCut", 6000, -3000, 3000);
128  m_PDGIDEff = new TH1F("m_PDGIDEff", "m_PDGIDEff", 6000, -3000, 3000);
129 
130  m_nCutHit = new TH1F("m_nCutHit", "m_nCutHit", 30, 0, 30);
131 
132 
133  m_noKickTree = new TTree("noKickTree", "noKickTree");
134  m_noKickTree->Branch("is_rejected", &m_isCutted);
135  m_noKickTree->Branch("p_mag", &m_pMag);
136  m_noKickTree->Branch("pt", &m_pt);
137  m_noKickTree->Branch("pdgID", &m_pdgID);
138  m_noKickTree->Branch("number_of_rejected_SP", &m_Ncuts);
139 
140  m_outputFlag = true;
141  }
142 
143  }
TFile * m_noKickOutputTFile
validartion output TFile
Definition: NoKickRTSel.h:44
int m_Ncuts
number of times the cut is applied on a particle
Definition: NoKickRTSel.h:56
double m_pdgID
pdg Code
Definition: NoKickRTSel.h:55
TH1F * m_momEff
histogram for efficiency
Definition: NoKickRTSel.h:47
TTree * m_noKickTree
TTree to which the information is written.
Definition: NoKickRTSel.h:57
TH1F * m_PDGIDCut
histogram for PDGID of cutted track
Definition: NoKickRTSel.h:48
TH1F * m_PDGIDSel
histogram for PDGID of selected track
Definition: NoKickRTSel.h:49
double m_pMag
momentum magnitut
Definition: NoKickRTSel.h:53
TH1F * m_PDGIDEff
histogram for efficiency for each PDGID
Definition: NoKickRTSel.h:50
TH1F * m_momCut
histrogram of cutted tracks
Definition: NoKickRTSel.h:46
bool m_isCutted
Indicator if cut is applied.
Definition: NoKickRTSel.h:52
bool m_outputFlag
true=produce validation output
Definition: NoKickRTSel.h:42
TH1F * m_nCutHit
histogram for number of cutted hist per track
Definition: NoKickRTSel.h:51
TH1F * m_momSel
histogram of selected tracks
Definition: NoKickRTSel.h:45
double m_pt
transverse momentum
Definition: NoKickRTSel.h:54

◆ 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 135 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 194 of file NoKickRTSel.cc.


The documentation for this class was generated from the following files: