]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/EVE/AliHLTEveCalo.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / HLT / EVE / AliHLTEveCalo.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   AliHLTEveCalo.cxx
18 /// @author Svein Lindal <slindal@fys.uio.no>
19 /// @brief  Calorimeter base class for the HLT EVE display
20
21 #include "TCollection.h"
22 #include "TObjArray.h"
23 #include "AliHLTEveCalo.h"
24 #include "AliHLTHOMERBlockDesc.h"
25 #include "TCanvas.h"
26 #include "AliHLTEveBase.h"
27 #include "TEveBoxSet.h"
28 #include "AliPHOSGeometry.h"
29 #include "TVector3.h"
30 #include "AliEveHLTEventManager.h"
31 #include "TEveManager.h"
32 #include "AliHLTCaloDigitDataStruct.h"
33 #include "AliHLTCaloClusterDataStruct.h"
34 #include "AliHLTCaloClusterReader.h"
35 #include "TEveTrans.h"
36 #include "TString.h"
37 #include "TH2F.h"
38 #include "TH1F.h"
39 #include "TRefArray.h"
40 #include "AliESDEvent.h"
41 #include "AliESDCaloCluster.h"
42
43 ClassImp(AliHLTEveCalo);
44
45 AliHLTEveCalo::AliHLTEveCalo(Int_t nm, TString name) : 
46   AliHLTEveBase(name), 
47   fBoxSetDigits(NULL),
48   fBoxSetClusters(NULL),
49   fNModules(nm),
50   fClustersArray(NULL),
51   fName(name), 
52   fPadTitles(NULL),
53   fInvMassCanvas(NULL),
54   fClusterReader(NULL)
55 {
56   // Constructor.
57
58   SetMaxHistograms(9);
59
60   fPadTitles = new TString[GetMaxHistograms()];
61  
62   for(int i = 0; i < GetMaxHistograms(); i++) {
63     fPadTitles[i] = "";
64   }
65
66   fClustersArray = new TRefArray();
67   fClusterReader = new AliHLTCaloClusterReader();
68
69
70
71
72 }
73
74 AliHLTEveCalo::~AliHLTEveCalo()
75 {
76   //Destructor
77   if(fBoxSetDigits)
78     delete fBoxSetDigits;
79   fBoxSetDigits = NULL;
80   
81   if(fBoxSetClusters)
82     delete fBoxSetClusters;
83   fBoxSetClusters = NULL;
84
85
86   if(fPadTitles)
87     delete [] fPadTitles;
88   fPadTitles = NULL;
89
90
91   if(fClusterReader)
92     delete fClusterReader;
93   fClusterReader = NULL;
94
95 }
96
97
98 void AliHLTEveCalo::ProcessBlock(AliHLTHOMERBlockDesc * block) {
99   //See header file for documentation
100
101   if ( block->GetDataType().CompareTo("ROOTHIST") == 0 ) { 
102     ProcessHistogram(block);
103    
104   } else {
105
106     
107     if ( block->GetDataType().CompareTo("CALOCLUS") == 0 ){
108       //cout <<"Skipping calo clusters"<<endl;
109       ProcessClusters( block );
110     } else if ( block->GetDataType().CompareTo("DIGITTYP") == 0 ) {
111       //ProcessDigits( block);
112       //
113     } else if ( block->GetDataType().CompareTo("CHANNELT") == 0 ) {
114       //ProcessClusters( block );
115     } else if (!block->GetDataType().CompareTo("ALIESDV0")) {
116       ProcessEsdBlock(block);
117     }
118   }
119 }
120
121 void AliHLTEveCalo::ProcessHistogram(AliHLTHOMERBlockDesc * block ) {
122   //See header file for documentation
123   
124   if(!fCanvas) {
125     fCanvas = CreateCanvas(Form("%s QA", fName.Data()), Form("%s QA", fName.Data()));
126     fCanvas->Divide(3, 3);
127   }
128
129   if(!fInvMassCanvas) {
130     fInvMassCanvas = CreateCanvas(Form("%s IM", fName.Data()), Form("%s IM", fName.Data()));
131     fInvMassCanvas->Divide(3, 2);
132   }
133
134
135   AddHistogramsToCanvas(block, fCanvas, fHistoCount);
136
137
138 }
139
140
141 // void AliHLTEveCalo::ProcessDigits(AliHLTHOMERBlockDesc* block) {
142 //   //See header file for documentation
143   
144 //   AliHLTCaloDigitDataStruct *ds = reinterpret_cast<AliHLTCaloDigitDataStruct*> (block->GetData());
145 //   UInt_t nDigits = block->GetSize()/sizeof(AliHLTCaloDigitDataStruct);
146     
147
148 //   for(UInt_t i = 0; i < nDigits; i++, ds++) {
149
150 //     Float_t x = (ds->fX - 32)* 2.2;
151 //       Float_t z = (ds->fZ - 28) * 2.2;
152
153
154 //     fBoxSetDigits[4-ds->fModule].AddBox(x, 0, z, 2.2, ds->fEnergy*200, 2.2);
155 //     fBoxSetDigits[4-ds->fModule].DigitValue(static_cast<Int_t>(ds->fEnergy*10));
156 //   }
157
158 // }
159
160 void AliHLTEveCalo::ProcessEsdBlock(AliHLTHOMERBlockDesc * block) {
161   AliESDEvent * event = dynamic_cast<AliESDEvent*>(block->GetTObject());
162   if (event) {
163     ProcessEvent(event);
164   } else {
165     cout << "problem getting the event!"<<endl;
166   }
167
168 }
169
170 void AliHLTEveCalo::ProcessEvent(AliESDEvent * event) {
171   //see header file for documentation
172
173
174
175   Int_t nClusters = GetClusters(event, fClustersArray);
176   for(int ic = 0; ic < nClusters; ic++) {
177     AliESDCaloCluster * cluster = dynamic_cast<AliESDCaloCluster*>(fClustersArray->At(ic));
178     ProcessESDCluster(cluster);
179   }
180   
181 }
182
183
184 void AliHLTEveCalo::ProcessClusters(AliHLTHOMERBlockDesc* block) {
185   //See header file for documentation
186
187   AliHLTCaloClusterHeaderStruct *dh = reinterpret_cast<AliHLTCaloClusterHeaderStruct*> (block->GetData());
188   fClusterReader->SetMemory(dh);  
189
190   AliHLTCaloClusterDataStruct * ds;
191
192   while( (ds = fClusterReader->NextCluster()) ){
193      AddClusters(ds->fGlobalPos, ds->fModule, ds->fEnergy);
194   }
195
196   AliHLTCaloDigitDataStruct *dg = fClusterReader->GetDigits();
197   UInt_t nDigits = fClusterReader->GetNDigits();;
198   for(UInt_t i = 0; i < nDigits; i++, dg++) {
199     AddDigits(dg->fX, dg->fZ, dg->fModule, dg->fEnergy);
200   }
201 }
202
203 void AliHLTEveCalo::UpdateElements() {
204   //See header file for documentation
205   if(fCanvas) fCanvas->Update();
206   if(fInvMassCanvas) fInvMassCanvas->Update();
207
208
209   if(fBoxSetDigits) {
210     for(int im = 0; im < fNModules; im++) {
211       fBoxSetDigits[im].ElementChanged();
212     }
213   }
214
215   if(fBoxSetClusters) {
216     for(int im = 0; im < fNModules; im++) {
217       fBoxSetClusters[im].ElementChanged();
218     }
219   }
220
221 }
222
223 void AliHLTEveCalo::ResetElements(){
224   //See header file for documentation
225   fHistoCount = 0;
226   
227   if ( fBoxSetDigits ){
228     for(int im = 0; im < fNModules; im++){
229       fBoxSetDigits[im].Reset();   
230     }
231   }
232
233   if ( fBoxSetClusters ){
234     for(int im = 0; im < fNModules; im++){
235       fBoxSetClusters[im].Reset();   
236     }
237   }
238
239 }
240
241 Int_t AliHLTEveCalo::GetPadNumber(TString name) {
242
243
244   for(int i = 0; i < GetMaxHistograms(); i++) {
245     if (!fPadTitles[i].CompareTo(name)){
246       return i+1;
247     }
248     else if (!fPadTitles[i].CompareTo("")) {
249       //cout <<"in empty title"<<endl;
250       fPadTitles[i] = name;
251       return i+1;
252     }
253   }
254   
255   cout << "AliHLTEveCalo()->GetPadNUmber():returning one"<<endl;
256   return 1;
257
258 }
259
260 void AliHLTEveCalo::DrawInvMassHistogram(TH1F * histo) {
261   
262   fInvMassCanvas->cd(++fHistoCount);
263
264   histo->SetAxisRange(histo->GetXaxis()->GetBinLowEdge(histo->FindFirstBinAbove(0, 1) - 3), histo->GetXaxis()->GetBinUpEdge(histo->FindLastBinAbove(0, 1) + 3), "X");
265   //histo->Fit("gaus", "", "", histo->GetXaxis()->GetBinLowEdge(histo->FindFirstBinAbove(0, 1)), histo->GetXaxis()->GetBinUpEdge(histo->FindLastBinAbove(0, 1)));
266
267   histo->Draw();
268
269 }
270
271 void AliHLTEveCalo::AddHistogramsToCanvas(AliHLTHOMERBlockDesc * block, TCanvas * canvas, Int_t &/*cdCount*/) {
272   //See header file for documentation
273
274   if ( ! block->GetClassName().CompareTo("TObjArray")) {
275     TIter next((TObjArray*)(block->GetTObject()));
276     TObject *object;
277
278
279    
280     while (( object = (TObject*) next())) {
281
282       TString name = static_cast<TH1*>(object)->GetName();
283       if(name.Contains("InvMass")){
284         DrawInvMassHistogram(static_cast<TH1F*>(object));
285         
286       } else {
287
288         
289         Int_t iPad = GetPadNumber(name);
290         canvas->cd(iPad);
291         
292         
293         //Check if histo is 2D histo
294         TH2F* histo2 = dynamic_cast<TH2F*>(object);
295         if(histo2){
296           
297           Int_t lb = histo2->FindLastBinAbove(0,1);
298           if(lb > -1) {
299             histo2->SetAxisRange(0, histo2->GetXaxis()->GetBinUpEdge(histo2->FindLastBinAbove(0, 1) + 3), "X");
300             histo2->SetAxisRange(0, histo2->GetYaxis()->GetBinUpEdge(histo2->FindLastBinAbove(0, 2) + 3), "Y");
301           }
302           
303           histo2->Draw("COLZ");
304         }
305         
306         
307         //Must be 1D histo
308         else {
309           TH1F* histo = dynamic_cast<TH1F*>(object);
310           if (histo) {
311             
312             TString name2 = histo->GetName();
313             
314             if(name2.Contains("Energy")) {
315               histo->SetAxisRange(0, histo->GetXaxis()->GetBinUpEdge(histo->FindLastBinAbove(0, 1) + 3), "X");
316             }
317             
318             
319             histo->Draw();
320           } else {
321             cout <<"AliHLTEveCaloBase::AddHistogramsTocCanvas: Histogram neither TH1F nor TH2F"<<endl;
322           }
323         }
324       }
325     }
326
327   } else if ( ! block->GetClassName().CompareTo("TH1F")) {
328
329     TH1F* histo = reinterpret_cast<TH1F*>(block->GetTObject());
330     if(histo) {
331       
332       Int_t iPad = GetPadNumber(histo->GetName());
333       canvas->cd(iPad);
334       histo->Draw();
335     }
336   
337   
338   } else if ( ! block->GetClassName().CompareTo("TH2F")) {
339     
340     TH2F *histo = reinterpret_cast<TH2F*>(block->GetTObject());
341     if(histo) {
342       
343       Int_t iPad = GetPadNumber(histo->GetName());
344       canvas->cd(iPad);
345       histo->Draw();
346     }
347
348
349   }
350   
351   canvas->cd();
352 }
353