Belle II Software development
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:

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.
 
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).
 
bool trackSelector (const RecoTrack &track)
 This method return true if every segment (see segmentSelector) of the input track respects the cuts contraints.
 
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.
 
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).
 
void initHistoNoKick (bool outHisto)
 This metod initialize some validation histograms of the Training Sample Selection.
 
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.

Constructor & Destructor Documentation

◆ NoKickRTSel() [1/2]

NoKickRTSel ( const std::string &  fileName,
bool  outputHisto 
)
inline

Constructor with input file for use specific cuts file and allows validation.

Definition at line 60 of file NoKickRTSel.h.

60 :
61 m_trackCuts(fileName)
62 {
63 m_outputFlag = false;
65 initHistoNoKick(outputHisto);
66 }
NoKickCuts m_trackCuts
auxiliary member to apply the cuts
Definition: NoKickRTSel.h:39
void initHistoNoKick(bool outHisto)
This metod initialize some validation histograms of the Training Sample Selection.
Definition: NoKickRTSel.h:118
void initNoKickRTSel()
Inizialize the class cleaning the member vectors.
Definition: NoKickRTSel.h:76
bool m_outputFlag
true=produce validation output
Definition: NoKickRTSel.h:42

◆ NoKickRTSel() [2/2]

NoKickRTSel ( )
inline

Empty Constructor that uses the defaults cuts file.

Definition at line 69 of file NoKickRTSel.h.

69 :
71 {
73 }

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.

97{
98 hitXPBuilder(track);
99 for (const hitXP& XP : m_hitXP) {
100 if (m_8hitTrack.size() < 1) {
101 m_8hitTrack.push_back(XP);
102 } else if (XP.m_sensorLayer != m_8hitTrack.back().m_sensorLayer ||
103 XP.m_sensorLadder != m_8hitTrack.back().m_sensorLadder ||
104 XP.m_sensorSensor != m_8hitTrack.back().m_sensorSensor) {
105 m_8hitTrack.push_back(XP);
106 }
107 }
108
109}
std::vector< hitXP > m_8hitTrack
vector of selected hit
Definition: NoKickRTSel.h:38
void hitXPBuilder(const RecoTrack &track)
this method build a vector of hitXP from a track.
Definition: NoKickRTSel.cc:21
std::vector< hitXP > m_hitXP
vector of hit, to convert the track
Definition: NoKickRTSel.h:36

◆ 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.

22{
23 StoreArray<SVDCluster> SVDClusters;
24 StoreArray<SVDTrueHit> SVDTrueHits;
25 StoreArray<MCParticle> MCParticles;
26 StoreArray<RecoTrack> recoTracks;
27
28 if (track.getRelationsTo<MCParticle>().size() > 0) {
29
30 const MCParticle* particle = track.getRelationsTo<MCParticle>()[0];
31
32 std::vector<Belle2::RecoHitInformation::UsedSVDHit*> clusterListSVD = track.getSVDHitList();
33 for (const SVDCluster* cluster : clusterListSVD) {
34 for (const SVDTrueHit& hit : cluster->getRelationsTo<SVDTrueHit>()) {
35 RelationVector<MCParticle> hitFromParticle = hit.getRelationsFrom<MCParticle>();
36 if (hitFromParticle.size() > 0) {
37 if (hitFromParticle[0]->getIndex() != particle->getIndex()) {
38 continue;
39 }
40 } else continue;
41 VxdID trueHitSensorID = hit.getSensorID();
42 const VXD::SensorInfoBase& sensorInfo = VXD::GeoCache::getInstance().getSensorInfo(trueHitSensorID);
43 hitXPDerivate entry(hit, *cluster, *particle, sensorInfo);
44 int NClusterU = 0;
45 int NClusterV = 0;
46 for (SVDCluster Ncluster : hit.getRelationsFrom<SVDCluster>()) {
47 if (Ncluster.isUCluster()) NClusterU++;
48 else NClusterV++;
49 }
50 entry.setClusterU(NClusterU);
51 entry.setClusterV(NClusterV);
52
53 bool isReconstructed(false);
54 for (const RecoTrack& aRecoTrack : particle->getRelationsFrom<RecoTrack>())
55 isReconstructed |= aRecoTrack.hasSVDHits();
56 entry.setReconstructed(isReconstructed);
57
58 m_setHitXP.insert(entry);
59 }
60 }
61
62
63 StoreArray<PXDCluster> PXDClusters;
64 StoreArray<PXDTrueHit> PXDTrueHits;
65
66 std::vector<Belle2::RecoHitInformation::UsedPXDHit*> clusterListPXD = track.getPXDHitList();
67 for (const PXDCluster* cluster : clusterListPXD) {
68 for (const PXDTrueHit& hit : cluster->getRelationsTo<PXDTrueHit>()) {
69 RelationVector<MCParticle> hitFromParticle = hit.getRelationsFrom<MCParticle>();
70 if (hitFromParticle.size() > 0) {
71 if (hitFromParticle[0]->getIndex() != particle->getIndex()) {
72 continue;
73 }
74 } else continue;
75 VxdID trueHitSensorID = hit.getSensorID();
76 const VXD::SensorInfoBase& sensorInfo = VXD::GeoCache::getInstance().getSensorInfo(trueHitSensorID);
77 hitXPDerivate entry(hit, *particle, sensorInfo);
78
79 bool isReconstructed(false);
80 for (const RecoTrack& aRecoTrack : particle->getRelationsFrom<RecoTrack>())
81 isReconstructed |= aRecoTrack.hasPXDHits();
82 entry.setReconstructed(isReconstructed);
83
84 m_setHitXP.insert(entry);
85 }
86 }
87
88 for (auto element : m_setHitXP) {
89 m_hitXP.push_back(element);
90 }
91
92 }
93
94}
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
std::set< hitXP, hitXP::timeCompare > m_setHitXP
set of hit to order the hit in time
Definition: NoKickRTSel.h:37
The PXD Cluster class This class stores all information about reconstructed PXD clusters The position...
Definition: PXDCluster.h:30
Class PXDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: PXDTrueHit.h:31
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
Class for type safe access to objects that are referred to in relations.
size_t size() const
Get number of relations.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:29
Class SVDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: SVDTrueHit.h:33
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
const SensorInfoBase & getSensorInfo(Belle2::VxdID id) const
Return a referecne to the SensorInfo of a given SensorID.
Definition: GeoCache.cc:67
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
Base class to provide Sensor Information for PXD and SVD.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
This class is the derivate of HitXP, and complete it with a constructor that use all other complex ty...
Definition: hitXPDerivate.h:27

