1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ////////////////////////////////////////////////////////////////////////////
20 // Produces the data needed to calculate the quality assurance. //
21 // All data must be mergeable objects. //
24 // Sylwester Radomski (radomski@physi.uni-heidelberg.de) //
26 ////////////////////////////////////////////////////////////////////////////
28 // --- ROOT system ---
29 #include <TClonesArray.h>
38 // --- AliRoot header files ---
39 #include "AliESDEvent.h"
41 #include "AliTRDcluster.h"
42 #include "AliTRDQADataMakerRec.h"
43 #include "AliTRDgeometry.h"
44 #include "AliTRDdataArrayI.h"
45 #include "AliTRDrawStreamTB.h"
47 #include "AliQAChecker.h"
49 ClassImp(AliTRDQADataMakerRec)
51 //____________________________________________________________________________
52 AliTRDQADataMakerRec::AliTRDQADataMakerRec() :
53 AliQADataMakerRec(AliQA::GetDetName(AliQA::kTRD), "TRD Quality Assurance Data Maker")
56 // Default constructor
59 //____________________________________________________________________________
60 AliTRDQADataMakerRec::AliTRDQADataMakerRec(const AliTRDQADataMakerRec& qadm) :
67 SetName((const char*)qadm.GetName()) ;
68 SetTitle((const char*)qadm.GetTitle());
72 //__________________________________________________________________
73 AliTRDQADataMakerRec& AliTRDQADataMakerRec::operator=(const AliTRDQADataMakerRec& qadm)
79 this->~AliTRDQADataMakerRec();
80 new(this) AliTRDQADataMakerRec(qadm);
85 //____________________________________________________________________________
86 void AliTRDQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX task, TObjArray * list)
89 // Detector specific actions at end of cycle
92 //AliInfo(Form("EndOfCycle", "Fitting RecPoints %d", task));
95 if (task == AliQA::kRECPOINTS) {
99 // Rec points full chambers
100 if (((TH2D*)list->At(1))->GetEntries() > 1e4) {
101 for (Int_t i=0; i<540; i++) {
103 TH1D *h = ((TH2D*)list->At(1))->ProjectionY(Form("qaTRD_recPoints_amp_%d",i), i+1, i+1);
104 if (h->GetSum() < 100) continue; // chamber not present
106 h->Fit("landau", "q0", "goff", 10, 180);
107 TF1 *fit = h->GetFunction("landau");
108 ((TH1D*)list->At(12))->Fill(fit->GetParameter(1));
109 ((TH1D*)list->At(13))->Fill(fit->GetParameter(2));
115 if (((TH2D*)list->At(10))->GetEntries() > 1e5) {
116 for (Int_t i=0; i<540; i++) {
118 TH1D *test = ((TH3D*)list->At(10))->ProjectionZ(Form("ampTime_%d",i), i+1, i+1, 0, 35);
119 if (test->GetSum() < 100) continue;
121 //AliInfo(Form("fitting det = %d", i));
123 for(Int_t j=0; j<35; j++) {
125 TH1D *h = ((TH3D*)list->At(10))->ProjectionZ(Form("ampTime_%d",i), i+1, i+1, j+1, j+1);
126 if (h->GetSum() < 50) continue;
128 h->Fit("landau", "q0", "goff", 10, 180);
129 TF1 *fit = h->GetFunction("landau");
133 TH2D *h2 = (TH2D*)list->At(14+sm);
134 Int_t bin = h2->FindBin(det,j);
135 // printf("%d %d %d\n", det, j, bin);
136 h2->SetBinContent(bin, fit->GetParameter(1));
143 AliQAChecker::Instance()->Run(AliQA::kTRD, task, list) ;
148 //____________________________________________________________________________
149 void AliTRDQADataMakerRec::InitESDs()
152 // Create ESDs histograms in ESDs subdir
155 const Int_t kNhist = 19;
157 Int_t histoCounter = -1 ;
159 hist[++histoCounter] = new TH1D("qaTRD_esd_ntracks", ":Number of tracks", 300, -0.5, 299.5);
160 hist[++histoCounter] = new TH1D("qaTRD_esd_sector", ":Sector", 18, -0.5, 17.7);
161 hist[++histoCounter] = new TH1D("qaTRD_esd_bits", ";Bits", 64, -0.5, 63.5);
163 const Int_t knbits = 6;
164 const char *suf[knbits] = {"TPCi", "TPCo", "TPCz", "TRDo", "TRDr", "TRDz"};
166 for(Int_t i=0; i<knbits; i++) {
167 hist[++histoCounter] = new TH1D(Form("qaTRD_esd_pt%s",suf[i]), ";p_{T} (GeV/c);", 50, 0, 10);
168 hist[++histoCounter] = new TH1D(Form("qaTRD_esd_trdz%s", suf[i]), ";z (cm)", 200, -400, 400);
171 hist[++histoCounter] = new TH1D("qaTRD_esd_clsTRDo", "TRDo;number of clusters", 130, -0.5, 129.5);;
172 hist[++histoCounter] = new TH1D("qaTRD_esd_clsTRDr", "TRDr;number of clusters", 130, -0.5, 129.5);;
173 hist[++histoCounter] = new TH1D("qaTRD_esd_clsTRDz", "TRDz;number of clusters", 130, -0.5, 129.5);;
174 //hist[++histoCounter] = new TH1D("qaTRD_esd_clsRatio", ";cluster ratio", 100, 0., 1.3);;
176 hist[++histoCounter] = new TH2D("qaTRD_esd_sigMom", ";momentum (GeV/c);signal", 100, 0, 5, 200, 0, 1e3);
178 for(Int_t i=0; i<=histoCounter; i++) {
180 Add2ESDsList(hist[i], i);
185 //____________________________________________________________________________
186 void AliTRDQADataMakerRec::InitRecPoints()
189 // Create Reconstructed Points histograms in RecPoints subdir
192 const Int_t kNhist = 14 + 18;
195 hist[0] = new TH1D("qaTRD_recPoints_det", ";Detector ID of the cluster", 540, -0.5, 539.5);
196 hist[1] = new TH2D("qaTRD_recPoints_amp", ";Amplitude", 540, -0.5, 539, 200, -0.5, 199.5);
197 hist[2] = new TH1D("qaTRD_recPoints_npad", ";Number of Pads", 12, -0.5, 11.5);
199 hist[3] = new TH1D("qaTRD_recPoints_dist2", ";residuals [2pad]", 100, -1, 1);
200 hist[4] = new TH1D("qaTRD_recPoints_dist3", ";residuals [3pad]", 100, -1, 1);
201 hist[5] = new TH1D("qaTRD_recPoints_dist4", ";residuals [4pad]", 100, -1, 1);
202 hist[6] = new TH1D("qaTRD_recPoints_dist5", ";residuals [5pad]", 100, -1, 1);
204 hist[7] = new TH2D("qaTRD_recPoints_rowCol", ";row;col", 16, -0.5, 15.5, 145, -0.5, 144.5);
205 hist[8] = new TH1D("qaTRD_recPoints_time", ";time bin", 35, -0.5, 34.5);
206 hist[9] = new TH1D("qaTRD_recPoints_nCls", ";number of clusters", 500, -0.5, 499.5);
208 hist[10] = new TH3D("qaTRD_recPoints_sigTime", ";chamber;time bin;signal",
209 540, -0.5, 539.5, 35, -0.5, 34.5, 100, 0, 200);
210 hist[11] = new TProfile("qaTRD_recPoints_prf", ";distance;center of gravity"
211 , 120, -0.6, 0.6, -1.2, 1.2, "");
213 hist[12] = new TH1D("qaTRD_recPoints_ampMPV", ";amplitude MPV", 100, 0, 100);
214 hist[13] = new TH1D("qaTRD_recPoints_ampSigma", ";amplitude Sigma", 100, 0, 100);
216 for(Int_t i=0; i<18; i++) {
217 hist[14+i] = new TH2D(Form("qaTRD_recPoints_sigTime_sm%d",i), Form("sm%d;det;time bin"),
218 30, -0.5, 29.5, 35, -0.5, 34.5);
219 hist[14+i]->SetMinimum(20);
220 hist[14+i]->SetMaximum(40);
223 for(Int_t i=0; i<kNhist; i++) {
225 Add2RecPointsList(hist[i], i);
230 //____________________________________________________________________________
231 void AliTRDQADataMakerRec::InitRaws()
234 // create Raws histograms in Raws subdir
237 const Int_t kSM = 18;
238 //const Int_t kNCh = 540;
239 const Int_t kNhist = 4+kSM;
242 // four histograms to be published
243 hist[0] = new TH1D("qaTRD_raws_det", ";detector", 540, -0.5, 539.5);
244 hist[1] = new TH1D("qaTRD_raws_sig", ";signal", 100, -0.5, 99.5);
245 hist[2] = new TH1D("qaTRD_raws_timeBin", ";time bin", 40, -0.5, 39.5);
246 hist[3] = new TH1D("qaTRD_raws_smId", ";supermodule", 18, -0.5, 17.5);
249 // one double per MCM (not published)
250 const Int_t kNMCM = 30 * 8 * 16;
251 for(Int_t i=0; i<kSM; i++)
252 hist[4+i] = new TH1D(Form("qaTRD_raws_sm%d",i),"",kNMCM, -0.5, kNMCM-0.5);
255 for(Int_t i=0; i<kNhist; i++) {
257 Add2RawsList(hist[i], i);
262 //____________________________________________________________________________
263 void AliTRDQADataMakerRec::MakeESDs(AliESDEvent * esd)
266 // Make QA data from ESDs
269 Int_t nTracks = esd->GetNumberOfTracks();
270 GetESDsData(0)->Fill(nTracks);
273 for (Int_t i=0; i<nTracks; i++) {
275 AliESDtrack *track = esd->GetTrack(i);
276 const AliExternalTrackParam *paramOut = track->GetOuterParam();
277 const AliExternalTrackParam *paramIn = track->GetInnerParam();
280 if (!paramIn) continue;
281 if (!paramOut) continue;
284 if (track->GetKinkIndex(0) > 0) continue;
286 Double_t extZ = GetExtZ(paramIn);
287 if (TMath::Abs(extZ) > 320) continue; // acceptance cut
289 // .. in the acceptance
290 Int_t sector = GetSector(paramOut->GetAlpha());
291 GetESDsData(1)->Fill(sector);
294 UInt_t status = track->GetStatus();
295 for(Int_t bit=0; bit<32; bit++)
296 if (u<<bit & status) GetESDsData(2)->Fill(bit);
298 const Int_t knbits = 6;
299 Int_t bit[6] = {0,0,0,0,0,0};
300 bit[0] = status & AliESDtrack::kTPCin;
301 bit[1] = status & AliESDtrack::kTPCout;
302 bit[2] = (status & AliESDtrack::kTPCout) && !(status & AliESDtrack::kTRDout);
303 bit[3] = status & AliESDtrack::kTRDout;
304 bit[4] = status & AliESDtrack::kTRDrefit;
305 bit[5] = (status & AliESDtrack::kTRDout) && !(status & AliESDtrack::kTRDrefit);
307 // transverse momentum
308 //const Double_t *val = paramOut->GetParameter(); // parameters at the Outer plane
309 Double_t pt = paramOut->Pt(); //1./TMath::Abs(val[4]);
311 for(Int_t b=0; b<knbits; b++) {
313 GetESDsData(2*b+3)->Fill(pt);
314 GetESDsData(2*b+4)->Fill(extZ);
319 for(Int_t b=0; b<3; b++)
320 if (bit[3+b]) GetESDsData(b+15)->Fill(track->GetTRDncls());
323 if (!bit[4]) continue;
325 //fQuality->Fill(track->GetTRDQuality());
326 //fBudget->Fill(track->GetTRDBudget());
327 //fSignal->Fill(track->GetTRDsignal());
329 GetESDsData(18)->Fill(track->GetP(), track->GetTRDsignal());
333 if (status & AliESDtrack::kTRDpid) {
335 for(Int_t l=0; l<6; l++) fTime->Fill(track->GetTRDTimBin(l));
337 // fill pid histograms
338 Double_t trdr0 = 0; //, tpcr0 = 0;
339 Int_t trdBestPid = 5; //, tpcBestPid = 5; // charged
340 const Double_t kminPidValue = 0.9;
343 //track->GetTPCpid(pp); // ESD inconsequence
345 for(Int_t pid=0; pid<5; pid++) {
347 trdr0 += track->GetTRDpid(pid);
350 fTrdPID[pid]->Fill(track->GetTRDpid(pid));
351 //fTpcPID[pid]->Fill(pp[pid]);
353 if (track->GetTRDpid(pid) > kminPidValue) trdBestPid = pid;
354 //if (pp[pid] > kminPidValue) tpcBestPid = pid;
357 fTrdPID[5]->Fill(trdr0); // check unitarity
358 fTrdSigMomPID[trdBestPid]->Fill(track->GetP(), track->GetTRDsignal());
360 //fTpcPID[5]->Fill(tpcr0); // check unitarity
361 //fTpcSigMomPID[tpcBestPid]->Fill(track->GetP(), track->GetTPCsignal());
369 //______________________________________________________________________________
370 Int_t AliTRDQADataMakerRec::GetSector(const Double_t alpha) const
373 // Gets the sector number
376 Double_t size = TMath::DegToRad() * 20.; // shall use TRDgeo
377 Int_t sector = (Int_t)((alpha + TMath::Pi())/size);
382 //______________________________________________________________________________
383 Double_t AliTRDQADataMakerRec::GetExtZ(const AliExternalTrackParam *in) const
386 // Returns the Z position at the entry to TRD
387 // using parameters from the TPC in
390 const Double_t kX0 = 300;
392 Double_t x = in->GetX();
393 const Double_t *par = in->GetParameter();
394 Double_t theta = par[3];
395 Double_t z = in->GetZ();
397 Double_t zz = z + (kX0-x) * TMath::Tan(theta);
402 //____________________________________________________________________________
403 void AliTRDQADataMakerRec::MakeRaws(AliRawReader* rawReader)
406 // Makes QA data from raw data
412 //const Int_t kSM = 18;
413 //const Int_t kROC = 30;
414 const Int_t kROB = 8;
415 //const Int_t kLayer = 6;
416 //const Int_t kStack = 5;
417 const Int_t kMCM = 16;
418 // const Int_t kADC = 22;
420 AliTRDrawStreamTB *raw = new AliTRDrawStreamTB(rawReader);
422 raw->SetRawVersion(3);
425 while (raw->Next()) {
427 GetRawsData(0)->Fill(raw->GetDet());
429 // possibly needs changes with the new reader !!
430 Int_t *sig = raw->GetSignals();
431 for(Int_t i=0; i<3; i++) GetRawsData(1)->Fill(sig[i]);
434 GetRawsData(2)->Fill(raw->GetTimeBin());
436 // calculate the index;
437 Int_t sm = raw->GetSM();
438 Int_t roc = raw->GetROC();
439 Int_t rob = raw->GetROB();
440 Int_t mcm = raw->GetMCM();
441 //Int_t adc = raw->GetADC();
443 //Int_t index = roc * (kROB*kMCM*kADC) + rob * (kMCM*kADC) + mcm * kADC + adc;
444 Int_t index = roc * (kROB*kMCM) + rob * kMCM + mcm;
445 GetRawsData(3)->Fill(sm);
446 GetRawsData(4+sm)->Fill(index);
450 //____________________________________________________________________________
451 void AliTRDQADataMakerRec::MakeRecPoints(TTree * clustersTree)
454 // Makes data from RecPoints
457 Int_t nsize = Int_t(clustersTree->GetTotBytes() / (sizeof(AliTRDcluster)));
458 TObjArray *clusterArray = new TObjArray(nsize+1000);
460 TBranch *branch = clustersTree->GetBranch("TRDcluster");
462 AliError("Can't get the branch !");
465 branch->SetAddress(&clusterArray);
467 // Loop through all entries in the tree
468 Int_t nEntries = (Int_t) clustersTree->GetEntries();
470 AliTRDcluster *c = 0;
472 for (Int_t i=0; i<540; i++) nDet[i] = 0;
474 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
477 nbytes += clustersTree->GetEvent(iEntry);
479 // Get the number of points in the detector
480 Int_t nCluster = clusterArray->GetEntriesFast();
482 // Loop through all TRD digits
483 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
484 c = (AliTRDcluster *) clusterArray->UncheckedAt(iCluster);
486 Int_t iDet = c->GetDetector();
488 GetRecPointsData(0)->Fill(iDet);
489 GetRecPointsData(1)->Fill(iDet, c->GetQ());
490 GetRecPointsData(2)->Fill(c->GetNPads());
491 if (c->GetNPads() < 6)
492 GetRecPointsData(1+c->GetNPads())->Fill(c->GetCenter());
494 //if (c->GetPadTime() < 5)
495 ((TH2D*)GetRecPointsData(7))->Fill(c->GetPadRow(), c->GetPadCol());
496 GetRecPointsData(8)->Fill(c->GetPadTime());
498 ((TH3D*)GetRecPointsData(10))->Fill(iDet, c->GetPadTime(), c->GetQ());
501 //if (c->GetNPads() == 2) {
502 Short_t *sig = c->GetSignals();
505 if (sig[0] == 0 && sig[1] == 0 && sig[2] == 0 && sig[5] == 0 && sig[6] == 0)
506 frac = 1. * sig[4] / (sig[3] + sig[4]);
508 if (sig[0] == 0 && sig[1] == 0 && sig[4] == 0 && sig[5] == 0 && sig[6] == 0)
509 frac = -1. * sig[2] / (sig[2] + sig[3]);
511 if (frac > -10) ((TProfile*)GetRecPointsData(11))->Fill(c->GetCenter(), frac);
517 for(Int_t i=0; i<540; i++)
518 if (nDet[i] > 0) GetRecPointsData(9)->Fill(nDet[i]);
524 //____________________________________________________________________________
525 void AliTRDQADataMakerRec::StartOfDetectorCycle()
528 // Detector specific actions at start of cycle
533 //__________________________________________________________________________
534 Int_t AliTRDQADataMakerRec::CheckPointer(TObject *obj, const char *name)
537 // Checks initialization of pointers
540 if (!obj) AliWarning(Form("null pointer: %s", name));