]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDQADataMakerRec.cxx
A small change in histogram names and a tiny bug fix
[u/mrichter/AliRoot.git] / FMD / AliFMDQADataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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 ---
16 #include <iostream>
17 #include <TClonesArray.h>
18 #include <TFile.h> 
19 #include <TH1F.h> 
20 #include <TH1I.h> 
21
22 // --- AliRoot header files ---
23 #include "AliESDEvent.h"
24 #include "AliLog.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"
31 #include "AliFMDRawReader.h"
32 #include "AliRawReader.h"
33 #include "AliFMDAltroMapping.h"
34 #include "AliFMDDebug.h"
35
36 //_____________________________________________________________________
37 // This is the class that collects the QA data for the FMD during
38 // reconstruction.  
39 //
40 // The following data types are picked up:
41 // - rec points
42 // - esd data
43 // - raws
44 // Author : Hans Hjersing Dalsgaard, hans.dalsgaard@cern.ch
45 //_____________________________________________________________________
46
47 ClassImp(AliFMDQADataMakerRec)
48 #if 0
49 ; // For Emacs - do not delete!
50 #endif
51            
52 //_____________________________________________________________________
53 AliFMDQADataMakerRec::AliFMDQADataMakerRec() : 
54   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kFMD), 
55                     "FMD Quality Assurance Data Maker"),
56   fRecPointsArray("AliFMDRecPoint", 1000)
57 {
58   // ctor
59  
60 }
61
62 //_____________________________________________________________________
63 AliFMDQADataMakerRec::AliFMDQADataMakerRec(const AliFMDQADataMakerRec& qadm) 
64   : AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kFMD), 
65                       "FMD Quality Assurance Data Maker"),
66     fRecPointsArray(qadm.fRecPointsArray)
67 {
68   // copy ctor 
69   // Parameters: 
70   //    qadm    Object to copy from
71   
72 }
73 //_____________________________________________________________________
74 AliFMDQADataMakerRec& AliFMDQADataMakerRec::operator = (const AliFMDQADataMakerRec& qadm ) 
75 {
76   fRecPointsArray = qadm.fRecPointsArray;
77   
78   return *this;
79 }
80 //_____________________________________________________________________
81 AliFMDQADataMakerRec::~AliFMDQADataMakerRec()
82 {
83  
84 }
85
86
87 //_____________________________________________________________________ 
88
89 void 
90 AliFMDQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, 
91                                          TObjArray ** list)
92 {
93   // Detector specific actions at end of cycle
94   // do the QA checking
95   AliLog::Message(5,"FMD: end of detector cycle",
96                   "AliFMDQADataMakerRec","AliFMDQADataMakerRec",
97                   "AliFMDQADataMakerRec::EndOfDetectorCycle",
98                   "AliFMDQADataMakerRec.cxx",95);
99   AliQAChecker::Instance()->Run(AliQAv1::kFMD, task, list);
100 }
101
102 //_____________________________________________________________________ 
103 void AliFMDQADataMakerRec::InitESDs()
104 {
105   // create Digits histograms in Digits subdir
106   const Bool_t expert   = kTRUE ; 
107   const Bool_t image    = kTRUE ; 
108   
109   TH1F* hEnergyOfESD = new TH1F("hEnergyOfESD","Energy distribution",100,0,3);
110   hEnergyOfESD->SetXTitle("Edep/Emip");
111   hEnergyOfESD->SetYTitle("Counts");
112   Add2ESDsList(hEnergyOfESD, 0, !expert, image);
113     
114 }
115
116 //_____________________________________________________________________
117 void AliFMDQADataMakerRec::InitDigits()
118 {
119   // create Digits histograms in Digits subdir
120   const Bool_t expert   = kTRUE ; 
121   const Bool_t image    = kTRUE ; 
122   
123   TH1I* hADCCounts = new TH1I("hADCCounts","Dist of ADC counts;ADC counts;Counts",1024,0,1024);
124   Add2DigitsList(hADCCounts, 0, !expert, image);
125 }
126
127
128 //_____________________________________________________________________ 
129 void AliFMDQADataMakerRec::InitRecPoints()
130 {
131   // create Reconstructed Points histograms in RecPoints subdir
132   const Bool_t expert   = kTRUE ; 
133   const Bool_t image    = kTRUE ; 
134
135   TH1F* hEnergyOfRecpoints = new TH1F("hEnergyOfRecpoints",
136                                       "Energy Distribution",100,0,3);
137   hEnergyOfRecpoints->SetXTitle("Edep/Emip");
138   hEnergyOfRecpoints->SetYTitle("Counts");
139   Add2RecPointsList(hEnergyOfRecpoints,0, !expert, image);
140 }
141
142 //_____________________________________________________________________ 
143 void AliFMDQADataMakerRec::InitRaws()
144 {
145   // create Raws histograms in Raws subdir  
146   const Bool_t expert   = kTRUE ; 
147   const Bool_t saveCorr = kTRUE ; 
148   const Bool_t image    = kTRUE ; 
149
150   TH1I* hADCCounts;
151   for(Int_t det = 1; det<=3; det++) {
152     Int_t firstring = (det==1 ? 1 : 0);
153     for(Int_t iring = firstring;iring<=1;iring++) {
154       Char_t ring = (iring == 1 ? 'I' : 'O');
155       hADCCounts      = new TH1I(Form("hADCCounts_FMD%d%c",
156                                       det, ring), "ADC counts;Amplitude [ADC counts];Counts",
157                                  1024,0,1023);
158       
159       Int_t index1 = GetHalfringIndex(det, ring, 0,1);
160       Add2RawsList(hADCCounts, index1, expert, !image, !saveCorr);
161       
162       for(Int_t b = 0; b<=1;b++) {
163         
164         //Hexadecimal board numbers 0x0, 0x1, 0x10, 0x11;
165         UInt_t board = (iring == 1 ? 0 : 1);
166         board = board + b*16;
167         
168         
169         hADCCounts      = new TH1I(Form("hADCCounts_FMD%d%c_board%d",
170                                         det, ring, board), "ADC counts;Amplitude [ADC counts];Counts",
171                                    1024,0,1023);
172         hADCCounts->SetXTitle("ADC counts");
173         hADCCounts->SetYTitle("");
174         Int_t index2 = GetHalfringIndex(det, ring, board/16,0);
175         Add2RawsList(hADCCounts, index2, !expert, image, !saveCorr);
176
177       }
178     }
179   }
180 }
181
182 #if 0
183 struct FillESDHist : public AliESDFMD::ForOne
184 {
185   FillESDHist(AliFMDQADataMakerRec* m) : fM(m) {}
186   FillESDHist(const FillESDHist& o) : fM(o.fM) {}
187   FillESDHist& operator=(const FillESDHist& o) { fM = o.fM; return *this;  }
188   Bool_t operator()(UShort_t, Char_t, UShort_t, UShort_t, Float_t m, Float_t) 
189   {
190     // Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
191     if(m == AliESDFMD::kInvalidMult) return true;
192     
193     fM->GetESDsData(0)->Fill(m);    
194     return true;
195   }
196   AliFMDQADataMakerRec* fM;
197 };
198 #endif
199
200 //_____________________________________________________________________
201 void AliFMDQADataMakerRec::MakeESDs(AliESDEvent * esd)
202 {
203   if(!esd) {
204     AliError("FMD ESD object not found!!") ; 
205     return;
206   }
207   AliFMDDebug(2, ("Will loop over ESD data and fill histogram"));
208
209   AliESDFMD* fmd = esd->GetFMDData();
210   if (!fmd) return;
211
212 #if 0
213   FillESDHist f(this);
214   fmd->ForEach(f);
215 #else
216
217
218
219   // FIXME - we should use AliESDFMD::ForOne subclass to do this!
220   for(UShort_t det=1;det<=3;det++) {
221     UShort_t nrng = (det == 1 ? 1 : 2);
222     for (UShort_t ir = 0; ir < nrng; ir++) {
223       Char_t   ring = (ir == 0 ? 'I' : 'O');
224       UShort_t nsec = (ir == 0 ? 20  : 40);
225       UShort_t nstr = (ir == 0 ? 512 : 256);
226       for(UShort_t sec =0; sec < nsec;  sec++)  {
227         for(UShort_t strip = 0; strip < nstr; strip++) {
228           Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
229           if(mult == AliESDFMD::kInvalidMult) continue;
230           
231           GetESDsData(0)->Fill(mult);
232         }
233       }
234     }
235   }
236 #endif
237 }
238
239
240 //_____________________________________________________________________
241 void AliFMDQADataMakerRec::MakeDigits()
242 {
243   // makes data from Digits  
244   if(!fDigitsArray)  {
245     AliError("FMD Digit object not found!!") ;
246     return;
247   }
248   
249   for(Int_t i=0;i<fDigitsArray->GetEntriesFast();i++) {
250     //Raw ADC counts
251     AliFMDDigit* digit = static_cast<AliFMDDigit*>(fDigitsArray->At(i));
252     GetDigitsData(0)->Fill(digit->Counts());
253   }
254 }
255
256 //_____________________________________________________________________
257 void AliFMDQADataMakerRec::MakeDigits(TTree * digitTree)
258 {
259   
260   if (fDigitsArray) 
261     fDigitsArray->Clear();
262   else 
263     fDigitsArray = new TClonesArray("AliFMDDigit", 1000);
264
265   TBranch*      branch = digitTree->GetBranch("FMD");
266   if (!branch) {
267     AliWarning("FMD branch in Digit Tree not found") ; 
268     return;
269   } 
270   branch->SetAddress(&fDigitsArray);
271   branch->GetEntry(0); 
272   MakeDigits();
273 }
274
275 //_____________________________________________________________________
276 void AliFMDQADataMakerRec::MakeRaws(AliRawReader* rawReader)
277 {
278  
279  AliFMDRawReader fmdReader(rawReader,0);
280   
281   if (fDigitsArray) 
282     fDigitsArray->Clear();
283   else 
284     fDigitsArray = new TClonesArray("AliFMDDigit", 1000);
285
286   TClonesArray* digitsAddress = fDigitsArray;
287   
288   rawReader->Reset();
289                 
290   digitsAddress->Clear();
291   fmdReader.ReadAdcs(digitsAddress);
292   for(Int_t i=0;i<digitsAddress->GetEntriesFast();i++) {
293     //Raw ADC counts
294     AliFMDDigit*      digit = static_cast<AliFMDDigit*>(digitsAddress->At(i));
295     UShort_t          det   = digit->Detector();
296     Char_t            ring  = digit->Ring();
297     UShort_t          sec   = digit->Sector();
298     // UShort_t strip = digit->Strip();
299     AliFMDParameters* pars  = AliFMDParameters::Instance();
300     Short_t           board = pars->GetAltroMap()->Sector2Board(ring, sec);
301     
302     Int_t index1 = GetHalfringIndex(det, ring, 0, 1);
303     GetRawsData(index1)->Fill(digit->Counts());
304     Int_t index2 = GetHalfringIndex(det, ring, board/16,0);
305     GetRawsData(index2)->Fill(digit->Counts());
306     
307   }
308 }
309
310 //_____________________________________________________________________
311 void AliFMDQADataMakerRec::MakeRecPoints(TTree* clustersTree)
312 {
313   // makes data from RecPoints
314   
315    AliFMDParameters* pars = AliFMDParameters::Instance();
316   fRecPointsArray.Clear();
317   TBranch *fmdbranch = clustersTree->GetBranch("FMD");
318   if (!fmdbranch) { 
319     AliError("can't get the branch with the FMD recpoints !");
320     return;
321   }
322   
323   TClonesArray* RecPointsAddress = &fRecPointsArray;
324   
325   fmdbranch->SetAddress(&RecPointsAddress);
326   fmdbranch->GetEntry(0);
327   TIter next(RecPointsAddress) ; 
328   AliFMDRecPoint * rp ; 
329   while ((rp = static_cast<AliFMDRecPoint*>(next()))) {
330     
331     GetRecPointsData(0)->Fill(rp->Edep()/pars->GetEdepMip()) ;
332   
333   }
334
335 }
336
337 //_____________________________________________________________________ 
338 void AliFMDQADataMakerRec::StartOfDetectorCycle()
339 {
340   // What 
341   // to 
342   // do?
343 }
344 //_____________________________________________________________________ 
345 Int_t AliFMDQADataMakerRec::GetHalfringIndex(UShort_t det, 
346                                              Char_t ring, 
347                                              UShort_t board, 
348                                              UShort_t monitor) {
349   
350   UShort_t iring  =  (ring == 'I' ? 1 : 0);
351   
352   Int_t index = ( ((det-1) << 3) | (iring << 2) | (board << 1) | (monitor << 0));
353   
354   return index-2;
355   
356 }
357
358 //_____________________________________________________________________ 
359
360
361 //
362 // EOF
363 //