]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSQADataMakerRec.cxx
A protection against a missing OCDB is introduced for bad channel map.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSQADataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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
16
17 /* $Id$ */
18
19 /*
20   Produces the data needed to calculate the quality assurance. 
21   All data must be mergeable objects.
22   Y. Schutz CERN July 2007
23 */
24
25 // --- ROOT system ---
26 #include <TClonesArray.h>
27 #include <TFile.h> 
28 #include <TH1F.h> 
29 #include <TH1I.h> 
30 #include <TH2F.h> 
31
32 // --- Standard library ---
33
34 // --- AliRoot header files ---
35 #include "AliESDCaloCluster.h"
36 #include "AliESDEvent.h"
37 #include "AliLog.h"
38 #include "AliPHOSQADataMakerRec.h"
39 #include "AliQAChecker.h"
40 #include "AliPHOSCpvRecPoint.h" 
41 #include "AliPHOSEmcRecPoint.h" 
42 #include "AliPHOSRecParticle.h" 
43 #include "AliPHOSTrackSegment.h" 
44 #include "AliPHOSRawDecoder.h"
45 #include "AliPHOSReconstructor.h"
46 #include "AliPHOSRecoParam.h"
47
48 ClassImp(AliPHOSQADataMakerRec)
49            
50 //____________________________________________________________________________ 
51   AliPHOSQADataMakerRec::AliPHOSQADataMakerRec() : 
52   AliQADataMakerRec(AliQA::GetDetName(AliQA::kPHOS), "PHOS Quality Assurance Data Maker")
53 {
54   // ctor
55 }
56
57 //____________________________________________________________________________ 
58 AliPHOSQADataMakerRec::AliPHOSQADataMakerRec(const AliPHOSQADataMakerRec& qadm) :
59   AliQADataMakerRec()
60 {
61   //copy ctor 
62   SetName((const char*)qadm.GetName()) ; 
63   SetTitle((const char*)qadm.GetTitle()); 
64 }
65
66 //__________________________________________________________________
67 AliPHOSQADataMakerRec& AliPHOSQADataMakerRec::operator = (const AliPHOSQADataMakerRec& qadm )
68 {
69   // Equal operator.
70   this->~AliPHOSQADataMakerRec();
71   new(this) AliPHOSQADataMakerRec(qadm);
72   return *this;
73 }
74  
75 //____________________________________________________________________________ 
76 void AliPHOSQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX task, TObjArray * list)
77 {
78   //Detector specific actions at end of cycle
79   // do the QA checking
80   AliQAChecker::Instance()->Run(AliQA::kPHOS, task, list) ;  
81 }
82
83 //____________________________________________________________________________ 
84 void AliPHOSQADataMakerRec::InitESDs()
85 {
86   //create ESDs histograms in ESDs subdir
87         
88   TH1F * h0 = new TH1F("hPhosESDs",    "ESDs energy distribution in PHOS",       100, 0., 100.) ;  
89   h0->Sumw2() ; 
90   Add2ESDsList(h0, 0) ;
91   TH1I * h1  = new TH1I("hPhosESDsMul", "ESDs multiplicity distribution in PHOS", 100, 0., 100) ; 
92   h1->Sumw2() ;
93   Add2ESDsList(h1, 1) ;
94 }
95
96 //____________________________________________________________________________ 
97 //void AliPHOSQADataMakerRec::InitRecParticles()
98 //{
99 //  // create Reconstructed particles histograms in RecParticles subdir
100 //  fhRecParticles     = new TH1F("hPhosRecParticles",    "RecParticles energy distribution in PHOS",       100, 0., 100.) ; 
101 //  fhRecParticles->Sumw2() ;
102 //  fhRecParticlesMul  = new TH1I("hPhosRecParticlesMul", "RecParticles multiplicity distribution in PHOS", 100, 0,  100) ; 
103 //  fhRecParticlesMul->Sumw2() ;
104 //}
105
106 //____________________________________________________________________________ 
107 void AliPHOSQADataMakerRec::InitRecPoints()
108 {
109   // create Reconstructed Points histograms in RecPoints subdir
110   TH1F * h0 = new TH1F("hEmcPhosRecPoints",    "EMCA RecPoints energy distribution in PHOS",       100, 0., 100.) ; 
111   h0->Sumw2() ;
112   Add2RecPointsList(h0, 0) ;
113   TH1I * h1 = new TH1I("hEmcPhosRecPointsMul", "EMCA RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
114   h1->Sumw2() ;
115   Add2RecPointsList(h1, 1) ;
116
117   TH1F * h2 = new TH1F("hCpvPhosRecPoints",    "CPV RecPoints energy distribution in PHOS",       100, 0., 100.) ; 
118   h2->Sumw2() ;
119   Add2RecPointsList(h2, 2) ;
120   TH1I * h3 = new TH1I("hCpvPhosRecPointsMul", "CPV RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
121   h3->Sumw2() ;
122   Add2RecPointsList(h3, 3) ;
123 }
124
125 //____________________________________________________________________________ 
126 void AliPHOSQADataMakerRec::InitRaws()
127 {
128   // create Raws histograms in Raws subdir
129   const Int_t modMax = 5 ; 
130   TH2I * h0[modMax*2] ; 
131   char name[32] ; 
132   char title[32] ; 
133   for (Int_t mod = 0; mod < modMax; mod++) {
134    sprintf(title, "Low Gain Rows x Columns for PHOS module %d", mod) ;  
135    sprintf(name, "hLowPHOSxyMod%d", mod) ; 
136    h0[mod] = new TH2I(name, title, 64, 1, 65, 56, 1, 57) ; 
137    Add2RawsList(h0[mod], mod) ;
138    sprintf(title, "High Gain Rows x Columns for PHOS module %d", mod) ;  
139    sprintf(name, "hHighPHOSxyMod%d", mod) ; 
140    h0[mod+modMax] = new TH2I(name, title, 64, 1, 65, 56, 1, 57) ; 
141    Add2RawsList(h0[mod+modMax], mod+modMax) ;
142   }
143   TH1I * h10 = new TH1I("hLowPhosModules",    "Low Gain Hits in EMCA PHOS modules",       6, 0, 6) ; 
144   h10->Sumw2() ;
145   Add2RawsList(h10, 10) ;
146   TH1I * h11 = new TH1I("hHighPhosModules",    "High Gain Hits in EMCA PHOS modules",       6, 0, 6) ; 
147   h11->Sumw2() ;
148   Add2RawsList(h11, 11) ;
149   TH1F * h12 = new TH1F("hLowPhosRawtime", "Low Gain Time of raw hits in PHOS", 100, 0, 100.) ; 
150   h12->Sumw2() ;
151   Add2RawsList(h12, 12) ;
152   TH1F * h13 = new TH1F("hHighPhosRawtime", "High Gain Time of raw hits in PHOS", 100, 0, 100.) ; 
153   h13->Sumw2() ;
154   Add2RawsList(h13, 13) ;
155   TH1F * h14 = new TH1F("hLowPhosRawEnergy", "Low Gain Energy of raw hits in PHOS", 100, 0., 100.) ; 
156   h14->Sumw2() ;
157   Add2RawsList(h14, 14) ;
158   TH1F * h15 = new TH1F("hHighPhosRawEnergy", "High Gain Energy of raw hits in PHOS", 100, 0., 100.) ; 
159   h15->Sumw2() ;
160   Add2RawsList(h15, 15) ;
161  
162 }
163
164 //____________________________________________________________________________ 
165 //void AliPHOSQADataMakerRec::InitTrackSegments()
166 //{
167 //  // create Track Segments histograms in TrackSegments subdir
168 //  fhTrackSegments     = new TH1F("hPhosTrackSegments",    "TrackSegments EMC-CPV distance in PHOS",       500, 0., 5000.) ; 
169 //  fhTrackSegments->Sumw2() ;
170 //  fhTrackSegmentsMul  = new TH1I("hPhosTrackSegmentsMul", "TrackSegments multiplicity distribution in PHOS", 100, 0,  100) ; 
171 //  fhTrackSegmentsMul->Sumw2() ;
172 //}
173
174 //____________________________________________________________________________
175 void AliPHOSQADataMakerRec::MakeESDs(AliESDEvent * esd)
176 {
177   // make QA data from ESDs
178
179   Int_t count = 0 ; 
180   for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
181         AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
182         if ( clu->IsPHOS() ) {
183                 GetESDsData(0)->Fill(clu->E()) ;
184                 count++ ;
185         } 
186   }
187   GetESDsData(1)->Fill(count) ;
188 }
189
190 //____________________________________________________________________________
191 // void AliPHOSQADataMakerRec::MakeRecParticles(TTree * recpar)
192 // {
193 //   // makes data from RecParticles
194
195 //   TClonesArray * recparticles = dynamic_cast<TClonesArray*>(fData) ; 
196 //   fhRecParticlesMul->Fill(recparticles->GetEntriesFast()) ; 
197 //   TIter next(recparticles) ; 
198 //   AliPHOSRecParticle * recparticle ; 
199 //   while ( (recparticle = dynamic_cast<AliPHOSRecParticle *>(next())) ) {
200 //     fhRecParticles->Fill( recparticle->Energy()) ;
201 //   }
202 // }
203
204 //____________________________________________________________________________
205 void AliPHOSQADataMakerRec::MakeRaws(AliRawReader* rawReader)
206 {
207   const Int_t modMax = 5 ; 
208   rawReader->Reset() ; 
209   AliPHOSRawDecoder decoder(rawReader);
210   decoder.SetOldRCUFormat(kFALSE);
211   decoder.SubtractPedestals(AliPHOSReconstructor::GetRecoParamEmc()->SubtractPedestals());
212   Int_t count = 0 ; 
213   while (decoder.NextDigit()) {
214    Int_t module  = decoder.GetModule() ;
215    Int_t row     = decoder.GetRow() ;
216    Int_t col     = decoder.GetColumn() ;
217    Double_t time = decoder.GetTime() ;
218    Double_t energy  = decoder.GetEnergy() ;     
219    Bool_t lowGain = decoder.IsLowGain();
220    if (lowGain) {
221      GetRawsData(module)->Fill(row, col) ; 
222          GetRawsData(10)->Fill(module) ; 
223      GetRawsData(12)->Fill(time) ; 
224      GetRawsData(14)->Fill(energy) ; 
225    } else {
226          GetRawsData(module+modMax)->Fill(row, col) ; 
227          GetRawsData(11)->Fill(module) ; 
228      GetRawsData(13)->Fill(time) ; 
229      GetRawsData(15)->Fill(energy) ; 
230    }
231    //AliInfo(Form(" %d %d %d %d %f %f\n", count, module, row, col, time, energy)) ;
232    count++ ; 
233   } 
234 }
235
236 //____________________________________________________________________________
237 void AliPHOSQADataMakerRec::MakeRecPoints(TTree * clustersTree)
238 {
239   {
240     // makes data from RecPoints
241     TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
242     if (!emcbranch) { 
243       AliError("can't get the branch with the PHOS EMC clusters !");
244       return;
245     }
246     TObjArray * emcrecpoints = new TObjArray(100) ;
247     emcbranch->SetAddress(&emcrecpoints);
248     emcbranch->GetEntry(0);
249     
250     GetRecPointsData(1)->Fill(emcrecpoints->GetEntriesFast()) ; 
251     TIter next(emcrecpoints) ; 
252     AliPHOSEmcRecPoint * rp ; 
253     while ( (rp = dynamic_cast<AliPHOSEmcRecPoint *>(next())) ) {
254       GetRecPointsData(0)->Fill( rp->GetEnergy()) ;
255     }
256     emcrecpoints->Delete();
257     delete emcrecpoints;
258   }
259   {
260     TBranch *cpvbranch = clustersTree->GetBranch("PHOSCpvRP");
261     if (!cpvbranch) { 
262       AliError("can't get the branch with the PHOS CPV clusters !");
263       return;
264     }
265     TObjArray *cpvrecpoints = new TObjArray(100) ;
266     cpvbranch->SetAddress(&cpvrecpoints);
267     cpvbranch->GetEntry(0);
268     
269     GetRecPointsData(1)->Fill(cpvrecpoints->GetEntriesFast()) ; 
270     TIter next(cpvrecpoints) ; 
271     AliPHOSCpvRecPoint * rp ; 
272     while ( (rp = dynamic_cast<AliPHOSCpvRecPoint *>(next())) ) {
273       GetRecPointsData(0)->Fill( rp->GetEnergy()) ;
274     }
275     cpvrecpoints->Delete();
276     delete cpvrecpoints;
277   }
278 }
279
280 //____________________________________________________________________________
281 // void AliPHOSQADataMakerRec::MakeTrackSegments(TTree * ts)
282 // {
283 //   // makes data from TrackSegments
284
285 //   TClonesArray * tracksegments = dynamic_cast<TClonesArray*>(fData) ;
286
287 //   fhTrackSegmentsMul->Fill(tracksegments->GetEntriesFast()) ; 
288 //   TIter next(tracksegments) ; 
289 //   AliPHOSTrackSegment * ts ; 
290 //   while ( (ts = dynamic_cast<AliPHOSTrackSegment *>(next())) ) {
291 //     fhTrackSegments->Fill( ts->GetCpvDistance()) ;
292 //   } 
293 // }
294
295 //____________________________________________________________________________ 
296 void AliPHOSQADataMakerRec::StartOfDetectorCycle()
297 {
298   //Detector specific actions at start of cycle
299   
300 }