]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSQADataMakerRec.cxx
Bad channels map and HG/LG are calculated in STANDALONE run only;
[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 #include <TParameter.h>
32
33 // --- Standard library ---
34
35 // --- AliRoot header files ---
36 #include "AliESDCaloCluster.h"
37 #include "AliESDEvent.h"
38 #include "AliLog.h"
39 #include "AliPHOSQADataMakerRec.h"
40 #include "AliQAChecker.h"
41 #include "AliPHOSCpvRecPoint.h" 
42 #include "AliPHOSEmcRecPoint.h" 
43 #include "AliPHOSRecParticle.h" 
44 #include "AliPHOSTrackSegment.h" 
45 #include "AliPHOSRawDecoder.h"
46 #include "AliPHOSRawDecoderv1.h"
47 #include "AliPHOSRawDecoderv2.h"
48 #include "AliPHOSReconstructor.h"
49
50 ClassImp(AliPHOSQADataMakerRec)
51            
52 //____________________________________________________________________________ 
53   AliPHOSQADataMakerRec::AliPHOSQADataMakerRec() : 
54   AliQADataMakerRec(AliQA::GetDetName(AliQA::kPHOS), "PHOS Quality Assurance Data Maker")
55 {
56   // ctor
57 }
58
59 //____________________________________________________________________________ 
60 AliPHOSQADataMakerRec::AliPHOSQADataMakerRec(const AliPHOSQADataMakerRec& qadm) :
61   AliQADataMakerRec()
62 {
63   //copy ctor 
64   SetName((const char*)qadm.GetName()) ; 
65   SetTitle((const char*)qadm.GetTitle()); 
66 }
67
68 //__________________________________________________________________
69 AliPHOSQADataMakerRec& AliPHOSQADataMakerRec::operator = (const AliPHOSQADataMakerRec& qadm )
70 {
71   // Equal operator.
72   this->~AliPHOSQADataMakerRec();
73   new(this) AliPHOSQADataMakerRec(qadm);
74   return *this;
75 }
76  
77 //____________________________________________________________________________ 
78 void AliPHOSQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
79 {
80   //Detector specific actions at end of cycle
81   // do the QA checking
82   AliQAChecker::Instance()->Run(AliQA::kPHOS, task, list) ;  
83 }
84
85 //____________________________________________________________________________ 
86 void AliPHOSQADataMakerRec::InitESDs()
87 {
88   //Create histograms to controll ESD
89  
90   TH1F * h1 = new TH1F("hESDPhosSpectrum",  "ESDs spectrum in PHOS"                ,  200, 0.,   20.) ;
91   h1->Sumw2() ;
92   Add2ESDsList(h1, kESDSpec)  ;
93
94   TH1I * h2 = new TH1I("hESDPhosMul",       "ESDs multiplicity distribution in PHOS", 100, 0,   100 ) ; 
95   h2->Sumw2() ;
96   Add2ESDsList(h2, kESDNtot) ;
97  
98   TH1F * h3 = new TH1F("hESDPhosEtot",      "ESDs total energy"                     , 100, 0,  1000.) ; 
99   h3->Sumw2() ;
100   Add2ESDsList(h3, kESDEtot) ;
101  
102   TH1F * h4 = new TH1F("hESDpid",           "ESDs PID distribution in PHOS"         , 100, 0.,    1.) ;
103   h4->Sumw2() ;
104   Add2ESDsList(h4, kESDpid) ;
105         
106 }
107
108 //____________________________________________________________________________ 
109 void AliPHOSQADataMakerRec::InitRecPoints()
110 {
111   // create Reconstructed Points histograms in RecPoints subdir
112   TH2I * h0 = new TH2I("hRpPHOSxyMod1","RecPoints Rows x Columns for PHOS module 1", 64, -72., 72., 56, -63., 63.) ;                             
113   Add2RecPointsList(h0,kRPmod1) ;
114   TH2I * h1 = new TH2I("hRpPHOSxyMod2","RecPoints Rows x Columns for PHOS module 2", 64, -72., 72., 56, -63., 63.) ;                             
115   Add2RecPointsList(h1,kRPmod2) ;
116   TH2I * h2 = new TH2I("hRpPHOSxyMod3","RecPoints Rows x Columns for PHOS module 3", 64, -72., 72., 56, -63., 63.) ;                             
117   Add2RecPointsList(h2,kRPmod3) ;
118   TH2I * h3 = new TH2I("hRpPHOSxyMod4","RecPoints Rows x Columns for PHOS module 4", 64, -72., 72., 56, -63., 63.) ;                             
119   Add2RecPointsList(h3,kRPmod4) ;
120   TH2I * h4 = new TH2I("hRpPHOSxyMod5","RecPoints Rows x Columns for PHOS module 5", 64, -72., 72., 56, -63., 63.) ;                             
121   Add2RecPointsList(h4,kRPmod5) ;
122  
123   TH1F * h5 = new TH1F("hEmcPhosRecPointsSpectrum",  "EMC RecPoints spectrum in PHOS",   2000, 0., 20.) ; 
124   h5->Sumw2() ;
125   Add2RecPointsList(h5, kRPSpec)  ;
126
127   TH1I * h6 = new TH1I("hEmcPhosRecPointsMul", "EMCA RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
128   h6->Sumw2() ;
129   Add2RecPointsList(h6, kRPNtot) ;
130
131   TH1I * h7 = new TH1I("hEmcPhosRecPointsEtot", "EMC RecPoints Etot", 200, 0,  200.) ; 
132   h7->Sumw2() ;
133   Add2RecPointsList(h7, kRPEtot) ;
134
135   TH1I * h8 = new TH1I("hCpvPhosRecPointsMul", "CPV RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
136   h8->Sumw2() ;
137   Add2RecPointsList(h8, kRPNcpv) ;
138 }
139
140 //____________________________________________________________________________ 
141 void AliPHOSQADataMakerRec::InitRaws()
142 {
143   // create Raws histograms in Raws subdir
144   TH2I * h0 = new TH2I("hHighPHOSxyMod1","High Gain Rows x Columns for PHOS module 1", 64, 0, 64, 56, 0, 56) ;
145   Add2RawsList(h0,kHGmod1) ;
146   TH2I * h1 = new TH2I("hHighPHOSxyMod2","High Gain Rows x Columns for PHOS module 2", 64, 0, 64, 56, 0, 56) ;
147   Add2RawsList(h1,kHGmod2) ;
148   TH2I * h2 = new TH2I("hHighPHOSxyMod3","High Gain Rows x Columns for PHOS module 3", 64, 0, 64, 56, 0, 56) ;
149   Add2RawsList(h2,kHGmod3) ;
150   TH2I * h3 = new TH2I("hHighPHOSxyMod4","High Gain Rows x Columns for PHOS module 4", 64, 0, 64, 56, 0, 56) ;
151   Add2RawsList(h3,kHGmod4) ;
152   TH2I * h4 = new TH2I("hHighPHOSxyMod5","High Gain Rows x Columns for PHOS module 5", 64, 0, 64, 56, 0, 56) ;
153   Add2RawsList(h4,kHGmod5) ;
154   TH2I * h5 = new TH2I("hLowPHOSxyMod1","Low Gain Rows x Columns for PHOS module 1", 64, 0, 64, 56, 0, 56) ;
155   Add2RawsList(h5,kLGmod1) ;
156   TH2I * h6 = new TH2I("hLowPHOSxyMod2","Low Gain Rows x Columns for PHOS module 2", 64, 0, 64, 56, 0, 56) ;
157   Add2RawsList(h6,kLGmod2) ;
158   TH2I * h7 = new TH2I("hLowPHOSxyMod3","Low Gain Rows x Columns for PHOS module 3", 64, 0, 64, 56, 0, 56) ;
159   Add2RawsList(h7,kLGmod3) ;
160   TH2I * h8 = new TH2I("hLowPHOSxyMod4","Low Gain Rows x Columns for PHOS module 4", 64, 0, 64, 56, 0, 56) ;
161   Add2RawsList(h8,kLGmod4) ;                                                                                                               
162   TH2I * h9 = new TH2I("hLowPHOSxyMod5","Low Gain Rows x Columns for PHOS module 5", 64, 0, 64, 56, 0, 56) ;                               
163   Add2RawsList(h9,kLGmod5) ;                                                                                                               
164                                                                                                                                            
165                                                                                                                                            
166   TH1I * h10 = new TH1I("hLowPhosModules",    "Low Gain Hits in EMCA PHOS modules",       6, 0, 6) ;                                       
167   h10->Sumw2() ;                                                                                                                           
168   Add2RawsList(h10, kNmodLG) ;                                                                                                             
169   TH1I * h11 = new TH1I("hHighPhosModules",   "High Gain Hits in EMCA PHOS modules",       6, 0, 6) ;                                      
170   h11->Sumw2() ;                                                                                                                           
171   Add2RawsList(h11, kNmodHG) ;                                                                                                             
172                                                                                                                                            
173   TH1F * h12 = new TH1F("hLowPhosRawtime", "Low Gain Time of raw hits in PHOS", 500, -50., 200.) ;                                            
174   h12->Sumw2() ;                                                                                                                           
175   Add2RawsList(h12, kLGtime) ;                                                                                                             
176   TH1F * h13 = new TH1F("hHighPhosRawtime", "High Gain Time of raw hits in PHOS", 500, -50., 200.) ;                                          
177   h13->Sumw2() ;                                                                                                                           
178   Add2RawsList(h13, kHGtime) ;                                                                                                             
179                                                                                                                                            
180   TH1F * h14 = new TH1F("hLowPhosRawEnergy", "Low Gain Energy of raw hits in PHOS", 500, 0., 1000.) ;                                      
181   h14->Sumw2() ;                                                                                                                           
182   Add2RawsList(h14, kSpecLG) ;                                                                                                             
183   TH1F * h15 = new TH1F("hHighPhosRawEnergy", "High Gain Energy of raw hits in PHOS",500,0., 1000.) ;                                      
184   h15->Sumw2() ;                                                                                                                           
185   Add2RawsList(h15, kSpecHG) ;                                                                                                             
186                                                                                                                                            
187   TH1F * h16 = new TH1F("hLowNtot", "Low Gain Total Number of raw hits in PHOS", 500, 0., 5000.) ;                                         
188   h16->Sumw2() ;                                                                                                                           
189   Add2RawsList(h16, kNtotLG, kTRUE) ;                                                                                                             
190   TH1F * h17 = new TH1F("hHighNtot", "High Gain Total Number of raw hits in PHOS",500,0., 5000.) ;                                         
191   h17->Sumw2() ;                                                                                                                           
192   Add2RawsList(h17, kNtotHG, kTRUE) ;                                                                                                             
193                                                                                                                                            
194   TH1F * h18 = new TH1F("hLowEtot", "Low Gain Total Energy of raw hits in PHOS", 500, 0., 5000.) ;                                       
195   h18->Sumw2() ;                                                                                                                           
196   Add2RawsList(h18, kEtotLG, kTRUE) ;                                                                                                             
197   TH1F * h19 = new TH1F("hHighEtot", "High Gain Total Energy of raw hits in PHOS",500,0., 100000.) ;                                       
198   h19->Sumw2() ;                                                                                                                           
199   Add2RawsList(h19, kEtotHG, kTRUE) ;                                                                                                             
200   
201 }
202
203 //____________________________________________________________________________
204 void AliPHOSQADataMakerRec::MakeESDs(AliESDEvent * esd)
205 {
206   // make QA data from ESDs
207
208   Int_t nTot = 0 ; 
209   Double_t eTot = 0 ; 
210   for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
211     AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
212     if( clu->IsPHOS() ) {
213       GetESDsData(kESDSpec)->Fill(clu->E()) ;
214       Double_t *pid=clu->GetPid() ;
215       GetESDsData(kESDpid)->Fill(pid[AliPID::kPhoton]) ;
216       eTot+=clu->E() ;
217       nTot++ ;
218     } 
219   }
220   GetESDsData(kESDNtot)->Fill(nTot) ;
221   GetESDsData(kESDEtot)->Fill(eTot) ;
222 }
223
224 //____________________________________________________________________________
225 void AliPHOSQADataMakerRec::MakeRaws(AliRawReader* rawReader)
226 {
227   //Fill prepared histograms with Raw digit properties
228   rawReader->Reset() ;
229   AliPHOSRawDecoder * decoder ;
230   if(strcmp(GetRecoParam()->EMCDecoderVersion(),"v1")==0)
231     decoder=new AliPHOSRawDecoderv1(rawReader);
232   else
233     if(strcmp(GetRecoParam()->EMCDecoderVersion(),"v2")==0)
234       decoder=new AliPHOSRawDecoderv2(rawReader);
235     else
236       decoder=new AliPHOSRawDecoder(rawReader);
237   decoder->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
238   Double_t lgEtot=0. ;
239   Double_t hgEtot=0. ;
240   Int_t lgNtot=0 ;
241   Int_t hgNtot=0 ;
242
243   while (decoder->NextDigit()) {
244    Int_t module  = decoder->GetModule() ;
245    Int_t row     = decoder->GetRow() ;
246    Int_t col     = decoder->GetColumn() ;
247    Double_t time = decoder->GetTime() ;
248    Double_t energy  = decoder->GetEnergy() ;
249    Bool_t lowGain = decoder->IsLowGain();
250    if (lowGain) {
251      if(energy<2.)
252        continue ;
253      switch(module){
254         case 1: GetRawsData(kLGmod1)->Fill(row-0.5,col-0.5) ; break ;
255         case 2: GetRawsData(kLGmod2)->Fill(row-0.5,col-0.5) ; break ;
256         case 3: GetRawsData(kLGmod3)->Fill(row-0.5,col-0.5) ; break ;
257         case 4: GetRawsData(kLGmod4)->Fill(row-0.5,col-0.5) ; break ;
258         case 5: GetRawsData(kLGmod5)->Fill(row-0.5,col-0.5) ; break ;
259      }                                  
260      GetRawsData(kNmodLG)->Fill(module) ;
261      GetRawsData(kLGtime)->Fill(time) ; 
262      GetRawsData(kSpecLG)->Fill(energy) ;    
263      lgEtot+=energy ;
264      lgNtot++ ;   
265    } else {        
266      if(energy<8.)
267        continue ;
268      switch (module){
269          case 1: GetRawsData(kHGmod1)->Fill(row-0.5,col-0.5) ; break ;
270          case 2: GetRawsData(kHGmod2)->Fill(row-0.5,col-0.5) ; break ;
271          case 3: GetRawsData(kHGmod3)->Fill(row-0.5,col-0.5) ; break ;
272          case 4: GetRawsData(kHGmod4)->Fill(row-0.5,col-0.5) ; break ;
273          case 5: GetRawsData(kHGmod5)->Fill(row-0.5,col-0.5) ; break ;
274      }              
275      GetRawsData(kNmodHG)->Fill(module) ; 
276      GetRawsData(kHGtime)->Fill(time) ;  
277      GetRawsData(kSpecHG)->Fill(energy) ;
278      hgEtot+=energy ; 
279      hgNtot++ ;  
280    }                 
281   }                    
282   GetRawsData(kEtotLG)->Fill(lgEtot) ; 
283   TParameter<double> * p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQA::GetTaskName(AliQA::kRAWS).Data(), GetRawsData(kEtotLG)->GetName()))) ; 
284   p->SetVal(lgEtot) ; 
285   GetRawsData(kEtotHG)->Fill(hgEtot) ;  
286   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQA::GetTaskName(AliQA::kRAWS).Data(), GetRawsData(kEtotHG)->GetName()))) ; 
287   p->SetVal(hgEtot) ; 
288   GetRawsData(kNtotLG)->Fill(lgNtot) ;
289   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQA::GetTaskName(AliQA::kRAWS).Data(), GetRawsData(kNtotLG)->GetName()))) ; 
290   p->SetVal(lgNtot) ; 
291   GetRawsData(kNtotHG)->Fill(hgNtot) ;
292   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQA::GetTaskName(AliQA::kRAWS).Data(), GetRawsData(kNtotHG)->GetName()))) ; 
293   p->SetVal(hgNtot) ; 
294 }
295 //____________________________________________________________________________
296 void AliPHOSQADataMakerRec::MakeRecPoints(TTree * clustersTree)
297 {
298   {
299     // makes data from RecPoints
300     TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
301     if (!emcbranch) { 
302       AliError("can't get the branch with the PHOS EMC clusters !");
303       return;
304     }
305     TObjArray * emcrecpoints = new TObjArray(100) ;
306     emcbranch->SetAddress(&emcrecpoints);
307     emcbranch->GetEntry(0);
308     
309     GetRecPointsData(kRPNtot)->Fill(emcrecpoints->GetEntriesFast()) ; 
310     TIter next(emcrecpoints) ; 
311     AliPHOSEmcRecPoint * rp ; 
312     Double_t eTot = 0. ; 
313     while ( (rp = dynamic_cast<AliPHOSEmcRecPoint *>(next())) ) {
314       GetRecPointsData(kRPSpec)->Fill( rp->GetEnergy()) ;
315       Int_t mod = rp->GetPHOSMod() ;
316       TVector3 pos ;
317       rp->GetLocalPosition(pos) ;
318       switch(mod){
319         case 1: GetRecPointsData(kRPmod1)->Fill(pos.X(),pos.Z()) ; break ;
320         case 2: GetRecPointsData(kRPmod2)->Fill(pos.X(),pos.Z()) ; break ;
321         case 3: GetRecPointsData(kRPmod3)->Fill(pos.X(),pos.Z()) ; break ;
322         case 4: GetRecPointsData(kRPmod4)->Fill(pos.X(),pos.Z()) ; break ;
323         case 5: GetRecPointsData(kRPmod5)->Fill(pos.X(),pos.Z()) ; break ;
324       }
325       eTot+= rp->GetEnergy() ;
326     }
327     GetRecPointsData(kRPEtot)->Fill(eTot) ;
328     emcrecpoints->Delete();
329     delete emcrecpoints;
330   }
331   {
332     TBranch *cpvbranch = clustersTree->GetBranch("PHOSCpvRP");
333     if (!cpvbranch) { 
334       AliError("can't get the branch with the PHOS CPV clusters !");
335       return;
336     }
337     TObjArray *cpvrecpoints = new TObjArray(100) ;
338     cpvbranch->SetAddress(&cpvrecpoints);
339     cpvbranch->GetEntry(0);
340     
341     GetRecPointsData(kRPNcpv)->Fill(cpvrecpoints->GetEntriesFast()) ; 
342     cpvrecpoints->Delete();
343     delete cpvrecpoints;
344   }
345 }
346
347 //____________________________________________________________________________ 
348 void AliPHOSQADataMakerRec::StartOfDetectorCycle()
349 {
350   //Detector specific actions at start of cycle
351   
352 }