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