FEDRA emulsion software from the OPERA Collaboration
tr_simple.C File Reference

Functions

void tr_simple (const char *rcut="nseg>2&&t.eFlag>-1")
 

Variables

float dpp = .5
 
EdbDataProcdproc =0
 
EdbPVRecgAli =0
 
EdbVertexRecgEVR =0
 
float momentum = 1.
 
float ProbMinP = 0.01
 
float ProbMinV = 0.01
 
bool ReFit = false
 
bool ReProp = false
 
bool UseSegPar = false
 

Function Documentation

◆ tr_simple()

void tr_simple ( const char *  rcut = "nseg>2&&t.eFlag>-1")
20{
21 dproc = new EdbDataProc("Carbonium_data_set.def");
22 dproc->InitVolume(100, rcut);
23
24 gAli = dproc->PVR();
25
26 int ntr = gAli->eTracks->GetEntries();
27 printf("\n%d tracks are in the list\n",ntr);
28
29 EdbTrackP *tra = 0;
30 int ntrgood = 0;
31 for(int i=0; i<ntr; i++)
32 {
33 tra = ((EdbTrackP*)(gAli->eTracks->At(i)));
34 if (tra->Flag() >= 0) ntrgood++;
35 }
36 printf("\n%d good tracks are in the list\n",ntrgood);
37
38 if (ReFit) // re-fit tracks
39 {
40 for(int i=0; i<ntr; i++)
41 {
42 tra = ((EdbTrackP*)(gAli->eTracks->At(i)));
43 if (tra->Flag() < 0) continue;
47 }
48 gAli->FitTracks(momentum,0.139);
49 }
50
51 gAli->FillCell(50,50,0.015,0.015);
52
53 if (ReProp) // apply additional propagation
54 {
56 ntrgood = 0;
57 for(int i=0; i<ntr; i++)
58 {
59 tra = ((EdbTrackP*)(gAli->eTracks->At(i)));
60 if (tra->Flag() >= 0) ntrgood++;
61 }
62 printf("\n%d good tracks are in the list after additional propagation\n",
63 ntrgood);
64 }
65
66 // tracks loop
67 ntr = gAli->eTracks->GetEntries();
68 ntrgood = 0;
69 int nseg = 0;
70 EdbSegP *seg = 0;
71 for(int i=0; i<ntr; i++)
72 {
73 tra = (EdbTrackP*)(gAli->eTracks->At(i));
74 if (!tra) continue;
75 if (tra->Flag() != -10) // skip "broken" tracks
76 {
77 if (tra->N() < 10) continue; // skip "short" tracks
78 // you can put here additional cuts
79
80 // you can put here track analysis
81 // number of segments is tra->N()
82 // track probability is tra->Prob()
83 // etc
84 ntrgood++;
85 nseg = tra->N(); // number of basetracks (segments) in this track
86 for (int j=0; j<nseg; j++) // j = 0 - segment with minimal Z
87 { // = nseg-1 - segment with maximal Z
88 seg = tra->GetSegment(j);
89 if (!seg) continue;
90 // you can put here basetracks analysis
91 }
92 }
93 }
94 printf("%d good tracks was found\n",ntrgood);
95
96 gEVR = new EdbVertexRec();
97 if (!(gAli->eTracks)) gAli->eTracks = new TObjArray(100);
99 if (!(gAli->eVTX)) gAli->eVTX = new TObjArray(100);
100 gEVR->eVTX = gAli->eVTX;
101 gEVR->ePVR = gAli;
102
103 gEVR->eDZmax = 3000.; // default is 3000. microns
104 gEVR->eProbMin = ProbMinV; // minimal vertex probability,
105 // default 0.01, i.e. 1%
106 gEVR->eImpMax = 55.; // maximal impact parameter-distance between
107 // tracks, for 2-tracks vertices, default is 25. microns
108 gEVR->eUseSegPar = UseSegPar; // default is FALSE - use fitted track parameters,
109 // not measured values for nearest basetrack
110 gEVR->eQualityMode= 0; // 0 - prob/(sx**2+sy**2), 1 = 1./aver_impact
111
112 int nvtx = gEVR->FindVertex(); // find 2-prongs vertices
113 printf("%d 2-track vertexes was found\n",nvtx);
114
115 if(nvtx != 0) {
116
117 int nadd = gEVR->ProbVertexN(); // find N-prongs vertices by joining 2-prongs
118 printf("%d vertexes with number of tracks >2 was found\n",nadd);
119 printf("----------------------------------------------\n");
120
121 EdbVertex *v = 0;
122 Vertex *vt = 0;
123 int ntv = 0;
124 for(int i=0; i<nvtx; i++)
125 {
126 v = (EdbVertex*)(gEVR->eVTX->At(i));
127 if (!v) continue;
128 if (v->Flag() < 0) continue; // skip used 2-tracks vertices
129 ntv=v->N(); // number of prongs (tracks) at vertex
130 // v->VX() is vertex X coordinate
131 // v->VY() is vertex Y coordinate
132 // v->VZ() is vertex Z coordinate
133 vt = v->V();
134 // vt->ndf() is number of degree of freedom
135 // vt->chi2() is vertex chi2
136 // vt->prob() is vertex probability
137 // vt->rmsDistAngle() is average distance from vertex to tracks
138 // vt->angle() is average angle between tracks
139
140 for(int j=0; j<ntv; j++)
141 {
142 tra = v->GetTrack(j);
143 if (!tra) continue;
144 // you can put here analysis of vertex tracks
145
146 nseg = tra->N(); // number of basetracks (segments) in this track
147 for (int k=0; k<nseg; k++)
148 {
149 seg = tra->GetSegment(k);
150 if (!seg) continue;
151 // you can put here basetracks analysis for vertex tracks
152 }
153 }
154 }
155 }
156}
emulsion data processing
Definition: EdbDataSet.h:181
int InitVolume(int datatype=0, const char *rcut="1")
Definition: EdbDataSet.cxx:2066
EdbPVRec * PVR() const
Definition: EdbDataSet.h:198
TObjArray * eVTX
array of vertex
Definition: EdbPVRec.h:162
TObjArray * eTracks
Definition: EdbPVRec.h:161
void FillCell(float stepx, float stepy, float steptx, float stepty)
Definition: EdbPVRec.cxx:1092
int PropagateTracks(int nplmax, int nplmin, float probMin=0.05, int ngapMax=3, int design=0)
Definition: EdbPVRec.cxx:2487
void FitTracks(float p=10., float mass=0.139, TObjArray *gener=0, int design=0)
Definition: EdbPVRec.cxx:1893
Definition: EdbSegP.h:21
void SetErrorP(float sp2)
Definition: EdbSegP.h:94
Int_t Flag() const
Definition: EdbSegP.h:149
Definition: EdbPattern.h:113
Int_t N() const
Definition: EdbPattern.h:177
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:195
EdbSegP * GetSegmentLast() const
Definition: EdbPattern.h:190
EdbSegP * GetSegmentFirst() const
Definition: EdbPattern.h:189
Int_t eQualityMode
vertex quality estimation method (0:=Prob/(sigVX^2+sigVY^2); 1:= inverse average track-vertex distanc...
Definition: EdbVertex.h:183
Bool_t eUseSegPar
use only the nearest measured segments for vertex fit (as Neuchatel)
Definition: EdbVertex.h:182
Float_t eImpMax
maximal acceptable impact parameter (preliminary check)
Definition: EdbVertex.h:179
Float_t eDZmax
maximum z-gap in the track-vertex group
Definition: EdbVertex.h:177
Float_t eProbMin
minimum acceptable probability for chi2-distance between tracks
Definition: EdbVertex.h:178
Definition: EdbVertex.h:194
Int_t ProbVertexN()
Definition: EdbVertex.cxx:1426
TObjArray * eVTX
array of vertex
Definition: EdbVertex.h:205
Int_t FindVertex()
Definition: EdbVertex.cxx:1065
EdbPVRec * ePVR
patterns volume (optional)
Definition: EdbVertex.h:206
TObjArray * eEdbTracks
Definition: EdbVertex.h:204
Definition: EdbVertex.h:69
VERTEX::Vertex * V() const
Definition: EdbVertex.h:154
Int_t N() const
Definition: EdbVertex.h:121
EdbTrackP * GetTrack(int i)
Definition: EdbVertex.h:141
Int_t Flag() const
Definition: EdbVertex.h:124
Definition: VtVertex.hh:88
bool ReFit
Definition: tr_simple.C:12
float momentum
Definition: tr_simple.C:9
float ProbMinV
Definition: tr_simple.C:7
bool UseSegPar
Definition: tr_simple.C:11
EdbDataProc * dproc
Definition: tr_simple.C:3
float dpp
Definition: tr_simple.C:10
bool ReProp
Definition: tr_simple.C:14
float ProbMinP
Definition: tr_simple.C:8
EdbVertexRec * gEVR
Definition: tr_simple.C:5
EdbPVRec * gAli
Definition: tr_simple.C:4

Variable Documentation

◆ dpp

float dpp = .5

◆ dproc

EdbDataProc* dproc =0

◆ gAli

EdbPVRec* gAli =0

◆ gEVR

EdbVertexRec* gEVR =0

◆ momentum

float momentum = 1.

◆ ProbMinP

float ProbMinP = 0.01

◆ ProbMinV

float ProbMinV = 0.01

◆ ReFit

bool ReFit = false

◆ ReProp

bool ReProp = false

◆ UseSegPar

bool UseSegPar = false