]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFQADataMakerRec.cxx
Fixes to remove compilation warnings
[u/mrichter/AliRoot.git] / TOF / AliTOFQADataMakerRec.cxx
CommitLineData
04236e67 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///////////////////////////////////////////////////////////////////////
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///////////////////////////////////////////////////////////////////////
23
24#include <TClonesArray.h>
5c7c93fa 25//#include <TFile.h>
26//#include <TH1I.h>
04236e67 27#include <TH1F.h>
28#include <TH2F.h>
5c7c93fa 29
30#include "AliLog.h"
04236e67 31#include "AliESDEvent.h"
32#include "AliESDtrack.h"
04236e67 33#include "AliQAChecker.h"
34#include "AliRawReader.h"
5c7c93fa 35
36#include "AliTOFcluster.h"
37#include "AliTOFQADataMakerRec.h"
04236e67 38#include "AliTOFRawStream.h"
39#include "AliTOFrawData.h"
5c7c93fa 40#include "AliTOFGeometry.h"
5192f264 41#include "AliTOFdigit.h"
04236e67 42
43ClassImp(AliTOFQADataMakerRec)
44
45//____________________________________________________________________________
46 AliTOFQADataMakerRec::AliTOFQADataMakerRec() :
4e25ac79 47 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker")
04236e67 48{
49 //
50 // ctor
51 //
52}
53
54//____________________________________________________________________________
55AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
56 AliQADataMakerRec()
57{
58 //
59 //copy ctor
60 //
61 SetName((const char*)qadm.GetName()) ;
62 SetTitle((const char*)qadm.GetTitle());
63}
64
65//__________________________________________________________________
66AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
67{
68 //
69 // assignment operator.
70 //
71 this->~AliTOFQADataMakerRec();
72 new(this) AliTOFQADataMakerRec(qadm);
73 return *this;
74}
75
76//____________________________________________________________________________
77void AliTOFQADataMakerRec::InitRaws()
78{
79 //
80 // create Raws histograms in Raws subdir
81 //
135306ea 82
7d297381 83 const Bool_t expert = kTRUE ;
84 const Bool_t saveCorr = kTRUE ;
85 const Bool_t image = kTRUE ;
86
db72ff3b 87 TH1F * h0 = new TH1F("hTOFRaws", "Number of TOF Raws;TOF raw number [10 power];Counts ",301, -1.02, 5.) ;
5561a41b 88 h0->Sumw2() ;
7d297381 89 Add2RawsList(h0, 0, !expert, image, !saveCorr) ;
04236e67 90
db72ff3b 91 TH1F * h1 = new TH1F("hTOFRawsTime", "Raws Time Spectrum in TOF (ns);Measured TOF time [ns];Counts", 2000, 0., 200) ;
04236e67 92 h1->Sumw2() ;
7d297381 93 Add2RawsList(h1, 1, !expert, image, !saveCorr) ;
04236e67 94
db72ff3b 95 TH1F * h2 = new TH1F("hTOFRawsToT", "Raws ToT Spectrum in TOF (ns);Measured TOT time [ns];Counts", 500, 0., 50) ;
04236e67 96 h2->Sumw2() ;
7d297381 97 Add2RawsList(h2, 2, !expert, image, !saveCorr) ;
04236e67 98
db72ff3b 99 TH2F * h3 = new TH2F("hTOFRawsClusMap","Raws vs TOF eta-phi;2*strip+padz (eta);48*sector+padx (phi)",183, -0.5, 182.5,865,-0.5,864.5) ;
04236e67 100 h3->Sumw2() ;
5561a41b 101 h3->GetYaxis()->SetTitleOffset(1.15);
7d297381 102 Add2RawsList(h3, 3, !expert, image, !saveCorr) ;
04236e67 103
104}
105
44ed7a66 106//____________________________________________________________________________
107void AliTOFQADataMakerRec::InitDigits()
108{
109 //
110 // create Digits histograms in Digits subdir
111 //
112
113 const Bool_t expert = kTRUE ;
114 const Bool_t image = kTRUE ;
115
db72ff3b 116 TH1F * h0 = new TH1F("hTOFDigits", "Number of TOF Digits;TOF digit number [10 power];Counts ",301, -1.02, 5.) ;
44ed7a66 117 h0->Sumw2() ;
44ed7a66 118 Add2DigitsList(h0, 0, !expert, image) ;
119
db72ff3b 120 TH1F * h1 = new TH1F("hTOFDigitsTime", "Digits Time Spectrum in TOF (ns);Digitized TOF time [ns];Counts", 2000, 0., 200) ;
44ed7a66 121 h1->Sumw2() ;
44ed7a66 122 Add2DigitsList(h1, 1, !expert, image) ;
123
db72ff3b 124 TH1F * h2 = new TH1F("hTOFDigitsToT", "Digits ToT Spectrum in TOF (ns);Digitized TOF time [ns];Counts", 500, 0., 50) ;
44ed7a66 125 h2->Sumw2() ;
44ed7a66 126 Add2DigitsList(h2, 2, !expert, image) ;
127
db72ff3b 128 TH2F * h3 = new TH2F("hTOFDigitsClusMap","Digits vs TOF eta-phi;2*strip+padz (eta);48*sector+padx (phi)",183, -0.5, 182.5,865,-0.5,864.5) ;
44ed7a66 129 h3->Sumw2() ;
44ed7a66 130 h3->GetYaxis()->SetTitleOffset(1.15);
131 Add2DigitsList(h3, 3, !expert, image) ;
44ed7a66 132}
133
04236e67 134//____________________________________________________________________________
135void AliTOFQADataMakerRec::InitRecPoints()
136{
137 //
138 // create RecPoints histograms in RecPoints subdir
139 //
135306ea 140
7d297381 141 const Bool_t expert = kTRUE ;
142 const Bool_t image = kTRUE ;
eca4fa66 143
db72ff3b 144 TH1F * h0 = new TH1F("hTOFRecPoints", "Number of TOF RecPoints;TOF recPoint number [10 power];Counts",301, -1.02, 5.) ;
5561a41b 145 h0->Sumw2() ;
7d297381 146 Add2RecPointsList(h0, 0, !expert, image) ;
04236e67 147
db72ff3b 148 TH1F * h1 = new TH1F("hTOFRecPointsTime", "RecPoints Time Spectrum in TOF (ns);Calibrated TOF time [ns];Counts", 2000, 0., 200) ;
04236e67 149 h1->Sumw2() ;
7d297381 150 Add2RecPointsList(h1, 1, !expert, image) ;
04236e67 151
db72ff3b 152 TH1F * h2 = new TH1F("hTOFRecPointsRawTime", "RecPoints raw Time Spectrum in TOF (ns);Measured TOF time [ns];Counts", 2000, 0., 200) ;
04236e67 153 h2->Sumw2() ;
7d297381 154 Add2RecPointsList(h2, 2, !expert, image) ;
04236e67 155
db72ff3b 156 TH1F * h3 = new TH1F("hTOFRecPointsToT", "RecPoints ToT Spectrum in TOF (ns);Measured TOT [ns];Counts", 500, 0., 50) ;
04236e67 157 h3->Sumw2() ;
7d297381 158 Add2RecPointsList(h3, 3, !expert, image) ;
04236e67 159
db72ff3b 160 TH2F * h4 = new TH2F("hTOFRecPointsClusMap","RecPoints vs TOF phi-eta;2*strip+padz (eta);48*sector+padx (phi)",183, -0.5, 182.5,865,-0.5,864.5) ;
04236e67 161 h4->Sumw2() ;
5561a41b 162 h4->GetYaxis()->SetTitleOffset(1.15);
7d297381 163 Add2RecPointsList(h4, 4, !expert, image) ;
04236e67 164}
165
166//____________________________________________________________________________
167void AliTOFQADataMakerRec::InitESDs()
168{
169 //
170 //create ESDs histograms in ESDs subdir
171 //
135306ea 172
7d297381 173 const Bool_t expert = kTRUE ;
174 const Bool_t image = kTRUE ;
db72ff3b 175 TH1F * h0 = new TH1F("hTOFESDs", "Number of matched TOF tracks over ESDs;Number of TOF matched ESD tracks [10 power];Counts", 250, -1., 4.) ;
04236e67 176 h0->Sumw2() ;
7d297381 177 Add2ESDsList(h0, 0, !expert, image) ;
04236e67 178
db72ff3b 179 TH1F * h1 = new TH1F("hTOFESDsTime", "Time Spectrum in TOF (ns);Calibrated TOF time [ns];Counts", 2000, 0., 200) ;
04236e67 180 h1->Sumw2() ;
7d297381 181 Add2ESDsList(h1, 1, !expert, image) ;
04236e67 182
db72ff3b 183 TH1F * h2 = new TH1F("hTOFESDsRawTime", "raw Time Spectrum in TOF (ns);Measured TOF time [ns];Counts", 2000, 0., 200) ;
04236e67 184 h2->Sumw2() ;
7d297381 185 Add2ESDsList(h2, 2, !expert, image) ;
04236e67 186
db72ff3b 187 TH1F * h3 = new TH1F("hTOFESDsToT", "ToT Spectrum in TOF (ns);Measured TOF time [ns];Counts", 500, 0., 50) ;
04236e67 188 h3->Sumw2() ;
7d297381 189 Add2ESDsList(h3, 3, !expert, image) ;
04236e67 190
db72ff3b 191 TH1F * h4 = new TH1F("hTOFESDsPID", "Fraction of matched TOF tracks with good PID flag (%);Fraction of TOF matched ESD tracks with good flag [%];Counts", 101, 0., 101.) ;
04236e67 192 h4->Sumw2() ;
7d297381 193 Add2ESDsList(h4, 4, !expert, image) ;
04236e67 194}
195
196//____________________________________________________________________________
197void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
198{
199 //
200 // makes data from Raws
201 //
eca4fa66 202
203 // Check id histograms already created for this Event Specie
204 if ( ! GetRawsData(0) )
205 InitRaws() ;
206
04236e67 207 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
208 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
209
210
211 Int_t ntof = 0 ;
212 Int_t in[5];
213 Int_t out[5];
214
215 TClonesArray * clonesRawData;
216 AliTOFRawStream tofInput(rawReader);
217 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
218 rawReader->Reset();
219 tofInput.LoadRawData(iDDL);
220 clonesRawData = (TClonesArray*)tofInput.GetRawData();
221 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
222 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
223 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
224 ntof++;
225 GetRawsData(1)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
226 GetRawsData(2)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
227
228 tofInput.EquipmentId2VolumeId(iDDL,
229 tofRawDatum->GetTRM(),
230 tofRawDatum->GetTRMchain(),
231 tofRawDatum->GetTDC(),
232 tofRawDatum->GetTDCchannel(),
233 in);
234
235 GetMapIndeces(in,out);
236 GetRawsData(3)->Fill( out[0],out[1]) ;//raw map
237
238 } // while loop
239
240 clonesRawData->Clear();
241
242 } // DDL Loop
243
244 Int_t nentries=ntof;
245 if(nentries<=0){
246 GetRawsData(0)->Fill(-1.) ;
247 }else{
248 GetRawsData(0)->Fill(TMath::Log10(nentries)) ;
249 }
250}
251
252//____________________________________________________________________________
44ed7a66 253void AliTOFQADataMakerRec::MakeDigits(TClonesArray * digits)
254{
255 //
256 // makes data from Digits
257 //
eca4fa66 258
44ed7a66 259 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
260 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
261 Int_t in[5];
262 Int_t out[5];
263
264 Int_t nentries=digits->GetEntriesFast();
265 if(nentries<=0){
266 GetDigitsData(0)->Fill(-1.) ;
267 }else{
268 GetDigitsData(0)->Fill(TMath::Log10(nentries)) ;
269 }
270
271 TIter next(digits) ;
272 AliTOFdigit * digit ;
273 while ( (digit = dynamic_cast<AliTOFdigit *>(next())) ) {
274
275 GetDigitsData(1)->Fill( digit->GetTdc()*tdc2ns) ;//in ns
276 GetDigitsData(2)->Fill( digit->GetToT()*tot2ns) ;//in ns
277
278 in[0] = digit->GetSector();
279 in[1] = digit->GetPlate();
280 in[2] = digit->GetStrip();
281 in[3] = digit->GetPadx();
282 in[4]= digit->GetPadz();
283 GetMapIndeces(in,out);
284 GetDigitsData(3)->Fill( out[0],out[1]) ;//digit map
285 }
44ed7a66 286}
287
288//____________________________________________________________________________
289void AliTOFQADataMakerRec::MakeDigits(TTree * digitTree)
290{
291 //
292 // makes data from Digit Tree
293 //
eca4fa66 294 // Check id histograms already created for this Event Specie
295 if ( ! GetDigitsData(0) )
296 InitDigits() ;
297
44ed7a66 298 TClonesArray * digits = new TClonesArray("AliTOFdigit", 1000) ;
299
300 TBranch * branch = digitTree->GetBranch("TOF") ;
301 if ( ! branch ) {
302 AliError("TOF branch in Digit Tree not found") ;
303 return;
304 }
305 branch->SetAddress(&digits) ;
306 branch->GetEntry(0) ;
307 MakeDigits(digits) ;
308}
309
310//____________________________________________________________________________
04236e67 311void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
312{
313 //
314 // Make data from Clusters
315 //
316
eca4fa66 317 // Check id histograms already created for this Event Specie
318 if ( ! GetRecPointsData(0) )
319 InitRecPoints() ;
320
04236e67 321 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
322 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
323
324 Int_t in[5];
325 Int_t out[5];
326
327 TBranch *branch=clustersTree->GetBranch("TOF");
328 if (!branch) {
329 AliError("can't get the branch with the TOF clusters !");
330 return;
331 }
332
b12e6fe4 333 static TClonesArray dummy("AliTOFcluster",10000);
334 dummy.Clear();
335 TClonesArray *clusters=&dummy;
04236e67 336 branch->SetAddress(&clusters);
337
338 // Import the tree
339 clustersTree->GetEvent(0);
340
341 Int_t nentries=clusters->GetEntriesFast();
342 if(nentries<=0){
343 GetRecPointsData(0)->Fill(-1.) ;
344 }else{
345 GetRecPointsData(0)->Fill(TMath::Log10(nentries)) ;
346 }
347
348 TIter next(clusters) ;
349 AliTOFcluster * c ;
350 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
351 GetRecPointsData(1)->Fill(c->GetTDC()*tdc2ns);
352 GetRecPointsData(2)->Fill(c->GetTDCRAW()*tdc2ns);
353 GetRecPointsData(3)->Fill(c->GetToT()*tot2ns);
354
355 in[0] = c->GetDetInd(0);
356 in[1] = c->GetDetInd(1);
357 in[2] = c->GetDetInd(2);
358 in[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
359 in[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
360
361 GetMapIndeces(in,out);
362 GetRecPointsData(4)->Fill(out[0],out[1]);
363
364 }
365}
366
367//____________________________________________________________________________
368void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
369{
370 //
371 // make QA data from ESDs
372 //
eca4fa66 373 // Check id histograms already created for this Event Specie
374 if ( ! GetESDsData(0) )
375 InitESDs() ;
376
04236e67 377 Int_t ntrk = esd->GetNumberOfTracks() ;
378 Int_t ntof=0;
379 Int_t ntofpid=0;
380 while (ntrk--) {
381 AliESDtrack *track=esd->GetTrack(ntrk);
382 Double_t tofTime=track->GetTOFsignal()*1E-3;//in ns
383 Double_t tofTimeRaw=track->GetTOFsignalRaw()*1E-3;//in ns
384 Double_t tofToT=track->GetTOFsignalToT(); //in ns
385 if(!(tofTime>0))continue;
386 ntof++;
387 GetESDsData(1)->Fill(tofTime);
388 GetESDsData(2)->Fill(tofTimeRaw);
389 GetESDsData(3)->Fill(tofToT);
390 //check how many tracks where ESD PID is ok
391 UInt_t status=track->GetStatus();
392 if (((status&AliESDtrack::kESDpid)==0) ||
393 ((status&AliESDtrack::kTOFpid)==0)) continue;
394 ntofpid++;
395 }
396
397 Int_t nentries=ntof;
398 if(nentries<=0){
399 GetESDsData(0)->Fill(-1.) ;
400 }else{
401 GetESDsData(0)->Fill(TMath::Log10(nentries)) ;
402 }
403
7d53ce61 404 if(ntof>0) {
405 Float_t ratio = (Float_t)ntofpid/(Float_t)ntof*100.;
406 GetESDsData(4)->Fill(ratio) ;
407 }
04236e67 408}
409
410//____________________________________________________________________________
411void AliTOFQADataMakerRec::StartOfDetectorCycle()
412{
413 //
414 //Detector specific actions at start of cycle
415 //to be implemented
416}
417
418//____________________________________________________________________________
4e25ac79 419void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
04236e67 420{
421 //Detector specific actions at end of cycle
422 // do the QA checking
4e25ac79 423 AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;
04236e67 424}
425//____________________________________________________________________________
426void AliTOFQADataMakerRec::GetMapIndeces(Int_t* in , Int_t* out)
427{
428 //
429 //return appropriate indeces for the theta-phi map
430 //
431
432 Int_t npadX = AliTOFGeometry::NpadX();
433 Int_t npadZ = AliTOFGeometry::NpadZ();
434 Int_t nStripA = AliTOFGeometry::NStripA();
435 Int_t nStripB = AliTOFGeometry::NStripB();
436 Int_t nStripC = AliTOFGeometry::NStripC();
437
438 Int_t isector = in[0];
439 Int_t iplate = in[1];
440 Int_t istrip = in[2];
441 Int_t ipadX = in[3];
442 Int_t ipadZ = in[4];
443
444 Int_t stripOffset = 0;
445 switch (iplate) {
446 case 0:
447 stripOffset = 0;
448 break;
449 case 1:
450 stripOffset = nStripC;
451 break;
452 case 2:
453 stripOffset = nStripC+nStripB;
454 break;
455 case 3:
456 stripOffset = nStripC+nStripB+nStripA;
457 break;
458 case 4:
459 stripOffset = nStripC+nStripB+nStripA+nStripB;
460 break;
461 default:
462 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
463 break;
464 };
465 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
466 Int_t phiindex=npadX*isector+ipadX+1;
467 out[0]=zindex;
468 out[1]=phiindex;
469
470}