◆ 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
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

◆ initNoKickRTSel()

void initNoKickRTSel ( )
inline

Inizialize the class cleaning the member vectors.

Definition at line 76 of file NoKickRTSel.h.

77 {
78 m_hitXP.clear();
79 m_setHitXP.clear();
80 m_8hitTrack.clear();
82 }
int m_numberOfCuts
number of catastrophic interaction for each track
Definition: NoKickRTSel.h:41

◆ produceHistoNoKick()

void produceHistoNoKick ( )

This method produce the validation histograms (to be used the endrun combined with the filling in trackSelector method)

Definition at line 277 of file NoKickRTSel.cc.

278{
279 if (m_outputFlag) {
281 m_momSel->Write();
282 m_momCut->Write();
283
284 m_momEff->Add(m_momSel, 1);
285 m_momEff->Add(m_momCut, 1);
286 m_momEff->Divide(m_momSel, m_momEff, 1, 1);
287 m_momEff->Write();
288
289 m_PDGIDSel->Write();
290 m_PDGIDCut->Write();
291
292 m_PDGIDEff->Add(m_PDGIDSel, 1);
293 m_PDGIDEff->Add(m_PDGIDCut, 1);
294 m_PDGIDEff->Divide(m_PDGIDSel, m_PDGIDEff, 1, 1);
295 m_PDGIDEff->Write();
296
297 m_nCutHit->Write();
298
299 m_noKickTree->Write();
300
301 delete m_noKickOutputTFile;
302 }
303}

◆ 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.

