]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFQADataMakerRec.cxx
Do not get main runloader from gAlice
[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"
04236e67 41
42ClassImp(AliTOFQADataMakerRec)
43
44//____________________________________________________________________________
45 AliTOFQADataMakerRec::AliTOFQADataMakerRec() :
46 AliQADataMakerRec(AliQA::GetDetName(AliQA::kTOF), "TOF Quality Assurance Data Maker")
47{
48 //
49 // ctor
50 //
51}
52
53//____________________________________________________________________________
54AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
55 AliQADataMakerRec()
56{
57 //
58 //copy ctor
59 //
60 SetName((const char*)qadm.GetName()) ;
61 SetTitle((const char*)qadm.GetTitle());
62}
63
64//__________________________________________________________________
65AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
66{
67 //
68 // assignment operator.
69 //
70 this->~AliTOFQADataMakerRec();
71 new(this) AliTOFQADataMakerRec(qadm);
72 return *this;
73}
74
75//____________________________________________________________________________
76void AliTOFQADataMakerRec::InitRaws()
77{
78 //
79 // create Raws histograms in Raws subdir
80 //
135306ea 81
82 Bool_t expert = kFALSE;
83
04236e67 84 TH1F * h0 = new TH1F("hTOFRaws", "Number of TOF Raws ",301, -1.02, 5.) ; h0->Sumw2() ;
135306ea 85 Add2RawsList(h0, 0, expert) ;
04236e67 86
87 TH1F * h1 = new TH1F("hTOFRawsTime", "Raws Time Spectrum in TOF (ns)", 2000, 0., 200) ;
88 h1->Sumw2() ;
135306ea 89 Add2RawsList(h1, 1, expert) ;
04236e67 90
91 TH1F * h2 = new TH1F("hTOFRawsToT", "Raws ToT Spectrum in TOF (ns)", 500, 0., 50) ;
92 h2->Sumw2() ;
135306ea 93 Add2RawsList(h2, 2, expert) ;
04236e67 94
95 TH2F * h3 = new TH2F("hTOFRawsClusMap","Raws vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
96 h3->Sumw2() ;
135306ea 97 Add2RawsList(h3, 3, expert) ;
04236e67 98
99}
100
101//____________________________________________________________________________
102void AliTOFQADataMakerRec::InitRecPoints()
103{
104 //
105 // create RecPoints histograms in RecPoints subdir
106 //
135306ea 107
108 Bool_t expert = kFALSE;
109
04236e67 110 TH1F * h0 = new TH1F("hTOFRecPoints", "Number of TOF RecPoints ",301, -1.02, 5.) ; h0->Sumw2() ;
135306ea 111 Add2RecPointsList(h0, 0, expert) ;
04236e67 112
113 TH1F * h1 = new TH1F("hTOFRecPointsTime", "RecPoints Time Spectrum in TOF (ns)", 2000, 0., 200) ;
114 h1->Sumw2() ;
135306ea 115 Add2RecPointsList(h1, 1, expert) ;
04236e67 116
117 TH1F * h2 = new TH1F("hTOFRecPointsRawTime", "RecPoints raw Time Spectrum in TOF (ns)", 2000, 0., 200) ;
118 h2->Sumw2() ;
135306ea 119 Add2RecPointsList(h2, 2, expert) ;
04236e67 120
121 TH1F * h3 = new TH1F("hTOFRecPointsToT", "RecPoints ToT Spectrum in TOF (ns)", 500, 0., 50) ;
122 h3->Sumw2() ;
135306ea 123 Add2RecPointsList(h3, 3, expert) ;
04236e67 124
125 TH2F * h4 = new TH2F("hTOFRecPointsClusMap","RecPoints vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
126 h4->Sumw2() ;
135306ea 127 Add2RecPointsList(h4, 4, expert) ;
04236e67 128
129}
130
131//____________________________________________________________________________
132void AliTOFQADataMakerRec::InitESDs()
133{
134 //
135 //create ESDs histograms in ESDs subdir
136 //
135306ea 137
138 Bool_t expert = kFALSE;
139
04236e67 140 TH1F * h0 = new TH1F("hTOFESDs", "Number of matched TOF tracks over ESDs", 250, -1., 4.) ;
141 h0->Sumw2() ;
135306ea 142 Add2ESDsList(h0, 0, expert) ;
04236e67 143
144 TH1F * h1 = new TH1F("hTOFESDsTime", "Time Spectrum in TOF (ns)", 2000, 0., 200) ;
145 h1->Sumw2() ;
135306ea 146 Add2ESDsList(h1, 1, expert) ;
04236e67 147
148 TH1F * h2 = new TH1F("hTOFESDsRawTime", "raw Time Spectrum in TOF (ns)", 2000, 0., 200) ;
149 h2->Sumw2() ;
135306ea 150 Add2ESDsList(h2, 2, expert) ;
04236e67 151
152 TH1F * h3 = new TH1F("hTOFESDsToT", "ToT Spectrum in TOF (ns)", 500, 0., 50) ;
153 h3->Sumw2() ;
135306ea 154 Add2ESDsList(h3, 3, expert) ;
04236e67 155
156 TH1F * h4 = new TH1F("hTOFESDsPID", "Fraction of matched TOF tracks with good PID glag", 100, 0., 1.) ;
157 h4->Sumw2() ;
135306ea 158 Add2ESDsList(h4, 4, expert) ;
04236e67 159}
160
161//____________________________________________________________________________
162void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
163{
164 //
165 // makes data from Raws
166 //
167
168 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
169 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
170
171
172 Int_t ntof = 0 ;
173 Int_t in[5];
174 Int_t out[5];
175
176 TClonesArray * clonesRawData;
177 AliTOFRawStream tofInput(rawReader);
178 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
179 rawReader->Reset();
180 tofInput.LoadRawData(iDDL);
181 clonesRawData = (TClonesArray*)tofInput.GetRawData();
182 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
183 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
184 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
185 ntof++;
186 GetRawsData(1)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
187 GetRawsData(2)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
188
189 tofInput.EquipmentId2VolumeId(iDDL,
190 tofRawDatum->GetTRM(),
191 tofRawDatum->GetTRMchain(),
192 tofRawDatum->GetTDC(),
193 tofRawDatum->GetTDCchannel(),
194 in);
195
196 GetMapIndeces(in,out);
197 GetRawsData(3)->Fill( out[0],out[1]) ;//raw map
198
199 } // while loop
200
201 clonesRawData->Clear();
202
203 } // DDL Loop
204
205 Int_t nentries=ntof;
206 if(nentries<=0){
207 GetRawsData(0)->Fill(-1.) ;
208 }else{
209 GetRawsData(0)->Fill(TMath::Log10(nentries)) ;
210 }
211}
212
213//____________________________________________________________________________
214void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
215{
216 //
217 // Make data from Clusters
218 //
219
220 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
221 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
222
223 Int_t in[5];
224 Int_t out[5];
225
226 TBranch *branch=clustersTree->GetBranch("TOF");
227 if (!branch) {
228 AliError("can't get the branch with the TOF clusters !");
229 return;
230 }
231
b12e6fe4 232 static TClonesArray dummy("AliTOFcluster",10000);
233 dummy.Clear();
234 TClonesArray *clusters=&dummy;
04236e67 235 branch->SetAddress(&clusters);
236
237 // Import the tree
238 clustersTree->GetEvent(0);
239
240 Int_t nentries=clusters->GetEntriesFast();
241 if(nentries<=0){
242 GetRecPointsData(0)->Fill(-1.) ;
243 }else{
244 GetRecPointsData(0)->Fill(TMath::Log10(nentries)) ;
245 }
246
247 TIter next(clusters) ;
248 AliTOFcluster * c ;
249 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
250 GetRecPointsData(1)->Fill(c->GetTDC()*tdc2ns);
251 GetRecPointsData(2)->Fill(c->GetTDCRAW()*tdc2ns);
252 GetRecPointsData(3)->Fill(c->GetToT()*tot2ns);
253
254 in[0] = c->GetDetInd(0);
255 in[1] = c->GetDetInd(1);
256 in[2] = c->GetDetInd(2);
257 in[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
258 in[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
259
260 GetMapIndeces(in,out);
261 GetRecPointsData(4)->Fill(out[0],out[1]);
262
263 }
264}
265
266//____________________________________________________________________________
267void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
268{
269 //
270 // make QA data from ESDs
271 //
272 Int_t ntrk = esd->GetNumberOfTracks() ;
273 Int_t ntof=0;
274 Int_t ntofpid=0;
275 while (ntrk--) {
276 AliESDtrack *track=esd->GetTrack(ntrk);
277 Double_t tofTime=track->GetTOFsignal()*1E-3;//in ns
278 Double_t tofTimeRaw=track->GetTOFsignalRaw()*1E-3;//in ns
279 Double_t tofToT=track->GetTOFsignalToT(); //in ns
280 if(!(tofTime>0))continue;
281 ntof++;
282 GetESDsData(1)->Fill(tofTime);
283 GetESDsData(2)->Fill(tofTimeRaw);
284 GetESDsData(3)->Fill(tofToT);
285 //check how many tracks where ESD PID is ok
286 UInt_t status=track->GetStatus();
287 if (((status&AliESDtrack::kESDpid)==0) ||
288 ((status&AliESDtrack::kTOFpid)==0)) continue;
289 ntofpid++;
290 }
291
292 Int_t nentries=ntof;
293 if(nentries<=0){
294 GetESDsData(0)->Fill(-1.) ;
295 }else{
296 GetESDsData(0)->Fill(TMath::Log10(nentries)) ;
297 }
298
299 if(ntof>0)GetESDsData(4)->Fill(ntofpid/ntof) ;
300
301}
302
303//____________________________________________________________________________
304void AliTOFQADataMakerRec::StartOfDetectorCycle()
305{
306 //
307 //Detector specific actions at start of cycle
308 //to be implemented
309}
310
311//____________________________________________________________________________
92a357bf 312void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
04236e67 313{
314 //Detector specific actions at end of cycle
315 // do the QA checking
316
317 AliQAChecker::Instance()->Run(AliQA::kTOF, task, list) ;
318}
319//____________________________________________________________________________
320void AliTOFQADataMakerRec::GetMapIndeces(Int_t* in , Int_t* out)
321{
322 //
323 //return appropriate indeces for the theta-phi map
324 //
325
326 Int_t npadX = AliTOFGeometry::NpadX();
327 Int_t npadZ = AliTOFGeometry::NpadZ();
328 Int_t nStripA = AliTOFGeometry::NStripA();
329 Int_t nStripB = AliTOFGeometry::NStripB();
330 Int_t nStripC = AliTOFGeometry::NStripC();
331
332 Int_t isector = in[0];
333 Int_t iplate = in[1];
334 Int_t istrip = in[2];
335 Int_t ipadX = in[3];
336 Int_t ipadZ = in[4];
337
338 Int_t stripOffset = 0;
339 switch (iplate) {
340 case 0:
341 stripOffset = 0;
342 break;
343 case 1:
344 stripOffset = nStripC;
345 break;
346 case 2:
347 stripOffset = nStripC+nStripB;
348 break;
349 case 3:
350 stripOffset = nStripC+nStripB+nStripA;
351 break;
352 case 4:
353 stripOffset = nStripC+nStripB+nStripA+nStripB;
354 break;
355 default:
356 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
357 break;
358 };
359 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
360 Int_t phiindex=npadX*isector+ipadX+1;
361 out[0]=zindex;
362 out[1]=phiindex;
363
364}