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) ;
93 //____________________________________________________________________________
94 void AliPHOSQADataMakerSim::InitDigits()
96 // create Digits histograms in Digits subdir
97 const Bool_t expert = kTRUE ;
98 const Bool_t image = kTRUE ;
99 TH1I * h0 = new TH1I("hPhosDigits", "Digits amplitude distribution in PHOS;Amplitude [ADC counts];Counts", 500, 0, 1000) ;
101 Add2DigitsList(h0, kDigits, !expert, image) ;
102 TH1I * h1 = new TH1I("hPhosDigitsMul", "Digits multiplicity distribution in PHOS;# of Digits;Entries", 2000, 0, 10000) ;
104 Add2DigitsList(h1, kDigitsMul, !expert, image) ;
107 //____________________________________________________________________________
108 void AliPHOSQADataMakerSim::InitSDigits()
110 // create SDigits histograms in SDigits subdir
111 const Bool_t expert = kTRUE ;
112 const Bool_t image = kTRUE ;
113 TH1F * h0 = new TH1F("hPhosSDigits", "SDigits energy distribution in PHOS; Energy [MeV];Counts", 500, 0., 1000.) ;
115 Add2SDigitsList(h0, kSDigits, !expert, image) ;
116 TH1I * h1 = new TH1I("hPhosSDigitsMul", "SDigits multiplicity distribution in PHOS;# of SDigits;Entries", 500, 0, 1000) ;
118 Add2SDigitsList(h1, kSDigitsMul, !expert, image) ;
121 //____________________________________________________________________________
122 void AliPHOSQADataMakerSim::MakeHits()
124 //make QA data from Hits
126 TIter next(fHitsArray) ;
128 while ( (hit = dynamic_cast<AliPHOSHit *>(next())) ) {
129 GetHitsData(kHits)->Fill( hit->GetEnergy()) ;
133 //____________________________________________________________________________
134 void AliPHOSQADataMakerSim::MakeHits(TTree * hitTree)
136 // make QA data from Hit Tree
139 fHitsArray->Clear() ;
141 fHitsArray = new TClonesArray("AliPHOSHit", 1000);
143 TBranch * branch = hitTree->GetBranch("PHOS") ;
145 AliWarning("PHOS branch in Hit Tree not found") ;
148 branch->SetAddress(&fHitsArray) ;
149 for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
150 branch->GetEntry(ientry) ;
151 nHits += fHitsArray->GetEntriesFast();
155 GetHitsData(1)->Fill(nHits) ;
159 //____________________________________________________________________________
160 void AliPHOSQADataMakerSim::MakeDigits()
162 // makes data from Digits
164 GetDigitsData(1)->Fill(fDigitsArray->GetEntriesFast()) ;
165 TIter next(fDigitsArray) ;
166 AliPHOSDigit * digit ;
167 while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
168 GetDigitsData(kDigits)->Fill( digit->GetEnergy()) ;
172 //____________________________________________________________________________
173 void AliPHOSQADataMakerSim::MakeDigits(TTree * digitTree)
175 // makes data from Digit Tree
177 fDigitsArray->Clear() ;
179 fDigitsArray = new TClonesArray("AliPHOSDigit", 1000) ;
181 TBranch * branch = digitTree->GetBranch("PHOS") ;
183 AliWarning("PHOS branch in Digit Tree not found") ;
185 branch->SetAddress(&fDigitsArray) ;
186 branch->GetEntry(0) ;
191 //____________________________________________________________________________
192 void AliPHOSQADataMakerSim::MakeSDigits()
194 // makes data from SDigits
196 GetSDigitsData(1)->Fill(fSDigitsArray->GetEntriesFast()) ;
197 TIter next(fSDigitsArray) ;
198 AliPHOSDigit * sdigit ;
199 while ( (sdigit = dynamic_cast<AliPHOSDigit *>(next())) ) {
200 GetSDigitsData(kSDigits)->Fill( sdigit->GetEnergy()) ;
204 //____________________________________________________________________________
205 void AliPHOSQADataMakerSim::MakeSDigits(TTree * sdigitTree)
207 // makes data from SDigit Tree
209 fSDigitsArray->Clear() ;
211 fSDigitsArray = new TClonesArray("AliPHOSDigit", 1000) ;
213 TBranch * branch = sdigitTree->GetBranch("PHOS") ;
215 AliWarning("PHOS branch in SDigit Tree not found") ;
217 branch->SetAddress(&fSDigitsArray) ;
218 branch->GetEntry(0) ;
223 //____________________________________________________________________________
224 void AliPHOSQADataMakerSim::StartOfDetectorCycle()
226 //Detector specific actions at start of cycle