136{
137 if (hit2.m_sensorLayer - hit1.m_sensorLayer > 1) return true;
138 else {
139 double deltaPar = 0;
140 switch (par) {
141 case NoKickCuts::c_Omega:
142 deltaPar = fabs(hit1.getOmegaEntry() - hit2.getOmegaEntry());
143 if (is0) deltaPar = fabs(hit1.getOmega0() - hit2.getOmegaEntry());
144 // selCutPXD =0.4;
145 break;
146
147 case NoKickCuts::c_D0:
148 deltaPar = hit1.getD0Entry() - hit2.getD0Entry();
149 if (is0) deltaPar = hit1.getD00() - hit2.getD0Entry();
150 // selCutPXD =1;
151 break;
152
153 case NoKickCuts::c_Phi0:
154 deltaPar = asin(sin(hit1.getPhi0Entry())) - asin(sin(hit2.getPhi0Entry()));
155 if (is0) deltaPar = asin(sin(hit1.getPhi00())) - asin(sin(hit2.getPhi0Entry()));
156 // selCutPXD =0.3;
157 break;
158
159 case NoKickCuts::c_Z0:
160 deltaPar = hit1.getZ0Entry() - hit2.getZ0Entry();
161 if (is0) deltaPar = hit1.getZ00() - hit2.getZ0Entry();
162 // selCutPXD =1;
163 break;
164
165 case NoKickCuts::c_Tanlambda:
166 deltaPar = hit1.getTanLambdaEntry() - hit2.getTanLambdaEntry();
167 if (is0) deltaPar = hit1.getTanLambda0() - hit2.getTanLambdaEntry();
168 // selCutPXD =0.3;
169 break;
170 }
171
172 double usedCut = 0;
173 if (fabs(selCut.at(0)) > fabs(selCut.at(1))) {
174 usedCut = fabs(selCut.at(0));
175 } else usedCut = fabs(selCut.at(1));
176
177 if (deltaPar > -usedCut && deltaPar < usedCut) return true;
178 else {
179 B2DEBUG(20, "--------------------------");
180 B2DEBUG(20, "lay1=" << hit1.m_sensorLayer);
181 B2DEBUG(20, "lay2=" << hit2.m_sensorLayer);
182 B2DEBUG(20, "parameter=" << par);
183 B2DEBUG(20, "Min=" << selCut.at(0));
184 B2DEBUG(20, "Max=" << selCut.at(1));
185 B2DEBUG(20, "deltaPar=" << deltaPar);
186 B2DEBUG(20, "momentum=" << hit1.m_momentum0.R());
187 return false;
188 }
189 }
190}
double getOmegaEntry() const
evaluate relative parameter using entrypoint position and momentum
Definition: hitXP.h:325
double getPhi00() const
evaluate relative parameter using IP position and momentum
Definition: hitXP.h:367
double getZ0Entry() const
evaluate relative parameter using entrypoint position and momentum
Definition: hitXP.h:373
double getOmega0() const
evaluate relative parameter using IP position and momentum
Definition: hitXP.h:331
int m_sensorLayer
layer of the hit
Definition: hitXP.h:52
double getD00() const
evaluate relative parameter using IP position and momentum
Definition: hitXP.h:355
double getPhi0Entry() const
evaluate relative parameter using entrypoint position and momentum
Definition: hitXP.h:361
double getZ00() const
evaluate relative parameter using IP position and momentum
Definition: hitXP.h:379
double getTanLambda0() const
evaluate relative parameter using IP position and momentum
Definition: hitXP.h:343
double getD0Entry() const
evaluate relative parameter using entrypoint position and momentum
Definition: hitXP.h:349
double getTanLambdaEntry() const
evaluate relative parameter using entrypoint position and momentum
Definition: hitXP.h:337
ROOT::Math::XYZVector m_momentum0
momentum at IP
Definition: hitXP.h:49

◆ 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.

195{
197 hit8TrackBuilder(track);
198
199 if (m_outputFlag) {
200 m_pdgID = m_8hitTrack[0].getPDGID();
201 m_pMag = track.getMomentumSeed().R();
202 m_pt = sqrt(track.getMomentumSeed().X() * track.getMomentumSeed().X() + track.getMomentumSeed().Y() * track.getMomentumSeed().Y());
203 }
204
205 bool good = globalCut(m_8hitTrack);
206 if (good == false) {
207 if (m_outputFlag) {
208 m_numberOfCuts = -1; //it means "no specific cuts applied, but rejected for global cuts"
210 m_isCutted = true;
212 m_momCut->Fill(track.getMomentumSeed().R());
213 m_PDGIDCut->Fill(track.getRelationsTo<MCParticle>()[0]->getPDG());
214 m_noKickTree->Fill();
215 }
216 return false;
217 }
218
219 if (track.getMomentumSeed().R() > m_pmax) {
220 if (m_outputFlag) {
222 m_isCutted = false;
224 m_momSel->Fill(track.getMomentumSeed().R());
225 m_PDGIDSel->Fill(track.getRelationsTo<MCParticle>()[0]->getPDG());
226 m_noKickTree->Fill();
227 }
228 return good;
229 }
230 for (int i = 0; i < (int)(m_8hitTrack.size() - 2); i++) {
231 double sinTheta = fabs(m_8hitTrack.at(i).m_momentum0.Y()) /
232 sqrt(pow(m_8hitTrack.at(i).m_momentum0.Y(), 2) +
233 pow(m_8hitTrack.at(i).m_momentum0.Z(), 2));
234
235 double momentum =
236 sqrt(pow(m_8hitTrack.at(i).m_momentum0.X(), 2) +
237 pow(m_8hitTrack.at(i).m_momentum0.Y(), 2) +
238 pow(m_8hitTrack.at(i).m_momentum0.Z(), 2));
239
240 for (int j = NoKickCuts::c_Omega; j <= NoKickCuts::c_Tanlambda; j++) { //track parameters loop
241 std::vector<double> selCut = m_trackCuts.cutSelector(sinTheta, momentum, m_8hitTrack.at(i).m_sensorLayer,
242 m_8hitTrack.at(i + 1).m_sensorLayer, (NoKickCuts::EParameters) j);
243 bool goodSeg = segmentSelector(m_8hitTrack.at(i), m_8hitTrack.at(i + 1), selCut, (NoKickCuts::EParameters) j);
244
245 if (!goodSeg) {
246 good = false;
248 }
249 if (i == 0) { //beampipe crossing
250 bool goodSeg0 = segmentSelector(m_8hitTrack.at(i), m_8hitTrack.at(i), selCut, (NoKickCuts::EParameters) j, true);
251 if (!goodSeg0) {
252 good = false;
254 }
255 }
256
257 }
258 }
259 if (m_outputFlag) {
262 if (good) {
263 m_isCutted = false;
264 m_momSel->Fill(track.getMomentumSeed().R());
265 m_PDGIDSel->Fill(track.getRelationsTo<MCParticle>()[0]->getPDG());
266 } else {
267 m_isCutted = true;
268 m_momCut->Fill(track.getMomentumSeed().R());
269 m_PDGIDCut->Fill(track.getRelationsTo<MCParticle>()[0]->getPDG());
270 }
271 m_noKickTree->Fill();
272 }
273 return good;
274}
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
EParameters
enum for parameters name
Definition: NoKickCuts.h:44
std::vector< double > cutSelector(double sintheta, double momentum, int layer1, int layer2, EParameters par)
This methods selects 2 cuts (minimum and maximum inside a vector) from the information of theta,...
Definition: NoKickCuts.cc:18
double m_pmax
range analyzed with cuts
Definition: NoKickRTSel.h:40
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).
Definition: NoKickRTSel.cc:111
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.
Definition: NoKickRTSel.cc:135
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 f...
Definition: NoKickRTSel.cc:96
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

