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