]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/testRawReaderFastDDL.C
Removing obsolete directory - the calibration classes directly in the TPC directory...
[u/mrichter/AliRoot.git] / TPC / testRawReaderFastDDL.C
1 #include <stdio.h>
2 #include <TString.h>
3 #include <TROOT.h>
4 #include <TStopwatch.h>
5 #include <TH2I.h>
6 #include <TFile.h>
7
8 #include "AliRawReader.h"
9 #include "AliRawReaderRoot.h"
10 #include "AliLog.h"
11
12 #include "AliAltroRawStreamFast.h"
13 #include "AliAltroRawStream.h"
14
15 /*
16  compare old and new AltroRawStream algorithms for each pad and timebin
17
18  check if bins are filled twice (which should not be the case!!!)
19 */
20
21 void testRawReaderFastDDL(const Char_t *file="/data.local/data/06000002142000.1A.root")
22 {
23     // set minimal screen output
24     AliLog::SetGlobalDebugLevel(0) ;
25     AliLog::SetGlobalLogLevel(AliLog::kFatal);
26
27 //    TString filename("/d/alice05/testtpc/raw/pulser/06000002142000.1A.root");
28     TString filename(file);
29
30     printf("File: %s\n", filename.Data());
31
32     AliRawReader *rawReader = new AliRawReaderRoot(filename);
33     if ( !rawReader ) return;
34     rawReader->RewindEvents();
35
36     AliAltroRawStreamFast *sf = new AliAltroRawStreamFast(rawReader);
37     AliAltroRawStream      *s = new AliAltroRawStream(rawReader);
38
39     s->SetOldRCUFormat(kTRUE);
40     s->SetNoAltroMapping(kFALSE);
41     s->SelectRawData("TPC");
42
43     Int_t ievent = 0;
44     Int_t count=0;
45
46     TH2I *h2ddlT1[216];
47     TH2I *h2ddlT2[216];
48
49     for ( Int_t i=0; i<216; i++ ){
50         h2ddlT1[i] = 0x0;
51         h2ddlT2[i] = 0x0;
52     }
53
54     TStopwatch timer1;
55     TStopwatch timer2;
56     TStopwatch timer3;
57
58     while (rawReader->NextEvent()){
59         printf("\nevent: %d\n",ievent);
60         Bool_t input=kFALSE;
61
62
63         //old algorithm
64         timer1.Start();timer2.Start(kFALSE);
65         while ( s->Next() ){
66             if ( !h2ddlT1[s->GetDDLNumber()] ) h2ddlT1[s->GetDDLNumber()] = new TH2I(Form("hddl1_%d",s->GetDDLNumber()),"h2c1",3584,0,3584,1024,0,1024);
67             TH2I *hist = h2ddlT1[s->GetDDLNumber()];
68
69             //fast filling, TH1::Fill takes awfully long
70             Int_t bin=(s->GetTime()+1)*(3584+2)+s->GetHWAddress()+1;
71             // check if this bin was allready filled
72             if ( hist->GetArray()[bin] > 0 )
73                 printf(" not 0:   |  %.3d : %.3d (%.3d)\n",
74                         s->GetHWAddress(), s->GetSignal(), hist->GetArray()[bin]);
75             else
76                 hist->GetArray()[bin]=s->GetSignal();
77
78
79             input=kTRUE;
80             count++;
81         }
82         timer1.Stop();timer2.Stop();
83         printf("old  --  Time: %.4f (%.4f)\n", timer1.RealTime(), timer1.CpuTime());
84         // END OLD
85
86         rawReader->Reset();
87
88         //new algorithm
89         timer1.Start();timer3.Start(kFALSE);
90         while ( sf->NextDDL() ){
91             if ( !h2ddlT2[sf->GetDDLNumber()] ) h2ddlT2[sf->GetDDLNumber()] = new TH2I(Form("hddl2_%d",s->GetDDLNumber()),"h2c1",3584,0,3584,1024,0,1024);
92             TH2I *hist = h2ddlT2[sf->GetDDLNumber()];
93             while ( sf->NextChannel() ){
94                 UInt_t signal=0;
95                 while ( sf->NextBunch() ){
96                     for (UInt_t timebin=sf->GetStartTimeBin(); timebin<sf->GetEndTimeBin(); timebin++){
97                         signal=sf->GetSignals()[timebin-sf->GetStartTimeBin()];
98
99                         //fast filling, TH1::Fill takes awfully long
100                         Int_t bin=(timebin+1+1)*(3584+2)+sf->GetHWAddress()+1; // timebins of old and new algorithm differ by 1!!!
101                         // check if this bin was allready filled
102                         if ( hist->GetArray()[bin] > 0 )
103                             printf(" not 0:   |  %.3d : %.3d (%.3d)\n",
104                                     sf->GetHWAddress(), signal, hist->GetArray()[bin]);
105                         else
106                             hist->GetArray()[bin]=signal;
107                     }
108                 }
109             }
110         }
111         timer1.Stop();timer3.Stop();
112         printf("new  --  Time: %.4f (%.4f)\n", timer1.RealTime(), timer1.CpuTime());
113         // END NEW
114
115         //check if all data are the same for both algorithms
116         for ( Int_t ddl=0; ddl<216; ddl++ ){
117             TH2I *hist = h2ddlT1[ddl];
118             if ( !hist ) continue;
119             TH2I *hist2 = h2ddlT2[ddl];
120             for ( Int_t hadd=0; hadd<3584; hadd++ )
121                 for ( Int_t time=0; time<1024; time++ ){
122                     Int_t bin=(time+1)*(3584+2)+hadd+1;
123                     Int_t val1 = hist->GetArray()[bin];
124                     Int_t val2 = hist2->GetArray()[bin];
125                     if ( val1 != val2 )
126                         printf("%.2d. %.3d %.4d %.4d: %d - %d = %d\n", ievent, ddl, hadd, time, val1, val2, val1-val2);
127
128                     //reset for the next event
129                     hist->GetArray()[bin]=0;
130                     hist2->GetArray()[bin]=0;
131                 }
132         }
133
134         if (input) ievent++;
135     }
136
137     printf("total old  --  Time: %.4f (%.4f)\n", timer2.RealTime(), timer2.CpuTime());
138     printf("total new  --  Time: %.4f (%.4f)\n", timer3.RealTime(), timer3.CpuTime());
139
140     delete rawReader;
141     delete [] h2ddlT1;
142     delete [] h2ddlT2;
143 }
144