3 #include "AliRunDigitizer.h"
4 #include "AliTPCDigitizer.h"
12 // test of the tpc merging using AliRunDigitizer and
13 // TPC Hits2Digits, Hits2SDigits and SDigits2Digits macros
16 // 0. make 2 directorys - ev1 and ev2
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
29 // 2. cp ev1/galice.root/galice.root.sdigit galice.root
30 // 3. load this macro and run testmerge()
32 // 4. run test function bellow to compare merged digits with original one
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
43 // merge two example events
45 //it merge two events -one from current directory -second from directory ev2
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);
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)
58 //if you think that there is memory leak -
59 //you are tru but othervise graphic doesn't work
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");
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);
81 TCanvas * c = new TCanvas(" Test merged digits", "test",600,900);
85 AliH2F * his = dig->DrawDigits("cont1",x1,x2,y1,y2);
86 his->SetTitle("MergedDigits");
87 his->SetName("MergedDigits");
91 AliH2F * his1 =dig1->DrawDigits("cont1",x1,x2,y1,y2);
92 his1->SetTitle("Event 1");
93 his1->SetName("Event1");
97 AliH2F * his2 =dig2->DrawDigits("cont1",x1,x2,y1,y2);
98 his2->SetTitle("Event 2");
99 his2->SetName("Event2");
104 void drawd(TFile * f, Int_t amp1, Int_t amp2)
106 TTree * tree = (TTree*)f->Get("TreeD_75x40_100x60_0;1");
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++){
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);
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");
135 AliSimDigits *dig1=0;
136 AliSimDigits *dig2=0;
138 tree->GetBranch("Segment")->SetAddress(&dig);
139 tree1->GetBranch("Segment")->SetAddress(&dig1);
140 tree2->GetBranch("Segment")->SetAddress(&dig2);
142 for (Int_t i=0;i<6000;i++){
147 dig1->ExpandBuffer();
148 dig2->ExpandBuffer();
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);
160 printf("%d\t%d\t%d\t%d\t%d\t%d\n",i,rows,col,d,d1,d2);
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");
176 AliSimDigits *dig1=0;
177 AliSimDigits *dig2=0;
179 tree->GetBranch("Segment")->SetAddress(&dig);
180 tree1->GetBranch("Segment")->SetAddress(&dig1);
181 tree2->GetBranch("Segment")->SetAddress(&dig2);
186 for (Int_t i=0;i<6000;i++){
191 dig1->ExpandBuffer();
192 dig2->ExpandBuffer();
194 Int_t nrows = dig->GetNRows();
195 Int_t ncols = dig->GetNCols();
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);
205 printf("%d\t%d\t%d\t%d\t%d\t%d\n",i,rows,col,d,d1,d2);
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");
220 AliSimDigits *dig1=0;
221 AliSimDigits *dig2=0;
223 tree->GetBranch("Segment")->SetAddress(&dig);
224 tree1->GetBranch("Segment")->SetAddress(&dig1);
225 tree2->GetBranch("Segment")->SetAddress(&dig2);
227 for (Int_t i=0;i<6000;i++){
231 if ( (dig1->GetID()!=i) || (dig2->GetID()!=i) ) {
232 printf("missed segments\n");
236 dig1->ExpandBuffer();
237 dig2->ExpandBuffer();
239 Int_t nrows = dig->GetNRows();
240 Int_t ncols = dig->GetNCols();
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;
249 printf("%d\t%d\t%d\t%d\t%d\t%d\n",i,rows,col,d,d1,d2);
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");
261 TTree * tree1 = (TTree*)f1.Get("TreeD_75x40_100x60_0;1");
262 TTree * tree2 = (TTree*)f2.Get("TreeD_75x40_100x60_0;1");
264 AliSimDigits *dig1=0;
265 AliSimDigits *dig2=0;
267 tree1->GetBranch("Segment")->SetAddress(&dig1);
268 tree2->GetBranch("Segment")->SetAddress(&dig2);
270 for (Int_t i=0;i<6000;i++){
274 //dig->ExpandBuffer();
275 if ( (dig1->GetID()!=i) || (dig2->GetID()!=i) ) {
276 printf("miised semgnets\n");
279 dig1->ExpandBuffer();
280 dig2->ExpandBuffer();
282 Int_t nrows = dig1->GetNRows();
283 Int_t ncols = dig1->GetNCols();
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);
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);
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");
305 AliSimDigits *dig1=0;
306 AliSimDigits *dig2=0;
308 tree1->GetBranch("Segment")->SetAddress(&dig1);
309 tree2->GetBranch("Segment")->SetAddress(&dig2);
311 for (Int_t i=0;i<6000;i++){
315 //dig->ExpandBuffer();
316 if ( (dig1->GetID()!=i) || (dig2->GetID()!=i) ) {
317 printf("miised semgnets\n");
320 dig1->ExpandBuffer();
321 dig2->ExpandBuffer();
323 Int_t nrows = dig1->GetNRows();
324 Int_t ncols = dig1->GetNCols();
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);