]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/ITStracks.C
New version with the Bari/Salerno model as default for SPD simulation
[u/mrichter/AliRoot.git] / ITS / ITStracks.C
CommitLineData
2e3c7808 1struct GoodTrack {
2 Int_t lab;
3 Int_t code;
4 Bool_t flag;
5 Float_t px,py,pz;
6 Float_t x,y,z;
7 Float_t pxg,pyg,pzg,ptg;
8};
9
10
11Int_t ITStracks(Int_t evNumber1=0,Int_t evNumber2=0,Int_t nclust=5) {
12
13 cerr<<"Select tracks which have nclust rec points in ITS...\n";
14
15 GoodTrack gt[15000];
16 Int_t ngood=0;
8ab89103 17 ifstream in("good_tracks_tpc");
2e3c7808 18
19 cerr<<"Reading good tracks...\n";
20 while (in>>gt[ngood].lab>>gt[ngood].code
21 >>gt[ngood].px >>gt[ngood].py>>gt[ngood].pz
22 >>gt[ngood].x >>gt[ngood].y >>gt[ngood].z
23 >>gt[ngood].pxg >>gt[ngood].pyg >>gt[ngood].pzg
24 >>gt[ngood].ptg >>gt[ngood].flag) {
25 ngood++;
26 cerr<<ngood<<'\r';
27 if (ngood==15000) {
28 cerr<<"Too many good tracks !\n";
29 break;
30 }
31 }
8ab89103 32 if (!in.eof()) cerr<<"Read error (good_tracks_tpc) !\n";
2e3c7808 33
34
35 if (gClassTable->GetID("AliRun") < 0) {
36 gROOT->LoadMacro("loadlibs.C");
37 loadlibs();
38 } else {
39 delete gAlice;
40 gAlice=0;
41 }
42
43// Connect the Root Galice file containing Geometry, Kine and Hits
44 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject("galice.root");
45 if (!file) file = new TFile("galice.root","UPDATE");
46 //if (!file) file = new TFile(filename);
47
48// Get AliRun object from file or create it if not on file
49 if (!gAlice) {
50 gAlice = (AliRun*)file->Get("gAlice");
51 if (gAlice) printf("AliRun object found on file\n");
52 if (!gAlice) gAlice = new AliRun("gAlice","Alice test program");
53 }
54
55
56 AliITS* ITS =(AliITS *)gAlice->GetDetector("ITS");
57 if (!ITS) return;
58
59 //TClonesArray *recPoints = ITS->RecPoints();
60 //cout<<" recPoints = "<<recPoints<<"\n";
61 Int_t numbpoints=0;
62 AliITSRecPoint *recp=0;
63 Int_t track=-3;
64
65
66 for (int nev=evNumber1; nev<= evNumber2; nev++) {
67 Int_t nparticles = gAlice->GetEvent(nev);
68 cout << "nev " << nev <<endl;
69 cout << "nparticles " << nparticles <<endl;
70 if (nev < evNumber1) continue;
71 if (nparticles <= 0) return;
72 TClonesArray *recPoints = ITS->RecPoints();
73 TObjArray *particles=gAlice->Particles();
74 // TClonesArray *particles=gAlice->Particles();
75 // Int_t np=particles->GetEntriesFast();
76 Int_t np=gAlice->GetNtrack(); //FCA correction
77 TMatrix goodITS(np,6);
78 Int_t modmin[]={0,80,240,324,500,1282};
79 Int_t modmax[]={79,239,323,499,1281,2269};
80 for(Int_t j=0; j<np; j++){
81 for(Int_t j1=0; j1<6; j1++) goodITS(j,j1)=0;
82 }
83
84
85 Int_t nent=gAlice->TreeR()->GetEntries();
86 printf("Found %d entries in the TreeR (must be one per module per event!)\n",nent); //
87
88 for (Int_t layer=1;layer<=6;layer++) {
89 for (Int_t mod=modmin[layer-1];mod<=modmax[layer-1];mod++) {
90 ITS->ResetRecPoints();
91 //gAlice->TreeR()->GetEvent(mod+1); //first entry in TreeR is empty
92 gAlice->TreeR()->GetEvent(mod); //modificato il 2-3-2001
93 numbpoints = recPoints->GetEntries();
94 if(!numbpoints) continue;
95 for (Int_t irec=0;irec<numbpoints;irec++) {
96 recp=(AliITSRecPoint*)recPoints->UncheckedAt(irec);
97 for (Int_t it=0;it<3;it++) {
98 track=recp->GetLabel(it);
99 if(track<0) continue;
100 if(track > np) {
101 cout<<" Error on track number \n";
102 getchar();
103 }
104 goodITS(track,layer-1)=1;
105 } //loop over tracks
106 } //loop over points
107 } //loop over modules
108 } //loop over layers
109
110 for(Int_t i=0; i<ngood; i++){
111 if(gt[i].lab>np) {cout<<" Error on TPC good tracks label \n"; getchar();}
112 Int_t ilab=gt[i].lab;
113 Int_t totclust=0;
114 for(Int_t j=0; j<6; j++) totclust+=goodITS(ilab,j);
115 if(totclust>=nclust) {
116 gt[i].flag=1;
117 printf("label nclusters %d %d\n",gt[i].lab,totclust); //
118 }
119 }
120
121 Int_t itsngood=0;
122 ofstream out("itsgood_tracks");
123 if (out) {
124 for (Int_t ngd=0; ngd<ngood; ngd++) {
125 if(gt[ngd].flag) {
126 out<<gt[ngd].lab<<' '<<gt[ngd].code<<' '
127 <<gt[ngd].px <<' '<<gt[ngd].py<<' '<<gt[ngd].pz<<' '
128 <<gt[ngd].x <<' '<<gt[ngd].y <<' '<<gt[ngd].z <<' '
129 <<gt[ngd].pxg <<' '<<gt[ngd].pyg <<' '<<gt[ngd].pzg <<' '
130 <<gt[ngd].ptg<<gt[ngd].flag<<endl;
131 itsngood++;
132 }
133 }
134 } else cerr<<"Can not open file (itsgood_tracks) !\n";
135 out.close();
136
137 cerr<<"Number of good tracks in TPC: "<<ngood<<endl;
138 cerr<<"Number of good tracks in ITS: "<<itsngood<<endl;
139
140 } // event loop
141
142 file->Close();
143
144 printf("after file close\n");
145
146 return 0;
147
148}
149