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 **************************************************************************/
18 Produces the data needed to calculate the quality assurance.
19 All data must be mergeable objects.
23 // --- ROOT system ---
24 #include <TClonesArray.h>
31 // --- Standard library ---
33 // --- AliRoot header files ---
36 #include "AliPMDhit.h"
37 #include "AliPMDsdigit.h"
38 #include "AliPMDdigit.h"
39 #include "AliPMDQADataMakerSim.h"
40 #include "AliQAChecker.h"
42 ClassImp(AliPMDQADataMakerSim)
44 //____________________________________________________________________________
45 AliPMDQADataMakerSim::AliPMDQADataMakerSim() :
46 AliQADataMakerSim(AliQAv1::GetDetName(AliQAv1::kPMD), "PMD Quality Assurance Data Maker")
51 //____________________________________________________________________________
52 AliPMDQADataMakerSim::AliPMDQADataMakerSim(const AliPMDQADataMakerSim& qadm) :
56 SetName((const char*)qadm.GetName()) ;
57 SetTitle((const char*)qadm.GetTitle());
60 //__________________________________________________________________
61 AliPMDQADataMakerSim& AliPMDQADataMakerSim::operator = (const AliPMDQADataMakerSim& qadm )
64 this->~AliPMDQADataMakerSim();
65 new(this) AliPMDQADataMakerSim(qadm);
69 //____________________________________________________________________________
70 void AliPMDQADataMakerSim::InitHits()
72 // create Hits histograms in Hits subdir
73 const Bool_t expert = kTRUE ;
74 const Bool_t image = kTRUE ;
76 TH1F *h0 = new TH1F("hPreHitsEdep","Hits energy distribution PRE(PMD);Energy [keV];Counts", 500, 0., 500.);
78 Add2HitsList(h0, 0, !expert, image) ;
80 TH1F *h1 = new TH1F("hCpvHitsEdep","Hits energy distribution CPV(PMD);Energy [keV];Counts", 500, 0., 500.);
82 Add2HitsList(h1, 1, !expert, image) ;
84 TH1I *h2 = new TH1I("hPreHitsMult","Hits multiplicity distribution in PRE(PMD);# of Hits;Entries", 500, 0, 3000) ;
86 Add2HitsList(h2, 2, !expert, image) ;
88 TH1I *h3 = new TH1I("hCpvHitsMult","Hits multiplicity distribution in PRE(PMD);# of Hits;Entries", 500, 0, 3000) ;
90 Add2HitsList(h3, 3, !expert, image) ;
93 //____________________________________________________________________________
94 void AliPMDQADataMakerSim::InitSDigits()
96 // create SDigits histograms in SDigits subdir
97 const Bool_t expert = kTRUE ;
98 const Bool_t image = kTRUE ;
100 TH1F *h0 = new TH1F("hPreSDigitsEdep","SDigits energy distribution in(keV) PRE(PMD);Energy [keV];Counts", 500, 0., 500.);
102 Add2SDigitsList(h0, 0, !expert, image);
104 TH1F *h1 = new TH1F("hCpvSDigitsEdep","SDigits energy distribution in (keV)CPV(PMD);Energy [keV];Counts", 500, 0., 500.);
106 Add2SDigitsList(h1, 1, !expert, image);
108 TH1I *h2 = new TH1I("hPreSDigitsMult","SDigits multiplicity distribution in PRE(PMD);# of SDigits;Entries", 500, 0., 1000.);
110 Add2SDigitsList(h2, 2, !expert, image);
112 TH1I *h3 = new TH1I("hCpvSDigitsMult","SDigits multiplicity distribution in CPV(PMD);# of SDigits;Entries", 500, 0., 1000.);
114 Add2SDigitsList(h3, 3, !expert, image);
118 //____________________________________________________________________________
119 void AliPMDQADataMakerSim::InitDigits()
121 // create Digits histograms in Digits subdir
122 const Bool_t expert = kTRUE ;
123 const Bool_t image = kTRUE ;
125 TH1F *h0 = new TH1F("hPreDigitsEdep","Digits energy distribution in PRE(PMD);Amplitude [ADC counts];Counts", 100, 0., 2000.);
127 Add2DigitsList(h0, 0, !expert, image);
129 TH1F *h1 = new TH1F("hCpvDigitsEdep","Digits energy distribution in CPV(PMD);Amplitude [ADC counts];Counts", 100, 0., 2000.);
131 Add2DigitsList(h1, 1, !expert, image);
133 TH1I *h2 = new TH1I("hPreDigitsMult","Digits multiplicity distribution in PRE(PMD);# of Digits;Entries", 500, 0, 1000) ;
135 Add2DigitsList(h2, 2, !expert, image);
137 TH1I *h3 = new TH1I("hCpvDigitsMult","Digits multiplicity distribution in CPV(PMD);# of Digits;Entries", 500, 0, 1000);
139 Add2DigitsList(h3, 3, !expert, image);
143 //____________________________________________________________________________
144 void AliPMDQADataMakerSim::MakeHits()
146 //make QA data from Hits
148 Int_t premul = 0, cpvmul = 0;
149 Float_t edepkev = 0.;
150 TIter next(fHitsArray);
153 while ( (hit = dynamic_cast<AliPMDhit *>(next())) )
155 if (hit->Z() > 361.5)
157 edepkev = (hit->GetEnergy())/1000.;
158 GetHitsData(0)->Fill(edepkev);
161 else if (hit->Z() < 361.5)
163 edepkev = (hit->GetEnergy())/1000.;
164 GetHitsData(1)->Fill(edepkev);
171 GetHitsData(2)->Fill(-1.);
175 GetHitsData(2)->Fill(premul);
180 GetHitsData(3)->Fill(-1.);
184 GetHitsData(3)->Fill(cpvmul);
189 //____________________________________________________________________________
190 void AliPMDQADataMakerSim::MakeHits(TTree * hitTree)
192 // make QA data from Hit Tree
194 TBranch * branch = hitTree->GetBranch("PMD") ;
197 AliWarning("PMD branch in Hit Tree not found") ;
202 fHitsArray->Clear() ;
204 fHitsArray = new TClonesArray("AliPMDhit", 1000);
206 branch->SetAddress(&fHitsArray);
208 for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
209 branch->GetEntry(ientry) ;
211 fHitsArray->Clear() ;
214 //____________________________________________________________________________
215 void AliPMDQADataMakerSim::MakeSDigits()
217 // makes data from SDigits
219 Int_t cpvmul = 0, premul = 0;
220 Float_t edepkev = 0.;
222 TIter next(fSDigitsArray) ;
223 AliPMDsdigit * sdigit ;
224 while ( (sdigit = dynamic_cast<AliPMDsdigit *>(next())) )
226 if(sdigit->GetDetector() == 0)
228 edepkev = (sdigit->GetCellEdep())/1000.;
229 GetSDigitsData(0)->Fill(edepkev);
232 if(sdigit->GetDetector() == 1)
234 edepkev = (sdigit->GetCellEdep())/1000.;
235 GetSDigitsData(1)->Fill(edepkev);
240 if (premul > 0) GetSDigitsData(2)->Fill(premul);
241 if (cpvmul > 0) GetSDigitsData(3)->Fill(cpvmul);
245 //____________________________________________________________________________
246 void AliPMDQADataMakerSim::MakeSDigits(TTree * sdigitTree)
248 // makes data from SDigit Tree
251 fSDigitsArray->Clear() ;
253 fSDigitsArray = new TClonesArray("AliPMDsdigit", 1000) ;
255 TBranch * branch = sdigitTree->GetBranch("PMDSDigit") ;
256 branch->SetAddress(&fSDigitsArray) ;
259 AliWarning("PMD branch in SDigit Tree not found") ;
261 branch->GetEntry(0) ;
266 //____________________________________________________________________________
267 void AliPMDQADataMakerSim::MakeDigits()
269 // makes data from Digits
271 Int_t cpvmul = 0, premul = 0;
273 TIter next(fDigitsArray) ;
274 AliPMDdigit * digit ;
275 while ( (digit = dynamic_cast<AliPMDdigit *>(next())) )
277 if(digit->GetDetector() == 0)
279 GetDigitsData(0)->Fill( digit->GetADC()) ;
282 if(digit->GetDetector() == 1)
284 GetDigitsData(1)->Fill( digit->GetADC());
289 if (premul > 0) GetDigitsData(2)->Fill(premul);
290 if (cpvmul > 0) GetDigitsData(3)->Fill(cpvmul);
295 //____________________________________________________________________________
296 void AliPMDQADataMakerSim::MakeDigits(TTree * digitTree)
298 // makes data from Digit Tree
301 fDigitsArray->Clear() ;
303 fDigitsArray = new TClonesArray("AliPMDdigit", 1000) ;
305 TBranch * branch = digitTree->GetBranch("PMDDigit") ;
306 branch->SetAddress(&fDigitsArray) ;
310 AliWarning("PMD branch in Digit Tree not found") ;
314 for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
316 branch->GetEntry(ient) ;
318 fDigitsArray->Clear() ;
326 //____________________________________________________________________________
328 void AliPMDQADataMakerSim::StartOfDetectorCycle()
330 //Detector specific actions at start of cycle
333 //____________________________________________________________________________
335 void AliPMDQADataMakerSim::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
337 //Detector specific actions at end of cycle
338 // do the QA checking
339 AliQAChecker::Instance()->Run(AliQAv1::kPMD, task, list) ;