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"),
51 fHits = new TClonesArray("AliPHOSHit", 1000);
54 //____________________________________________________________________________
55 AliPHOSQADataMakerSim::AliPHOSQADataMakerSim(const AliPHOSQADataMakerSim& qadm) :
60 SetName((const char*)qadm.GetName()) ;
61 SetTitle((const char*)qadm.GetTitle());
62 fHits = new TClonesArray("AliPHOSHit", 1000);
65 //__________________________________________________________________
66 AliPHOSQADataMakerSim& AliPHOSQADataMakerSim::operator = (const AliPHOSQADataMakerSim& qadm )
69 this->~AliPHOSQADataMakerSim();
70 new(this) AliPHOSQADataMakerSim(qadm);
74 //____________________________________________________________________________
75 void AliPHOSQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
77 //Detector specific actions at end of cycle
79 AliQAChecker::Instance()->Run(AliQAv1::kPHOS, task, list) ;
82 //____________________________________________________________________________
83 void AliPHOSQADataMakerSim::InitHits()
85 // create Hits histograms in Hits subdir
86 const Bool_t expert = kTRUE ;
87 const Bool_t image = kTRUE ;
88 TH1F * h0 = new TH1F("hPhosHits", "Hits energy distribution in PHOS;Energy [MeV];Counts", 100, 0., 100.) ;
90 Add2HitsList(h0, kHits, !expert, image) ;
91 TH1I * h1 = new TH1I("hPhosHitsMul", "Hits multiplicity distribution in PHOS;# of Hits;Entries", 500, 0., 10000) ;
93 Add2HitsList(h1, kHitsMul, !expert, image) ;
97 //____________________________________________________________________________
98 void AliPHOSQADataMakerSim::InitDigits()
100 // create Digits histograms in Digits subdir
101 const Bool_t expert = kTRUE ;
102 const Bool_t image = kTRUE ;
103 TH1I * h0 = new TH1I("hPhosDigits", "Digits amplitude distribution in PHOS;Amplitude [ADC counts];Counts", 500, 0, 1000) ;
105 Add2DigitsList(h0, kDigits, !expert, image) ;
106 TH1I * h1 = new TH1I("hPhosDigitsMul", "Digits multiplicity distribution in PHOS;# of Digits;Entries", 2000, 0, 10000) ;
108 Add2DigitsList(h1, kDigitsMul, !expert, image) ;
111 //____________________________________________________________________________
112 void AliPHOSQADataMakerSim::InitSDigits()
114 // create SDigits histograms in SDigits subdir
115 const Bool_t expert = kTRUE ;
116 const Bool_t image = kTRUE ;
117 TH1F * h0 = new TH1F("hPhosSDigits", "SDigits energy distribution in PHOS; Energy [MeV];Counts", 500, 0., 1000.) ;
119 Add2SDigitsList(h0, kSDigits, !expert, image) ;
120 TH1I * h1 = new TH1I("hPhosSDigitsMul", "SDigits multiplicity distribution in PHOS;# of SDigits;Entries", 500, 0, 1000) ;
122 Add2SDigitsList(h1, kSDigitsMul, !expert, image) ;
125 //____________________________________________________________________________
126 void AliPHOSQADataMakerSim::MakeHits()
128 //make QA data from Hits
130 // Check id histograms already created for this Event Specie
131 if ( ! GetHitsData(kHits) )
137 while ( (hit = dynamic_cast<AliPHOSHit *>(next())) ) {
138 GetHitsData(kHits)->Fill( hit->GetEnergy()) ;
142 //____________________________________________________________________________
143 void AliPHOSQADataMakerSim::MakeHits(TTree * hitTree)
145 // make QA data from Hit Tree
147 TBranch * branch = hitTree->GetBranch("PHOS") ;
149 AliWarning("PHOS branch in Hit Tree not found") ;
152 branch->SetAddress(&fHits) ;
153 for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
154 branch->GetEntry(ientry) ;
155 nHits += fHits->GetEntriesFast();
159 GetHitsData(1)->Fill(nHits) ;
163 //____________________________________________________________________________
164 void AliPHOSQADataMakerSim::MakeDigits(TClonesArray * digits)
166 // makes data from Digits
168 // Check id histograms already created for this Event Specie
169 if ( ! GetDigitsData(kDigits) )
172 GetDigitsData(1)->Fill(digits->GetEntriesFast()) ;
174 AliPHOSDigit * digit ;
175 while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
176 GetDigitsData(kDigits)->Fill( digit->GetEnergy()) ;
180 //____________________________________________________________________________
181 void AliPHOSQADataMakerSim::MakeDigits(TTree * digitTree)
183 // makes data from Digit Tree
184 TClonesArray * digits = new TClonesArray("AliPHOSDigit", 1000) ;
186 TBranch * branch = digitTree->GetBranch("PHOS") ;
188 AliWarning("PHOS branch in Digit Tree not found") ;
190 branch->SetAddress(&digits) ;
191 branch->GetEntry(0) ;
196 //____________________________________________________________________________
197 void AliPHOSQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
199 // makes data from SDigits
202 // Check id histograms already created for this Event Specie
203 if ( ! GetSDigitsData(kSDigits) )
206 GetSDigitsData(1)->Fill(sdigits->GetEntriesFast()) ;
207 TIter next(sdigits) ;
208 AliPHOSDigit * sdigit ;
209 while ( (sdigit = dynamic_cast<AliPHOSDigit *>(next())) ) {
210 GetSDigitsData(kSDigits)->Fill( sdigit->GetEnergy()) ;
214 //____________________________________________________________________________
215 void AliPHOSQADataMakerSim::MakeSDigits(TTree * sdigitTree)
217 // makes data from SDigit Tree
218 TClonesArray * sdigits = new TClonesArray("AliPHOSDigit", 1000) ;
220 TBranch * branch = sdigitTree->GetBranch("PHOS") ;
222 AliWarning("PHOS branch in SDigit Tree not found") ;
224 branch->SetAddress(&sdigits) ;
225 branch->GetEntry(0) ;
226 MakeSDigits(sdigits) ;
230 //____________________________________________________________________________
231 void AliPHOSQADataMakerSim::StartOfDetectorCycle()
233 //Detector specific actions at start of cycle