]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/EVE/AliHLTEveTRD.cxx
4c0c94963ad9ed0802350debd876393c667048f8
[u/mrichter/AliRoot.git] / HLT / EVE / AliHLTEveTRD.cxx
1 /**************************************************************************
2  * This file is property of and copyright by the ALICE HLT Project        *
3  * ALICE Experiment at CERN, All rights reserved.                         *
4  *                                                                        *
5  * Primary Authors: Svein Lindal <slindal@fys.uio.no   >                  *
6  *                  for The ALICE HLT Project.                            *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /// @file   AliHLTEvePhos.cxx
18 /// @author Svein Lindal <slindal@fys.uio.no>
19 /// @brief  TPC processor for the HLT EVE display
20
21 #include "AliHLTEveTRD.h"
22 #include "AliHLTHOMERBlockDesc.h"
23 #include "TCanvas.h"
24 #include "AliHLTEveBase.h"
25 #include "AliEveHOMERManager.h"
26 #include "TEveManager.h"
27 #include "TEvePointSet.h"
28 #include "TColor.h"
29 #include "TMath.h"
30 #include "TH1F.h"
31 #include "AliHLTTRDCluster.h"
32 #include "AliTRDcluster.h"
33
34 ClassImp(AliHLTEveTRD)
35
36 AliHLTEveTRD::AliHLTEveTRD() : 
37   AliHLTEveBase(), 
38   fEveClusters(NULL),
39   fEveColClusters(NULL),
40   fNColorBins(15)
41 {
42   // Constructor.
43
44 }
45
46 AliHLTEveTRD::~AliHLTEveTRD()
47 {
48   //Destructor, not implemented
49   if(fEveColClusters)
50     delete fEveColClusters;
51   fEveColClusters = NULL;
52   
53   if(fEveClusters)
54     delete fEveClusters;
55   fEveClusters = NULL;
56 }
57
58
59 void AliHLTEveTRD::ProcessBlock(AliHLTHOMERBlockDesc * block) {
60   //See header file for documentation
61
62   if ( ! block->GetDataType().CompareTo("CLUSTERS") ) {
63     
64     if(!fEveColClusters){
65       fEveColClusters = CreatePointSetArray();
66       fEventManager->GetEveManager()->AddElement(fEveColClusters);
67     } 
68     
69     ProcessClusters(block, fEveColClusters);
70    
71   } else if ( ! block->GetDataType().CompareTo("ROOTHIST") ) {
72   
73     if(!fCanvas) {
74       fCanvas = CreateCanvas("TRD QA", "TRD QA");
75       fCanvas->Divide(3, 2);
76     }
77
78     AddHistogramsToCanvas(block, fCanvas, fHistoCount);
79                    
80   }
81   
82 }
83
84
85 void AliHLTEveTRD::AddHistogramsToCanvas(AliHLTHOMERBlockDesc* block, TCanvas * canvas, Int_t &cdCount ) {
86   //See header file for documentation
87   if ( ! block->GetClassName().CompareTo("TH1F")) {
88     TH1F* histo = reinterpret_cast<TH1F*>(block->GetTObject());
89     ++cdCount;
90   
91     TVirtualPad* pad = canvas->cd(cdCount);
92     histo->Draw();
93     pad->SetGridy();
94     pad->SetGridx();
95
96     if ( ! strcmp(histo->GetName(), "nscls") ) {
97       histo->GetXaxis()->SetRangeUser(0.,15.);
98     }
99
100     if ( ! strcmp(histo->GetName(),"sclsdist") ||
101          ! strcmp(histo->GetName(),"evSize") )
102       pad->SetLogy();
103   }
104
105 }
106
107
108
109 TEvePointSet * AliHLTEveTRD::CreatePointSet() {
110   //See header file for documentation
111   TEvePointSet * ps = new TEvePointSet("TRD Clusters");
112   ps->SetMainColor(kBlue);
113   ps->SetMarkerStyle((Style_t)kFullDotSmall);
114  
115   return ps;
116
117 }
118
119 TEvePointSetArray * AliHLTEveTRD::CreatePointSetArray(){
120   //See header file for documentation
121   TEvePointSetArray * cc = new TEvePointSetArray("TRD Clusters Colorized");
122   cc->SetMainColor(kRed);
123   cc->SetMarkerStyle(4); // antialiased circle
124   cc->SetMarkerSize(0.4);
125   cc->InitBins("Cluster Charge", fNColorBins, 0., fNColorBins*100.);
126   
127   const Int_t nCol = TColor::GetNumberOfColors();
128   for (Int_t ii = 0; ii < fNColorBins + 1; ++ii) {
129     cc->GetBin(ii)->SetMainColor(TColor::GetColorPalette(ii * nCol / (fNColorBins+2)));
130   }
131
132   return cc;
133      
134 }
135
136
137 void AliHLTEveTRD::UpdateElements() {
138   //See header file for documentation
139   if(fCanvas) fCanvas->Update();
140   // if(fEveClusters) fEveClusters->ResetBBox();
141
142   for (Int_t ib = 0; ib <= fNColorBins + 1; ++ib) {
143     fEveColClusters->GetBin(ib)->ResetBBox();
144   }
145
146 }
147
148 void AliHLTEveTRD::ResetElements(){
149   //See header file for documentation
150   // if(fEveClusters) fEveClusters->Reset();
151  
152   if(fEveColClusters){
153     for (Int_t ib = 0; ib <= fNColorBins + 1; ++ib) {
154       fEveColClusters->GetBin(ib)->Reset();
155     }
156   }
157
158   fHistoCount = 0;
159 }
160
161 Int_t AliHLTEveTRD::ProcessClusters( AliHLTHOMERBlockDesc * block, TEvePointSetArray * contCol ){
162   //See header file for documentation
163
164   Int_t iResult = 0;
165
166   Int_t sm = block->GetSubDetector();
167   if ( sm == 6 ) sm = 7;
168   
169   Float_t phi   = ( sm + 0.5 ) * TMath::Pi() / 9.0;  
170   Float_t cos   = TMath::Cos( phi );
171   Float_t sin   = TMath::Sin( phi );
172   
173   Byte_t* ptrData = reinterpret_cast<Byte_t*>(block->GetData());
174   UInt_t ptrSize = block->GetSize();
175
176   for (UInt_t size = 0; size+sizeof(AliHLTTRDCluster) <= ptrSize; size+=sizeof(AliHLTTRDCluster) ) {
177     AliHLTTRDCluster *cluster = reinterpret_cast<AliHLTTRDCluster*>(&(ptrData[size]));
178    
179     AliTRDcluster *trdCluster = new AliTRDcluster;
180     cluster->ExportTRDCluster( trdCluster );
181    
182     contCol->Fill(cos*trdCluster->GetX() - sin*trdCluster->GetY(), 
183                    sin*trdCluster->GetX() + cos*trdCluster->GetY(), 
184                    trdCluster->GetZ(),
185                    trdCluster->GetQ() );    
186     
187     //cont->SetNextPoint(cos*trdCluster->GetX() - sin*trdCluster->GetY(), 
188     //         sin*trdCluster->GetX() + cos*trdCluster->GetY(), trdCluster->GetZ());
189   }
190   
191   return iResult;
192
193
194
195
196
197  
198   return 0;  
199
200
201 }
202
203