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: AliVZEROQADataMakerSim.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 "AliVZEROdigit.h"
36 #include "AliVZEROhit.h"
37 #include "AliVZEROQADataMakerSim.h"
38 #include "AliQAChecker.h"
40 ClassImp(AliVZEROQADataMakerSim)
42 //____________________________________________________________________________
43 AliVZEROQADataMakerSim::AliVZEROQADataMakerSim() :
44 AliQADataMakerSim(AliQA::GetDetName(AliQA::kVZERO), "VZERO Quality Assurance Data Maker")
52 //____________________________________________________________________________
53 AliVZEROQADataMakerSim::AliVZEROQADataMakerSim(const AliVZEROQADataMakerSim& qadm) :
58 SetName((const char*)qadm.GetName()) ;
59 SetTitle((const char*)qadm.GetTitle());
62 //__________________________________________________________________
63 AliVZEROQADataMakerSim& AliVZEROQADataMakerSim::operator = (const AliVZEROQADataMakerSim& qadm )
66 this->~AliVZEROQADataMakerSim();
67 new(this) AliVZEROQADataMakerSim(qadm);
70 //____________________________________________________________________________
71 void AliVZEROQADataMakerSim::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray ** list)
73 //Detector specific actions at end of cycle
75 AliQAChecker::Instance()->Run(AliQA::kVZERO, task, list) ;
79 //____________________________________________________________________________
80 void AliVZEROQADataMakerSim::InitHits()
83 // create Hits histograms in Hits subdir
85 TH1I * h0 = new TH1I("hHitMultiplicity", "Hit multiplicity distribution in VZERO", 300, 0, 299) ;
89 TH1I * h1 = new TH1I("hHitCellNumber", "Hit cell distribution in VZERO", 80, 0, 79) ;
95 //____________________________________________________________________________
96 void AliVZEROQADataMakerSim::InitDigits()
98 // create Digits histograms in Digits subdir
107 // create Digits histograms in Digits subdir
108 TH1I * h0 = new TH1I("hDigitMultiplicity", "Digits multiplicity distribution in VZERO", 100, 0, 99) ;
110 Add2DigitsList(h0, 0) ;
112 for (Int_t i=0; i<64; i++)
114 sprintf(TDCname, "hDigitTDC%d", i);
115 sprintf(texte,"Digit TDC in cell %d",i);
116 fhDigTDC[i] = new TH1I(TDCname,texte,300,0.,149.);
118 sprintf(ADCname,"hDigitADC%d",i);
119 sprintf(texte,"Digit ADC in cell %d",i);
120 fhDigADC[i]= new TH1I(ADCname,texte,1024,0.,1023.);
122 Add2DigitsList(fhDigTDC[i],i+1);
123 Add2DigitsList(fhDigADC[i],i+1+64);
128 //____________________________________________________________________________
129 void AliVZEROQADataMakerSim::MakeHits(TClonesArray * hits)
131 //make QA data from Hits
133 GetHitsData(0)->Fill(hits->GetEntriesFast()) ; // fills Hit multiplicity
134 Int_t nhits = hits->GetEntriesFast();
135 for (Int_t ihit=0;ihit<nhits;ihit++)
137 AliVZEROhit * VZEROHit = (AliVZEROhit*) hits->UncheckedAt(ihit);
139 AliError("The unchecked hit doesn't exist");
142 GetHitsData(1)->Fill(VZEROHit->Cell());
147 //____________________________________________________________________________
149 void AliVZEROQADataMakerSim::MakeHits(TTree *hitTree)
151 //fills QA histos for Hits
152 TClonesArray * hits = new TClonesArray("AliVZEROhit", 1000);
154 TBranch * branch = hitTree->GetBranch("VZERO") ;
156 AliWarning("VZERO branch in Hit Tree not found") ;
160 branch->SetAddress(&hits);
162 AliError("Branch VZERO hit not found");
165 Int_t ntracks = (Int_t) hitTree->GetEntries();
167 if (ntracks<=0) return;
168 // Start loop on tracks in the hits containers
169 for (Int_t track=0; track<ntracks;track++) {
170 branch->GetEntry(track);
171 GetHitsData(0)->Fill(hits->GetEntriesFast()) ; // fills Hit multiplicity
172 Int_t nhits = hits->GetEntriesFast();
173 for (Int_t ihit=0;ihit<nhits;ihit++)
175 AliVZEROhit * VZEROHit = (AliVZEROhit*) hits->UncheckedAt(ihit);
177 AliError("The unchecked hit doesn't exist");
180 GetHitsData(1)->Fill(VZEROHit->Cell());
187 //____________________________________________________________________________
188 void AliVZEROQADataMakerSim::MakeDigits(TClonesArray * digits)
190 // makes data from Digits
192 GetDigitsData(0)->Fill(digits->GetEntriesFast()) ;
194 AliVZEROdigit *VZERODigit ;
195 while ( (VZERODigit = dynamic_cast<AliVZEROdigit *>(next())) ) {
196 Int_t PMNumber = VZERODigit->PMNumber();
197 GetDigitsData(PMNumber +1)->Fill( VZERODigit->Time()) ; // in 100 of picoseconds
198 GetDigitsData(PMNumber +1+64)->Fill( VZERODigit->ADC()) ;
203 //____________________________________________________________________________
204 void AliVZEROQADataMakerSim::MakeDigits(TTree *digitTree)
206 // makes data from Digit Tree
208 TClonesArray * digits = new TClonesArray("AliVZEROdigit", 1000) ;
210 TBranch * branch = digitTree->GetBranch("VZERODigit") ;
212 AliWarning("VZERO branch in Digit Tree not found") ;
214 branch->SetAddress(&digits) ;
215 branch->GetEntry(0) ;
221 //____________________________________________________________________________
222 void AliVZEROQADataMakerSim::StartOfDetectorCycle()
224 //Detector specific actions at start of cycle