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

Functions

void VTExample (char *data_set="july2003_data_set.def")
 

Variables

float AngleAcceptance = 1.0
 
EdbDataProcdset =0
 
EdbPVRecgAli =0
 
int maxgaps [6] = {0,3,3,6,3,6}
 
int nsegMin = 3
 
float ProbMinP = 0.01
 
float ProbMinT = 0.01
 
float ProbMinV = 0.01
 
bool usemom = false
 

Function Documentation

◆ VTExample()

void VTExample ( char *  data_set = "july2003_data_set.def")
17{
18 dset=new EdbDataProc(data_set);
19 dset->InitVolume(0);
20 gAli = dset->PVR();
21
22 gAli->Link();
23
26 if (!gAli->eTracks) return;
27
28 int ntr = gAli->eTracks->GetEntries();
29 EdbTrackP *tr = 0;
30 float p = 0.;
31 // loop on tracks for momentum error setting
32 for(int i=0; i<ntr; i++)
33 {
34 tr = (EdbTrackP*)(gAli->eTracks->At(i));
35 p = 4.; // use 20% error for 4 GeV/c
36 tr->SetErrorP(0.2*0.2*p*p);
37 }
38 // fit all tracks assume momentum 4 GeV/c and pion mass
39 gAli->FitTracks(4.);
40
41 gAli->FillCell(50,50,0.015,0.015);
42 // merge consistient tracks
44 // find 2-tracks vertexes
45 int nvtx = gAli->ProbVertex(maxgaps, AngleAcceptance, ProbMinV, ProbMinT, nsegMin, usemom);
46 printf("%d 2-track vertexes was found\n",nvtx);
47
48 if(nvtx != 0) {
49 // merge 2-tracks vertexes into N-tracks (if common track exist)
50 int nadd = gAli->ProbVertexN(ProbMinV);
51 printf("%d vertexes with number of tracks >2 was found\n",nadd);
52 }
53
54 if (!(gAli->eVTX)) return;
55 EdbVertex *v = 0;
56 EdbVertex *vc = 0;
57 Vertex *vt = 0;
58 int nv = gAli->eVTX->GetEntries();
59 int nvg = 0;
60 float xv,yv,zv,xve,yve,zve,probv, rmsdist,rmsang,chindfv;
61 int ndfv, ntrv;
62 for(int i=0; i<nv; i++)
63 {
64 v = (EdbVertex*)(gAli->eVTX->At(i));
65 if (v->Flag() != -10) // skip "broken" vertexes
66 {
67 nvg++;
68 xv = v->VX(); // X vertex coordinate
69 yv = v->VY(); // Y vertex coordinate
70 zv = v->VZ(); // Z vertex coordinate
71 vt = v->V(); // Vt++ vertex pointer
72 xve = vt->vxerr(); // error on x-coordinate
73 yve = vt->vyerr(); // error on y-coordinate
74 zve = vt->vzerr(); // error on z-coordinate
75 chindfv = vt->ndf() > 0 ? vt->chi2()/vt->ndf() : 0;
76 probv = vt->prob(); // vertex chi square probability
77 rmsdist = vt->rmsDistAngle(); // RMS of track-vertex distances
78 rmsang = vt->angle(); // RMS of angles between tracks
79 ntrv = vt->ntracks(); // number of tracks in vertex
80 v->Print();
81 }
82 }
83 printf("Total %d vertexes was found\n",nvg);
84
85 if (nvg)
86 {
87 int nlv = gAli->LinkedVertexes();
88 printf("%d linked vertexes found\n", nlv);
89 if (nlv)
90 {
91 for(int i=0; i<nv; i++)
92 {
93 v = (EdbVertex*)(gAli->eVTX->At(i));
94 if (v->Flag() > 2)
95 {
96 vc = v->GetConnectedVertex(0);
97 if (vc->ID() > v->ID())
98 {
99 v->Print();
100 vc->Print();
101 }
102 }
103 }
104 }
105 }
106
107 int nn = gAli->VertexNeighboor(1000.,2);
108 printf("%d neighbooring tracks found\n", nn);
109}
TTree * tr
Definition: Shower_E_FromShowerRoot.C:5
int maxgaps[6]
Definition: VTExample.C:7
bool usemom
Definition: VTExample.C:13
float ProbMinV
Definition: VTExample.C:9
float AngleAcceptance
Definition: VTExample.C:8
EdbDataProc * dset
Definition: VTExample.C:4
float ProbMinT
Definition: VTExample.C:11
int nsegMin
Definition: VTExample.C:12
float ProbMinP
Definition: VTExample.C:10
EdbPVRec * gAli
Definition: VTExample.C:5
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
int MakeTracks(int nsegments=2, int flag=0)
Definition: EdbPVRec.cxx:1989
TObjArray * eVTX
array of vertex
Definition: EdbPVRec.h:162
void FillTracksCell()
Definition: EdbPVRec.cxx:1319
TObjArray * eTracks
Definition: EdbPVRec.h:161
void FillCell(float stepx, float stepy, float steptx, float stepty)
Definition: EdbPVRec.cxx:1092
int Link()
Definition: EdbPVRec.cxx:1187
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: EdbPattern.h:113
Definition: EdbVertex.h:69
void Print()
Definition: EdbVertex.cxx:328
Int_t ID() const
Definition: EdbVertex.h:126
Float_t VX() const
Definition: EdbVertex.h:133
VERTEX::Vertex * V() const
Definition: EdbVertex.h:154
Float_t VY() const
Definition: EdbVertex.h:134
EdbVertex * GetConnectedVertex(int nv)
Definition: EdbVertex.cxx:405
Float_t VZ() const
Definition: EdbVertex.h:135
Int_t Flag() const
Definition: EdbVertex.h:124
Definition: VtVertex.hh:88
double vyerr() const
$\sqrt{\sigma_{vy}^2}$ vertex $y$-error
unsigned short int ndf() const
degrees of freedom of vertex fit
Definition: VtVertex.C:238
float chi2() const
$\chi^2$ of vertex fit
Definition: VtVertex.C:236
double vzerr() const
$\sqrt{\sigma_{vz}^2}$ vertex $z$-error
double angle() const
unsigned short int ntracks() const
no of tracks in vertex $n$
Definition: VtVertex.C:240
double rmsDistAngle() const
calc rms dist and rms angle
Definition: VtVertex.C:1073
float prob() const
upper tail $\chi^2$ probability
Definition: VtVertex.C:237
double vxerr() const
$\sqrt{\sigma_{vx}^2}$ vertex $x$-error
p
Definition: testBGReduction_AllMethods.C:8

Variable Documentation

◆ AngleAcceptance

float AngleAcceptance = 1.0

◆ dset

EdbDataProc* dset =0

◆ gAli

EdbPVRec* gAli =0

◆ maxgaps

int maxgaps[6] = {0,3,3,6,3,6}

◆ nsegMin

int nsegMin = 3

◆ ProbMinP

float ProbMinP = 0.01

◆ ProbMinT

float ProbMinT = 0.01

◆ ProbMinV

float ProbMinV = 0.01

◆ usemom

bool usemom = false