]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSQADataMakerRec.cxx
Expert/shifter separation added, new raw histos added
[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   if(GetRawsData(kHGqualMod1) && GetRawsData(kHGmod1))
82     GetRawsData(kHGqualMod1)->Divide( GetRawsData(kHGmod1) ) ;
83   if(GetRawsData(kHGqualMod2) && GetRawsData(kHGmod2))
84     GetRawsData(kHGqualMod2)->Divide( GetRawsData(kHGmod2) ) ;
85   if(GetRawsData(kHGqualMod3) && GetRawsData(kHGmod3))
86     GetRawsData(kHGqualMod3)->Divide( GetRawsData(kHGmod3) ) ;
87   if(GetRawsData(kHGqualMod4) && GetRawsData(kHGmod4))
88     GetRawsData(kHGqualMod4)->Divide( GetRawsData(kHGmod4) ) ;
89   if(GetRawsData(kHGqualMod5) && GetRawsData(kHGmod5))
90     GetRawsData(kHGqualMod5)->Divide( GetRawsData(kHGmod5) ) ;
91   // do the QA checking
92   AliQAChecker::Instance()->Run(AliQA::kPHOS, task, list) ;  
93 }
94
95 //____________________________________________________________________________ 
96 void AliPHOSQADataMakerRec::InitESDs()
97 {
98   //Create histograms to controll ESD
99  
100   TH1F * h1 = new TH1F("hESDPhosSpectrum",  "ESDs spectrum in PHOS"                ,  200, 0.,   20.) ;
101   h1->Sumw2() ;
102   Add2ESDsList(h1, kESDSpec)  ;
103
104   TH1I * h2 = new TH1I("hESDPhosMul",       "ESDs multiplicity distribution in PHOS", 100, 0,   100 ) ; 
105   h2->Sumw2() ;
106   Add2ESDsList(h2, kESDNtot) ;
107  
108   TH1F * h3 = new TH1F("hESDPhosEtot",      "ESDs total energy"                     , 100, 0,  1000.) ; 
109   h3->Sumw2() ;
110   Add2ESDsList(h3, kESDEtot,1) ;  //Expert histo
111  
112   TH1F * h4 = new TH1F("hESDpid",           "ESDs PID distribution in PHOS"         , 100, 0.,    1.) ;
113   h4->Sumw2() ;
114   Add2ESDsList(h4, kESDpid,1) ; //Expert histo
115         
116 }
117
118 //____________________________________________________________________________ 
119 void AliPHOSQADataMakerRec::InitRecPoints()
120 {
121   // create Reconstructed Points histograms in RecPoints subdir
122   TH2I * h0 = new TH2I("hRpPHOSxyMod1","RecPoints Rows x Columns for PHOS module 1", 64, -72., 72., 56, -63., 63.) ;                             
123   Add2RecPointsList(h0,kRPmod1,1) ;
124   TH2I * h1 = new TH2I("hRpPHOSxyMod2","RecPoints Rows x Columns for PHOS module 2", 64, -72., 72., 56, -63., 63.) ;                             
125   Add2RecPointsList(h1,kRPmod2,1) ;
126   TH2I * h2 = new TH2I("hRpPHOSxyMod3","RecPoints Rows x Columns for PHOS module 3", 64, -72., 72., 56, -63., 63.) ;                             
127   Add2RecPointsList(h2,kRPmod3) ;
128   TH2I * h3 = new TH2I("hRpPHOSxyMod4","RecPoints Rows x Columns for PHOS module 4", 64, -72., 72., 56, -63., 63.) ;                             
129   Add2RecPointsList(h3,kRPmod4,1) ;
130   TH2I * h4 = new TH2I("hRpPHOSxyMod5","RecPoints Rows x Columns for PHOS module 5", 64, -72., 72., 56, -63., 63.) ;                             
131   Add2RecPointsList(h4,kRPmod5,1) ;
132  
133   TH1F * h5 = new TH1F("hEmcPhosRecPointsSpectrum",  "EMC RecPoints spectrum in PHOS",   2000, 0., 20.) ; 
134   h5->Sumw2() ;
135   Add2RecPointsList(h5, kRPSpec)  ;
136
137   TH1I * h6 = new TH1I("hEmcPhosRecPointsMul", "EMCA RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
138   h6->Sumw2() ;
139   Add2RecPointsList(h6, kRPNtot) ;
140
141   TH1I * h7 = new TH1I("hEmcPhosRecPointsEtot", "EMC RecPoints Etot", 200, 0,  200.) ; 
142   h7->Sumw2() ;
143   Add2RecPointsList(h7, kRPEtot,1) ;
144
145   TH1I * h8 = new TH1I("hCpvPhosRecPointsMul", "CPV RecPoints multiplicity distribution in PHOS", 100, 0,  100) ; 
146   h8->Sumw2() ;
147   Add2RecPointsList(h8, kRPNcpv,1) ;
148 }
149
150 //____________________________________________________________________________ 
151 void AliPHOSQADataMakerRec::InitRaws()
152 {
153   // create Raws histograms in Raws subdir
154   TH2I * h0 = new TH2I("hHighPHOSxyMod1","High Gain Rows x Columns for PHOS module 1", 64, 0, 64, 56, 0, 56) ;
155   Add2RawsList(h0,kHGmod1,1) ;
156   TH2I * h1 = new TH2I("hHighPHOSxyMod2","High Gain Rows x Columns for PHOS module 2", 64, 0, 64, 56, 0, 56) ;
157   Add2RawsList(h1,kHGmod2,1) ;
158   TH2I * h2 = new TH2I("hHighPHOSxyMod3","High Gain Rows x Columns for PHOS module 3", 64, 0, 64, 56, 0, 56) ;
159   Add2RawsList(h2,kHGmod3) ;
160   TH2I * h3 = new TH2I("hHighPHOSxyMod4","High Gain Rows x Columns for PHOS module 4", 64, 0, 64, 56, 0, 56) ;
161   Add2RawsList(h3,kHGmod4,1) ;
162   TH2I * h4 = new TH2I("hHighPHOSxyMod5","High Gain Rows x Columns for PHOS module 5", 64, 0, 64, 56, 0, 56) ;
163   Add2RawsList(h4,kHGmod5,1) ;
164   TH2I * h5 = new TH2I("hLowPHOSxyMod1","Low Gain Rows x Columns for PHOS module 1", 64, 0, 64, 56, 0, 56) ;
165   Add2RawsList(h5,kLGmod1,1) ;
166   TH2I * h6 = new TH2I("hLowPHOSxyMod2","Low Gain Rows x Columns for PHOS module 2", 64, 0, 64, 56, 0, 56) ;
167   Add2RawsList(h6,kLGmod2,1) ;
168   TH2I * h7 = new TH2I("hLowPHOSxyMod3","Low Gain Rows x Columns for PHOS module 3", 64, 0, 64, 56, 0, 56) ;
169   Add2RawsList(h7,kLGmod3) ;
170   TH2I * h8 = new TH2I("hLowPHOSxyMod4","Low Gain Rows x Columns for PHOS module 4", 64, 0, 64, 56, 0, 56) ;
171   Add2RawsList(h8,kLGmod4,1) ;                                                                                                               
172   TH2I * h9 = new TH2I("hLowPHOSxyMod5","Low Gain Rows x Columns for PHOS module 5", 64, 0, 64, 56, 0, 56) ;                               
173   Add2RawsList(h9,kLGmod5,1) ;                                                                                                               
174                                                                                                                                            
175                                                                                                                                            
176   TH1I * h10 = new TH1I("hLowPhosModules",    "Low Gain Hits in EMCA PHOS modules",       6, 0, 6) ;                                       
177   h10->Sumw2() ;                                                                                                                           
178   Add2RawsList(h10, kNmodLG,1) ;                                                                                                             
179   TH1I * h11 = new TH1I("hHighPhosModules",   "High Gain Hits in EMCA PHOS modules",       6, 0, 6) ;                                      
180   h11->Sumw2() ;                                                                                                                           
181   Add2RawsList(h11, kNmodHG,1) ;                                                                                                             
182                                                                                                                                            
183   TH1F * h12 = new TH1F("hLowPhosRawtime", "Low Gain Time of raw hits in PHOS", 500, -50., 200.) ;                                            
184   h12->Sumw2() ;                                                                                                                           
185   Add2RawsList(h12, kLGtime,1) ;                                                                                                             
186   TH1F * h13 = new TH1F("hHighPhosRawtime", "High Gain Time of raw hits in PHOS", 500, -50., 200.) ;                                          
187   h13->Sumw2() ;                                                                                                                           
188   Add2RawsList(h13, kHGtime,1) ;                                                                                                             
189                                                                                                                                            
190   TH1F * h14 = new TH1F("hLowPhosRawEnergy", "Low Gain Energy of raw hits in PHOS", 500, 0., 1000.) ;                                      
191   h14->Sumw2() ;                                                                                                                           
192   Add2RawsList(h14, kSpecLG,1) ;                                                                                                             
193   TH1F * h15 = new TH1F("hHighPhosRawEnergy", "High Gain Energy of raw hits in PHOS",500,0., 1000.) ;                                      
194   h15->Sumw2() ;                                                                                                                           
195   Add2RawsList(h15, kSpecHG,1) ;                                                                                                             
196                                                                                                                                            
197   TH1F * h16 = new TH1F("hLowNtot", "Low Gain Total Number of raw hits in PHOS", 500, 0., 5000.) ;                                         
198   h16->Sumw2() ;                                                                                                                           
199   Add2RawsList(h16, kNtotLG, kTRUE) ;                                                                                                             
200   TH1F * h17 = new TH1F("hHighNtot", "High Gain Total Number of raw hits in PHOS",500,0., 5000.) ;                                         
201   h17->Sumw2() ;                                                                                                                           
202   Add2RawsList(h17, kNtotHG, kTRUE) ;                                                                                                             
203                                                                                                                                            
204   TH1F * h18 = new TH1F("hLowEtot", "Low Gain Total Energy of raw hits in PHOS", 500, 0., 5000.) ;                                       
205   h18->Sumw2() ;                                                                                                                           
206   Add2RawsList(h18, kEtotLG, kTRUE) ;                                                                                                             
207   TH1F * h19 = new TH1F("hHighEtot", "High Gain Total Energy of raw hits in PHOS",500,0., 100000.) ;                                       
208   h19->Sumw2() ;                                                                                                                           
209   Add2RawsList(h19, kEtotHG, kTRUE) ;                                                                                                             
210   
211   TH2F * h20 = new TH2F("hQualHGxyMod1","High Gain signal quality Rows x Columns module 1", 64, 0, 64, 56, 0, 56) ;
212   Add2RawsList(h20,kHGqualMod1,1) ;
213   TH2F * h21 = new TH2F("hQualHGxyMod2","High Gain signal quality Rows x Columns module 2", 64, 0, 64, 56, 0, 56) ;
214   Add2RawsList(h21,kHGqualMod2,1) ;
215   TH2F * h22 = new TH2F("hQualHGxyMod3","High Gain signal quality Rows x Columns module 3", 64, 0, 64, 56, 0, 56) ;
216   Add2RawsList(h22,kHGqualMod3) ;
217   TH2F * h23 = new TH2F("hQualHGxyMod4","High Gain signal quality Rows x Columns module 4", 64, 0, 64, 56, 0, 56) ;
218   Add2RawsList(h23,kHGqualMod4,1) ;
219   TH2F * h24 = new TH2F("hQualHGxyMod5","High Gain signal quality Rows x Columns module 5", 64, 0, 64, 56, 0, 56) ;
220   Add2RawsList(h24,kHGqualMod5,1) ;
221  
222   TH1F * h25 = new TH1F("hHGpedRMS","High Gain pedestal RMS",200,0.,20.) ;
223   Add2RawsList(h25,kHGpedRMS) ;
224 }
225
226 //____________________________________________________________________________
227 void AliPHOSQADataMakerRec::MakeESDs(AliESDEvent * esd)
228 {
229   // make QA data from ESDs
230
231   Int_t nTot = 0 ; 
232   Double_t eTot = 0 ; 
233   for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
234     AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
235     if( clu->IsPHOS() ) {
236       GetESDsData(kESDSpec)->Fill(clu->E()) ;
237       Double_t *pid=clu->GetPid() ;
238       GetESDsData(kESDpid)->Fill(pid[AliPID::kPhoton]) ;
239       eTot+=clu->E() ;
240       nTot++ ;
241     } 
242   }
243   GetESDsData(kESDNtot)->Fill(nTot) ;
244   GetESDsData(kESDEtot)->Fill(eTot) ;
245 }
246
247 //____________________________________________________________________________
248 void AliPHOSQADataMakerRec::MakeRaws(AliRawReader* rawReader)
249 {
250   //Fill prepared histograms with Raw digit properties
251   rawReader->Reset() ;
252   AliPHOSRawDecoder * decoder ;
253   if(strcmp(GetRecoParam()->EMCDecoderVersion(),"v1")==0)
254     decoder=new AliPHOSRawDecoderv1(rawReader);
255   else
256     if(strcmp(GetRecoParam()->EMCDecoderVersion(),"v2")==0)
257       decoder=new AliPHOSRawDecoderv2(rawReader);
258     else
259       decoder=new AliPHOSRawDecoder(rawReader);
260   decoder->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
261   Double_t lgEtot=0. ;
262   Double_t hgEtot=0. ;
263   Int_t lgNtot=0 ;
264   Int_t hgNtot=0 ;
265
266   while (decoder->NextDigit()) {
267    Int_t module  = decoder->GetModule() ;
268    Int_t row     = decoder->GetRow() ;
269    Int_t col     = decoder->GetColumn() ;
270    Double_t time = decoder->GetTime() ;
271    Double_t energy  = decoder->GetEnergy() ;
272    Bool_t lowGain = decoder->IsLowGain();
273    if (lowGain) {
274      if(energy<2.)
275        continue ;
276      switch(module){
277         case 1: GetRawsData(kLGmod1)->Fill(row-0.5,col-0.5) ; break ;
278         case 2: GetRawsData(kLGmod2)->Fill(row-0.5,col-0.5) ; break ;
279         case 3: GetRawsData(kLGmod3)->Fill(row-0.5,col-0.5) ; break ;
280         case 4: GetRawsData(kLGmod4)->Fill(row-0.5,col-0.5) ; break ;
281         case 5: GetRawsData(kLGmod5)->Fill(row-0.5,col-0.5) ; break ;
282      }                                  
283      GetRawsData(kNmodLG)->Fill(module) ;
284      GetRawsData(kLGtime)->Fill(time) ; 
285      GetRawsData(kSpecLG)->Fill(energy) ;    
286      lgEtot+=energy ;
287      lgNtot++ ;   
288    } else {        
289      //if this isnon-ZS run - fill pedestal RMS
290      if(GetRecoParam()->EMCSubtractPedestals())
291        GetRawsData(kHGpedRMS)->Fill(decoder->GetPedestalRMS()) ;
292      if(energy<8.)
293        continue ;
294      switch (module){
295          case 1: GetRawsData(kHGmod1)->Fill(row-0.5,col-0.5) ; break ;
296          case 2: GetRawsData(kHGmod2)->Fill(row-0.5,col-0.5) ; break ;
297          case 3: GetRawsData(kHGmod3)->Fill(row-0.5,col-0.5) ; break ;
298          case 4: GetRawsData(kHGmod4)->Fill(row-0.5,col-0.5) ; break ;
299          case 5: GetRawsData(kHGmod5)->Fill(row-0.5,col-0.5) ; break ;
300      }              
301      //if quality was evaluated, fill histo
302      if(strcmp(GetRecoParam()->EMCDecoderVersion(),"v1")==0){
303        switch (module){
304          case 1: ((TH2F*)GetRawsData(kHGqualMod1))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
305          case 2: ((TH2F*)GetRawsData(kHGqualMod2))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
306          case 3: ((TH2F*)GetRawsData(kHGqualMod3))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
307          case 4: ((TH2F*)GetRawsData(kHGqualMod4))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
308          case 5: ((TH2F*)GetRawsData(kHGqualMod5))->Fill(row-0.5,col-0.5,decoder->GetSampleQuality()) ; break ;
309        }
310      }
311      GetRawsData(kNmodHG)->Fill(module) ; 
312      GetRawsData(kHGtime)->Fill(time) ;  
313      GetRawsData(kSpecHG)->Fill(energy) ;
314      hgEtot+=energy ; 
315      hgNtot++ ;  
316    }                 
317   }                    
318   GetRawsData(kEtotLG)->Fill(lgEtot) ; 
319   TParameter<double> * p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQA::GetTaskName(AliQA::kRAWS).Data(), GetRawsData(kEtotLG)->GetName()))) ; 
320   p->SetVal(lgEtot) ; 
321   GetRawsData(kEtotHG)->Fill(hgEtot) ;  
322   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQA::GetTaskName(AliQA::kRAWS).Data(), GetRawsData(kEtotHG)->GetName()))) ; 
323   p->SetVal(hgEtot) ; 
324   GetRawsData(kNtotLG)->Fill(lgNtot) ;
325   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQA::GetTaskName(AliQA::kRAWS).Data(), GetRawsData(kNtotLG)->GetName()))) ; 
326   p->SetVal(lgNtot) ; 
327   GetRawsData(kNtotHG)->Fill(hgNtot) ;
328   p = dynamic_cast<TParameter<double>*>(GetParameterList()->FindObject(Form("%s_%s_%s", GetName(), AliQA::GetTaskName(AliQA::kRAWS).Data(), GetRawsData(kNtotHG)->GetName()))) ; 
329   p->SetVal(hgNtot) ; 
330 }
331 //____________________________________________________________________________
332 void AliPHOSQADataMakerRec::MakeRecPoints(TTree * clustersTree)
333 {
334   {
335     // makes data from RecPoints
336     TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
337     if (!emcbranch) { 
338       AliError("can't get the branch with the PHOS EMC clusters !");
339       return;
340     }
341     TObjArray * emcrecpoints = new TObjArray(100) ;
342     emcbranch->SetAddress(&emcrecpoints);
343     emcbranch->GetEntry(0);
344     
345     GetRecPointsData(kRPNtot)->Fill(emcrecpoints->GetEntriesFast()) ; 
346     TIter next(emcrecpoints) ; 
347     AliPHOSEmcRecPoint * rp ; 
348     Double_t eTot = 0. ; 
349     while ( (rp = dynamic_cast<AliPHOSEmcRecPoint *>(next())) ) {
350       GetRecPointsData(kRPSpec)->Fill( rp->GetEnergy()) ;
351       Int_t mod = rp->GetPHOSMod() ;
352       TVector3 pos ;
353       rp->GetLocalPosition(pos) ;
354       switch(mod){
355         case 1: GetRecPointsData(kRPmod1)->Fill(pos.X(),pos.Z()) ; break ;
356         case 2: GetRecPointsData(kRPmod2)->Fill(pos.X(),pos.Z()) ; break ;
357         case 3: GetRecPointsData(kRPmod3)->Fill(pos.X(),pos.Z()) ; break ;
358         case 4: GetRecPointsData(kRPmod4)->Fill(pos.X(),pos.Z()) ; break ;
359         case 5: GetRecPointsData(kRPmod5)->Fill(pos.X(),pos.Z()) ; break ;
360       }
361       eTot+= rp->GetEnergy() ;
362     }
363     GetRecPointsData(kRPEtot)->Fill(eTot) ;
364     emcrecpoints->Delete();
365     delete emcrecpoints;
366   }
367   {
368     TBranch *cpvbranch = clustersTree->GetBranch("PHOSCpvRP");
369     if (!cpvbranch) { 
370       AliError("can't get the branch with the PHOS CPV clusters !");
371       return;
372     }
373     TObjArray *cpvrecpoints = new TObjArray(100) ;
374     cpvbranch->SetAddress(&cpvrecpoints);
375     cpvbranch->GetEntry(0);
376     
377     GetRecPointsData(kRPNcpv)->Fill(cpvrecpoints->GetEntriesFast()) ; 
378     cpvrecpoints->Delete();
379     delete cpvrecpoints;
380   }
381 }
382
383 //____________________________________________________________________________ 
384 void AliPHOSQADataMakerRec::StartOfDetectorCycle()
385 {
386   //Detector specific actions at start of cycle
387   
388 }