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(AliQA::GetDetName(AliQA::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
74 TH1F *h0 = new TH1F("hPreHitsEdep","Hits energy distribution in (keV)PRE(PMD)", 500, 0., 500.);
78 TH1F *h1 = new TH1F("hCpvHitsEdep","Hits energy distribution in (keV)CPV(PMD)", 500, 0., 500.);
82 TH1I *h2 = new TH1I("hPreHitsMult","Hits multiplicity distribution in PRE(PMD)", 500, 0, 3000) ;
86 TH1I *h3 = new TH1I("hCpvHitsMult","Hits multiplicity distribution in PRE(PMD)", 500, 0, 3000) ;
91 //____________________________________________________________________________
92 void AliPMDQADataMakerSim::InitSDigits()
94 // create SDigits histograms in SDigits subdir
96 TH1F *h0 = new TH1F("hPreSDigitsEdep","SDigits energy distribution in(keV) PRE(PMD)", 500, 0., 500.);
98 Add2SDigitsList(h0, 0);
100 TH1F *h1 = new TH1F("hCpvSDigitsEdep","SDigits energy distribution in (keV)CPV(PMD)", 500, 0., 500.);
102 Add2SDigitsList(h1, 1);
104 TH1I *h2 = new TH1I("hPreSDigitsMult","SDigits multiplicity distribution in PRE(PMD)", 500, 0., 1000.);
106 Add2SDigitsList(h2, 2);
108 TH1I *h3 = new TH1I("hCpvSDigitsMult","SDigits multiplicity distribution in CPV(PMD)", 500, 0., 1000.);
110 Add2SDigitsList(h3, 3);
114 //____________________________________________________________________________
115 void AliPMDQADataMakerSim::InitDigits()
117 // create Digits histograms in Digits subdir
119 TH1F *h0 = new TH1F("hPreDigitsEdep","Digits energy distribution in PRE(PMD)", 100, 0., 2000.);
121 Add2DigitsList(h0, 0);
123 TH1F *h1 = new TH1F("hCpvDigitsEdep","Digits energy distribution in CPV(PMD)", 100, 0., 2000.);
125 Add2DigitsList(h1, 1);
127 TH1I *h2 = new TH1I("hPreDigitsMult","Digits multiplicity distribution in PRE(PMD)", 500, 0, 1000) ;
129 Add2DigitsList(h2, 2);
131 TH1I *h3 = new TH1I("hCpvDigitsMult","Digits multiplicity distribution in CPV(PMD)", 500, 0, 1000);
133 Add2DigitsList(h3, 3);
137 //____________________________________________________________________________
138 void AliPMDQADataMakerSim::MakeHits(TClonesArray *hits)
140 //make QA data from Hits
142 Int_t premul = 0, cpvmul = 0;
143 Float_t edepkev = 0.;
147 while ( (hit = dynamic_cast<AliPMDhit *>(next())) )
149 if (hit->Z() > 361.5)
151 edepkev = (hit->GetEnergy())/1000.;
152 GetHitsData(0)->Fill(edepkev);
155 else if (hit->Z() < 361.5)
157 edepkev = (hit->GetEnergy())/1000.;
158 GetHitsData(1)->Fill(edepkev);
165 GetHitsData(2)->Fill(-1.);
169 GetHitsData(2)->Fill(premul);
174 GetHitsData(3)->Fill(-1.);
178 GetHitsData(3)->Fill(cpvmul);
183 //____________________________________________________________________________
184 void AliPMDQADataMakerSim::MakeHits(TTree * hitTree)
186 // make QA data from Hit Tree
188 TBranch * branch = hitTree->GetBranch("PMD") ;
191 AliWarning("PMD branch in Hit Tree not found") ;
195 static TClonesArray statichits("AliPMDhit", 1000);
197 TClonesArray *hits = &statichits;
198 static TClonesArray staticdummy("AliPMDhit", 1000);
200 TClonesArray *dummy = &staticdummy;
201 branch->SetAddress(&dummy);
204 for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
205 branch->GetEntry(ientry) ;
206 for (Int_t ihit = 0 ; ihit < dummy->GetEntries() ; ihit++) {
207 AliPMDhit * hit = dynamic_cast<AliPMDhit *> (dummy->At(ihit)) ;
208 new((*hits)[index]) AliPMDhit(*hit) ;
217 //____________________________________________________________________________
218 void AliPMDQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
220 // makes data from SDigits
222 Int_t cpvmul = 0, premul = 0;
223 Float_t edepkev = 0.;
225 TIter next(sdigits) ;
226 AliPMDsdigit * sdigit ;
227 while ( (sdigit = dynamic_cast<AliPMDsdigit *>(next())) )
229 if(sdigit->GetDetector() == 0)
231 edepkev = (sdigit->GetCellEdep())/1000.;
232 GetSDigitsData(0)->Fill(edepkev);
235 if(sdigit->GetDetector() == 1)
237 edepkev = (sdigit->GetCellEdep())/1000.;
238 GetSDigitsData(1)->Fill(edepkev);
243 if (premul > 0) GetSDigitsData(2)->Fill(premul);
244 if (cpvmul > 0) GetSDigitsData(3)->Fill(cpvmul);
248 //____________________________________________________________________________
249 void AliPMDQADataMakerSim::MakeSDigits(TTree * sdigitTree)
251 // makes data from SDigit Tree
253 TClonesArray * sdigits = new TClonesArray("AliPMDsdigit", 1000) ;
255 TBranch * branch = sdigitTree->GetBranch("PMDSDigit") ;
256 branch->SetAddress(&sdigits) ;
260 AliWarning("PMD branch in SDigit Tree not found") ;
264 for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
266 branch->GetEntry(ient) ;
267 MakeSDigits(sdigits) ;
272 //____________________________________________________________________________
273 void AliPMDQADataMakerSim::MakeDigits(TClonesArray * digits)
275 // makes data from Digits
277 Int_t cpvmul = 0, premul = 0;
280 AliPMDdigit * digit ;
281 while ( (digit = dynamic_cast<AliPMDdigit *>(next())) )
283 if(digit->GetDetector() == 0)
285 GetDigitsData(0)->Fill( digit->GetADC()) ;
288 if(digit->GetDetector() == 1)
290 GetDigitsData(1)->Fill( digit->GetADC());
295 if (premul > 0) GetDigitsData(2)->Fill(premul);
296 if (cpvmul > 0) GetDigitsData(3)->Fill(cpvmul);
301 //____________________________________________________________________________
302 void AliPMDQADataMakerSim::MakeDigits(TTree * digitTree)
304 // makes data from Digit Tree
306 TClonesArray * digits = new TClonesArray("AliPMDdigit", 1000) ;
308 TBranch * branch = digitTree->GetBranch("PMDDigit") ;
309 branch->SetAddress(&digits) ;
313 AliWarning("PMD branch in Digit Tree not found") ;
317 for (Int_t ient = 0; ient < branch->GetEntries(); ient++)
319 branch->GetEntry(ient) ;
327 //____________________________________________________________________________
329 void AliPMDQADataMakerSim::StartOfDetectorCycle()
331 //Detector specific actions at start of cycle
334 //____________________________________________________________________________
336 void AliPMDQADataMakerSim::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
338 //Detector specific actions at end of cycle
339 // do the QA checking
340 AliQAChecker::Instance()->Run(AliQA::kPMD, task, list) ;