]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TPCLib/OnlineDisplay/AliHLTTPCDisplayPadRow.cxx
- configure adapted to the new directory structure of the HOMER module in PubSub
[u/mrichter/AliRoot.git] / HLT / TPCLib / OnlineDisplay / AliHLTTPCDisplayPadRow.cxx
1 // $Id$
2
3 /** \class AliHLTTPCDisplayPadRow
4 <pre>
5 //_____________________________________________________________
6 // AliHLTTPCDisplayPadRow
7 //
8 // Display class for the HLT TPC-PadRow events.
9 </pre>
10 */
11 // Author: Jochen Thaeder <mailto:thaeder@kip.uni-heidelberg.de>
12 //*-- Copyright &copy ALICE HLT Group 
13
14 #include <TH2.h>
15 #include <TFile.h>
16 #include <TStyle.h>
17 #include <TAttText.h>
18 #include <TAxis.h>
19 #include <TCanvas.h>
20 #include <TPolyMarker3D.h>
21
22 #ifdef use_aliroot
23 #include <TClonesArray.h>
24 #include <AliRun.h>
25 #include <AliSimDigits.h>
26 #include <AliTPCParam.h>
27 #endif
28
29 #include "AliHLTStdIncludes.h"
30 #include "AliHLTTPCDefinitions.h"
31 #include "AliHLTDataTypes.h"
32 #include "AliHLTTPCSpacePointData.h"
33 #include "AliHLTTPCClusterDataFormat.h"
34 #include "AliHLTTPCLogging.h"
35 #include "AliHLTTPCTransform.h"
36 #include "AliHLTTPCDigitReaderPacked.h"
37 #include "AliHLTTPCDigitReaderRaw.h"
38
39
40 #include "AliHLTTPCDisplayMain.h"
41 #include "AliHLTTPCDisplayPadRow.h"
42 #include "AliHLTTPCDisplayPad.h"
43
44 #if __GNUC__ >= 3
45 using namespace std;
46 #endif
47
48 ClassImp(AliHLTTPCDisplayPadRow)
49
50 //____________________________________________________________________________________________________
51 AliHLTTPCDisplayPadRow::AliHLTTPCDisplayPadRow(AliHLTTPCDisplayMain* display) {
52     // constructor
53     fDisplay = display;
54
55     fNTimes = AliHLTTPCTransform::GetNTimeBins();
56     
57     fBinX[0] = 0;      
58     fBinX[1] = fDisplay->GetNPads(); 
59     fBinY[0] = 0;
60     fBinY[1] = fNTimes -1;
61     fTmpEvent = 0;         
62
63     for (Int_t ii = 0;ii < 20; ii++)
64       fpmarr[ii] = 0;
65
66     Int_t maxpads = 150;
67     Int_t padbinning = maxpads*10;
68
69     fHistraw = new TH2F("fHistraw","Selected PadRow with found Clusters;Pad #;Timebin #",maxpads,0,maxpads-1,fNTimes,0,fNTimes-1);
70     fHistrawcl = new TH1F("fHistrawcl","cvcv;ddd;kkk",padbinning,0,maxpads-1);
71
72     fHistraw->SetOption("COLZ");    
73     fHistraw->SetTitleSize(0.03);
74     fHistraw->GetXaxis()->SetLabelSize(0.03);
75     fHistraw->GetXaxis()->SetTitleSize(0.03);
76     fHistraw->GetYaxis()->SetLabelSize(0.03);
77     fHistraw->GetYaxis()->SetTitleSize(0.03);
78
79     gStyle->SetPalette(1);
80 }
81
82 //____________________________________________________________________________________________________
83 AliHLTTPCDisplayPadRow::~AliHLTTPCDisplayPadRow() {
84     // destructor     
85     if ( fHistraw ){
86         delete fHistraw;
87         fHistraw = NULL;
88     }
89     if ( fHistrawcl ){
90         delete fHistrawcl;
91         fHistrawcl = NULL;
92     }
93     for (Int_t ii=19; ii <= 0; ii--){ // try fwd
94       if ( fpmarr[ii] ){
95         delete [] fpmarr[ii];
96         fpmarr[ii] = NULL;
97       }
98     }
99 }
100
101 //____________________________________________________________________________________________________
102 void AliHLTTPCDisplayPadRow::Reset(){  
103   fHistraw->Reset();   
104   fHistrawcl->Reset(); 
105 }
106
107 //____________________________________________________________________________________________________
108 void AliHLTTPCDisplayPadRow::Save(){  
109   fDisplay->GetCanvasCharge()->SaveAs("HLT-PadRowView.eps"); 
110 }
111
112 //____________________________________________________________________________________________________
113 void AliHLTTPCDisplayPadRow::Fill(Int_t patch, ULong_t dataBlock, ULong_t dataLen){
114     // Fill PadRow Histogram
115     
116 #if defined(HAVE_TPC_MAPPING)
117     AliHLTTPCDigitReaderRaw digitReader(0);
118
119     bool readValue = true;
120     Int_t rowOffset = 0;
121     
122     Int_t padRow = fDisplay->GetPadRow();
123     Int_t slice = fDisplay->GetSlicePadRow();
124
125     // Initialize RAW DATA
126     Int_t firstRow = AliHLTTPCTransform::GetFirstRow(patch);
127     Int_t lastRow = AliHLTTPCTransform::GetLastRow(patch);
128
129     // Outer sector, patches 2, 3, 4, 5 -  start counting in patch 2 with row 0
130     if ( patch >= 2 ) rowOffset = AliHLTTPCTransform::GetFirstRow( 2 );
131
132     // Initialize block for reading packed data
133     void* tmpdataBlock = (void*) dataBlock;
134     digitReader.InitBlock(tmpdataBlock,dataLen,firstRow,lastRow,patch,slice);
135
136     readValue = digitReader.Next();
137
138     if (!readValue){    
139         LOG(AliHLTTPCLog::kError,"AliHLTTPCDisplayPadRow::Fill","Read first value") << "No value in data block" << ENDLOG;
140         return;
141     }
142
143     if ( fDisplay->GetZeroSuppression() ){
144       for (Int_t pad=0; pad < AliHLTTPCTransform::GetNPads(padRow); pad++){
145         //      for (Int_t timeBin=fDisplay->GetTimeBinMin(); timeBin <= fDisplay->GetTimeBinMax(); timeBin++){
146         for (Int_t timeBin=0; timeBin <= fDisplay->GetNTimeBins(); timeBin++){
147
148           fHistraw->Fill(pad,timeBin,fDisplay->fRawDataZeroSuppressed[padRow][pad][timeBin]);
149         } // end time
150       }  // end pad
151       
152     }  // end - if ( fDisplay->GetZeroSuppression() ){
153
154     else {
155       for (Int_t pad=0; pad < AliHLTTPCTransform::GetNPads(padRow); pad++){
156         //      for (Int_t timeBin=fDisplay->GetTimeBinMin(); timeBin <= fDisplay->GetTimeBinMax(); timeBin++){
157         for (Int_t timeBin=0; timeBin <= fDisplay->GetNTimeBins(); timeBin++){
158           fHistraw->Fill(pad,timeBin,fDisplay->fRawData[padRow][pad][timeBin]);
159         } // end time
160       }  // end pad
161     }  // end - else of if ( fDisplay->GetZeroSuppression() ){
162
163
164     // FILL PADROW 3D --- Initialize the colorbins
165     if ( fDisplay->Get3DSwitchPadRow() ){
166         for (UInt_t ii=0;ii < 20;ii++){
167             fbinct[ii] = 0;
168             fcolorbin[ii] = 0;
169         }
170
171         // read number of entries in colorbin
172         while ( readValue ){ 
173             
174             Int_t row = digitReader.GetRow() + rowOffset;
175             
176             if (row == padRow){    
177                 UInt_t charge = digitReader.GetSignal();
178                 
179                 for (UInt_t ii=0;ii < 19;ii++){
180                     if ( charge > (ii*15) && charge <= ((ii*15) + 15) ) fcolorbin[ii]++;
181                 }
182                 // larger than 19 * 15  
183                 if (charge > 285 ) fcolorbin[19]++;
184             }
185
186             // read next value
187             readValue = digitReader.Next();
188       
189             if(!readValue) break; //No more value
190         } 
191
192         // Initialize fpmarr[color][3*colorbin[ii]] 
193         for (Int_t ii = 0; ii < 20; ii++) {
194           if ( fpmarr[ii] ) delete[] fpmarr[ii];
195           fpmarr[ii] = new Float_t[fcolorbin[ii]*3]; 
196         }
197
198         // Rewind the raw reader and fill the polymarker3D
199         digitReader.InitBlock(tmpdataBlock,dataLen,firstRow,lastRow,patch,slice);
200         
201         readValue = digitReader.Next();
202     } // END if (fDisplay->Get3DSwitchPadRow())
203
204     // -- Fill Raw Data
205     while ( readValue ){ 
206
207         Int_t row = digitReader.GetRow() + rowOffset;
208
209         // select padrow to fill in histogramm
210         if (row == padRow){    
211             UChar_t pad = digitReader.GetPad();
212             UShort_t time = digitReader.GetTime();
213             UInt_t charge = digitReader.GetSignal();
214             Float_t xyz[3];
215             //fHistraw->Fill(pad,time,charge);
216
217             if ( fDisplay->Get3DSwitchPadRow() ) {
218                 // Transform raw coordinates to local coordinates
219                 AliHLTTPCTransform::RawHLT2Global(xyz, slice, padRow, pad, time);
220
221                 for (UInt_t ii=0;ii < 19;ii++){
222                     if ( charge > (ii*15) && charge <= ((ii*15) + 15) ){
223                         fpmarr[ii][fbinct[ii]] = xyz[0];
224                         fpmarr[ii][fbinct[ii]+1] = xyz[1];
225                         fpmarr[ii][fbinct[ii]+2] = xyz[2];
226                         fbinct[ii] += 3;
227                     }
228                 }
229                 // larger than 19 * 15
230                 if (charge > 285 ) {
231                     fpmarr[19][fbinct[19]] = xyz[0];
232                     fpmarr[19][fbinct[19]+1] = xyz[1];
233                     fpmarr[19][fbinct[19]+2] = xyz[2];
234                     fbinct[19] += 3;
235                 }
236             } // END if (fSwitch3DPadRow)
237         
238         } // END if (row == padRow){ 
239         
240         // read next value
241         readValue = digitReader.Next();
242       
243         //Check where to stop:
244         if(!readValue) break; //No more value
245     } 
246     
247      if (fDisplay->ExistsClusterData()){
248         AliHLTTPCSpacePointData *points = fDisplay->GetSpacePointDataPointer(slice,patch);
249         if(!points) return;
250         
251         Float_t xyz[3];
252         for(Int_t i=0; i< fDisplay->GetNumberSpacePoints(slice,patch); i++){
253             xyz[0] = points[i].fX;
254             xyz[1] = points[i].fY;
255             xyz[2] = points[i].fZ;
256         
257             // select padrow to fill in histogramm
258             if (padRow == AliHLTTPCTransform::GetPadRow(xyz[0])){
259                 AliHLTTPCTransform::LocHLT2Raw(xyz, slice, padRow);
260                 fHistrawcl->Fill(xyz[1],xyz[2]);
261             }
262         }
263      } // END if (fDisplay->ExistsClusterData()){
264 #else //! defined(HAVE_TPC_MAPPING)
265       HLTFatal("DigitReaderRaw not available - check your build");
266 #endif //defined(HAVE_TPC_MAPPING)
267  
268 }
269
270 //____________________________________________________________________________________________________
271 void AliHLTTPCDisplayPadRow::Draw(){
272
273     fDisplay->GetCanvasPadRow()->cd();
274     fDisplay->GetCanvasPadRow()->Clear();
275
276     if (fDisplay->GetSplitPadRow()){
277         fDisplay->GetCanvasPadRow()->Divide(1,2);
278         fDisplay->GetCanvasPadRow()->cd(1);
279     }
280   
281     Char_t title[256];
282     sprintf(title,"Selected PadRow %d with found Clusters",fDisplay->GetPadRow());
283   
284     // Keep Zoom
285     if (!fDisplay->GetKeepView() || (fBinX[1]>fDisplay->GetNPads()) || (fBinY[1]>fNTimes)){
286         fBinX[0] = 0;
287         fBinX[1] = fDisplay->GetNPads();
288         fBinY[0] = 0;
289         fBinY[1] = fNTimes - 1;
290     }
291     
292     fHistraw->SetAxisRange(fBinX[0],fBinX[1]); 
293     fHistraw->SetAxisRange(fBinY[0],fBinY[1],"Y");
294     fHistraw->SetTitle(title);
295     fHistraw->SetStats(kFALSE);
296     fHistraw->Draw("COLZ");
297     
298     if ( fDisplay->ExistsClusterData() ){
299       //cout << "XX" <<  endl;
300         fHistrawcl->SetAxisRange(fBinX[0],fBinX[1]);
301         fHistrawcl->SetAxisRange(fBinY[0],fBinY[1],"Y");
302         fHistrawcl->SetStats(kFALSE);
303         fHistrawcl->SetMarkerStyle(28);
304         fHistrawcl->SetMarkerSize(2);
305         fHistrawcl->SetMarkerColor(1);
306         fHistrawcl->Draw("psame");
307     }
308     
309     if (fDisplay->GetSplitPadRow()){
310         fDisplay->GetCanvasPadRow()->cd(2);
311         fDisplay->GetPadPointer()->fHistpad2->Draw();
312     }
313     
314     fDisplay->GetCanvasPadRow()->Modified();
315     fDisplay->GetCanvasPadRow()->Update();
316     
317     // Select Pad
318     fDisplay->GetCanvasPadRow()->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)","AliHLTTPCDisplayMain",(void*) fDisplay,"ExecPadEvent(Int_t,Int_t,Int_t,TObject*)");
319     // Keep Zoom
320     fDisplay->GetCanvasPadRow()->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)","AliHLTTPCDisplayPadRow",(void*) this,"ExecEvent(Int_t,Int_t,Int_t,TObject*)");
321 }
322
323 //____________________________________________________________________________________________________
324 void AliHLTTPCDisplayPadRow::Draw3D(){
325     Int_t markercolor = 51;
326     
327     for (UInt_t ii=0;ii < 20;ii++){
328         if (fcolorbin[ii]> 0){
329             
330             TPolyMarker3D* pm = new  TPolyMarker3D(fcolorbin[ii], fpmarr[ii], 7 );
331             pm->SetBit(kCanDelete);
332
333             pm->SetMarkerColor(markercolor); 
334             pm->Draw(""); 
335         }
336
337         // in order to have the SetPalette(1), so called "pretty"
338         if (ii % 2 == 0 ) markercolor += 2;
339         else  markercolor += 3;
340     }
341 }
342
343 //____________________________________________________________________________________________________
344 void AliHLTTPCDisplayPadRow::ExecEvent(Int_t event, Int_t px, Int_t py, TObject *selected){
345    // Saves the Zoom Position of the Histogram 
346
347    // - Mouse down on Axis : StartPoint of Range
348    if (event == 1 && selected->InheritsFrom("TAxis"))
349        fTmpEvent = 1;
350
351    // - Mouse pressed on Axis : Real Zoom process not only click
352    else if (event == 21 && selected->InheritsFrom("TAxis") && fTmpEvent == 1)
353        fTmpEvent = 21;
354
355    // - Mouse pressed on Axis : Still pressed
356    else if (event == 21 && selected->InheritsFrom("TAxis") && fTmpEvent == 21) 
357        return;
358
359     // - Mouse up on Axis : End Point of Rangex
360    else if (event == 11 && selected->InheritsFrom("TAxis") && fTmpEvent == 21){
361        TAxis * axis = (TAxis*) selected;
362        
363        if (selected == fHistraw->GetXaxis() || selected == fHistrawcl->GetXaxis()){ 
364            fBinX[0] = axis->GetFirst() -1;
365            fBinX[1] = axis->GetLast() -1;
366        }
367        else {
368            fBinY[0] = axis->GetFirst() -1;
369            fBinY[1] = axis->GetLast() -1;
370        }
371
372        
373        //    Reset();
374 //     Fill();
375        //   Draw();
376
377        fTmpEvent = 0; 
378    }
379    else fTmpEvent = 0;
380 }