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", 17,-0.5,16.5) ;
126 Add2SDigitsList(h0, 0, !expert, image) ;
128 TH2D * h1 = new TH2D("hSDigitChargePerPM", "SDigits total amplified charge per PM in AD;PM number;Charge [C]", 16, 0, 16, 1000,1e-13, 1e-9);
130 Add2SDigitsList(h1, 1, !expert, image) ;
133 ClonePerTrigClass(AliQAv1::kSDIGITS); // this should be the last line
138 //____________________________________________________________________________
139 void AliADQADataMakerSim::InitDigits()
141 // create Digits histograms in Digits subdir
142 const Bool_t expert = kTRUE ;
143 const Bool_t image = kTRUE ;
145 // create Digits histograms in Digits subdir
146 TH1I * h0 = new TH1I("hDigitMultiplicity", "Digits multiplicity distribution in AD;# of Digits;Entries", 17,-0.5,16.5) ;
148 Add2DigitsList(h0, 0, !expert, image) ;
150 TH2D * h1 = new TH2D("hDigitLeadingTimePerPM", "Leading time distribution per PM in AD;PM number;Leading Time [ns]",16,0,16, 1000, 200, 300);
152 Add2DigitsList(h1, 1, !expert, image) ;
154 TH2D * h2 = new TH2D("hDigitTimeWidthPerPM", "Time width distribution per PM in AD;PM number;Time width [ns]",16,0,16, 1000, 0, 100);
156 Add2DigitsList(h2, 2, !expert, image) ;
158 TH2I * h3 = new TH2I("hDigitChargePerClockPerPM", "Charge array per PM in AD;PM number; Clock",16,0,16,21, -10.5, 10.5);
160 Add2DigitsList(h3, 3, !expert, image) ;
162 TH1I * h4 = new TH1I("hDigitBBflagsAD","Number of BB flags in AD; # of BB flags; Entries",17,-0.5,16.5);
164 Add2DigitsList(h4, 4, !expert, image) ;
166 TH1I * h5 = new TH1I("hDigitBBflagsADA","Number of BB flags in ADA; # of BB flags; Entries",9,-0.5,8.5);
168 Add2DigitsList(h5, 5, !expert, image) ;
170 TH1I * h6 = new TH1I("hDigitBBflagsADC","Number of BB flags in ADC; # of BB flags; Entries",9,-0.5,8.5);
172 Add2DigitsList(h6, 6, !expert, image) ;
174 TH2D * h7 = new TH2D("hDigitTotalChargePerPM", "Total Charge per PM in AD;PM number; Charge [ADC counts]",16,0,16,10000,0,10000);
176 Add2DigitsList(h7, 7, !expert, image) ;
178 TH2I * h8 = new TH2I("hDigitMaxChargeClockPerPM", "Clock with maximum charge per PM in AD;PM number; Clock ",16,0,16,21, -10.5, 10.5);
180 Add2DigitsList(h8, 8, !expert, image) ;
183 ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
187 //____________________________________________________________________________
188 void AliADQADataMakerSim::MakeHits()
190 //make QA data from Hits
192 Int_t nhits = fHitsArray->GetEntriesFast();
193 FillHitsData(0,nhits) ; // fills Hit multiplicity
194 for (Int_t ihit=0;ihit<nhits;ihit++)
196 AliADhit * ADHit = (AliADhit*) fHitsArray->UncheckedAt(ihit);
198 AliError("The unchecked hit doesn't exist");
201 FillHitsData(1,ADHit->GetCell());
202 FillHitsData(2,ADHit->GetNphot());
203 FillHitsData(3,ADHit->GetCell(),ADHit->GetNphot());
204 FillHitsData(4,ADHit->GetCell(),ADHit->GetTof());
209 //____________________________________________________________________________
211 void AliADQADataMakerSim::MakeHits(TTree *hitTree)
213 //fills QA histos for Hits
215 fHitsArray->Clear() ;
217 fHitsArray = new TClonesArray("AliADhit", 1000);
219 TBranch * branch = hitTree->GetBranch("AD") ;
221 AliWarning("AD branch in Hit Tree not found") ;
225 branch->SetAddress(&fHitsArray);
227 AliError("Branch AD hit not found");
230 // Check id histograms already created for this Event Specie
231 if ( ! GetHitsData(0) )
234 Int_t ntracks = (Int_t) hitTree->GetEntries();
236 if (ntracks<=0) return;
237 // Start loop on tracks in the hits containers
238 for (Int_t track=0; track<ntracks;track++) {
239 branch->GetEntry(track);
240 Int_t nhits = fHitsArray->GetEntriesFast();
241 FillHitsData(0,nhits) ; // fills Hit multiplicity
242 for (Int_t ihit=0;ihit<nhits;ihit++)
244 AliADhit * ADHit = (AliADhit*) fHitsArray->UncheckedAt(ihit);
246 AliError("The unchecked hit doesn't exist");
249 FillHitsData(1,ADHit->GetCell());
250 FillHitsData(2,ADHit->GetNphot());
251 FillHitsData(3,ADHit->GetCell(),ADHit->GetNphot());
252 FillHitsData(4,ADHit->GetCell(),ADHit->GetTof());
257 IncEvCountCycleHits();
258 IncEvCountTotalHits();
264 //____________________________________________________________________________
265 void AliADQADataMakerSim::MakeSDigits(TTree *sdigitTree)
267 // makes data from Digit Tree
270 fSDigitsArray->Clear() ;
272 fSDigitsArray = new TClonesArray("AliADSDigit", 1000) ;
274 TBranch * branch = sdigitTree->GetBranch("ADSDigit") ;
276 AliWarning("AD branch in SDigit Tree not found") ;
278 branch->SetAddress(&fSDigitsArray) ;
279 branch->GetEntry(0) ;
283 IncEvCountCycleDigits();
284 IncEvCountTotalDigits();
288 //____________________________________________________________________________
289 void AliADQADataMakerSim::MakeSDigits()
291 // makes data from SDigits
293 FillSDigitsData(0,fSDigitsArray->GetEntriesFast()) ;
294 TIter next(fSDigitsArray) ;
295 AliADSDigit *ADSDigit ;
296 while ( (ADSDigit = dynamic_cast<AliADSDigit *>(next())) ) {
297 Int_t PMNumber = ADSDigit->PMNumber();
298 Int_t Nbins = ADSDigit->GetNBins();
299 Double_t totCharge = 0;
301 for(Int_t i = 0; i<Nbins; i++)totCharge += ADSDigit->GetCharges()[i];
302 FillSDigitsData(1, PMNumber, totCharge) ;
306 //____________________________________________________________________________
307 void AliADQADataMakerSim::MakeDigits()
309 // makes data from Digits
311 FillDigitsData(0,fDigitsArray->GetEntriesFast()) ;
312 TIter next(fDigitsArray) ;
313 AliADdigit *ADDigit ;
314 Int_t nBBflagsADA = 0;
315 Int_t nBBflagsADC = 0;
317 while ( (ADDigit = dynamic_cast<AliADdigit *>(next())) ) {
319 Int_t PMNumber = ADDigit->PMNumber();
321 if(PMNumber<8 && ADDigit->GetBBflag()) nBBflagsADC++;
322 if(PMNumber>7 && ADDigit->GetBBflag()) nBBflagsADA++;
325 for(Int_t iClock=0; iClock<21; iClock++) {
326 adc[iClock]= ADDigit->ChargeADC(iClock);
327 FillDigitsData(3, PMNumber,(float)iClock-10,(float)adc[iClock]);
328 totCharge += adc[iClock];
331 FillDigitsData(1,PMNumber,ADDigit->Time());
332 FillDigitsData(2,PMNumber,ADDigit->Width());
333 FillDigitsData(7,PMNumber,totCharge);
334 FillDigitsData(8,PMNumber,TMath::LocMax(21,adc)-10);
337 FillDigitsData(4,nBBflagsADA+nBBflagsADC);
338 FillDigitsData(5,nBBflagsADA);
339 FillDigitsData(6,nBBflagsADC);
342 //____________________________________________________________________________
343 void AliADQADataMakerSim::MakeDigits(TTree *digitTree)
345 // makes data from Digit Tree
348 fDigitsArray->Clear() ;
350 fDigitsArray = new TClonesArray("AliADdigit", 1000) ;
352 TBranch * branch = digitTree->GetBranch("ADDigit") ;
354 AliWarning("AD branch in Digit Tree not found") ;
356 branch->SetAddress(&fDigitsArray) ;
357 branch->GetEntry(0) ;
361 IncEvCountCycleDigits();
362 IncEvCountTotalDigits();
367 //____________________________________________________________________________
368 void AliADQADataMakerSim::StartOfDetectorCycle()
370 //Detector specific actions at start of cycle