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 **************************************************************************/
17 /* $Id: AliADQADataMakerSim.cxx 23123 2007-12-18 09:08:18Z hristov $ */
20 // Produces the data needed to calculate the quality assurance.
21 // All data must be mergeable objects.
25 // --- ROOT system ---
26 #include <TClonesArray.h>
29 #include <TDirectory.h>
30 // --- Standard library ---
32 // --- AliRoot header files ---
33 #include "AliESDEvent.h"
35 #include "AliADdigit.h"
36 #include "AliADSDigit.h"
38 #include "AliADQADataMakerSim.h"
39 #include "AliQAChecker.h"
41 ClassImp(AliADQADataMakerSim)
43 //____________________________________________________________________________
44 AliADQADataMakerSim::AliADQADataMakerSim() :
45 AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kAD), "AD Quality Assurance Data Maker")
53 //____________________________________________________________________________
54 AliADQADataMakerSim::AliADQADataMakerSim(const AliADQADataMakerSim& qadm) :
59 SetName((const char*)qadm.GetName()) ;
60 SetTitle((const char*)qadm.GetTitle());
63 //__________________________________________________________________
64 AliADQADataMakerSim& AliADQADataMakerSim::operator = (const AliADQADataMakerSim& qadm )
67 this->~AliADQADataMakerSim();
68 new(this) AliADQADataMakerSim(qadm);
71 //____________________________________________________________________________
72 void AliADQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
74 //Detector specific actions at end of cycle
76 ResetEventTrigClasses();
77 AliQAChecker::Instance()->Run(AliQAv1::kAD, task, list) ;
81 //____________________________________________________________________________
82 void AliADQADataMakerSim::InitHits()
85 // create Hits histograms in Hits subdir
86 const Bool_t expert = kTRUE ;
87 const Bool_t image = kTRUE ;
89 TH1I * h0 = new TH1I("hHitMultiplicity", "Hit multiplicity distribution in AD;# of Hits;Entries", 300, 0, 299) ;
91 Add2HitsList(h0, 0, !expert, image) ;
93 TH1I * h1 = new TH1I("hHitCellNumber", "Hit cell distribution in AD;Cell;# of Hits", 16, 0, 16) ;
95 Add2HitsList(h1, 1, !expert, image) ;
97 TH1I * h2 = new TH1I("hHitNPhotons", "Number of photons per hit in AD;# of Photons;Entries", 100000, 0, 100000) ;
99 Add2HitsList(h2, 2, !expert, image) ;
101 TH2I * h3 = new TH2I("hCellNPhotons", "Number of photons per cell in AD;Cell;# of Photons", 16, 0, 16, 100000, 0, 100000) ;
103 Add2HitsList(h3, 3, !expert, image) ;
105 TH2D * h4 = new TH2D("hCellTof", "Time of flight per cell in AD;Cell;Time of Flight [ns]", 16, 0, 16, 6000,40,100) ;
107 Add2HitsList(h4, 4, !expert, image) ;
110 ClonePerTrigClass(AliQAv1::kHITS); // this should be the last line
114 //____________________________________________________________________________
115 void AliADQADataMakerSim::InitSDigits()
117 // create Digits histograms in Digits subdir
118 const Bool_t expert = kTRUE ;
119 const Bool_t image = kTRUE ;
121 TH1I *fhSDigCharge[16];
123 // create SDigits histograms in SDigits subdir
124 TH1I * h0 = new TH1I("hSDigitMultiplicity", "SDigits multiplicity distribution in AD;# of Digits;Entries", 100, 0, 99) ;
126 Add2SDigitsList(h0, 0, !expert, image) ;
128 for (Int_t i=0; i<16; i++)
130 fhSDigCharge[i] = new TH1I(Form("hSDigitCharge%d", i),Form("SDigit charges in cell %d; Time;Entries",i),1700,0.,1700);
132 Add2SDigitsList(fhSDigCharge[i],i+1, !expert, image);
136 ClonePerTrigClass(AliQAv1::kSDIGITS); // this should be the last line
141 //____________________________________________________________________________
142 void AliADQADataMakerSim::InitDigits()
144 // create Digits histograms in Digits subdir
145 const Bool_t expert = kTRUE ;
146 const Bool_t image = kTRUE ;
148 // create Digits histograms in Digits subdir
149 TH1I * h0 = new TH1I("hDigitMultiplicity", "Digits multiplicity distribution in AD;# of Digits;Entries", 100, 0, 99) ;
151 Add2DigitsList(h0, 0, !expert, image) ;
153 TH2D * h1 = new TH2D("hDigitLeadingTimePerPM", "Leading time distribution per PM in AD;PM number;Leading Time [ns]",16,0,16, 1000, 200, 300);
155 Add2DigitsList(h1, 1, !expert, image) ;
157 TH2D * h2 = new TH2D("hDigitTimeWidthPerPM", "Time width distribution per PM in AD;PM number;Time width [ns]",16,0,16, 1000, 0, 100);
159 Add2DigitsList(h2, 2, !expert, image) ;
161 TH2I * h3 = new TH2I("hDigitChargePerClockPerPM", "Charge array per PM in AD;PM number; Clock",16,0,16,21, -10.5, 10.5);
163 Add2DigitsList(h3, 3, !expert, image) ;
165 TH1I * h4 = new TH1I("hDigitBBflagsAD","Number of BB flags in AD; # of BB flags; Entries",17,-0.5,16.5);
167 Add2DigitsList(h4, 4, !expert, image) ;
169 TH1I * h5 = new TH1I("hDigitBBflagsADA","Number of BB flags in ADA; # of BB flags; Entries",9,-0.5,8.5);
171 Add2DigitsList(h5, 5, !expert, image) ;
173 TH1I * h6 = new TH1I("hDigitBBflagsADC","Number of BB flags in ADC; # of BB flags; Entries",9,-0.5,8.5);
175 Add2DigitsList(h6, 6, !expert, image) ;
177 TH2D * h7 = new TH2D("hDigitTotalChargePerPM", "Total Charge per PM in AD;PM number; Charge [ADC counts]",16,0,16,10000,0,10000);
179 Add2DigitsList(h7, 7, !expert, image) ;
181 TH2I * h8 = new TH2I("hDigitMaxChargeClockPerPM", "Clock with maximum charge per PM in AD;PM number; Clock ",16,0,16,21, -10.5, 10.5);
183 Add2DigitsList(h8, 8, !expert, image) ;
186 ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
190 //____________________________________________________________________________
191 void AliADQADataMakerSim::MakeHits()
193 //make QA data from Hits
195 Int_t nhits = fHitsArray->GetEntriesFast();
196 FillHitsData(0,nhits) ; // fills Hit multiplicity
197 for (Int_t ihit=0;ihit<nhits;ihit++)
199 AliADhit * ADHit = (AliADhit*) fHitsArray->UncheckedAt(ihit);
201 AliError("The unchecked hit doesn't exist");
204 FillHitsData(1,ADHit->GetCell());
205 FillHitsData(2,ADHit->GetNphot());
206 FillHitsData(3,ADHit->GetCell(),ADHit->GetNphot());
207 FillHitsData(4,ADHit->GetCell(),ADHit->GetTof());
212 //____________________________________________________________________________
214 void AliADQADataMakerSim::MakeHits(TTree *hitTree)
216 //fills QA histos for Hits
218 fHitsArray->Clear() ;
220 fHitsArray = new TClonesArray("AliADhit", 1000);
222 TBranch * branch = hitTree->GetBranch("AD") ;
224 AliWarning("AD branch in Hit Tree not found") ;
228 branch->SetAddress(&fHitsArray);
230 AliError("Branch AD hit not found");
233 // Check id histograms already created for this Event Specie
234 if ( ! GetHitsData(0) )
237 Int_t ntracks = (Int_t) hitTree->GetEntries();
239 if (ntracks<=0) return;
240 // Start loop on tracks in the hits containers
241 for (Int_t track=0; track<ntracks;track++) {
242 branch->GetEntry(track);
243 Int_t nhits = fHitsArray->GetEntriesFast();
244 FillHitsData(0,nhits) ; // fills Hit multiplicity
245 for (Int_t ihit=0;ihit<nhits;ihit++)
247 AliADhit * ADHit = (AliADhit*) fHitsArray->UncheckedAt(ihit);
249 AliError("The unchecked hit doesn't exist");
252 FillHitsData(1,ADHit->GetCell());
253 FillHitsData(2,ADHit->GetNphot());
254 FillHitsData(3,ADHit->GetCell(),ADHit->GetNphot());
255 FillHitsData(4,ADHit->GetCell(),ADHit->GetTof());
260 IncEvCountCycleHits();
261 IncEvCountTotalHits();
267 //____________________________________________________________________________
268 void AliADQADataMakerSim::MakeSDigits(TTree *sdigitTree)
270 // makes data from Digit Tree
273 fSDigitsArray->Clear() ;
275 fSDigitsArray = new TClonesArray("AliADSDigit", 1000) ;
277 TBranch * branch = sdigitTree->GetBranch("ADSDigit") ;
279 AliWarning("AD branch in SDigit Tree not found") ;
281 branch->SetAddress(&fSDigitsArray) ;
282 branch->GetEntry(0) ;
286 IncEvCountCycleDigits();
287 IncEvCountTotalDigits();
291 //____________________________________________________________________________
292 void AliADQADataMakerSim::MakeSDigits()
294 // makes data from SDigits
296 FillSDigitsData(0,fSDigitsArray->GetEntriesFast()) ;
297 TIter next(fSDigitsArray) ;
298 AliADSDigit *ADSDigit ;
299 while ( (ADSDigit = dynamic_cast<AliADSDigit *>(next())) ) {
300 Int_t PMNumber = ADSDigit->PMNumber();
301 FillSDigitsData(PMNumber +1, ADSDigit->GetNBins()) ;
305 //____________________________________________________________________________
306 void AliADQADataMakerSim::MakeDigits()
308 // makes data from Digits
310 FillDigitsData(0,fDigitsArray->GetEntriesFast()) ;
311 TIter next(fDigitsArray) ;
312 AliADdigit *ADDigit ;
313 Int_t nBBflagsADA = 0;
314 Int_t nBBflagsADC = 0;
316 while ( (ADDigit = dynamic_cast<AliADdigit *>(next())) ) {
318 Int_t PMNumber = ADDigit->PMNumber();
320 if(PMNumber<8 && ADDigit->GetBBflag()) nBBflagsADA++;
321 if(PMNumber>7 && ADDigit->GetBBflag()) nBBflagsADC++;
324 for(Int_t iClock=0; iClock<21; iClock++) {
325 adc[iClock]= ADDigit->ChargeADC(iClock);
326 FillDigitsData(3, PMNumber,(float)iClock-10,(float)adc[iClock]);
327 totCharge += adc[iClock];
330 FillDigitsData(1,PMNumber,ADDigit->Time());
331 FillDigitsData(2,PMNumber,ADDigit->Width());
332 FillDigitsData(7,PMNumber,totCharge);
333 FillDigitsData(8,PMNumber,TMath::LocMax(21,adc)-10);
336 FillDigitsData(4,nBBflagsADA+nBBflagsADC);
337 FillDigitsData(5,nBBflagsADA);
338 FillDigitsData(6,nBBflagsADC);
341 //____________________________________________________________________________
342 void AliADQADataMakerSim::MakeDigits(TTree *digitTree)
344 // makes data from Digit Tree
347 fDigitsArray->Clear() ;
349 fDigitsArray = new TClonesArray("AliADdigit", 1000) ;
351 TBranch * branch = digitTree->GetBranch("ADDigit") ;
353 AliWarning("AD branch in Digit Tree not found") ;
355 branch->SetAddress(&fDigitsArray) ;
356 branch->GetEntry(0) ;
360 IncEvCountCycleDigits();
361 IncEvCountTotalDigits();
366 //____________________________________________________________________________
367 void AliADQADataMakerSim::StartOfDetectorCycle()
369 //Detector specific actions at start of cycle