]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDQADataMakerRec.cxx
Fixes for bug #49914: Compilation breaks in trunk, and bug #48629: Trunk cannot read...
[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   TH1F* hEnergyOfESD = new TH1F("hEnergyOfESD","Energy distribution",100,0,3);
109   hEnergyOfESD->SetXTitle("Edep/Emip");
110   hEnergyOfESD->SetYTitle("Counts");
111   Add2ESDsList(hEnergyOfESD, 0);
112     
113 }
114
115 //_____________________________________________________________________ 
116 /*void AliFMDQADataMakerRec::InitDigits()
117 {
118   // create Digits histograms in Digits subdir
119   TH1I* hADCCounts      = new TH1I("hADCCounts","Dist of ADC counts",
120                                    1024,0,1024);
121   hADCCounts->SetXTitle("ADC counts");
122   hADCCounts->SetYTitle("");
123   Add2DigitsList(hADCCounts, 0);
124   
125 }
126 */
127 //_____________________________________________________________________ 
128 void AliFMDQADataMakerRec::InitRecPoints()
129 {
130   TH1F* hEnergyOfRecpoints = new TH1F("hEnergyOfRecpoints",
131                                       "Energy Distribution",100,0,3);
132   hEnergyOfRecpoints->SetXTitle("Edep/Emip");
133   hEnergyOfRecpoints->SetYTitle("");
134   Add2RecPointsList(hEnergyOfRecpoints,0);
135 }
136
137 //_____________________________________________________________________ 
138 void AliFMDQADataMakerRec::InitRaws()
139 {
140   
141   TH1I* hADCCounts;
142   for(Int_t det = 1; det<=3; det++) {
143     Int_t firstring = (det==1 ? 1 : 0);
144     for(Int_t iring = firstring;iring<=1;iring++) {
145       Char_t ring = (iring == 1 ? 'I' : 'O');
146       hADCCounts      = new TH1I(Form("hADCCounts_FMD%d%c",
147                                       det, ring), "ADC counts",
148                                  1024,0,1023);
149       
150       Int_t index1 = GetHalfringIndex(det, ring, 0,1);
151       Add2RawsList(hADCCounts, index1,kTRUE);
152       
153       for(Int_t b = 0; b<=1;b++) {
154         
155         //Hexadecimal board numbers 0x0, 0x1, 0x10, 0x11;
156         UInt_t board = (iring == 1 ? 0 : 1);
157         board = board + b*16;
158         
159         
160         hADCCounts      = new TH1I(Form("hADCCounts_FMD%d%c_board%d",
161                                         det, ring, board), "ADC counts",
162                                    1024,0,1023);
163         hADCCounts->SetXTitle("ADC counts");
164         hADCCounts->SetYTitle("");
165         Int_t index2 = GetHalfringIndex(det, ring, board/16,0);
166         Add2RawsList(hADCCounts, index2);
167
168       }
169     }
170   }
171 }
172
173 //_____________________________________________________________________
174 void AliFMDQADataMakerRec::MakeESDs(AliESDEvent * esd)
175 {
176   if(!esd) {
177     AliError("FMD ESD object not found!!") ; 
178     return;
179   }
180   AliESDFMD* fmd = esd->GetFMDData();
181   if (!fmd) return;
182   
183   for(UShort_t det=1;det<=3;det++) {
184     for (UShort_t ir = 0; ir < 2; ir++) {
185       Char_t   ring = (ir == 0 ? 'I' : 'O');
186       UShort_t nsec = (ir == 0 ? 20  : 40);
187       UShort_t nstr = (ir == 0 ? 512 : 256);
188       for(UShort_t sec =0; sec < nsec;  sec++)  {
189         for(UShort_t strip = 0; strip < nstr; strip++) {
190           Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
191           if(mult == AliESDFMD::kInvalidMult) continue;
192           
193           GetESDsData(0)->Fill(mult);
194         }
195       }
196     }
197   }
198 }
199
200 /*
201 //_____________________________________________________________________
202 void AliFMDQADataMakerRec::MakeDigits(TClonesArray * digits)
203 {
204   // makes data from Digits  
205   if(!digits)  {
206     AliError("FMD Digit object not found!!") ;
207     return;
208   }
209   for(Int_t i=0;i<digits->GetEntriesFast();i++) {
210     //Raw ADC counts
211     AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
212     GetDigitsData(0)->Fill(digit->Counts());
213   }
214 }
215
216 //_____________________________________________________________________
217 void AliFMDQADataMakerRec::MakeDigits(TTree * digitTree)
218 {
219   
220   fDigitsArray.Clear();
221   TBranch*      branch = digitTree->GetBranch("FMD");
222   if (!branch) {
223     AliWarning("FMD branch in Digit Tree not found") ; 
224     return;
225   } 
226   TClonesArray* digitsAddress = &fDigitsArray;
227   branch->SetAddress(&digitsAddress);
228   branch->GetEntry(0); 
229   MakeDigits(digitsAddress);
230 }
231 */
232 //_____________________________________________________________________
233 void AliFMDQADataMakerRec::MakeRaws(AliRawReader* rawReader)
234 {
235  
236   AliFMDRawReader fmdReader(rawReader,0);
237   TClonesArray* digitsAddress = &fDigitsArray;
238   
239   rawReader->Reset();
240                 
241   digitsAddress->Clear();
242   fmdReader.ReadAdcs(digitsAddress);
243   for(Int_t i=0;i<digitsAddress->GetEntriesFast();i++) {
244     //Raw ADC counts
245     AliFMDDigit*      digit = static_cast<AliFMDDigit*>(digitsAddress->At(i));
246     UShort_t          det   = digit->Detector();
247     Char_t            ring  = digit->Ring();
248     UShort_t          sec   = digit->Sector();
249     // UShort_t strip = digit->Strip();
250     AliFMDParameters* pars  = AliFMDParameters::Instance();
251     Short_t           board = pars->GetAltroMap()->Sector2Board(ring, sec);
252     
253     Int_t index1 = GetHalfringIndex(det, ring, 0, 1);
254     GetRawsData(index1)->Fill(digit->Counts());
255     Int_t index2 = GetHalfringIndex(det, ring, board/16,0);
256     GetRawsData(index2)->Fill(digit->Counts());
257     
258   }
259 }
260
261 //_____________________________________________________________________
262 void AliFMDQADataMakerRec::MakeRecPoints(TTree* clustersTree)
263 {
264   // makes data from RecPoints
265   AliFMDParameters* pars = AliFMDParameters::Instance();
266   fRecPointsArray.Clear();
267   TBranch *fmdbranch = clustersTree->GetBranch("FMD");
268   if (!fmdbranch) { 
269     AliError("can't get the branch with the FMD recpoints !");
270     return;
271   }
272   
273   TClonesArray* RecPointsAddress = &fRecPointsArray;
274   
275   fmdbranch->SetAddress(&RecPointsAddress);
276   fmdbranch->GetEntry(0);
277   TIter next(RecPointsAddress) ; 
278   AliFMDRecPoint * rp ; 
279   while ((rp = static_cast<AliFMDRecPoint*>(next()))) {
280     
281     GetRecPointsData(0)->Fill(rp->Edep()/pars->GetEdepMip()) ;
282   
283   }
284
285 }
286
287 //_____________________________________________________________________ 
288 void AliFMDQADataMakerRec::StartOfDetectorCycle()
289 {
290   // What 
291   // to 
292   // do?
293 }
294 //_____________________________________________________________________ 
295 Int_t AliFMDQADataMakerRec::GetHalfringIndex(UShort_t det, 
296                                              Char_t ring, 
297                                              UShort_t board, 
298                                              UShort_t monitor) {
299   
300   UShort_t iring  =  (ring == 'I' ? 1 : 0);
301   
302   Int_t index = ( ((det-1) << 3) | (iring << 2) | (board << 1) | (monitor << 0));
303   
304   return index-2;
305   
306 }
307
308 //_____________________________________________________________________ 
309
310
311 //
312 // EOF
313 //