]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDQADataMakerRec.cxx
1.The QA data created on demand according to the event species at filling time. 2...
[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
35 //_____________________________________________________________________
36 // This is the class that collects the QA data for the FMD during
37 // reconstruction.  
38 //
39 // The following data types are picked up:
40 // - rec points
41 // - esd data
42 // - raws
43 // Author : Hans Hjersing Dalsgaard, hans.dalsgaard@cern.ch
44 //_____________________________________________________________________
45
46 ClassImp(AliFMDQADataMakerRec)
47 #if 0
48 ; // For Emacs - do not delete!
49 #endif
50            
51 //_____________________________________________________________________
52 AliFMDQADataMakerRec::AliFMDQADataMakerRec() : 
53   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kFMD), 
54                     "FMD Quality Assurance Data Maker"),
55   fDigitsArray("AliFMDDigit", 0),
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     fDigitsArray(qadm.fDigitsArray),
67     fRecPointsArray(qadm.fRecPointsArray)
68 {
69   // copy ctor 
70   // Parameters: 
71   //    qadm    Object to copy from
72   
73 }
74 //_____________________________________________________________________
75 AliFMDQADataMakerRec& AliFMDQADataMakerRec::operator = (const AliFMDQADataMakerRec& qadm ) 
76 {
77   fDigitsArray = qadm.fDigitsArray;
78   fRecPointsArray = qadm.fRecPointsArray;
79   
80   return *this;
81 }
82 //_____________________________________________________________________
83 AliFMDQADataMakerRec::~AliFMDQADataMakerRec()
84 {
85  
86 }
87
88
89 //_____________________________________________________________________ 
90
91 void 
92 AliFMDQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, 
93                                          TObjArray ** list)
94 {
95   // Detector specific actions at end of cycle
96   // do the QA checking
97   AliLog::Message(5,"FMD: end of detector cycle",
98                   "AliFMDQADataMakerRec","AliFMDQADataMakerRec",
99                   "AliFMDQADataMakerRec::EndOfDetectorCycle",
100                   "AliFMDQADataMakerRec.cxx",95);
101   AliQAChecker::Instance()->Run(AliQAv1::kFMD, task, list);
102 }
103
104 //_____________________________________________________________________ 
105 void AliFMDQADataMakerRec::InitESDs()
106 {
107   // create Digits histograms in Digits subdir
108   const Bool_t expert   = kTRUE ; 
109   const Bool_t image    = kTRUE ; 
110   
111   TH1F* hEnergyOfESD = new TH1F("hEnergyOfESD","Energy distribution",100,0,3);
112   hEnergyOfESD->SetXTitle("Edep/Emip");
113   hEnergyOfESD->SetYTitle("Counts");
114   Add2ESDsList(hEnergyOfESD, 0, !expert, image);
115     
116 }
117
118 //_____________________________________________________________________
119 void AliFMDQADataMakerRec::InitDigits()
120 {
121   // create Digits histograms in Digits subdir
122   const Bool_t expert   = kTRUE ; 
123   const Bool_t image    = kTRUE ; 
124   
125   TH1I* hADCCounts = new TH1I("hADCCounts","Dist of ADC counts;ADC counts;Counts",1024,0,1024);
126   Add2DigitsList(hADCCounts, 0, !expert, image);
127 }
128
129
130 //_____________________________________________________________________ 
131 void AliFMDQADataMakerRec::InitRecPoints()
132 {
133   // create Reconstructed Points histograms in RecPoints subdir
134   const Bool_t expert   = kTRUE ; 
135   const Bool_t image    = kTRUE ; 
136
137   TH1F* hEnergyOfRecpoints = new TH1F("hEnergyOfRecpoints",
138                                       "Energy Distribution",100,0,3);
139   hEnergyOfRecpoints->SetXTitle("Edep/Emip");
140   hEnergyOfRecpoints->SetYTitle("Counts");
141   Add2RecPointsList(hEnergyOfRecpoints,0, !expert, image);
142 }
143
144 //_____________________________________________________________________ 
145 void AliFMDQADataMakerRec::InitRaws()
146 {
147   // create Raws histograms in Raws subdir  
148   const Bool_t expert   = kTRUE ; 
149   const Bool_t saveCorr = kTRUE ; 
150   const Bool_t image    = kTRUE ; 
151
152   TH1I* hADCCounts;
153   for(Int_t det = 1; det<=3; det++) {
154     Int_t firstring = (det==1 ? 1 : 0);
155     for(Int_t iring = firstring;iring<=1;iring++) {
156       Char_t ring = (iring == 1 ? 'I' : 'O');
157       hADCCounts      = new TH1I(Form("hADCCounts_FMD%d%c",
158                                       det, ring), "ADC counts;Amplitude [ADC counts];Counts",
159                                  1024,0,1023);
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 //_____________________________________________________________________
185 void AliFMDQADataMakerRec::MakeESDs(AliESDEvent * esd)
186 {
187   if(!esd) {
188     AliError("FMD ESD object not found!!") ; 
189     return;
190   }
191   AliESDFMD* fmd = esd->GetFMDData();
192   if (!fmd) return;
193   
194   // Check id histograms already created for this Event Specie
195   if ( ! GetESDsData(0) )
196     InitESDs() ;
197
198   for(UShort_t det=1;det<=3;det++) {
199     for (UShort_t ir = 0; ir < 2; ir++) {
200       Char_t   ring = (ir == 0 ? 'I' : 'O');
201       UShort_t nsec = (ir == 0 ? 20  : 40);
202       UShort_t nstr = (ir == 0 ? 512 : 256);
203       for(UShort_t sec =0; sec < nsec;  sec++)  {
204         for(UShort_t strip = 0; strip < nstr; strip++) {
205           Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
206           if(mult == AliESDFMD::kInvalidMult) continue;
207           
208           GetESDsData(0)->Fill(mult);
209         }
210       }
211     }
212   }
213 }
214
215
216 //_____________________________________________________________________
217 void AliFMDQADataMakerRec::MakeDigits(TClonesArray * digits)
218 {
219   // makes data from Digits  
220   if(!digits)  {
221     AliError("FMD Digit object not found!!") ;
222     return;
223   }
224   
225   // Check id histograms already created for this Event Specie
226   if ( ! GetDigitsData(0) )
227     InitDigits() ;
228
229   for(Int_t i=0;i<digits->GetEntriesFast();i++) {
230     //Raw ADC counts
231     AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
232     GetDigitsData(0)->Fill(digit->Counts());
233   }
234 }
235
236 //_____________________________________________________________________
237 void AliFMDQADataMakerRec::MakeDigits(TTree * digitTree)
238 {
239   
240   fDigitsArray.Clear();
241   TBranch*      branch = digitTree->GetBranch("FMD");
242   if (!branch) {
243     AliWarning("FMD branch in Digit Tree not found") ; 
244     return;
245   } 
246   TClonesArray* digitsAddress = &fDigitsArray;
247   branch->SetAddress(&digitsAddress);
248   branch->GetEntry(0); 
249   MakeDigits(digitsAddress);
250 }
251
252 //_____________________________________________________________________
253 void AliFMDQADataMakerRec::MakeRaws(AliRawReader* rawReader)
254 {
255  
256   // Check id histograms already created for this Event Specie
257   if ( ! GetRawsData(0) )
258     InitRaws() ;
259
260   AliFMDRawReader fmdReader(rawReader,0);
261   TClonesArray* digitsAddress = &fDigitsArray;
262   
263   rawReader->Reset();
264                 
265   digitsAddress->Clear();
266   fmdReader.ReadAdcs(digitsAddress);
267   for(Int_t i=0;i<digitsAddress->GetEntriesFast();i++) {
268     //Raw ADC counts
269     AliFMDDigit*      digit = static_cast<AliFMDDigit*>(digitsAddress->At(i));
270     UShort_t          det   = digit->Detector();
271     Char_t            ring  = digit->Ring();
272     UShort_t          sec   = digit->Sector();
273     // UShort_t strip = digit->Strip();
274     AliFMDParameters* pars  = AliFMDParameters::Instance();
275     Short_t           board = pars->GetAltroMap()->Sector2Board(ring, sec);
276     
277     Int_t index1 = GetHalfringIndex(det, ring, 0, 1);
278     GetRawsData(index1)->Fill(digit->Counts());
279     Int_t index2 = GetHalfringIndex(det, ring, board/16,0);
280     GetRawsData(index2)->Fill(digit->Counts());
281     
282   }
283 }
284
285 //_____________________________________________________________________
286 void AliFMDQADataMakerRec::MakeRecPoints(TTree* clustersTree)
287 {
288   // makes data from RecPoints
289   
290   // Check id histograms already created for this Event Specie
291   if ( ! GetRecPointsData(0) )
292     InitRecPoints() ;
293
294   AliFMDParameters* pars = AliFMDParameters::Instance();
295   fRecPointsArray.Clear();
296   TBranch *fmdbranch = clustersTree->GetBranch("FMD");
297   if (!fmdbranch) { 
298     AliError("can't get the branch with the FMD recpoints !");
299     return;
300   }
301   
302   TClonesArray* RecPointsAddress = &fRecPointsArray;
303   
304   fmdbranch->SetAddress(&RecPointsAddress);
305   fmdbranch->GetEntry(0);
306   TIter next(RecPointsAddress) ; 
307   AliFMDRecPoint * rp ; 
308   while ((rp = static_cast<AliFMDRecPoint*>(next()))) {
309     
310     GetRecPointsData(0)->Fill(rp->Edep()/pars->GetEdepMip()) ;
311   
312   }
313
314 }
315
316 //_____________________________________________________________________ 
317 void AliFMDQADataMakerRec::StartOfDetectorCycle()
318 {
319   // What 
320   // to 
321   // do?
322 }
323 //_____________________________________________________________________ 
324 Int_t AliFMDQADataMakerRec::GetHalfringIndex(UShort_t det, 
325                                              Char_t ring, 
326                                              UShort_t board, 
327                                              UShort_t monitor) {
328   
329   UShort_t iring  =  (ring == 'I' ? 1 : 0);
330   
331   Int_t index = ( ((det-1) << 3) | (iring << 2) | (board << 1) | (monitor << 0));
332   
333   return index-2;
334   
335 }
336
337 //_____________________________________________________________________ 
338
339
340 //
341 // EOF
342 //