]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDQADataMakerRec.cxx
Adding more intelligent names to the raw histos shown at P2
[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   Int_t colors[3] = {kRed,kGreen,kBlue};
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), Form("FMD%d%c ADC counts;Amplitude [ADC counts];Counts",det,ring),
157                                  1024,0,1023);
158       hADCCounts->SetFillStyle(3001);
159       hADCCounts->SetFillColor(colors[det-1]+2-2*iring);
160       
161       Int_t index1 = GetHalfringIndex(det, ring, 0,1);
162       Add2RawsList(hADCCounts, index1, !expert, image, !saveCorr);
163       
164       for(Int_t b = 0; b<=1;b++) {
165         
166         //Hexadecimal board numbers 0x0, 0x1, 0x10, 0x11;
167         UInt_t board = (iring == 1 ? 0 : 1);
168         board = board + b*16;
169         
170         
171         hADCCounts      = new TH1I(Form("hADCCounts_FMD%d%c_board%d",
172                                         det, ring, board), "ADC counts;Amplitude [ADC counts];Counts",
173                                    1024,0,1023);
174         hADCCounts->SetXTitle("ADC counts");
175         hADCCounts->SetYTitle("");
176         Int_t index2 = GetHalfringIndex(det, ring, board/16,0);
177         Add2RawsList(hADCCounts, index2, expert, !image, !saveCorr);
178
179       }
180     }
181   }
182 }
183
184 #if 0
185 struct FillESDHist : public AliESDFMD::ForOne
186 {
187   FillESDHist(AliFMDQADataMakerRec* m) : fM(m) {}
188   FillESDHist(const FillESDHist& o) : fM(o.fM) {}
189   FillESDHist& operator=(const FillESDHist& o) { fM = o.fM; return *this;  }
190   Bool_t operator()(UShort_t, Char_t, UShort_t, UShort_t, Float_t m, Float_t) 
191   {
192     // Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
193     if(m == AliESDFMD::kInvalidMult) return true;
194     
195     fM->GetESDsData(0)->Fill(m);    
196     return true;
197   }
198   AliFMDQADataMakerRec* fM;
199 };
200 #endif
201
202 //_____________________________________________________________________
203 void AliFMDQADataMakerRec::MakeESDs(AliESDEvent * esd)
204 {
205   if(!esd) {
206     AliError("FMD ESD object not found!!") ; 
207     return;
208   }
209   AliFMDDebug(2, ("Will loop over ESD data and fill histogram"));
210
211   AliESDFMD* fmd = esd->GetFMDData();
212   if (!fmd) return;
213
214 #if 0
215   FillESDHist f(this);
216   fmd->ForEach(f);
217 #else
218
219
220
221   // FIXME - we should use AliESDFMD::ForOne subclass to do this!
222   for(UShort_t det=1;det<=3;det++) {
223     UShort_t nrng = (det == 1 ? 1 : 2);
224     for (UShort_t ir = 0; ir < nrng; ir++) {
225       Char_t   ring = (ir == 0 ? 'I' : 'O');
226       UShort_t nsec = (ir == 0 ? 20  : 40);
227       UShort_t nstr = (ir == 0 ? 512 : 256);
228       for(UShort_t sec =0; sec < nsec;  sec++)  {
229         for(UShort_t strip = 0; strip < nstr; strip++) {
230           Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
231           if(mult == AliESDFMD::kInvalidMult) continue;
232           
233           GetESDsData(0)->Fill(mult);
234         }
235       }
236     }
237   }
238 #endif
239 }
240
241
242 //_____________________________________________________________________
243 void AliFMDQADataMakerRec::MakeDigits()
244 {
245   // makes data from Digits  
246   if(!fDigitsArray)  {
247     AliError("FMD Digit object not found!!") ;
248     return;
249   }
250   
251   for(Int_t i=0;i<fDigitsArray->GetEntriesFast();i++) {
252     //Raw ADC counts
253     AliFMDDigit* digit = static_cast<AliFMDDigit*>(fDigitsArray->At(i));
254     GetDigitsData(0)->Fill(digit->Counts());
255   }
256 }
257
258 //_____________________________________________________________________
259 void AliFMDQADataMakerRec::MakeDigits(TTree * digitTree)
260 {
261   
262   if (fDigitsArray) 
263     fDigitsArray->Clear();
264   else 
265     fDigitsArray = new TClonesArray("AliFMDDigit", 1000);
266
267   TBranch*      branch = digitTree->GetBranch("FMD");
268   if (!branch) {
269     AliWarning("FMD branch in Digit Tree not found") ; 
270     return;
271   } 
272   branch->SetAddress(&fDigitsArray);
273   branch->GetEntry(0); 
274   MakeDigits();
275 }
276
277 //_____________________________________________________________________
278 void AliFMDQADataMakerRec::MakeRaws(AliRawReader* rawReader)
279 {
280  
281  AliFMDRawReader fmdReader(rawReader,0);
282   
283   if (fDigitsArray) 
284     fDigitsArray->Clear();
285   else 
286     fDigitsArray = new TClonesArray("AliFMDDigit", 1000);
287
288   TClonesArray* digitsAddress = fDigitsArray;
289   
290   rawReader->Reset();
291                 
292   digitsAddress->Clear();
293   fmdReader.ReadAdcs(digitsAddress);
294   for(Int_t i=0;i<digitsAddress->GetEntriesFast();i++) {
295     //Raw ADC counts
296     AliFMDDigit*      digit = static_cast<AliFMDDigit*>(digitsAddress->At(i));
297     UShort_t          det   = digit->Detector();
298     Char_t            ring  = digit->Ring();
299     UShort_t          sec   = digit->Sector();
300     // UShort_t strip = digit->Strip();
301     AliFMDParameters* pars  = AliFMDParameters::Instance();
302     Short_t           board = pars->GetAltroMap()->Sector2Board(ring, sec);
303     
304     Int_t index1 = GetHalfringIndex(det, ring, 0, 1);
305     GetRawsData(index1)->Fill(digit->Counts());
306     Int_t index2 = GetHalfringIndex(det, ring, board/16,0);
307     GetRawsData(index2)->Fill(digit->Counts());
308     
309   }
310 }
311
312 //_____________________________________________________________________
313 void AliFMDQADataMakerRec::MakeRecPoints(TTree* clustersTree)
314 {
315   // makes data from RecPoints
316   
317    AliFMDParameters* pars = AliFMDParameters::Instance();
318   fRecPointsArray.Clear();
319   TBranch *fmdbranch = clustersTree->GetBranch("FMD");
320   if (!fmdbranch) { 
321     AliError("can't get the branch with the FMD recpoints !");
322     return;
323   }
324   
325   TClonesArray* RecPointsAddress = &fRecPointsArray;
326   
327   fmdbranch->SetAddress(&RecPointsAddress);
328   fmdbranch->GetEntry(0);
329   TIter next(RecPointsAddress) ; 
330   AliFMDRecPoint * rp ; 
331   while ((rp = static_cast<AliFMDRecPoint*>(next()))) {
332     
333     GetRecPointsData(0)->Fill(rp->Edep()/pars->GetEdepMip()) ;
334   
335   }
336
337 }
338
339 //_____________________________________________________________________ 
340 void AliFMDQADataMakerRec::StartOfDetectorCycle()
341 {
342   // What 
343   // to 
344   // do?
345 }
346 //_____________________________________________________________________ 
347 Int_t AliFMDQADataMakerRec::GetHalfringIndex(UShort_t det, 
348                                              Char_t ring, 
349                                              UShort_t board, 
350                                              UShort_t monitor) {
351   
352   UShort_t iring  =  (ring == 'I' ? 1 : 0);
353   
354   Int_t index = ( ((det-1) << 3) | (iring << 2) | (board << 1) | (monitor << 0));
355   
356   return index-2;
357   
358 }
359
360 //_____________________________________________________________________ 
361
362
363 //
364 // EOF
365 //