]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDQADataMakerRec.cxx
Fix for bug #78633 (Diego)
[u/mrichter/AliRoot.git] / FMD / AliFMDQADataMakerRec.cxx
CommitLineData
c9dd1c4d 1/**************************************************************************
2 * Copyright(c) 2004, 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// --- ROOT system ---
16#include <iostream>
17#include <TClonesArray.h>
18#include <TFile.h>
19#include <TH1F.h>
20#include <TH1I.h>
21
22// --- AliRoot header files ---
23#include "AliESDEvent.h"
24#include "AliLog.h"
25#include "AliFMDQADataMakerRec.h"
26#include "AliFMDDigit.h"
27#include "AliFMDRecPoint.h"
28#include "AliQAChecker.h"
29#include "AliESDFMD.h"
30#include "AliFMDParameters.h"
86c6cec8 31#include "AliFMDRawReader.h"
3ceaa9ad 32#include "AliRawReader.h"
b995fc28 33#include "AliFMDAltroMapping.h"
f4b5062e 34#include "AliFMDDebug.h"
c9dd1c4d 35
36//_____________________________________________________________________
37// This is the class that collects the QA data for the FMD during
38// reconstruction.
39//
40// The following data types are picked up:
c9dd1c4d 41// - rec points
42// - esd data
c9dd1c4d 43// - raws
44// Author : Hans Hjersing Dalsgaard, hans.dalsgaard@cern.ch
45//_____________________________________________________________________
46
47ClassImp(AliFMDQADataMakerRec)
48#if 0
49; // For Emacs - do not delete!
50#endif
51
52//_____________________________________________________________________
53AliFMDQADataMakerRec::AliFMDQADataMakerRec() :
4e25ac79 54 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kFMD),
a7e41e8d 55 "FMD Quality Assurance Data Maker"),
198942c7 56 fRecPointsArray("AliFMDRecPoint", 1000)
c9dd1c4d 57{
58 // ctor
56236ce9 59
c9dd1c4d 60}
61
62//_____________________________________________________________________
a7e41e8d 63AliFMDQADataMakerRec::AliFMDQADataMakerRec(const AliFMDQADataMakerRec& qadm)
4e25ac79 64 : AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kFMD),
a7e41e8d 65 "FMD Quality Assurance Data Maker"),
a7e41e8d 66 fRecPointsArray(qadm.fRecPointsArray)
c9dd1c4d 67{
ffa78f64 68 // copy ctor
c9dd1c4d 69 // Parameters:
70 // qadm Object to copy from
71
72}
9bd2ccc2 73//_____________________________________________________________________
63720979 74AliFMDQADataMakerRec&
75AliFMDQADataMakerRec::operator = (const AliFMDQADataMakerRec& qadm )
a7e41e8d 76{
a7e41e8d 77 fRecPointsArray = qadm.fRecPointsArray;
78
79 return *this;
80}
81//_____________________________________________________________________
9bd2ccc2 82AliFMDQADataMakerRec::~AliFMDQADataMakerRec()
83{
56236ce9 84
9bd2ccc2 85}
c9dd1c4d 86
87
88//_____________________________________________________________________
89
90void
4e25ac79 91AliFMDQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task,
57acd2d2 92 TObjArray ** list)
c9dd1c4d 93{
94 // Detector specific actions at end of cycle
95 // do the QA checking
ffa78f64 96 AliLog::Message(5,"FMD: end of detector cycle",
97 "AliFMDQADataMakerRec","AliFMDQADataMakerRec",
98 "AliFMDQADataMakerRec::EndOfDetectorCycle",
99 "AliFMDQADataMakerRec.cxx",95);
4e25ac79 100 AliQAChecker::Instance()->Run(AliQAv1::kFMD, task, list);
c9dd1c4d 101}
102
103//_____________________________________________________________________
104void AliFMDQADataMakerRec::InitESDs()
105{
106 // create Digits histograms in Digits subdir
7d297381 107 const Bool_t expert = kTRUE ;
108 const Bool_t image = kTRUE ;
109
c9dd1c4d 110 TH1F* hEnergyOfESD = new TH1F("hEnergyOfESD","Energy distribution",100,0,3);
111 hEnergyOfESD->SetXTitle("Edep/Emip");
112 hEnergyOfESD->SetYTitle("Counts");
7d297381 113 Add2ESDsList(hEnergyOfESD, 0, !expert, image);
c9dd1c4d 114
115}
116
44ed7a66 117//_____________________________________________________________________
118void AliFMDQADataMakerRec::InitDigits()
119{
120 // create Digits histograms in Digits subdir
121 const Bool_t expert = kTRUE ;
122 const Bool_t image = kTRUE ;
123
63720979 124 TH1I* hADCCounts = new TH1I("hADCCounts",
125 "Dist of ADC counts;ADC counts;Counts",
126 1024,0,1024);
44ed7a66 127 Add2DigitsList(hADCCounts, 0, !expert, image);
128}
129
130
c9dd1c4d 131//_____________________________________________________________________
132void AliFMDQADataMakerRec::InitRecPoints()
133{
7d297381 134 // create Reconstructed Points histograms in RecPoints subdir
135 const Bool_t expert = kTRUE ;
136 const Bool_t image = kTRUE ;
137
c9dd1c4d 138 TH1F* hEnergyOfRecpoints = new TH1F("hEnergyOfRecpoints",
139 "Energy Distribution",100,0,3);
140 hEnergyOfRecpoints->SetXTitle("Edep/Emip");
db72ff3b 141 hEnergyOfRecpoints->SetYTitle("Counts");
7d297381 142 Add2RecPointsList(hEnergyOfRecpoints,0, !expert, image);
c9dd1c4d 143}
144
145//_____________________________________________________________________
146void AliFMDQADataMakerRec::InitRaws()
147{
7d297381 148 // create Raws histograms in Raws subdir
149 const Bool_t expert = kTRUE ;
150 const Bool_t saveCorr = kTRUE ;
151 const Bool_t image = kTRUE ;
17a67721 152 Int_t colors[3] = {kRed,kGreen,kBlue};
b4da452d 153 TH1I* hADCCounts;
154 for(Int_t det = 1; det<=3; det++) {
155 Int_t firstring = (det==1 ? 1 : 0);
156 for(Int_t iring = firstring;iring<=1;iring++) {
c18fe410 157 Char_t ring = (iring == 1 ? 'I' : 'O');
63720979 158 hADCCounts =
159 new TH1I(Form("hADCCounts_FMD%d%c", det, ring),
160 Form("FMD%d%c ADC counts;Amplitude [ADC counts];Counts",
161 det,ring),1024,0,1023);
17a67721 162 hADCCounts->SetFillStyle(3001);
63720979 163 hADCCounts->SetFillColor(colors[det-1]+2+iring);
164 hADCCounts->SetDrawOption("LOGY");
165
c18fe410 166 Int_t index1 = GetHalfringIndex(det, ring, 0,1);
6ceca4ef 167 Add2RawsList(hADCCounts, index1, !expert, image, !saveCorr);
c18fe410 168
b4da452d 169 for(Int_t b = 0; b<=1;b++) {
170
c18fe410 171 //Hexadecimal board numbers 0x0, 0x1, 0x10, 0x11;
b4da452d 172 UInt_t board = (iring == 1 ? 0 : 1);
173 board = board + b*16;
c18fe410 174
b4da452d 175
63720979 176 hADCCounts =
177 new TH1I(Form("hADCCounts_FMD%d%c_0x%02x", det, ring, board),
178 Form("FMD%d%c_0x%x ADC counts;Amplitude [ADC counts];Counts",
179 det,ring, board),1024,0,1023);
b4da452d 180 hADCCounts->SetXTitle("ADC counts");
181 hADCCounts->SetYTitle("");
63720979 182 hADCCounts->SetFillStyle(3001);
183 hADCCounts->SetFillColor(colors[det-1]+2+iring);
184 hADCCounts->SetDrawOption("LOGY");
c18fe410 185 Int_t index2 = GetHalfringIndex(det, ring, board/16,0);
6ceca4ef 186 Add2RawsList(hADCCounts, index2, expert, !image, !saveCorr);
c18fe410 187
b4da452d 188 }
189 }
190 }
c9dd1c4d 191}
192
408bf2b4 193#if 0
194struct FillESDHist : public AliESDFMD::ForOne
195{
196 FillESDHist(AliFMDQADataMakerRec* m) : fM(m) {}
197 FillESDHist(const FillESDHist& o) : fM(o.fM) {}
198 FillESDHist& operator=(const FillESDHist& o) { fM = o.fM; return *this; }
199 Bool_t operator()(UShort_t, Char_t, UShort_t, UShort_t, Float_t m, Float_t)
200 {
201 // Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
202 if(m == AliESDFMD::kInvalidMult) return true;
203
204 fM->GetESDsData(0)->Fill(m);
205 return true;
206 }
207 AliFMDQADataMakerRec* fM;
208};
209#endif
210
c9dd1c4d 211//_____________________________________________________________________
212void AliFMDQADataMakerRec::MakeESDs(AliESDEvent * esd)
213{
214 if(!esd) {
215 AliError("FMD ESD object not found!!") ;
216 return;
217 }
f4b5062e 218 AliFMDDebug(2, ("Will loop over ESD data and fill histogram"));
219
c9dd1c4d 220 AliESDFMD* fmd = esd->GetFMDData();
221 if (!fmd) return;
f4b5062e 222
408bf2b4 223#if 0
224 FillESDHist f(this);
225 fmd->ForEach(f);
226#else
227
228
229
f4b5062e 230 // FIXME - we should use AliESDFMD::ForOne subclass to do this!
c9dd1c4d 231 for(UShort_t det=1;det<=3;det++) {
f4b5062e 232 UShort_t nrng = (det == 1 ? 1 : 2);
233 for (UShort_t ir = 0; ir < nrng; ir++) {
c9dd1c4d 234 Char_t ring = (ir == 0 ? 'I' : 'O');
235 UShort_t nsec = (ir == 0 ? 20 : 40);
236 UShort_t nstr = (ir == 0 ? 512 : 256);
237 for(UShort_t sec =0; sec < nsec; sec++) {
238 for(UShort_t strip = 0; strip < nstr; strip++) {
239 Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
240 if(mult == AliESDFMD::kInvalidMult) continue;
241
242 GetESDsData(0)->Fill(mult);
243 }
244 }
245 }
246 }
408bf2b4 247#endif
c9dd1c4d 248}
249
44ed7a66 250
c9dd1c4d 251//_____________________________________________________________________
6252ceeb 252void AliFMDQADataMakerRec::MakeDigits()
c9dd1c4d 253{
254 // makes data from Digits
6252ceeb 255 if(!fDigitsArray) {
ffa78f64 256 AliError("FMD Digit object not found!!") ;
257 return;
c9dd1c4d 258 }
eca4fa66 259
6252ceeb 260 for(Int_t i=0;i<fDigitsArray->GetEntriesFast();i++) {
c9dd1c4d 261 //Raw ADC counts
6252ceeb 262 AliFMDDigit* digit = static_cast<AliFMDDigit*>(fDigitsArray->At(i));
c9dd1c4d 263 GetDigitsData(0)->Fill(digit->Counts());
264 }
265}
266
267//_____________________________________________________________________
268void AliFMDQADataMakerRec::MakeDigits(TTree * digitTree)
269{
270
6252ceeb 271 if (fDigitsArray)
272 fDigitsArray->Clear();
273 else
274 fDigitsArray = new TClonesArray("AliFMDDigit", 1000);
275
ffa78f64 276 TBranch* branch = digitTree->GetBranch("FMD");
277 if (!branch) {
c9dd1c4d 278 AliWarning("FMD branch in Digit Tree not found") ;
279 return;
280 }
6252ceeb 281 branch->SetAddress(&fDigitsArray);
ffa78f64 282 branch->GetEntry(0);
6252ceeb 283 MakeDigits();
c9dd1c4d 284}
44ed7a66 285
c9dd1c4d 286//_____________________________________________________________________
86c6cec8 287void AliFMDQADataMakerRec::MakeRaws(AliRawReader* rawReader)
c9dd1c4d 288{
198942c7 289
6252ceeb 290 AliFMDRawReader fmdReader(rawReader,0);
291
292 if (fDigitsArray)
293 fDigitsArray->Clear();
294 else
295 fDigitsArray = new TClonesArray("AliFMDDigit", 1000);
eca4fa66 296
6252ceeb 297 TClonesArray* digitsAddress = fDigitsArray;
3ceaa9ad 298
299 rawReader->Reset();
78328afd 300
c18fe410 301 digitsAddress->Clear();
302 fmdReader.ReadAdcs(digitsAddress);
303 for(Int_t i=0;i<digitsAddress->GetEntriesFast();i++) {
304 //Raw ADC counts
b995fc28 305 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digitsAddress->At(i));
306 UShort_t det = digit->Detector();
307 Char_t ring = digit->Ring();
308 UShort_t sec = digit->Sector();
309 // UShort_t strip = digit->Strip();
310 AliFMDParameters* pars = AliFMDParameters::Instance();
311 Short_t board = pars->GetAltroMap()->Sector2Board(ring, sec);
c18fe410 312
313 Int_t index1 = GetHalfringIndex(det, ring, 0, 1);
314 GetRawsData(index1)->Fill(digit->Counts());
315 Int_t index2 = GetHalfringIndex(det, ring, board/16,0);
316 GetRawsData(index2)->Fill(digit->Counts());
317
318 }
c9dd1c4d 319}
320
321//_____________________________________________________________________
322void AliFMDQADataMakerRec::MakeRecPoints(TTree* clustersTree)
323{
324 // makes data from RecPoints
eca4fa66 325
6252ceeb 326 AliFMDParameters* pars = AliFMDParameters::Instance();
56236ce9 327 fRecPointsArray.Clear();
c9dd1c4d 328 TBranch *fmdbranch = clustersTree->GetBranch("FMD");
329 if (!fmdbranch) {
330 AliError("can't get the branch with the FMD recpoints !");
331 return;
332 }
333
56236ce9 334 TClonesArray* RecPointsAddress = &fRecPointsArray;
9bd2ccc2 335
56236ce9 336 fmdbranch->SetAddress(&RecPointsAddress);
c9dd1c4d 337 fmdbranch->GetEntry(0);
56236ce9 338 TIter next(RecPointsAddress) ;
c9dd1c4d 339 AliFMDRecPoint * rp ;
340 while ((rp = static_cast<AliFMDRecPoint*>(next()))) {
198942c7 341
342 GetRecPointsData(0)->Fill(rp->Edep()/pars->GetEdepMip()) ;
343
c9dd1c4d 344 }
c9dd1c4d 345
346}
347
348//_____________________________________________________________________
349void AliFMDQADataMakerRec::StartOfDetectorCycle()
350{
351 // What
352 // to
353 // do?
354}
355//_____________________________________________________________________
b4da452d 356Int_t AliFMDQADataMakerRec::GetHalfringIndex(UShort_t det,
357 Char_t ring,
c18fe410 358 UShort_t board,
359 UShort_t monitor) {
b4da452d 360
361 UShort_t iring = (ring == 'I' ? 1 : 0);
362
c18fe410 363 Int_t index = ( ((det-1) << 3) | (iring << 2) | (board << 1) | (monitor << 0));
b4da452d 364
365 return index-2;
366
367}
368
369//_____________________________________________________________________
370
371
c9dd1c4d 372//
373// EOF
374//