]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFQADataMaker.cxx
Update timestamp for new data points simulation
[u/mrichter/AliRoot.git] / TOF / AliTOFQADataMaker.cxx
CommitLineData
380c04bf 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 **************************************************************************/
8fd6fd6c 15
16///////////////////////////////////////////////////////////////////////
17// //
18// Produces the data needed to calculate the TOF quality assurance. //
19// QA objects are 1 & 2 Dimensional histograms. //
20// author: S.Arcelli //
21// //
22///////////////////////////////////////////////////////////////////////
380c04bf 23
24#include <TClonesArray.h>
25#include <TFile.h>
26#include <TH1I.h>
27#include <TH1F.h>
28#include <TH2F.h>
29#include "AliESDEvent.h"
30#include "AliESDtrack.h"
31#include "AliTOFcluster.h"
32#include "AliTOFdigit.h"
33#include "AliTOFSDigit.h"
34#include "AliTOFhitT0.h"
35#include "AliTOFQADataMaker.h"
693a790d 36#include "AliQAChecker.h"
380c04bf 37#include "AliRawReader.h"
38#include "AliTOFRawStream.h"
39#include "AliTOFrawData.h"
40#include "AliLog.h"
41
42ClassImp(AliTOFQADataMaker)
43
44//____________________________________________________________________________
45 AliTOFQADataMaker::AliTOFQADataMaker() :
97430dfe 46 AliQADataMaker(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker"),
47 fTOFRawStream(AliTOFRawStream())
380c04bf 48{
49 //
50 // ctor
51 //
52}
53
54//____________________________________________________________________________
55AliTOFQADataMaker::AliTOFQADataMaker(const AliTOFQADataMaker& qadm) :
97430dfe 56 AliQADataMaker(),
57 fTOFRawStream(AliTOFRawStream(qadm.fTOFRawStream))
380c04bf 58{
59 //
60 //copy ctor
61 //
62 SetName((const char*)qadm.GetName()) ;
63 SetTitle((const char*)qadm.GetTitle());
64}
65
66//__________________________________________________________________
67AliTOFQADataMaker& AliTOFQADataMaker::operator = (const AliTOFQADataMaker& qadm )
68{
69 //
70 // assignment operator.
71 //
72 this->~AliTOFQADataMaker();
73 new(this) AliTOFQADataMaker(qadm);
74 return *this;
75}
76
77//____________________________________________________________________________
78void AliTOFQADataMaker::InitHits()
79{
80 //
81 // create Hits histograms in Hits subdir
82 //
135306ea 83
84 Bool_t expert = kFALSE;
85
5561a41b 86 TH1F * h0 = new TH1F("hTOFHits", "Number of TOF Hits ",301, -1.02, 5.) ;
380c04bf 87 h0->Sumw2() ;
5561a41b 88 h0->GetXaxis()->SetTitle("TOF hit number [10 power]");
135306ea 89 Add2HitsList(h0, 0, expert) ;
380c04bf 90
5561a41b 91 TH1F * h1 = new TH1F("hTOFHitsTime", "Hits Time Spectrum in TOF (ns)", 2000, 0., 200) ;
380c04bf 92 h1->Sumw2() ;
5561a41b 93 h1->GetXaxis()->SetTitle("Simulated TOF time [ns]");
135306ea 94 Add2HitsList(h1, 1, expert) ;
380c04bf 95
5561a41b 96 TH1F * h2 = new TH1F("hTOFHitsLength", "Length Spectrum in TOF (cm)", 500, 0., 500) ;
380c04bf 97 h2->Sumw2() ;
5561a41b 98 h2->GetXaxis()->SetTitle("Track length from primary vertex till hit TOF pad [cm]");
135306ea 99 Add2HitsList(h2, 2, expert) ;
380c04bf 100
5561a41b 101 TH2F * h3 = new TH2F("hTOFHitsClusMap","Hits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
380c04bf 102 h3->Sumw2() ;
5561a41b 103 h3->GetXaxis()->SetTitle("2*strip+padz (eta)");
104 h3->GetYaxis()->SetTitle("48*sector+padx (phi)");
105 h3->GetYaxis()->SetTitleOffset(1.15);
135306ea 106 Add2HitsList(h3, 3, expert) ;
92664bc8 107 //
108 ClonePerTrigClass(AliQAv1::kHITS); // this should be the last line
380c04bf 109}
110
111//____________________________________________________________________________
112void AliTOFQADataMaker::InitDigits()
113{
114 //
115 // create Digits histograms in Digits subdir
116 //
135306ea 117
118 Bool_t expert = kFALSE;
119
5561a41b 120 TH1F * h0 = new TH1F("hTOFDigits", "Number of TOF Digits ",301, -1.02, 5.) ;
121 h0->Sumw2() ;
122 h0->GetXaxis()->SetTitle("TOF digit number [10 power]");
135306ea 123 Add2DigitsList(h0, 0, expert) ;
380c04bf 124
5561a41b 125 TH1F * h1 = new TH1F("hTOFDigitsTime", "Digits Time Spectrum in TOF (ns)", 2000, 0., 200) ;
380c04bf 126 h1->Sumw2() ;
5561a41b 127 h1->GetXaxis()->SetTitle("Digitized TOF time [ns]");
135306ea 128 Add2DigitsList(h1, 1, expert) ;
380c04bf 129
5561a41b 130 TH1F * h2 = new TH1F("hTOFDigitsToT", "Digits ToT Spectrum in TOF (ns)", 500, 0., 50) ;
380c04bf 131 h2->Sumw2() ;
5561a41b 132 h2->GetXaxis()->SetTitle("Digitized TOT [ns]");
133 h2->GetYaxis()->SetTitle("Digitized TOF time [ns]");
135306ea 134 Add2DigitsList(h2, 2, expert) ;
380c04bf 135
5561a41b 136 TH2F * h3 = new TH2F("hTOFDigitsClusMap","Digits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
380c04bf 137 h3->Sumw2() ;
5561a41b 138 h3->GetXaxis()->SetTitle("2*strip+padz (eta)");
139 h3->GetYaxis()->SetTitle("48*sector+padx (phi)");
140 h3->GetYaxis()->SetTitleOffset(1.15);
135306ea 141 Add2DigitsList(h3, 3, expert) ;
92664bc8 142 //
143 ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
380c04bf 144}
145
146//____________________________________________________________________________
147void AliTOFQADataMaker::InitSDigits()
148{
149 //
150 // create SDigits histograms in SDigits subdir
151 //
135306ea 152
153 Bool_t expert = kFALSE;
154
5561a41b 155 TH1F * h0 = new TH1F("hTOFSDigits", "Number of TOF SDigits ",301, -1.02, 5.) ;
156 h0->Sumw2() ;
157 h0->GetXaxis()->SetTitle("TOF sdigit number [10 power]");
135306ea 158 Add2SDigitsList(h0, 0, expert) ;
380c04bf 159
5561a41b 160 TH1F * h1 = new TH1F("hTOFSDigitsTime", "SDigits Time Spectrum in TOF (ns)", 2000, 0., 200) ;
380c04bf 161 h1->Sumw2() ;
5561a41b 162 h1->GetXaxis()->SetTitle("SDigitized TOF time [ns]");
135306ea 163 Add2SDigitsList(h1, 1, expert) ;
380c04bf 164
5561a41b 165 TH2F * h2 = new TH2F("hTOFSDigitsClusMap","SDigits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
380c04bf 166 h2->Sumw2() ;
5561a41b 167 h2->GetXaxis()->SetTitle("2*strip+padz (eta)");
168 h2->GetYaxis()->SetTitle("48*sector+padx (phi)");
169 h2->GetYaxis()->SetTitleOffset(1.15);
135306ea 170 Add2SDigitsList(h2, 2, expert) ;
92664bc8 171 //
172 ClonePerTrigClass(AliQAv1::kSDIGITS); // this should be the last line
380c04bf 173}
174
175//____________________________________________________________________________
176void AliTOFQADataMaker::InitRaws()
177{
178 //
179 // create Raws histograms in Raws subdir
180 //
135306ea 181
182 Bool_t expert = kFALSE;
183
5561a41b 184 TH1F * h0 = new TH1F("hTOFRaws", "Number of TOF Raws ",301, -1.02, 5.) ;
185 h0->Sumw2() ;
186 h0->GetXaxis()->SetTitle("TOF raw number [10 power]");
135306ea 187 Add2RawsList(h0, 0, expert) ;
380c04bf 188
0a2a478c 189 TH1F * h1 = new TH1F("hTOFRawsTime", "Raws Time Spectrum in TOF (ns)", 20000, 0., 2000) ;
380c04bf 190 h1->Sumw2() ;
5561a41b 191 h1->GetXaxis()->SetTitle("Measured TOF time [ns]");
135306ea 192 Add2RawsList(h1, 1, expert) ;
380c04bf 193
5561a41b 194 TH1F * h2 = new TH1F("hTOFRawsToT", "Raws ToT Spectrum in TOF (ns)", 500, 0., 50) ;
380c04bf 195 h2->Sumw2() ;
5561a41b 196 h2->GetXaxis()->SetTitle("Measured TOT [ns]");
197 h2->GetYaxis()->SetTitle("Measured TOF time [ns]");
135306ea 198 Add2RawsList(h2, 2, expert) ;
380c04bf 199
5561a41b 200 TH2F * h3 = new TH2F("hTOFRawsClusMap","Raws vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
380c04bf 201 h3->Sumw2() ;
5561a41b 202 h3->GetXaxis()->SetTitle("2*strip+padz (eta)");
203 h3->GetYaxis()->SetTitle("48*sector+padx (phi)");
204 h3->GetYaxis()->SetTitleOffset(1.15);
135306ea 205 Add2RawsList(h3, 3, expert) ;
92664bc8 206 //
207 ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
380c04bf 208}
209
210//____________________________________________________________________________
211void AliTOFQADataMaker::InitRecPoints()
212{
213 //
214 // create RecPoints histograms in RecPoints subdir
215 //
135306ea 216
217 Bool_t expert = kFALSE;
218
5561a41b 219 TH1F * h0 = new TH1F("hTOFRecPoints", "Number of TOF RecPoints ",301, -1.02, 5.) ;
220 h0->Sumw2() ;
221 h0->GetXaxis()->SetTitle("TOF recPoint number [10 power]");
135306ea 222 Add2RecPointsList(h0, 0, expert) ;
380c04bf 223
0a2a478c 224 TH1F * h1 = new TH1F("hTOFRecPointsTime", "RecPoints Time Spectrum in TOF (ns)", 20000, 0., 2000) ;
380c04bf 225 h1->Sumw2() ;
5561a41b 226 h1->GetXaxis()->SetTitle("Calibrated TOF time [ns]");
135306ea 227 Add2RecPointsList(h1, 1, expert) ;
380c04bf 228
0a2a478c 229 TH1F * h2 = new TH1F("hTOFRecPointsRawTime", "RecPoints raw Time Spectrum in TOF (ns)", 20000, 0., 2000) ;
380c04bf 230 h2->Sumw2() ;
5561a41b 231 h2->GetXaxis()->SetTitle("Measured TOT [ns]");
232 h2->GetYaxis()->SetTitle("Measured TOF time [ns]");
135306ea 233 Add2RecPointsList(h2, 2, expert) ;
380c04bf 234
5561a41b 235 TH1F * h3 = new TH1F("hTOFRecPointsToT", "RecPoints ToT Spectrum in TOF (ns)", 500, 0., 50) ;
380c04bf 236 h3->Sumw2() ;
5561a41b 237 h3->GetXaxis()->SetTitle("Measured TOT [ns]");
135306ea 238 Add2RecPointsList(h3, 3, expert) ;
380c04bf 239
5561a41b 240 TH2F * h4 = new TH2F("hTOFRecPointsClusMap","RecPoints vs TOF phi-eta",183, -0.5, 182.5,865,-0.5,864.5) ;
380c04bf 241 h4->Sumw2() ;
0ea3c44e 242 h4->GetXaxis()->SetTitle("2*strip+padz (eta)");
243 h4->GetYaxis()->SetTitle("48*sector+padx (phi)");
5561a41b 244 h4->GetYaxis()->SetTitleOffset(1.15);
135306ea 245 Add2RecPointsList(h4, 4, expert) ;
92664bc8 246 //
247 ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
380c04bf 248}
249
250//____________________________________________________________________________
251void AliTOFQADataMaker::InitESDs()
252{
253 //
254 //create ESDs histograms in ESDs subdir
255 //
135306ea 256
257 Bool_t expert = kFALSE;
258
5561a41b 259 TH1F * h0 = new TH1F("hTOFESDs", "Number of matched TOF tracks over ESDs", 250, -1., 4.) ;
380c04bf 260 h0->Sumw2() ;
5561a41b 261 h0->GetXaxis()->SetTitle("Number of TOF matched ESD tracks [10 power]");
135306ea 262 Add2ESDsList(h0, 0, expert) ;
380c04bf 263
0a2a478c 264 TH1F * h1 = new TH1F("hTOFESDsTime", "Time Spectrum in TOF (ns)", 20000, 0., 2000) ;
380c04bf 265 h1->Sumw2() ;
5561a41b 266 h1->GetXaxis()->SetTitle("Calibrated TOF time [ns]");
135306ea 267 Add2ESDsList(h1, 1, expert) ;
380c04bf 268
0a2a478c 269 TH1F * h2 = new TH1F("hTOFESDsRawTime", "raw Time Spectrum in TOF (ns)", 20000, 0., 2000) ;
380c04bf 270 h2->Sumw2() ;
5561a41b 271 h2->GetXaxis()->SetTitle("Measured TOF time [ns]");
135306ea 272 Add2ESDsList(h2, 2, expert) ;
380c04bf 273
5561a41b 274 TH1F * h3 = new TH1F("hTOFESDsToT", "ToT Spectrum in TOF (ns)", 500, 0., 50) ;
380c04bf 275 h3->Sumw2() ;
5561a41b 276 h3->GetXaxis()->SetTitle("Measured TOT [ns]");
135306ea 277 Add2ESDsList(h3, 3, expert) ;
380c04bf 278
5561a41b 279 TH1F * h4 = new TH1F("hTOFESDsPID", "Fraction of matched TOF tracks with good PID flag", 101, 0., 101.) ;
280 h4->Sumw2() ;
281 h4->GetXaxis()->SetTitle("Fraction of TOF matched ESD tracks with good flag [%]");
135306ea 282 Add2ESDsList(h4, 4, expert) ;
92664bc8 283 //
284 ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line
380c04bf 285}
286
287
288//____________________________________________________________________________
289void AliTOFQADataMaker::MakeHits(TClonesArray * hits)
290{
291 //
292 //make QA data from Hits
293 //
294
295 Int_t in[5];
296 Int_t out[5];
297
298 Int_t nentries=hits->GetEntriesFast();
c332d0d1 299 if(nentries<=0) {
92664bc8 300 FillHitsData(0,-1.) ;
ed77256d 301 } else{
92664bc8 302 FillHitsData(0,TMath::Log10(nentries)) ;
ed77256d 303 }
380c04bf 304 TIter next(hits) ;
305 AliTOFhitT0 * hit ;
306 while ( (hit = dynamic_cast<AliTOFhitT0 *>(next())) ) {
307
92664bc8 308 FillHitsData(1, hit->GetTof()*1.E9) ;//in ns
309 FillHitsData(2, hit->GetLen()) ;//in cm
380c04bf 310
311 in[0] = hit->GetSector();
312 in[1] = hit->GetPlate();
313 in[2]= hit->GetStrip();
314 in[3]= hit->GetPadx();
315 in[4]= hit->GetPadz();
316 GetMapIndeces(in,out);
92664bc8 317 FillHitsData(3, out[0],out[1]) ;//hit map
380c04bf 318 }
319
320}
321
322
323//____________________________________________________________________________
324void AliTOFQADataMaker::MakeHits(TTree * hitTree)
325{
326 //
327 // make QA data from Hit Tree
328 //
329 if(!hitTree){
330 AliError("can't get the tree with TOF hits !");
331 return;
332 }
333
334 TBranch * branch = hitTree->GetBranch("TOF") ;
335
336 if (!branch ) {
337 AliError("TOF branch in Hit Tree not found") ;
338 return;
339 }
340
341 TClonesArray * hits = new TClonesArray("AliTOFhitT0", 1000);
342 TClonesArray * dummy = new TClonesArray("AliTOFhitT0", 1000);
343 branch->SetAddress(&dummy);
344 Int_t index = 0 ;
345 for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
346 branch->GetEntry(ientry) ;
347 for (Int_t ihit = 0 ; ihit < dummy->GetEntries() ; ihit++) {
348 AliTOFhitT0 * hit = dynamic_cast<AliTOFhitT0 *> (dummy->At(ihit)) ;
349 new((*hits)[index]) AliTOFhitT0(*hit) ;
350 index++ ;
351 }
352 }
353
354 dummy->Delete();
355 delete dummy;
356 MakeHits(hits) ;
92664bc8 357 //
358 IncEvCountCycleHits();
359 IncEvCountTotalHits();
360 //
380c04bf 361}
362
363//____________________________________________________________________________
364void AliTOFQADataMaker::MakeDigits(TClonesArray * digits)
365{
366 //
367 // makes data from Digits
368 //
369 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
370 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
371 Int_t in[5];
372 Int_t out[5];
373
374 Int_t nentries=digits->GetEntriesFast();
c332d0d1 375 if(nentries<=0){
92664bc8 376 FillDigitsData(0,-1.) ;
ed77256d 377 }else{
92664bc8 378 FillDigitsData(0, TMath::Log10(nentries)) ;
ed77256d 379 }
380
380c04bf 381 TIter next(digits) ;
382 AliTOFdigit * digit ;
383 while ( (digit = dynamic_cast<AliTOFdigit *>(next())) ) {
384
92664bc8 385 FillDigitsData(1, digit->GetTdc()*tdc2ns) ;//in ns
386 FillDigitsData(2, digit->GetToT()*tot2ns) ;//in ns
380c04bf 387
388 in[0] = digit->GetSector();
389 in[1] = digit->GetPlate();
390 in[2] = digit->GetStrip();
391 in[3] = digit->GetPadx();
392 in[4]= digit->GetPadz();
393 GetMapIndeces(in,out);
92664bc8 394 FillDigitsData(3, out[0],out[1]) ;//digit map
380c04bf 395 }
396
397}
398
399
400//____________________________________________________________________________
401void AliTOFQADataMaker::MakeDigits(TTree * digitTree)
402{
403 //
404 // makes data from Digit Tree
405 //
406 TClonesArray * digits = new TClonesArray("AliTOFdigit", 1000) ;
407
408 TBranch * branch = digitTree->GetBranch("TOF") ;
409 if ( ! branch ) {
410 AliError("TOF branch in Digit Tree not found") ;
411 return;
412 }
413 branch->SetAddress(&digits) ;
414 branch->GetEntry(0) ;
415 MakeDigits(digits) ;
92664bc8 416 //
417 IncEvCountCycleDigits();
418 IncEvCountTotalDigits();
419 //
380c04bf 420}
421
422//____________________________________________________________________________
423void AliTOFQADataMaker::MakeSDigits(TClonesArray * sdigits)
424{
425 //
426 // makes data from SDigits
427 //
428
429 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
430 Int_t in[5];
431 Int_t out[5];
432
433 Int_t nentries=sdigits->GetEntriesFast();
c332d0d1 434 if(nentries<=0){
92664bc8 435 FillSDigitsData(0,-1.) ;
ed77256d 436 }else{
92664bc8 437 FillSDigitsData(0,TMath::Log10(nentries)) ;
ed77256d 438 }
380c04bf 439
440 TIter next(sdigits) ;
441 AliTOFSDigit * sdigit ;
442 while ( (sdigit = dynamic_cast<AliTOFSDigit *>(next())) ) {
443
444 for(Int_t i=0;i<sdigit->GetNDigits();i++){
92664bc8 445 FillSDigitsData(1, sdigit->GetTdc(i)*tdc2ns) ;//in ns
380c04bf 446 }
447
448 in[0] = sdigit->GetSector();
449 in[1] = sdigit->GetPlate();
450 in[2] = sdigit->GetStrip();
451 in[3] = sdigit->GetPadx();
452 in[4]= sdigit->GetPadz();
453 GetMapIndeces(in,out);
92664bc8 454 FillSDigitsData(2, out[0],out[1]) ;//sdigit map
380c04bf 455 }
456}
457
458//____________________________________________________________________________
459void AliTOFQADataMaker::MakeSDigits(TTree * sdigitTree)
460{
461 //
462 // makes data from SDigit Tree
463 //
464 TClonesArray * sdigits = new TClonesArray("AliTOFSDigit", 1000) ;
465
466 TBranch * branch = sdigitTree->GetBranch("TOF") ;
467 if ( ! branch ) {
468 AliError("TOF branch in SDigit Tree not found") ;
469 return;
470 }
471 branch->SetAddress(&sdigits) ;
472 branch->GetEntry(0) ;
473 MakeSDigits(sdigits) ;
92664bc8 474 //
475 IncEvCountCycleSDigits();
476 IncEvCountTotalSDigits();
477 //
380c04bf 478}
479
480//____________________________________________________________________________
481void AliTOFQADataMaker::MakeRaws(AliRawReader* rawReader)
482{
483 //
484 // makes data from Raws
485 //
486
487 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
488 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
489
490
491 Int_t ntof = 0 ;
492 Int_t in[5];
493 Int_t out[5];
494
495 TClonesArray * clonesRawData;
97430dfe 496 //AliTOFRawStream tofInput(rawReader);
497 fTOFRawStream.SetRawReader(rawReader);
380c04bf 498 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
499 rawReader->Reset();
97430dfe 500 fTOFRawStream.LoadRawData(iDDL);
501 clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
380c04bf 502 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
503 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
5561a41b 504 //if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
505 if (!tofRawDatum->GetTOF()) continue;
380c04bf 506 ntof++;
92664bc8 507 FillRawsData(1, tofRawDatum->GetTOF()*tdc2ns) ;//in ns
508 FillRawsData(2, tofRawDatum->GetTOT()*tot2ns) ;//in ns
380c04bf 509
97430dfe 510 fTOFRawStream.EquipmentId2VolumeId(iDDL,
511 tofRawDatum->GetTRM(),
512 tofRawDatum->GetTRMchain(),
513 tofRawDatum->GetTDC(),
514 tofRawDatum->GetTDCchannel(),
515 in);
380c04bf 516
517 GetMapIndeces(in,out);
92664bc8 518 FillRawsData(3, out[0],out[1]) ;//raw map
380c04bf 519
520 } // while loop
521
522 clonesRawData->Clear();
523
524 } // DDL Loop
525
97430dfe 526 fTOFRawStream.Clear();
527
380c04bf 528 Int_t nentries=ntof;
c332d0d1 529 if(nentries<=0){
92664bc8 530 FillRawsData(0,-1.) ;
ed77256d 531 }else{
92664bc8 532 FillRawsData(0,TMath::Log10(nentries)) ;
ed77256d 533 }
92664bc8 534 //
535 IncEvCountCycleRaws();
536 IncEvCountTotalRaws();
537 //
380c04bf 538}
539
540//____________________________________________________________________________
541void AliTOFQADataMaker::MakeRecPoints(TTree * clustersTree)
542{
543 //
544 // Make data from Clusters
545 //
546
547 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
548 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
549
550 Int_t in[5];
551 Int_t out[5];
552
553 TBranch *branch=clustersTree->GetBranch("TOF");
554 if (!branch) {
555 AliError("can't get the branch with the TOF clusters !");
556 return;
557 }
558
ed77256d 559 TClonesArray dummy("AliTOFcluster",10000), *clusters=&dummy;
380c04bf 560 branch->SetAddress(&clusters);
380c04bf 561
ed77256d 562 // Import the tree
563 clustersTree->GetEvent(0);
564
380c04bf 565 Int_t nentries=clusters->GetEntriesFast();
c332d0d1 566 if(nentries<=0){
92664bc8 567 FillRecPointsData(0,-1.) ;
ed77256d 568 }else{
92664bc8 569 FillRecPointsData(0,TMath::Log10(nentries)) ;
ed77256d 570 }
571
380c04bf 572 TIter next(clusters) ;
573 AliTOFcluster * c ;
574 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
92664bc8 575 FillRecPointsData(1,c->GetTDC()*tdc2ns);
576 FillRecPointsData(2,c->GetTDCRAW()*tdc2ns);
577 FillRecPointsData(3,c->GetToT()*tot2ns);
ed77256d 578
380c04bf 579 in[0] = c->GetDetInd(0);
580 in[1] = c->GetDetInd(1);
581 in[2] = c->GetDetInd(2);
582 in[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
583 in[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
ed77256d 584
380c04bf 585 GetMapIndeces(in,out);
92664bc8 586 FillRecPointsData(4,out[0],out[1]);
380c04bf 587 }
92664bc8 588 //
589 IncEvCountCycleRecPoints();
590 IncEvCountTotalRecPoints();
591 //
380c04bf 592}
593
594//____________________________________________________________________________
595void AliTOFQADataMaker::MakeESDs(AliESDEvent * esd)
596{
597 //
598 // make QA data from ESDs
599 //
600 Int_t ntrk = esd->GetNumberOfTracks() ;
601 Int_t ntof=0;
602 Int_t ntofpid=0;
603 while (ntrk--) {
604 AliESDtrack *track=esd->GetTrack(ntrk);
605 Double_t tofTime=track->GetTOFsignal()*1E-3;//in ns
606 Double_t tofTimeRaw=track->GetTOFsignalRaw()*1E-3;//in ns
607 Double_t tofToT=track->GetTOFsignalToT(); //in ns
608 if(!(tofTime>0))continue;
609 ntof++;
92664bc8 610 FillESDsData(1,tofTime);
611 FillESDsData(2,tofTimeRaw);
612 FillESDsData(3,tofToT);
380c04bf 613 //check how many tracks where ESD PID is ok
614 UInt_t status=track->GetStatus();
615 if (((status&AliESDtrack::kESDpid)==0) ||
616 ((status&AliESDtrack::kTOFpid)==0)) continue;
617 ntofpid++;
618 }
619
620 Int_t nentries=ntof;
c332d0d1 621 if(nentries<=0){
92664bc8 622 FillESDsData(0,-1.) ;
ed77256d 623 }else{
92664bc8 624 FillESDsData(0,TMath::Log10(nentries)) ;
ed77256d 625 }
626
5561a41b 627 Float_t ratio = 0.;
628 if(ntof>0) ratio=(Float_t)ntofpid/(Float_t)ntof*100.;
92664bc8 629 FillESDsData(4,ratio) ;
630 //
631 IncEvCountCycleESDs();
632 IncEvCountTotalESDs();
633 //
380c04bf 634}
635
636//____________________________________________________________________________
637void AliTOFQADataMaker::StartOfDetectorCycle()
638{
639 //
640 //Detector specific actions at start of cycle
641 //to be implemented
642}
643
693a790d 644//____________________________________________________________________________
4e25ac79 645void AliTOFQADataMaker::EndOfDetectorCycle(AliQAv1::TASKINDEX task, TObjArray * list)
693a790d 646{
647 //Detector specific actions at end of cycle
648 // do the QA checking
649
4e25ac79 650 AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;
693a790d 651}
380c04bf 652//____________________________________________________________________________
653void AliTOFQADataMaker::GetMapIndeces(Int_t* in , Int_t* out)
654{
655 //
656 //return appropriate indeces for the theta-phi map
657 //
658
659 Int_t npadX = AliTOFGeometry::NpadX();
660 Int_t npadZ = AliTOFGeometry::NpadZ();
661 Int_t nStripA = AliTOFGeometry::NStripA();
662 Int_t nStripB = AliTOFGeometry::NStripB();
663 Int_t nStripC = AliTOFGeometry::NStripC();
664
665 Int_t isector = in[0];
666 Int_t iplate = in[1];
667 Int_t istrip = in[2];
668 Int_t ipadX = in[3];
669 Int_t ipadZ = in[4];
670
671 Int_t stripOffset = 0;
672 switch (iplate) {
673 case 0:
674 stripOffset = 0;
675 break;
676 case 1:
677 stripOffset = nStripC;
678 break;
679 case 2:
680 stripOffset = nStripC+nStripB;
681 break;
682 case 3:
683 stripOffset = nStripC+nStripB+nStripA;
684 break;
685 case 4:
686 stripOffset = nStripC+nStripB+nStripA+nStripB;
687 break;
688 default:
689 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
690 break;
691 };
692 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
693 Int_t phiindex=npadX*isector+ipadX+1;
694 out[0]=zindex;
695 out[1]=phiindex;
696
697}