]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCTestMerge.C
Updated
[u/mrichter/AliRoot.git] / TPC / AliTPCTestMerge.C
1 #ifndef __CINT__
2 #include "alles.h"
3 #include "AliRunDigitizer.h"
4 #include "AliTPCDigitizer.h"
5 #include "AliH2F.h"
6 #include "TTree.h"
7
8
9 #endif
10
11 // Marian Ivanov
12 // test of the tpc merging using AliRunDigitizer and 
13 // TPC  Hits2Digits, Hits2SDigits and SDigits2Digits macros
14
15 // preparation  
16 // 0. make 2 directorys - ev1 and ev2
17
18 // 1.make hits, digits,sdigits and sdigits-digits in both directories 
19 //  1.a  galice -b -q grun.C and produce hits
20 //  1.b. cp galice.root galice.root.hits
21 //  1.c  run AliTPCHits2Digits.C
22 //  1.d  cp galice.root galice.root.digits
23 //  1.e  copy back cp galice.root.hits galice.root
24 //  1.f  run AliTPCSDigits2Digits.C
25 //  1.g  cp galice.root  galice.root.sdigits
26 //  1.h  run AliTPCSDigits2Digit.C 
27 //  1.i  cp galice.root galice.root.dig2
28
29 // 2. cp ev1/galice.root/galice.root.sdigit galice.root
30 // 3. load this macro and run testmerge()
31
32 // 4. run test function bellow to compare merged digits with original one
33
34 // 5. to be  noticed output ftom testx() function - should be bigger than
35 //    noise because in proces of digitisation we use different seed
36 //    of random numbers for diffusion and gas gain
37 //    -anyway in place where we have the signal should occur the signal in both casses
38 //    - only the amplitude should be different - by factor of sqrt
39
40       
41 void testmerge()
42 {
43   // merge two example events
44   //
45   //it merge two events -one from current directory -second from directory ev2
46   
47   AliRunDigitizer * manager = new AliRunDigitizer(2,1);
48   manager->SetInputStream(0,"galice.root");
49   manager->SetInputStream(1,"ev2/galice.root.sdigits");
50   AliTPCDigitizer dTPC(manager);
51   manager->SetNrOfEventsToWrite(1);
52   manager->Exec("");
53
54 }
55
56 void drawmerged(Int_t sec, Int_t row, Int_t x1=-1, Int_t x2=-1, Int_t y1=-1, Int_t y2=-1)
57 {
58   //if you think that there is memory leak -
59   //you are tru but othervise graphic doesn't work
60   // sec=0; row =0;
61   TFile * f = new TFile("galice.root");
62   TFile * f1= new TFile("ev1/galice.root.dig2");
63   TFile * f2= new TFile("ev2/galice.root.dig2");
64   TTree * tree = (TTree*)f->Get("TreeD_75x40_100x60_0;1");
65   TTree * tree1 = (TTree*)f1->Get("TreeD_75x40_100x60_0;1");
66   TTree * tree2 = (TTree*)f2->Get("TreeD_75x40_100x60_0;1");
67   //
68   AliSimDigits *dig=0;
69   AliSimDigits *dig1=0;
70   AliSimDigits *dig2=0;
71   //
72   tree->GetBranch("Segment")->SetAddress(&dig);
73   tree1->GetBranch("Segment")->SetAddress(&dig1);
74   tree2->GetBranch("Segment")->SetAddress(&dig2);
75   AliTPCParam * param =(AliTPCParam*) f->Get("75x40_100x60");
76   Int_t index = param->GetIndex(sec,row);
77   tree->GetEvent(index);
78   tree1->GetEvent(index);
79   tree2->GetEvent(index);
80
81   TCanvas * c = new TCanvas(" Test merged digits", "test",600,900);
82   c->Divide(1,3);
83   //
84   c->cd(1);
85   AliH2F * his = dig->DrawDigits("cont1",x1,x2,y1,y2);
86   his->SetTitle("MergedDigits");
87   his->SetName("MergedDigits");
88   his->DrawClone();
89   //
90   c->cd(2);
91   AliH2F * his1 =dig1->DrawDigits("cont1",x1,x2,y1,y2); 
92   his1->SetTitle("Event 1");
93   his1->SetName("Event1");
94   his1->DrawClone();
95   //
96   c->cd(3);
97   AliH2F * his2 =dig2->DrawDigits("cont1",x1,x2,y1,y2); 
98   his2->SetTitle("Event 2");
99   his2->SetName("Event2");
100   his2->DrawClone();  
101 }
102
103
104 void drawd(TFile * f, Int_t amp1, Int_t amp2)
105 {
106   TTree * tree = (TTree*)f->Get("TreeD_75x40_100x60_0;1");
107   AliSimDigits *dig=0;
108   tree->GetBranch("Segment")->SetAddress(&dig); 
109   TH1F * his = new TH1F("his","his",amp2-amp1,amp1,amp2);
110   for (Int_t i=0;i<60;i++){  
111     tree->GetEvent(i);   
112     dig->ExpandBuffer();
113     Int_t nrows = dig->GetNRows();
114     Int_t ncols = dig->GetNCols();
115     for (Int_t rows=0;rows<nrows; rows++)
116       for (Int_t col=0;col<ncols; col++){    
117         Int_t d  = dig->GetDigitFast(rows,col);
118         his->Fill(d);
119       }
120   }
121   his->Draw();
122 }
123
124 void test1(){
125   //test of the merged digits
126   //compare merged digits with standard digits
127   TFile f("galice.root");
128   TFile f1("ev1/galice.root.digits");
129   TFile f2("ev2/galice.root.digits");
130   TTree * tree = (TTree*)f.Get("TreeD_75x40_100x60_0;1");
131   TTree * tree1 = (TTree*)f1.Get("TreeD_75x40_100x60_0;1");
132   TTree * tree2 = (TTree*)f2.Get("TreeD_75x40_100x60_0;1");
133   //
134   AliSimDigits *dig=0;
135   AliSimDigits *dig1=0;
136   AliSimDigits *dig2=0;
137   //
138   tree->GetBranch("Segment")->SetAddress(&dig);
139   tree1->GetBranch("Segment")->SetAddress(&dig1);
140   tree2->GetBranch("Segment")->SetAddress(&dig2);
141   //
142   for (Int_t i=0;i<6000;i++){
143     tree->GetEvent(i);
144     tree1->GetEvent(i);
145     tree2->GetEvent(i);
146     dig->ExpandBuffer();
147     dig1->ExpandBuffer();
148     dig2->ExpandBuffer();
149     //
150     Int_t nrows = dig->GetNRows();
151     Int_t ncols = dig->GetNCols();
152     for (Int_t rows=0;rows<nrows; rows++)
153       for (Int_t col=0;col<ncols; col++){
154         Int_t d  = dig->GetDigitFast(rows,col);
155         Int_t d1 = dig1->GetDigitFast(rows,col);
156         Int_t d2 = dig2->GetDigitFast(rows,col);
157         
158         if (abs(d-d1-d2)>4)
159         //if (d2>5)
160           printf("%d\t%d\t%d\t%d\t%d\t%d\n",i,rows,col,d,d1,d2);
161       }
162   }
163 }
164
165 void test5(){
166   //
167   //compare merged digits with digits obtained hits2sdig->sdigtodig
168   TFile f("galice.root");
169   TFile f1("ev1/galice.root.dig2");
170   TFile f2("ev2/galice.root.dig2");
171   TTree * tree = (TTree*)f.Get("TreeD_75x40_100x60_0;1");
172   TTree * tree1 = (TTree*)f1.Get("TreeD_75x40_100x60_0;1");
173   TTree * tree2 = (TTree*)f2.Get("TreeD_75x40_100x60_0;1");
174
175   AliSimDigits *dig=0;
176   AliSimDigits *dig1=0;
177   AliSimDigits *dig2=0;
178   
179   tree->GetBranch("Segment")->SetAddress(&dig);
180   tree1->GetBranch("Segment")->SetAddress(&dig1);
181   tree2->GetBranch("Segment")->SetAddress(&dig2);
182
183
184
185
186   for (Int_t i=0;i<6000;i++){
187     tree->GetEvent(i);
188     tree1->GetEvent(i);
189     tree2->GetEvent(i);
190     dig->ExpandBuffer();
191     dig1->ExpandBuffer();
192     dig2->ExpandBuffer();
193     
194     Int_t nrows = dig->GetNRows();
195     Int_t ncols = dig->GetNCols();
196
197     for (Int_t rows=0;rows<nrows; rows++)
198       for (Int_t col=0;col<ncols; col++){
199         Int_t d  = dig->GetDigitFast(rows,col);
200         Int_t d1 = dig1->GetDigitFast(rows,col);
201         Int_t d2 = dig2->GetDigitFast(rows,col);
202         
203         if (abs(d-d1-d2)>4)
204         //if (d2>5)
205           printf("%d\t%d\t%d\t%d\t%d\t%d\n",i,rows,col,d,d1,d2);
206       }
207   }
208 }
209
210 void test3(){
211   //test of the merged digits
212   TFile f("galice.root");
213   TFile f1("ev1/galice.root.sdigits");
214   TFile f2("ev2/galice.root.sdigits");
215   TTree * tree = (TTree*)f.Get("TreeD_75x40_100x60_0;1");
216   TTree * tree1 = (TTree*)f1.Get("TreeS_75x40_100x60_0;1");
217   TTree * tree2 = (TTree*)f2.Get("TreeS_75x40_100x60_0;1");
218   //
219   AliSimDigits *dig=0;
220   AliSimDigits *dig1=0;
221   AliSimDigits *dig2=0;
222   //  
223   tree->GetBranch("Segment")->SetAddress(&dig);
224   tree1->GetBranch("Segment")->SetAddress(&dig1);
225   tree2->GetBranch("Segment")->SetAddress(&dig2);
226   //
227   for (Int_t i=0;i<6000;i++){
228     tree->GetEvent(i);
229     tree1->GetEvent(i);
230     tree2->GetEvent(i);
231     if ( (dig1->GetID()!=i) || (dig2->GetID()!=i) ) {
232       printf("missed segments\n");
233     }
234     //
235     dig->ExpandBuffer();
236     dig1->ExpandBuffer();
237     dig2->ExpandBuffer();
238     //
239     Int_t nrows = dig->GetNRows();
240     Int_t ncols = dig->GetNCols();
241     //
242     for (Int_t rows=0;rows<nrows; rows++)
243       for (Int_t col=0;col<ncols; col++){
244         Int_t d  = dig->GetDigitFast(rows,col);
245         Int_t d1 = dig1->GetDigitFast(rows,col)/16;
246         Int_t d2 = dig2->GetDigitFast(rows,col)/16;     
247         if (abs(d-d1-d2)>4)
248         //if (d2>5)
249           printf("%d\t%d\t%d\t%d\t%d\t%d\n",i,rows,col,d,d1,d2);
250       }
251   }
252 }
253
254
255 void test2(){
256   //test of the digits produced by the Hits2Digits 
257   //and Hits2SDigits - SDigits2Digits chain
258   TFile f1("galice.root.digits");
259   TFile f2("galice.root.dig2");
260   //
261   TTree * tree1 = (TTree*)f1.Get("TreeD_75x40_100x60_0;1");
262   TTree * tree2 = (TTree*)f2.Get("TreeD_75x40_100x60_0;1");
263   //
264   AliSimDigits *dig1=0;
265   AliSimDigits *dig2=0;
266   //  
267   tree1->GetBranch("Segment")->SetAddress(&dig1);
268   tree2->GetBranch("Segment")->SetAddress(&dig2);
269   //
270   for (Int_t i=0;i<6000;i++){
271     //tree->GetEvent(i);
272     tree1->GetEvent(i);
273     tree2->GetEvent(i);
274     //dig->ExpandBuffer();
275     if ( (dig1->GetID()!=i) || (dig2->GetID()!=i) ) {
276       printf("miised semgnets\n");
277     }
278     //        
279     dig1->ExpandBuffer();
280     dig2->ExpandBuffer();
281     //
282     Int_t nrows = dig1->GetNRows();
283     Int_t ncols = dig1->GetNCols();
284     //
285     for (Int_t rows=0;rows<nrows; rows++)
286       for (Int_t col=0;col<ncols; col++){
287         //      Int_t d  = dig->GetDigitFast(rows,col);
288         Int_t d1 = dig1->GetDigitFast(rows,col);
289         Int_t d2 = dig2->GetDigitFast(rows,col);
290         //      
291         if ((d2>5) &&abs(d1-d2)>2)        printf("%d\t%d\t%d\t\t%d\t%d\n",i,rows,col,d1,d2);
292       }
293   }
294 }
295
296
297
298 void test4(){
299   //TPC internal test
300   TFile f1("galice.root.sdigits");
301   TFile f2("galice.root.digits");
302   TTree * tree1 = (TTree*)f1.Get("TreeS_75x40_100x60_0;1");
303   TTree * tree2 = (TTree*)f2.Get("TreeD_75x40_100x60_0;1");
304   //
305   AliSimDigits *dig1=0;
306   AliSimDigits *dig2=0;
307   //
308   tree1->GetBranch("Segment")->SetAddress(&dig1);
309   tree2->GetBranch("Segment")->SetAddress(&dig2);
310   //
311   for (Int_t i=0;i<6000;i++){
312     //tree->GetEvent(i);
313     tree1->GetEvent(i);
314     tree2->GetEvent(i);
315     //dig->ExpandBuffer();
316     if ( (dig1->GetID()!=i) || (dig2->GetID()!=i) ) {
317       printf("miised semgnets\n");
318     }
319     //        
320     dig1->ExpandBuffer();
321     dig2->ExpandBuffer();
322     //
323     Int_t nrows = dig1->GetNRows();
324     Int_t ncols = dig1->GetNCols();
325     //
326     for (Int_t rows=0;rows<nrows; rows++)
327       for (Int_t col=0;col<ncols; col++){
328         Int_t d1 = dig1->GetDigitFast(rows,col)/16.;
329         Int_t d2 = dig2->GetDigitFast(rows,col);                
330         if ((d2>5) &&abs(d1-d2)>2)        printf("%d\t%d\t%d\t\t%d\t%d\n",i,rows,col,d1,d2);
331       }
332   }
333 }
334
335
336