]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDQADataMakerRec.cxx
Restoration of earlier checked in version of the DataMaker
[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
34 //_____________________________________________________________________
35 // This is the class that collects the QA data for the FMD during
36 // reconstruction.  
37 //
38 // The following data types are picked up:
39 // - rec points
40 // - esd data
41 // - raws
42 // Author : Hans Hjersing Dalsgaard, hans.dalsgaard@cern.ch
43 //_____________________________________________________________________
44
45 ClassImp(AliFMDQADataMakerRec)
46 #if 0
47 ; // For Emacs - do not delete!
48 #endif
49            
50 //_____________________________________________________________________
51 AliFMDQADataMakerRec::AliFMDQADataMakerRec() : 
52   AliQADataMakerRec(AliQA::GetDetName(AliQA::kFMD), 
53                     "FMD Quality Assurance Data Maker"),
54   fDigitsArray("AliFMDDigit", 0),
55   fRecPointsArray("AliFMDRecPoint", 1000)
56 {
57   // ctor
58  
59 }
60
61 //_____________________________________________________________________
62 AliFMDQADataMakerRec::AliFMDQADataMakerRec(const AliFMDQADataMakerRec& qadm) 
63   : AliQADataMakerRec(AliQA::GetDetName(AliQA::kFMD), 
64                       "FMD Quality Assurance Data Maker"),
65     fDigitsArray(qadm.fDigitsArray),
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   fDigitsArray = qadm.fDigitsArray;
77   fRecPointsArray = qadm.fRecPointsArray;
78   
79   return *this;
80 }
81 //_____________________________________________________________________
82 AliFMDQADataMakerRec::~AliFMDQADataMakerRec()
83 {
84  
85 }
86
87
88 //_____________________________________________________________________ 
89
90 void 
91 AliFMDQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, 
92                                          TObjArray * list)
93 {
94   // Detector specific actions at end of cycle
95   // do the QA checking
96   AliLog::Message(5,"FMD: end of detector cycle",
97                   "AliFMDQADataMakerRec","AliFMDQADataMakerRec",
98                   "AliFMDQADataMakerRec::EndOfDetectorCycle",
99                   "AliFMDQADataMakerRec.cxx",95);
100   AliQAChecker::Instance()->Run(AliQA::kFMD, task, list);
101 }
102
103 //_____________________________________________________________________ 
104 void AliFMDQADataMakerRec::InitESDs()
105 {
106   // create Digits histograms in Digits subdir
107   TH1F* hEnergyOfESD = new TH1F("hEnergyOfESD","Energy distribution",100,0,3);
108   hEnergyOfESD->SetXTitle("Edep/Emip");
109   hEnergyOfESD->SetYTitle("Counts");
110   Add2ESDsList(hEnergyOfESD, 0);
111     
112 }
113
114 //_____________________________________________________________________ 
115 /*void AliFMDQADataMakerRec::InitDigits()
116 {
117   // create Digits histograms in Digits subdir
118   TH1I* hADCCounts      = new TH1I("hADCCounts","Dist of ADC counts",
119                                    1024,0,1024);
120   hADCCounts->SetXTitle("ADC counts");
121   hADCCounts->SetYTitle("");
122   Add2DigitsList(hADCCounts, 0);
123   
124 }
125 */
126 //_____________________________________________________________________ 
127 void AliFMDQADataMakerRec::InitRecPoints()
128 {
129   TH1F* hEnergyOfRecpoints = new TH1F("hEnergyOfRecpoints",
130                                       "Energy Distribution",100,0,3);
131   hEnergyOfRecpoints->SetXTitle("Edep/Emip");
132   hEnergyOfRecpoints->SetYTitle("");
133   Add2RecPointsList(hEnergyOfRecpoints,0);
134 }
135
136 //_____________________________________________________________________ 
137 void AliFMDQADataMakerRec::InitRaws()
138 {
139   
140   TH1I* hADCCounts;
141   for(Int_t det = 1; det<=3; det++) {
142     Int_t firstring = (det==1 ? 1 : 0);
143     for(Int_t iring = firstring;iring<=1;iring++) {
144       for(Int_t b = 0; b<=1;b++) {
145         
146         //Hexadecimal board numbers 0x0, 0x1, 0x10, 0x10;
147         UInt_t board = (iring == 1 ? 0 : 1);
148         board = board + b*16;
149         Char_t ring = (iring == 1 ? 'I' : 'O');
150         
151         hADCCounts      = new TH1I(Form("hADCCounts_FMD%d%c_board%d",
152                                         det, ring, board), "ADC counts",
153                                    1024,0,1023);
154         hADCCounts->SetXTitle("ADC counts");
155         hADCCounts->SetYTitle("");
156         Int_t index = GetHalfringIndex(det, ring, board/16);
157         Add2RawsList(hADCCounts, index);
158         
159       }
160     }
161   }
162 }
163
164 //_____________________________________________________________________
165 void AliFMDQADataMakerRec::MakeESDs(AliESDEvent * esd)
166 {
167   if(!esd) {
168     AliError("FMD ESD object not found!!") ; 
169     return;
170   }
171   AliESDFMD* fmd = esd->GetFMDData();
172   if (!fmd) return;
173   
174   for(UShort_t det=1;det<=3;det++) {
175     for (UShort_t ir = 0; ir < 2; ir++) {
176       Char_t   ring = (ir == 0 ? 'I' : 'O');
177       UShort_t nsec = (ir == 0 ? 20  : 40);
178       UShort_t nstr = (ir == 0 ? 512 : 256);
179       for(UShort_t sec =0; sec < nsec;  sec++)  {
180         for(UShort_t strip = 0; strip < nstr; strip++) {
181           Float_t mult = fmd->Multiplicity(det,ring,sec,strip);
182           if(mult == AliESDFMD::kInvalidMult) continue;
183           
184           GetESDsData(0)->Fill(mult);
185         }
186       }
187     }
188   }
189 }
190
191 /*
192 //_____________________________________________________________________
193 void AliFMDQADataMakerRec::MakeDigits(TClonesArray * digits)
194 {
195   // makes data from Digits  
196   if(!digits)  {
197     AliError("FMD Digit object not found!!") ;
198     return;
199   }
200   for(Int_t i=0;i<digits->GetEntriesFast();i++) {
201     //Raw ADC counts
202     AliFMDDigit* digit = static_cast<AliFMDDigit*>(digits->At(i));
203     GetDigitsData(0)->Fill(digit->Counts());
204   }
205 }
206
207 //_____________________________________________________________________
208 void AliFMDQADataMakerRec::MakeDigits(TTree * digitTree)
209 {
210   
211   fDigitsArray.Clear();
212   TBranch*      branch = digitTree->GetBranch("FMD");
213   if (!branch) {
214     AliWarning("FMD branch in Digit Tree not found") ; 
215     return;
216   } 
217   TClonesArray* digitsAddress = &fDigitsArray;
218   branch->SetAddress(&digitsAddress);
219   branch->GetEntry(0); 
220   MakeDigits(digitsAddress);
221 }
222 */
223 //_____________________________________________________________________
224 void AliFMDQADataMakerRec::MakeRaws(AliRawReader* rawReader)
225 {
226  
227   AliFMDRawReader fmdReader(rawReader,0);
228   TClonesArray* digitsAddress = &fDigitsArray;
229   
230   rawReader->Reset();
231                 
232         digitsAddress->Clear();
233         fmdReader.ReadAdcs(digitsAddress);
234         for(Int_t i=0;i<digitsAddress->GetEntriesFast();i++) {
235           //Raw ADC counts
236           AliFMDDigit* digit = static_cast<AliFMDDigit*>(digitsAddress->At(i));
237           UShort_t det = digit->Detector();
238           Char_t ring  = digit->Ring();
239           UShort_t sec = digit->Sector();
240           UShort_t strip = digit->Strip();
241           UInt_t ddl, board, chip, channel;
242           AliFMDParameters* pars = AliFMDParameters::Instance();
243           pars->Detector2Hardware(det,ring,sec,strip,ddl,board,chip,channel);
244           Int_t index = GetHalfringIndex(det, ring, board/16);
245           
246           GetRawsData(index)->Fill(digit->Counts());
247           
248         }
249 }
250
251 //_____________________________________________________________________
252 void AliFMDQADataMakerRec::MakeRecPoints(TTree* clustersTree)
253 {
254   // makes data from RecPoints
255   AliFMDParameters* pars = AliFMDParameters::Instance();
256   fRecPointsArray.Clear();
257   TBranch *fmdbranch = clustersTree->GetBranch("FMD");
258   if (!fmdbranch) { 
259     AliError("can't get the branch with the FMD recpoints !");
260     return;
261   }
262   
263   TClonesArray* RecPointsAddress = &fRecPointsArray;
264   
265   fmdbranch->SetAddress(&RecPointsAddress);
266   fmdbranch->GetEntry(0);
267   TIter next(RecPointsAddress) ; 
268   AliFMDRecPoint * rp ; 
269   while ((rp = static_cast<AliFMDRecPoint*>(next()))) {
270     
271     GetRecPointsData(0)->Fill(rp->Edep()/pars->GetEdepMip()) ;
272   
273   }
274
275 }
276
277 //_____________________________________________________________________ 
278 void AliFMDQADataMakerRec::StartOfDetectorCycle()
279 {
280   // What 
281   // to 
282   // do?
283 }
284 //_____________________________________________________________________ 
285 Int_t AliFMDQADataMakerRec::GetHalfringIndex(UShort_t det, 
286                                              Char_t ring, 
287                                              UShort_t board) {
288   
289   UShort_t iring  =  (ring == 'I' ? 1 : 0);
290   
291   Int_t index = (((det-1) << 2) | (iring << 1) | (board << 0));
292   
293   return index-2;
294   
295 }
296
297 //_____________________________________________________________________ 
298
299
300 //
301 // EOF
302 //