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 **************************************************************************/
20 Produces the data needed to calculate the quality assurance.
21 All data must be mergeable objects.
22 Y. Schutz CERN July 2007
25 // --- ROOT system ---
26 #include <TClonesArray.h>
33 // --- Standard library ---
35 // --- AliRoot header files ---
36 #include "AliESDCaloCluster.h"
38 #include "AliPHOSDigit.h"
39 #include "AliPHOSHit.h"
40 #include "AliPHOSQADataMakerSim.h"
41 #include "AliQAChecker.h"
43 ClassImp(AliPHOSQADataMakerSim)
45 //____________________________________________________________________________
46 AliPHOSQADataMakerSim::AliPHOSQADataMakerSim() :
47 AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kPHOS), "PHOS Quality Assurance Data Maker")
52 //____________________________________________________________________________
53 AliPHOSQADataMakerSim::AliPHOSQADataMakerSim(const AliPHOSQADataMakerSim& qadm) :
57 SetName((const char*)qadm.GetName()) ;
58 SetTitle((const char*)qadm.GetTitle());
61 //__________________________________________________________________
62 AliPHOSQADataMakerSim& AliPHOSQADataMakerSim::operator = (const AliPHOSQADataMakerSim& qadm )
65 this->~AliPHOSQADataMakerSim();
66 new(this) AliPHOSQADataMakerSim(qadm);
70 //____________________________________________________________________________
71 void AliPHOSQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
73 //Detector specific actions at end of cycle
75 AliQAChecker::Instance()->Run(AliQAv1::kPHOS, task, list) ;
78 //____________________________________________________________________________
79 void AliPHOSQADataMakerSim::InitHits()
81 // create Hits histograms in Hits subdir
82 const Bool_t expert = kTRUE ;
83 const Bool_t image = kTRUE ;
84 TH1F * h0 = new TH1F("hPhosHits", "Hits energy distribution in PHOS;Energy [MeV];Counts", 100, 0., 100.) ;
86 Add2HitsList(h0, kHits, !expert, image) ;
87 TH1I * h1 = new TH1I("hPhosHitsMul", "Hits multiplicity distribution in PHOS;# of Hits;Entries", 500, 0., 10000) ;
89 Add2HitsList(h1, kHitsMul, !expert, image) ;
91 ClonePerTrigClass(AliQAv1::kHITS); // this should be the last line
94 //____________________________________________________________________________
95 void AliPHOSQADataMakerSim::InitDigits()
97 // create Digits histograms in Digits subdir
98 const Bool_t expert = kTRUE ;
99 const Bool_t image = kTRUE ;
100 TH1I * h0 = new TH1I("hPhosDigits", "Digits amplitude distribution in PHOS;Amplitude [ADC counts];Counts", 500, 0, 1000) ;
102 Add2DigitsList(h0, kDigits, !expert, image) ;
103 TH1I * h1 = new TH1I("hPhosDigitsMul", "Digits multiplicity distribution in PHOS;# of Digits;Entries", 2000, 0, 10000) ;
105 Add2DigitsList(h1, kDigitsMul, !expert, image) ;
107 ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
110 //____________________________________________________________________________
111 void AliPHOSQADataMakerSim::InitSDigits()
113 // create SDigits histograms in SDigits subdir
114 const Bool_t expert = kTRUE ;
115 const Bool_t image = kTRUE ;
116 TH1F * h0 = new TH1F("hPhosSDigits", "SDigits energy distribution in PHOS; Energy [MeV];Counts", 500, 0., 1000.) ;
118 Add2SDigitsList(h0, kSDigits, !expert, image) ;
119 TH1I * h1 = new TH1I("hPhosSDigitsMul", "SDigits multiplicity distribution in PHOS;# of SDigits;Entries", 500, 0, 1000) ;
121 Add2SDigitsList(h1, kSDigitsMul, !expert, image) ;
123 ClonePerTrigClass(AliQAv1::kSDIGITS); // this should be the last line
126 //____________________________________________________________________________
127 void AliPHOSQADataMakerSim::MakeHits()
129 //make QA data from Hits
131 TIter next(fHitsArray) ;
133 while ( (hit = dynamic_cast<AliPHOSHit *>(next())) ) {
134 FillHitsData(kHits, hit->GetEnergy()) ;
138 //____________________________________________________________________________
139 void AliPHOSQADataMakerSim::MakeHits(TTree * hitTree)
141 // make QA data from Hit Tree
144 fHitsArray->Clear() ;
146 fHitsArray = new TClonesArray("AliPHOSHit", 1000);
148 TBranch * branch = hitTree->GetBranch("PHOS") ;
149 if ( ! branch ) { AliWarning("PHOS branch in Hit Tree not found"); return;}
152 branch->SetAddress(&fHitsArray) ;
153 for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
154 branch->GetEntry(ientry) ;
155 nHits += fHitsArray->GetEntriesFast();
159 FillHitsData(1,nHits) ;
161 IncEvCountCycleHits();
162 IncEvCountTotalHits();
166 //____________________________________________________________________________
167 void AliPHOSQADataMakerSim::MakeDigits()
169 // makes data from Digits
171 FillDigitsData(1,fDigitsArray->GetEntriesFast()) ;
172 TIter next(fDigitsArray) ;
173 AliPHOSDigit * digit ;
174 while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
175 FillDigitsData(kDigits, digit->GetEnergy()) ;
179 //____________________________________________________________________________
180 void AliPHOSQADataMakerSim::MakeDigits(TTree * digitTree)
182 // makes data from Digit Tree
184 fDigitsArray->Clear() ;
186 fDigitsArray = new TClonesArray("AliPHOSDigit", 1000) ;
188 TBranch * branch = digitTree->GetBranch("PHOS") ;
189 if ( ! branch ) {AliWarning("PHOS branch in Digit Tree not found"); return;}
190 branch->SetAddress(&fDigitsArray) ;
191 branch->GetEntry(0) ;
194 IncEvCountCycleDigits();
195 IncEvCountTotalDigits();
199 //____________________________________________________________________________
200 void AliPHOSQADataMakerSim::MakeSDigits()
202 // makes data from SDigits
204 FillSDigitsData(1,fSDigitsArray->GetEntriesFast()) ;
205 TIter next(fSDigitsArray) ;
206 AliPHOSDigit * sdigit ;
207 while ( (sdigit = dynamic_cast<AliPHOSDigit *>(next())) ) {
208 FillSDigitsData(kSDigits, sdigit->GetEnergy()) ;
212 //____________________________________________________________________________
213 void AliPHOSQADataMakerSim::MakeSDigits(TTree * sdigitTree)
215 // makes data from SDigit Tree
217 fSDigitsArray->Clear() ;
219 fSDigitsArray = new TClonesArray("AliPHOSDigit", 1000) ;
221 TBranch * branch = sdigitTree->GetBranch("PHOS") ;
222 if ( ! branch ) {AliWarning("PHOS branch in SDigit Tree not found"); return;}
223 branch->SetAddress(&fSDigitsArray) ;
224 branch->GetEntry(0) ;
227 IncEvCountCycleSDigits();
228 IncEvCountTotalSDigits();
232 //____________________________________________________________________________
233 void AliPHOSQADataMakerSim::StartOfDetectorCycle()
235 //Detector specific actions at start of cycle