Member Data Documentation

◆ m_8hitTrack

std::vector<hitXP> m_8hitTrack

vector of selected hit

Definition at line 38 of file NoKickRTSel.h.

◆ m_hitXP

std::vector<hitXP> m_hitXP

vector of hit, to convert the track

Definition at line 36 of file NoKickRTSel.h.

◆ m_isCutted

bool m_isCutted

Indicator if cut is applied.

Definition at line 52 of file NoKickRTSel.h.

◆ m_momCut

TH1F* m_momCut

histrogram of cutted tracks

Definition at line 46 of file NoKickRTSel.h.

◆ m_momEff

TH1F* m_momEff

histogram for efficiency

Definition at line 47 of file NoKickRTSel.h.

◆ m_momSel

TH1F* m_momSel

histogram of selected tracks

Definition at line 45 of file NoKickRTSel.h.

◆ m_nCutHit

TH1F* m_nCutHit

histogram for number of cutted hist per track

Definition at line 51 of file NoKickRTSel.h.

◆ m_Ncuts

int m_Ncuts

number of times the cut is applied on a particle

Definition at line 56 of file NoKickRTSel.h.

◆ m_noKickOutputTFile

TFile* m_noKickOutputTFile

validartion output TFile

Definition at line 44 of file NoKickRTSel.h.

◆ m_noKickTree

TTree* m_noKickTree

TTree to which the information is written.

Definition at line 57 of file NoKickRTSel.h.

◆ m_numberOfCuts

int m_numberOfCuts

number of catastrophic interaction for each track

Definition at line 41 of file NoKickRTSel.h.

◆ m_outputFlag

bool m_outputFlag

true=produce validation output

Definition at line 42 of file NoKickRTSel.h.

◆ m_pdgID

double m_pdgID

pdg Code

Definition at line 55 of file NoKickRTSel.h.

◆ m_PDGIDCut

TH1F* m_PDGIDCut

histogram for PDGID of cutted track

Definition at line 48 of file NoKickRTSel.h.

◆ m_PDGIDEff

TH1F* m_PDGIDEff

histogram for efficiency for each PDGID

Definition at line 50 of file NoKickRTSel.h.

◆ m_PDGIDSel

TH1F* m_PDGIDSel

histogram for PDGID of selected track

Definition at line 49 of file NoKickRTSel.h.

◆ m_pMag

double m_pMag

momentum magnitut

Definition at line 53 of file NoKickRTSel.h.

◆ m_pmax

double m_pmax = 10.

range analyzed with cuts

Definition at line 40 of file NoKickRTSel.h.

◆ m_pt

double m_pt

transverse momentum

Definition at line 54 of file NoKickRTSel.h.

◆ m_setHitXP

std::set<hitXP, hitXP::timeCompare> m_setHitXP

set of hit to order the hit in time

Definition at line 37 of file NoKickRTSel.h.

◆ m_trackCuts

NoKickCuts m_trackCuts

auxiliary member to apply the cuts

Definition at line 39 of file NoKickRTSel.h.


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