1 /**************************************************************************
2 * Copyright(c) 2004, 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 **************************************************************************/
15 // --- ROOT system ---
17 #include <TClonesArray.h>
22 // --- AliRoot header files ---
23 #include "AliESDEvent.h"
25 #include "AliFMDQADataMakerRec.h"
26 #include "AliFMDDigit.h"
27 #include "AliFMDRecPoint.h"
28 #include "AliQAChecker.h"
29 #include "AliESDFMD.h"
30 #include "AliFMDParameters.h"
32 //_____________________________________________________________________
33 // This is the class that collects the QA data for the FMD during
36 // The following data types are picked up:
40 // The following data types are not supported (yet):
42 // Author : Hans Hjersing Dalsgaard, hans.dalsgaard@cern.ch
43 //_____________________________________________________________________
45 ClassImp(AliFMDQADataMakerRec)
47 ; // For Emacs - do not delete!
50 //_____________________________________________________________________
51 AliFMDQADataMakerRec::AliFMDQADataMakerRec() :
52 AliQADataMakerRec(AliQA::GetDetName(AliQA::kFMD),
53 "FMD Quality Assurance Data Maker")
58 //_____________________________________________________________________
59 AliFMDQADataMakerRec::AliFMDQADataMakerRec(const AliFMDQADataMakerRec& /*qadm*/)
64 // qadm Object to copy from
69 //_____________________________________________________________________
72 AliFMDQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task,
75 // Detector specific actions at end of cycle
77 AliLog::Message(5,"FMD: end of detector cycle",
78 "AliFMDQADataMakerRec","AliFMDQADataMakerRec",
79 "AliFMDQADataMakerRec::EndOfDetectorCycle",
80 "AliFMDQADataMakerRec.cxx",95);
81 AliQAChecker::Instance()->Run(AliQA::kFMD, task, list);
84 //_____________________________________________________________________
85 void AliFMDQADataMakerRec::InitESDs()
87 // create Digits histograms in Digits subdir
88 TH1F* hEnergyOfESD = new TH1F("hEnergyOfESD","Energy distribution",100,0,3);
89 hEnergyOfESD->SetXTitle("Edep/Emip");
90 hEnergyOfESD->SetYTitle("Counts");
91 Add2ESDsList(hEnergyOfESD, 0);
95 //_____________________________________________________________________
96 void AliFMDQADataMakerRec::InitDigits()
98 // create Digits histograms in Digits subdir
99 TH1I* hADCCounts = new TH1I("hADCCounts","Dist of ADC counts",
101 hADCCounts->SetXTitle("ADC counts");
102 hADCCounts->SetYTitle("");
103 Add2DigitsList(hADCCounts, 0);
107 //_____________________________________________________________________
108 void AliFMDQADataMakerRec::InitRecPoints()
110 TH1F* hEnergyOfRecpoints = new TH1F("hEnergyOfRecpoints",
111 "Energy Distribution",100,0,3);
112 hEnergyOfRecpoints->SetXTitle("Edep/Emip");
113 hEnergyOfRecpoints->SetYTitle("");
114 Add2RecPointsList(hEnergyOfRecpoints,0);
117 //_____________________________________________________________________
118 void AliFMDQADataMakerRec::InitRaws()
123 //_____________________________________________________________________
124 void AliFMDQADataMakerRec::MakeESDs(AliESDEvent * esd)
127 AliError("FMD ESD object not found!!") ;
130 AliESDFMD* fmd = esd->GetFMDData();
133 for(UShort_t det=1;det<=3;det++) {
134 for (UShort_t ir = 0; ir < 2; ir++) {
135 Char_t ring = (ir == 0 ? 'I' : 'O');
136 UShort_t nsec = (ir == 0 ? 20 : 40);
137 UShort_t nstr = (ir == 0 ? 512 : 256);
138 for(UShort_t sec =0; sec < nsec; sec++) {
139 for(UShort_t strip = 0; strip < nstr; strip++) {
140 Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
141 if(mult == AliESDFMD::kInvalidMult) continue;
143 GetESDsData(0)->Fill(mult);
150 //_____________________________________________________________________
151 void AliFMDQADataMakerRec::MakeDigits(TClonesArray * digits)
153 // makes data from Digits
155 AliError("FMD Digit object not found!!") ;
158 for(Int_t i=0;i<digits->GetEntries();i++) {
160 AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
161 GetDigitsData(0)->Fill(digit->Counts());
165 //_____________________________________________________________________
166 void AliFMDQADataMakerRec::MakeDigits(TTree * digitTree)
169 TClonesArray* digits = new TClonesArray("AliFMDDigit", 1000);
170 TBranch* branch = digitTree->GetBranch("FMD");
172 AliWarning("FMD branch in Digit Tree not found") ;
176 branch->SetAddress(&digits);
181 //_____________________________________________________________________
182 void AliFMDQADataMakerRec::MakeRaws(AliRawReader* /*rawReader*/)
187 //_____________________________________________________________________
188 void AliFMDQADataMakerRec::MakeRecPoints(TTree* clustersTree)
190 // makes data from RecPoints
191 TBranch *fmdbranch = clustersTree->GetBranch("FMD");
193 AliError("can't get the branch with the FMD recpoints !");
197 TClonesArray * fmdrecpoints = new TClonesArray("AliFMDRecPoint", 1000);
198 fmdbranch->SetAddress(&fmdrecpoints);
199 fmdbranch->GetEntry(0);
201 TIter next(fmdrecpoints) ;
202 AliFMDRecPoint * rp ;
203 while ((rp = static_cast<AliFMDRecPoint*>(next()))) {
204 GetRecPointsData(0)->Fill(rp->Particles()) ;
206 fmdrecpoints->Delete();
211 //_____________________________________________________________________
212 void AliFMDQADataMakerRec::StartOfDetectorCycle()
218 //_____________________________________________________________________