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