]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONTrackerDataHistogrammer.cxx
correction of trivial typo preventing compilation
[u/mrichter/AliRoot.git] / MUON / AliMUONTrackerDataHistogrammer.cxx
CommitLineData
10eb3d17 1/**************************************************************************
2* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3* *
4* Author: The ALICE Off-line Project. *
5* Contributors are mentioned in the code where appropriate. *
6* *
7* Permission to use, copy, modify and distribute this software and its *
8* documentation strictly for non-commercial purposes is hereby granted *
9* without fee, provided that the above copyright notice appears in all *
10* copies and that both the copyright notice and this permission notice *
11* appear in the supporting documentation. The authors make no claims *
12* about the suitability of this software for any purpose. It is *
13* provided "as is" without express or implied warranty. *
14**************************************************************************/
15
16// $Id$
17
18#include "AliMUONTrackerDataHistogrammer.h"
19
49419555 20#include "AliLog.h"
10eb3d17 21#include "AliMUONPainterGroup.h"
22#include "AliMUONSparseHisto.h"
23#include "AliMUONVPainter.h"
24#include "AliMUONVTrackerData.h"
25#include "AliMpBusPatch.h"
26#include "AliMpConstants.h"
27#include "AliMpDDLStore.h"
28#include "AliMpDEIterator.h"
29#include "AliMpDetElement.h"
30#include "AliMpManuUID.h"
49419555 31#include <TClass.h>
10eb3d17 32#include <TH1.h>
33#include <TObjArray.h>
49419555 34#include <TROOT.h>
10eb3d17 35
36///\class AliMUONTrackerDataHistogrammer
37///
38/// Class to generate histograms from AliMUONVTrackerData
39/// (and AliMUONVPainter) objects
40///
41/// \author Laurent Aphecetche, Subatech
42///
43
44///\cond CLASSIMP
45ClassImp(AliMUONTrackerDataHistogrammer)
46///\endcond CLASSIMP
47
48//_____________________________________________________________________________
49AliMUONTrackerDataHistogrammer::AliMUONTrackerDataHistogrammer(const AliMUONVTrackerData& data,
49419555 50 Int_t externalDim,
51 Int_t internalDim)
10eb3d17 52: TObject(),
53fData(data),
49419555 54fExternalDim(externalDim),
55fInternalDim(internalDim)
10eb3d17 56{
57 /// ctor
58}
59
60//_____________________________________________________________________________
61AliMUONTrackerDataHistogrammer::~AliMUONTrackerDataHistogrammer()
62{
63 /// dtor
64}
65
66//_____________________________________________________________________________
67void
68AliMUONTrackerDataHistogrammer::Add(TH1& h, const AliMUONSparseHisto& sh) const
69{
70 /// Add sparse histo content to histogram.
71
72 Double_t entries(h.GetEntries());
73
74 for ( Int_t i = 0; i < sh.GetNbins(); ++i )
75 {
76 Int_t count = sh.GetBinContent(i);
77
78 h.Fill(sh.GetBinCenter(i),count);
79
80 entries += count;
81 }
82
83 h.SetEntries(entries);
84
85 if (sh.HasUnderflow()) h.SetBinContent(0,1);
86 if (sh.HasOverflow()) h.SetBinContent(h.GetNbinsX()+1,1);
87}
88
89//_____________________________________________________________________________
90void
91AliMUONTrackerDataHistogrammer::AddBusPatchHisto(TH1& h, Int_t busPatchId) const
92{
93 /// Add data from one bus patch to the histogram
94
95 if ( fData.HasBusPatch(busPatchId ) )
96 {
97 AliMpBusPatch* busPatch = AliMpDDLStore::Instance()->GetBusPatch(busPatchId);
98 for ( Int_t i = 0; i < busPatch->GetNofManus(); ++i )
99 {
100 Int_t manuId = busPatch->GetManuId(i);
101 AddManuHisto(h,busPatch->GetDEId(),manuId);
102 }
103 }
104}
105//_____________________________________________________________________________
106void
107AliMUONTrackerDataHistogrammer::AddDEHisto(TH1& h, Int_t detElemId) const
108{
109 /// Add data from one detection element to the histogram
110
111 if ( fData.HasDetectionElement(detElemId) )
112 {
113 AliMpDetElement* de = AliMpDDLStore::Instance()->GetDetElement(detElemId);
114 for ( Int_t i = 0; i < de->GetNofBusPatches(); ++ i )
115 {
116 Int_t busPatchId = de->GetBusPatchId(i);
117 AddBusPatchHisto(h,busPatchId);
118 }
119 }
120}
121
122//_____________________________________________________________________________
123void
124AliMUONTrackerDataHistogrammer::AddManuHisto(TH1& h, Int_t detElemId, Int_t manuId) const
125{
126 /// Add data from a given manu to histogram
127
128 if ( fData.HasManu(detElemId,manuId) )
129 {
130 for ( Int_t i = 0; i < AliMpConstants::ManuNofChannels(); ++i )
131 {
132 if ( fData.HasChannel(detElemId,manuId,i) )
133 {
49419555 134 if ( IsInternalMode() )
135 {
136 h.Fill(fData.Channel(detElemId,manuId,i,fInternalDim));
137 }
138 else
139 {
140 AliMUONSparseHisto* sh = fData.GetChannelSparseHisto(detElemId,manuId,i);
10eb3d17 141
49419555 142 if ( sh )
143 {
144 Add(h,*sh);
145 }
10eb3d17 146 }
147 }
148 }
149 }
150}
151
10eb3d17 152//_____________________________________________________________________________
153TH1*
49419555 154AliMUONTrackerDataHistogrammer::CreateChannelHisto(Int_t detElemId,
155 Int_t manuId,
10eb3d17 156 Int_t manuChannel) const
157{
158 /// Create histogram of a given channel. Note that in order
159 /// to keep memory footprint as low as possible, you should delete
160 /// the returned pointer as soon as possible...
161
49419555 162 if ( fData.HasChannel(detElemId, manuId, manuChannel) && fData.IsHistogrammed(fExternalDim) )
10eb3d17 163 {
164 AliMUONSparseHisto* sh = fData.GetChannelSparseHisto(detElemId,manuId,manuChannel);
165
166 if ( sh )
167 {
49419555 168 Int_t nbins((1<<sh->Nbits()));
169 Double_t xmin,xmax;
170 fData.HistogramRange(xmin,xmax);
171
172 TH1* h = CreateHisto(Form("DE%04dMANU%04dCH%02d",detElemId,manuId,manuChannel),
173 nbins,xmin,xmax);
10eb3d17 174 if (h )
175 {
176 Add(*h,*sh);
177 }
178 return h;
179 }
180 }
181 return 0x0;
182}
183
184//_____________________________________________________________________________
185TH1*
49419555 186AliMUONTrackerDataHistogrammer::CreateHisto(const char* name,
187 Int_t nbins,
188 Double_t xmin,
189 Double_t xmax) const
10eb3d17 190{
49419555 191 /// Create a single histogram
10eb3d17 192
49419555 193 TH1* h(0);
10eb3d17 194
49419555 195 if ( xmin < xmax )
10eb3d17 196 {
49419555 197 h = new TH1F(name,name,nbins,xmin,xmax);
198 h->SetDirectory(gROOT);
10eb3d17 199 }
10eb3d17 200 return h;
201}
202
203//_____________________________________________________________________________
204TH1*
49419555 205AliMUONTrackerDataHistogrammer::CreateHisto(const AliMUONVPainter& painter,
206 Int_t externalDim,
207 Int_t internalDim)
10eb3d17 208{
209 /// Create an histogram, from given dim of given data,
210 /// for all the channels handled by painter
211
212 AliMUONPainterGroup* group = painter.Master()->PlotterGroup();
213
214 if ( !group ) return 0x0; // no data to histogram in this painter
215
216 AliMUONVTrackerData* data = group->Data();
10eb3d17 217
49419555 218 if ( externalDim >= data->ExternalDimension() )
219 {
220 AliErrorClass(Form("externalDim %d is out of bounds",externalDim));
221 return 0x0;
222 }
223
224 if ( internalDim >= data->NumberOfDimensions() )
225 {
226 AliErrorClass(Form("internalDim %d is out of bounds",internalDim));
227 return 0x0;
228 }
229
230 if ( internalDim < 0 && externalDim < 0 )
231 {
232 AliErrorClass("Both internal and external dim are < 0 !!!");
233 return 0x0;
234 }
235
236 AliMUONTrackerDataHistogrammer tdh(*data,externalDim,internalDim);
10eb3d17 237
238 TObjArray manuArray;
239
240 painter.FillManuList(manuArray);
49419555 241
10eb3d17 242 AliMpManuUID* mid;
243 TIter next(&manuArray);
49419555 244
245 TString basename(Form("%s-%s",painter.PathName().Data(),painter.Attributes().GetName()));
246 TString ext;
247 Int_t nbins((1<<12));
248 Double_t xmin(0.0);
249 Double_t xmax(0.0);
10eb3d17 250
49419555 251 if ( !tdh.IsInternalMode() )
252 {
253 data->HistogramRange(xmin,xmax);
254
255 xmin -= 0.5;
256 xmax -= 0.5;
257
258 ext = data->ExternalDimensionName(externalDim).Data();
259 }
260 else
261 {
262 tdh.GetDataRange(manuArray,xmin,xmax);
263 ext = data->DimensionName(internalDim).Data();
264 nbins = 100;
265 }
10eb3d17 266
49419555 267 TString name(Form("%s-%s",basename.Data(),ext.Data()));
268
269 TH1* histo = tdh.CreateHisto(name.Data(),nbins,xmin,xmax);
270
271 if ( histo )
10eb3d17 272 {
49419555 273 while ( ( mid = static_cast<AliMpManuUID*>(next()) ) )
10eb3d17 274 {
49419555 275 TH1* h = tdh.CreateManuHisto(mid->DetElemId(),mid->ManuId(),nbins,xmin,xmax);
276 if ( h )
277 {
278 histo->Add(h);
279 }
280 delete h;
10eb3d17 281 }
49419555 282 }
283 else
284 {
285 AliErrorClass(Form("Could not create histo for painter %s external dim %d internal dim %d",
286 painter.PathName().Data(),externalDim,internalDim));
10eb3d17 287 }
288
289 return histo;
290}
291
292//_____________________________________________________________________________
293TH1*
49419555 294AliMUONTrackerDataHistogrammer::CreateManuHisto(Int_t detElemId, Int_t manuId,
295 Int_t nbins,
296 Double_t xmin,
297 Double_t xmax) const
10eb3d17 298{
49419555 299 /// Create histogram of a given manu. Note that in order
300 /// to keep memory footprint as low as possible, you should delete
301 /// the returned pointer as soon as possible...
10eb3d17 302
303 TH1* h(0x0);
304
49419555 305 if ( !fData.HasManu(detElemId,manuId) ) return 0x0;
306
307 if ( ( fExternalDim >= 0 && fData.IsHistogrammed(fExternalDim) ) ||
308 ( fInternalDim >= 0 && fInternalDim < fData.NumberOfDimensions() ) )
10eb3d17 309 {
49419555 310 h = CreateHisto(Form("DE%04dMANU%04d",detElemId,manuId),
311 nbins,xmin,xmax);
312 if ( h ) AddManuHisto(*h,detElemId,manuId);
10eb3d17 313 }
49419555 314
10eb3d17 315 return h;
316}
317
10eb3d17 318//_____________________________________________________________________________
49419555 319void
320AliMUONTrackerDataHistogrammer::GetDataRange(const TObjArray& manuArray,
321 Double_t& xmin, Double_t& xmax) const
10eb3d17 322{
49419555 323 /// Get data range (in case of InternalMode() only) spanned by the manus in
324 /// manuArray
10eb3d17 325
49419555 326 xmin = FLT_MAX;
327 xmax = -FLT_MAX;
10eb3d17 328
49419555 329 if (!IsInternalMode())
10eb3d17 330 {
49419555 331 AliError("Cannot use this method for external mode !");
10eb3d17 332 }
49419555 333
334 AliMpManuUID* mid;
335 TIter next(&manuArray);
10eb3d17 336
49419555 337 while ( ( mid = static_cast<AliMpManuUID*>(next()) ) )
338 {
339 Int_t detElemId = mid->DetElemId();
340 Int_t manuId = mid->ManuId();
341
342 for ( Int_t i = 0; i < AliMpConstants::ManuNofChannels(); ++i )
343 {
344 if ( fData.HasChannel(detElemId,manuId,i) )
345 {
346 Double_t value = fData.Channel(detElemId,manuId,i,fInternalDim);
347 xmin = TMath::Min(xmin,value);
348 xmax = TMath::Max(xmax,value);
349 }
350 }
351 }
352
10eb3d17 353}
354