]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/testRawReaderFastTPC.C
time bin limits introduced
[u/mrichter/AliRoot.git] / TPC / testRawReaderFastTPC.C
1 #include <stdio.h>
2 #include <TString.h>
3 #include <TROOT.h>
4 #include <TStopwatch.h>
5 #include <TCanvas.h>
6 #include <TFile.h>
7 #include <TH2F.h>
8
9 #include "AliLog.h"
10
11 #include "AliRawReader.h"
12 #include "AliRawReaderRoot.h"
13
14 #include "AliTPCRawStream.h"
15 #include "AliTPCRawStreamFast.h"
16
17 #include "AliAltroRawStream.h"
18 #include "AliAltroRawStreamFast.h"
19
20 #include "AliTPCCalPad.h"
21 #include "AliTPCCalROC.h"
22
23 /*
24     .L testRawReaderFastTPC.C+
25     testRawReaderFastTPC()
26 */
27
28
29 AliTPCCalPad *testRawReaderFastTPC(const Char_t *file="/data.local/data/06000002142000.1A.root")
30 {
31
32     AliLog::SetGlobalDebugLevel(0) ;
33     AliLog::SetGlobalLogLevel(AliLog::kFatal);
34 //    TString filename("/d/alice05/testtpc/raw/pulser/06000002142000.1A.root");  //nfs
35 //    TString filename("root://lxfs35.gsi.de:1094//alice/testtpc/raw2006/06000001537001.001.root");
36     TString filename(file);  //local
37     // on castor: /castor/cern.ch/alice/data/2006/09/18/15/06000002142000.1A.root
38
39
40     printf("File: %s\n", filename.Data());
41
42     AliRawReader *rawReader = new AliRawReaderRoot(filename);
43     if ( !rawReader ) return 0x0;
44     rawReader->RewindEvents();
45
46     AliTPCRawStreamFast *sf = new AliTPCRawStreamFast(rawReader);
47     AliTPCRawStream      *s = new AliTPCRawStream(rawReader);
48
49     s->SetOldRCUFormat(kTRUE);
50     s->SetNoAltroMapping(kFALSE);
51
52     Int_t ievent = 0;
53     Int_t ievent2 = 0;
54     Int_t count=0;
55
56     AliTPCCalPad *padOld = new AliTPCCalPad("old","old");
57     AliTPCCalPad *padNew = new AliTPCCalPad("new","new");
58     AliTPCCalPad *padSum = new AliTPCCalPad("sum","sum");
59     AliTPCCalPad *padhadd = new AliTPCCalPad("sum","sum");
60
61     for (Int_t sec=0; sec<72; sec++){
62         AliTPCCalROC *roc = padhadd->GetCalROC(sec);
63         for (UInt_t ch=0; ch<roc->GetNchannels(); ch++)
64             roc->SetValue(ch,-1);
65     }
66
67     TStopwatch timer1;
68     TStopwatch timer2;
69     TStopwatch timer3;
70
71     TCanvas *c1 = (TCanvas*)gROOT->FindObject("c1");
72     if ( !c1 ) c1 = new TCanvas("c1","c1");
73     c1->Clear();
74     c1->Divide(2,2);
75
76
77     while (rawReader->NextEvent()){
78         //old algorithm
79         printf("\nevent: %d (%d)\n",ievent, ievent2);
80
81         Bool_t input=kFALSE;
82
83 //      if ( ievent != 19 ) input=kTRUE;
84 //      else{
85             timer1.Start();timer2.Start(kFALSE);
86
87             Int_t sum1 = 0;
88             Int_t sum2 = 0;
89         while ( s->Next() ){
90             AliTPCCalROC *roc = padhadd->GetCalROC(s->GetSector());
91
92             //check if hwaddress gets overwritten
93             if ( roc->GetValue(s->GetRow(), s->GetPad()) == -1 )
94                 roc->SetValue(s->GetRow(), s->GetPad(), s->GetHWAddress());
95             else
96                 if ( roc->GetValue(s->GetRow(), s->GetPad()) != s->GetHWAddress()){
97                     printf("#%.2d: %.2d.%.3d.%.2d,%.4d - old [%.4d] now[%.4d]",
98                            ievent, s->GetSector(), s->GetRow(),s->GetPad(), s->GetTime(),
99                            (Int_t)roc->GetValue(s->GetRow(), s->GetPad()), s->GetHWAddress());
100                 }
101
102
103             Float_t val=padOld->GetCalROC(s->GetSector())->GetValue(s->GetRow(), s->GetPad());
104             padOld->GetCalROC(s->GetSector())->SetValue(s->GetRow(), s->GetPad(), s->GetSignal()+val);
105
106 /*          if ( ievent == 19 && s->GetSector() == 25 && s->GetRow() == 00 && s->GetPad()==9 ){
107                 printf("old         | (%.2d.%.3d.%.2d.%.4d | %.3d ): %.3d\n",
108                        s->GetSector(), s->GetRow(), s->GetPad(), s->GetTime(), s->GetHWAddress(), s->GetSignal());
109                 sum1+=s->GetSignal();
110             }
111 */
112
113             input=kTRUE;
114             count++;
115         }
116         timer1.Stop();timer2.Stop();
117         printf("old  --  Time: %.4f (%.4f)\n", timer1.RealTime(), timer1.CpuTime());
118
119         rawReader->Reset();
120
121         //new algorithm
122         timer1.Start();timer3.Start(kFALSE);
123         while ( sf->NextDDL() ){
124             while ( sf->NextChannel() ){
125                 UInt_t signal=0;
126                 while ( sf->NextBunch() ){
127                     for (UInt_t timebin=sf->GetStartTimeBin(); timebin<sf->GetEndTimeBin(); timebin++){
128
129                         AliTPCCalROC *roc = padhadd->GetCalROC(sf->GetSector());
130 //                      if ( roc->GetValue(sf->GetRow(), sf->GetPad()) == -1 )
131 //                          roc->SetValue(sf->GetRow(), sf->GetPad(), sf->GetHWAddress());
132 //                      else
133                             if ( roc->GetValue(sf->GetRow(), sf->GetPad()) != sf->GetHWAddress()){
134                                 printf("#%.2d: %.2d.%.3d.%.2d,%.4d - old [%.4d] now[%.4d]",
135                                        ievent, sf->GetSector(), sf->GetRow(),sf->GetPad(), timebin+1,
136                                        (Int_t)roc->GetValue(sf->GetRow(), sf->GetPad()), sf->GetHWAddress());
137                             }
138                             /*
139                             Int_t sig = sf->GetSignals()[timebin-sf->GetStartTimeBin()];
140                             if ( ievent == 19 && sf->GetSector() == 25 && sf->GetRow() == 00 && sf->GetPad()==9  ){
141                                 printf("new         | (%.2d.%.3d.%.2d.%.4d | %.3d ): %.3d\n",
142                                        sf->GetSector(), sf->GetRow(), sf->GetPad(), timebin+1, sf->GetHWAddress(), sig );
143                                 sum2+=sig;
144                             }
145                             */
146                         signal+=sf->GetSignals()[timebin-sf->GetStartTimeBin()];
147                     }
148                     padNew->GetCalROC(sf->GetSector())->SetValue(sf->GetRow(),sf->GetPad(),signal);
149                 }
150             }
151         }
152         timer1.Stop();timer3.Stop();
153         printf("new  --  Time: %.4f (%.4f)\n", timer1.RealTime(), timer1.CpuTime());
154
155         printf("sum1: %d, sum2: %d, diff: %d\n",sum1,sum2,sum1-sum2);
156
157         AliTPCCalPad com(*padOld);
158         com.Add(padNew, -1);
159
160         c1->cd(1);
161         padOld->MakeHisto2D(1)->Draw("colz");
162         c1->cd(2);
163         padNew->MakeHisto2D(1)->Draw("colz");
164         c1->cd(3);
165         com.MakeHisto2D(1)->Draw("colz");
166
167         //loop over all sectors, rows, pads
168         for ( UInt_t iSec=0; iSec<72; iSec++ )
169             for ( UInt_t iRow=0; iRow<com.GetCalROC(iSec)->GetNrows(); iRow++ )
170                 for ( UInt_t iPad=0; iPad<com.GetCalROC(iSec)->GetNPads(iRow); iPad++){
171                     Float_t val = com.GetCalROC(iSec)->GetValue(iRow, iPad);
172                     Float_t valo = padNew->GetCalROC(iSec)->GetValue(iRow, iPad);
173                     Float_t valn = padOld->GetCalROC(iSec)->GetValue(iRow, iPad);
174                     //check if values for old and new algorithm differ
175 //                  if ( val != 0 ){
176                     if ( (Int_t)valo != (Int_t)valn ){
177                         Float_t hadd = padhadd->GetCalROC(iSec)->GetValue(iRow,iPad);
178                         printf("Event: %.2d | (%.2d.%.3d.%.2d | %f ): %.3f\n", ievent, iSec, iRow, iPad, hadd, val);
179                         padSum->GetCalROC(iSec)->SetValue(iRow, iPad, padSum->GetCalROC(iSec)->GetValue(iRow, iPad)+1);
180                     }
181                     com.GetCalROC(iSec)->SetValue(iRow, iPad, 0);
182                     padOld->GetCalROC(iSec)->SetValue(iRow, iPad, 0);
183                     padNew->GetCalROC(iSec)->SetValue(iRow, iPad, 0);
184                 }
185
186
187         c1->cd(4);
188         padSum->MakeHisto2D(1)->Draw("colz");
189
190         c1->Modified();
191         c1->Update();
192   //      }// end event sel
193         if (input) ievent++;
194         ievent2++;
195     }
196     TFile f("output.root","recreate");
197     padSum->Write();
198     f.Save();
199     f.Close();
200
201     printf("total old  --  Time: %.4f (%.4f)\n", timer2.RealTime(), timer2.CpuTime());
202     printf("total new  --  Time: %.4f (%.4f)\n", timer3.RealTime(), timer3.CpuTime());
203
204     delete rawReader;
205     return padSum;
206 }
207