]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDQADataMakerRec.cxx
Fixes for coverity
[u/mrichter/AliRoot.git] / PMD / AliPMDQADataMakerRec.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 /*
18   Produces the data needed to calculate the quality assurance. 
19   All data must be mergeable objects.
20   B.K. Nandi
21 */
22
23 // --- ROOT system ---
24 #include <TClonesArray.h>
25 #include <TFile.h>
26 #include <TH1D.h>  
27 #include <TH1F.h> 
28 #include <TH1I.h> 
29 #include <TH2F.h> 
30
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34
35 #include "AliESDEvent.h"
36 #include "AliLog.h"
37 #include "AliPMDQADataMakerRec.h"
38 #include "AliQAChecker.h"
39 #include "AliPMDdigit.h" 
40 #include "AliPMDrecpoint1.h" 
41 #include "AliPMDRawStream.h"
42 #include "AliPMDddldata.h"
43 #include "AliPMDUtility.h"
44 #include "AliESDPmdTrack.h"
45 //#include "AliPMDRecoParam.h"
46
47 ClassImp(AliPMDQADataMakerRec)
48            
49 //____________________________________________________________________________ 
50   AliPMDQADataMakerRec::AliPMDQADataMakerRec() : 
51   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kPMD), "PMD Quality Assurance Data Maker")
52 {
53   // ctor
54 }
55
56 //____________________________________________________________________________ 
57 AliPMDQADataMakerRec::AliPMDQADataMakerRec(const AliPMDQADataMakerRec& qadm) :
58   AliQADataMakerRec()
59 {
60   //copy ctor 
61   SetName((const char*)qadm.GetName()) ; 
62   SetTitle((const char*)qadm.GetTitle()); 
63 }
64
65 //__________________________________________________________________
66 AliPMDQADataMakerRec& AliPMDQADataMakerRec::operator = (const AliPMDQADataMakerRec& qadm )
67 {
68   // Equal operator.
69   this->~AliPMDQADataMakerRec();
70   new(this) AliPMDQADataMakerRec(qadm);
71   return *this;
72 }
73  
74
75 //____________________________________________________________________________ 
76 void AliPMDQADataMakerRec::InitRaws() {
77   // create Raws histograms in Raws subdir
78
79   const Bool_t expert   = kTRUE ; 
80   const Bool_t saveCorr = kTRUE ; 
81   const Bool_t image    = kTRUE ; 
82   
83
84   TH1F *h0 = new TH1F("hDdl2304","Cell Wise Frequency for PRE Plane DDL 2304;DDL2304;Frequency",50000 ,0,50000.0);
85   TH1F *h1 = new TH1F("hDdl2305","Cell Wise Frequency for PRE Plane DDL 2305;DDL2305;Frequency",50000 ,0,50000.0);
86   TH1F *h2 = new TH1F("hDdl2308","Cell Wise Frequency for CPV Plane DDL 2308;DDL2308;Frequency",50000 ,0,50000.0);
87   TH1F *h3 = new TH1F("hDdl2309","Cell Wise Frequency for CPV Plane DDL 2309;DDL2309;Frequency",50000 ,0,50000.0);
88     
89   for(int i = 0; i < 50; i++ )  {
90     h0->GetXaxis()->SetBinLabel(i*1000+300,Form("%2d",i+1));
91     h1->GetXaxis()->SetBinLabel(i*1000+300,Form("%2d",i+1));
92     h2->GetXaxis()->SetBinLabel(i*1000+300,Form("%2d",i+1));
93     h3->GetXaxis()->SetBinLabel(i*1000+300,Form("%2d",i+1));
94   }
95
96
97   Add2RawsList(h0, 0, expert, image, saveCorr);  
98   Add2RawsList(h1, 1, expert, image, saveCorr);  
99   Add2RawsList(h2, 2, expert, image, saveCorr);  
100   Add2RawsList(h3, 3, expert, image, saveCorr);  
101
102
103   TH2F *hPreXY = new TH2F("hPreXY","PRE plane;X [cm];Y [cm]",200,-100.,100.,200,-100.,100.);
104   Add2RawsList(hPreXY, 4, expert, !image, saveCorr);
105   
106   TH2F *hCpvXY = new TH2F("hCpvXY","CPV plane;X [cm];Y [cm]",200,-100.,100.,200,-100.,100.);
107   Add2RawsList(hCpvXY, 5, expert, !image, saveCorr);
108   
109   TH1F *hQuality1 = new TH1F("hPmdQualityWAdc","Quality Plot SumWAdc/SumWAdc;#eta Bins (PRE   <--->    CPV);Quality",25,0,25);
110   TH1F *hQuality2 = new TH1F("hPmdQualityHit","Quality Plot Hit/Hit;#eta Bins (PRE  <--->   CPV);Quality",25,0,25);
111
112   for(int i = 0; i < 10; i++ )  {
113     hQuality1->GetXaxis()->SetBinLabel(i+3,Form("%2d",i+1));
114     hQuality2->GetXaxis()->SetBinLabel(i+3,Form("%2d",i+1));
115
116     hQuality1->GetXaxis()->SetBinLabel(i+15,Form("%2d",i+1));
117     hQuality2->GetXaxis()->SetBinLabel(i+15,Form("%2d",i+1));
118   }
119
120   Add2RawsList(hQuality1, 6, !expert, image, saveCorr);
121   Add2RawsList(hQuality2, 7, !expert, image, saveCorr);
122
123   TH1F *hMP = new TH1F("hPmdMultiplicityP"," Cell Hit Distribution with adc>40; Multiplicity; Entries",1000,0,7000); 
124   TH1F *hMC = new TH1F("hPmdMultiplicityC"," Cell Hit Distribution with adc>40; Multiplicity; Entries",1000,0,7000); 
125   Add2RawsList(hMP, 8, expert, !image, saveCorr);
126   Add2RawsList(hMC, 9, expert, !image, saveCorr);
127   
128
129
130
131   ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
132 }
133
134 //____________________________________________________________________________
135 void AliPMDQADataMakerRec::InitDigits()
136 {
137   // create Digits histograms in Digits subdir
138   const Bool_t expert   = kTRUE ; 
139   const Bool_t image    = kTRUE ; 
140   
141   TH1F *h0 = new TH1F("hPreDigitsEdep","Digits energy distribution in PRE(PMD);Amplitude [ADC counts];Counts", 100, 0., 2000.);
142   h0->Sumw2();
143   Add2DigitsList(h0, 0, !expert, image);
144   
145   TH1F *h1 = new TH1F("hCpvDigitsEdep","Digits energy distribution in CPV(PMD);Amplitude [ADC counts];Counts", 100, 0., 2000.); 
146   h1->Sumw2();
147   Add2DigitsList(h1, 1, !expert, image);
148   
149   TH1I *h2 = new TH1I("hPreDigitsMult","Digits multiplicity distribution in PRE(PMD);# of Digits;Entries", 500, 0, 1000) ; 
150   h2->Sumw2();
151   Add2DigitsList(h2, 2, !expert, image);
152   
153   TH1I *h3 = new TH1I("hCpvDigitsMult","Digits multiplicity distribution in CPV(PMD);# of Digits;Entries", 500, 0, 1000);
154   h3->Sumw2();
155   Add2DigitsList(h3, 3, !expert, image);  
156   //
157   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
158 }
159
160 //____________________________________________________________________________ 
161 void AliPMDQADataMakerRec::InitRecPoints()
162 {
163   // create Reconstructed Points histograms in RecPoints subdir
164   
165   /*
166     TH2F * h0 = new TH2F("hPreXY","RecPoints Y vs X PRE(PMD)", 100,-100.,100.,100,-100.,100.);
167     Add2RecPointsList(h0,0);
168     
169     TH2F * h1 = new TH2F("hCpvXY","RecPoints Y vs X CPV(PMD)", 100,-100.,100.,100,-100.,100.);
170     Add2RecPointsList(h1,1);
171   */
172
173     //  Ncell distribution in a cluster
174   
175     // PRE plane
176
177   const Bool_t expert   = kTRUE ; 
178   const Bool_t image    = kTRUE ; 
179
180   TH1F * h0 = new TH1F("hPreDdl0Ncell","PRE: Ddl0 Ncell in a cluster;# of cells;Counts",50,0.,50.);
181   h0->Sumw2();
182   Add2RecPointsList(h0, 0, !expert, image);
183
184
185   TH1F * h1 = new TH1F("hPreDdl1Ncell","PRE: Ddl1 Ncell in a cluste;# of cells;Countsr",50,0.,50.);
186   h1->Sumw2();
187   Add2RecPointsList(h1, 1, !expert, image);
188
189
190   TH1F * h2 = new TH1F("hPreDdl2Ncell","PRE: Ddl2 Ncell in a cluster;# of cells;Counts",50,0.,50.);
191   h2->Sumw2();
192   Add2RecPointsList(h2, 2, !expert, image);
193
194
195   TH1F * h3 = new TH1F("hPreDdl3Ncell","PRE: Ddl3 Ncell in a cluster;# of cells;Counts",50,0.,50.);
196   h3->Sumw2();
197   Add2RecPointsList(h3, 3, !expert, image);
198
199   // CPV plane
200
201   TH1F * h4 = new TH1F("hCpvDdl4Ncell","CPV: Ddl4 Ncell in a cluster;# of cells;Counts",50,0.,50.);
202   h4->Sumw2();
203   Add2RecPointsList(h4, 4, !expert, image);
204
205   TH1F * h5 = new TH1F("hCpvDdl5Ncell","CPV: Ddl5 Ncell in a cluster;# of cells;Counts",50,0.,50.);
206   h5->Sumw2();
207   Add2RecPointsList(h5, 5, !expert, image);
208
209   // Correlation plot
210
211   TH2I *h6 = new TH2I("hPre10","Cluster - DDL1 vs DDL0;DDL0;DDL1", 100,0,200,100,0,200);
212   Add2RecPointsList(h6,6, !expert, image);
213
214   TH2I *h7 = new TH2I("hPre32","Cluster - DDL3 vs DDL2;DDL2;DDL3", 100,0,200,100,0,200);
215   Add2RecPointsList(h7,7, !expert, image);
216
217
218   TH2I *h8 = new TH2I("hCpv54","Cluster - DDL5 vs DDL4;DDL4;DDL5", 100,0,200,100,0,200);
219   Add2RecPointsList(h8,8, !expert, image);
220   //
221   ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
222 }
223
224 //____________________________________________________________________________ 
225
226 void AliPMDQADataMakerRec::InitESDs()
227 {
228   //Create histograms to controll ESD
229
230   const Bool_t expert   = kTRUE ; 
231   const Bool_t image    = kTRUE ; 
232
233   TH1F *h0 = new TH1F("hPreClADC","Cluster ADC of PRE plane;# of clusters;Counts",500,0.,10000.);
234   h0->Sumw2();
235   Add2ESDsList(h0, 0, !expert, image)  ;
236
237   TH1F *h1 = new TH1F("hCpvClADC","Cluster ADC of CPV plane;# of clusters;Counts",500,0.,10000.);
238   h1->Sumw2();
239   Add2ESDsList(h1, 1, !expert, image)  ;
240
241   TH2I *h2 = new TH2I("hPmdClMult","Cluster Multiplicity: PRE vs. CPVplane;CPV multiplicity;PRE Multiplicity",100,0,1000,100,0,1000);
242   h2->Sumw2();
243   Add2ESDsList(h2, 2, !expert, image)  ;
244
245   /*
246     TH1I * h3 = new TH1I("hCpvClMult","Cluster Multiplicity of CPV plane",100,0.,1000.);
247     h3->Sumw2();
248     Add2ESDsList(h3, 3, !expert, image)  ;
249   */
250   //
251   ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line
252 }
253
254 //____________________________________________________________________________
255 void AliPMDQADataMakerRec::MakeRaws(AliRawReader* rawReader)
256 {
257   //Fill prepared histograms with Raw digit properties
258
259   TObjArray *pmdddlcont = 0x0;
260     pmdddlcont = new TObjArray();
261     AliPMDRawStream stream(rawReader);
262     
263     AliPMDddldata *pmdddl = 0x0;
264
265     Int_t   iddl =  1;
266     Int_t   sddl = -1;
267     Float_t xx = 0., yy = 0.;
268
269     Double_t eta;
270     Int_t etaindex = 0;
271
272     Float_t nDdl1[10] = {0,0,0,0,0,0,0,0,0,0};
273     Float_t nDdl2[10] = {1,1,1,1,1,1,1,1,1,1};
274     Float_t nDdl3[10] = {0,0,0,0,0,0,0,0,0,0};
275     Float_t nDdl4[10] = {1,1,1,1,1,1,1,1,1,1};
276
277     Float_t nDdl1a[10] = {1,1,1,1,1,1,1,1,1,1};
278     Float_t nDdl2a[10] = {1,1,1,1,1,1,1,1,1,1};
279     Float_t nDdl3a[10] = {1,1,1,1,1,1,1,1,1,1};
280     Float_t nDdl4a[10] = {1,1,1,1,1,1,1,1,1,1};
281     
282     Float_t nDdl11[10] = {1,1,1,1,1,1,1,1,1,1};
283     Float_t nDdl22[10] = {1,1,1,1,1,1,1,1,1,1};
284     Float_t nDdl33[10] = {1,1,1,1,1,1,1,1,1,1};
285     Float_t nDdl44[10] = {1,1,1,1,1,1,1,1,1,1};
286
287     Float_t nddlp = 0;
288     Float_t nddlc  = 0;
289
290     AliPMDUtility cc;
291
292        
293
294     while ((iddl = stream.DdlData(pmdddlcont)) >=0)
295     {
296         Int_t ientries = pmdddlcont->GetEntries();
297         //printf(" ======= DDLNO = %d ientries = %d \n", iddl, ientries);
298         if(iddl == 0 || iddl == 1) sddl = iddl;
299         if(iddl == 4 ) sddl = 2;
300         if(iddl == 5 ) sddl = 3;
301         
302         for (Int_t ient = 0; ient < ientries; ient++)
303           {
304             pmdddl = (AliPMDddldata*)pmdddlcont->UncheckedAt(ient);
305
306             Int_t det  = pmdddl->GetDetector();
307             Int_t smn  = pmdddl->GetSMN();
308             Int_t mcm  = pmdddl->GetMCM();
309             Int_t chno = pmdddl->GetChannel();
310             Int_t pbus = pmdddl->GetPatchBusId();
311             Int_t row  = pmdddl->GetRow();
312             Int_t col  = pmdddl->GetColumn();
313             Int_t sig  = pmdddl->GetSignal();
314
315             if (mcm == 0) continue;
316
317             if (det < 0 || det > 1)  continue;
318             if (smn < 0 || smn > 23) continue;
319             if (row < 0 || row > 47) continue;
320             if (col < 0 || col > 95) continue;
321             
322             Int_t rc = pbus*1000 + mcm*64 + chno;
323
324             if(sddl == 0) FillRawsData(0,rc);
325             else if(sddl == 1) FillRawsData(1,rc);
326             else if(sddl == 2) FillRawsData(2,rc);
327             else if(sddl == 3) FillRawsData(3,rc);
328
329             cc.GetEtaIndexXY(smn,row,col,xx,yy,eta,etaindex);
330             
331             if(etaindex >= 0 && etaindex <10) {
332
333               if(sddl == 0) { 
334                 nDdl1a[etaindex] += 1; 
335                 nDdl1[etaindex] += sig; 
336                 if(sig >40) {
337                   nDdl11[etaindex] += 1;
338                   nddlp += 1;
339                 }
340                   
341               }
342               else if(sddl == 1) { 
343                 nDdl2a[etaindex] += 1; 
344                 nDdl2[etaindex] += sig; 
345                 if(sig >40) {
346                   nDdl22[etaindex] += 1;
347                   nddlp += 1;
348                 }
349  
350
351               }
352               else if(sddl == 2) { 
353                 nDdl3a[etaindex] += 1; 
354                 nDdl3[etaindex] += sig; 
355                 if(sig >40) {
356                   nDdl33[etaindex] += 1;
357                   nddlc += 1;
358                 } 
359               }
360               else if(sddl == 3) { 
361                 nDdl4a[etaindex] += 1; 
362                 nDdl4[etaindex] += sig; 
363                 if(sig >40) {
364                   nDdl44[etaindex] += 1; 
365                   nddlc += 1;
366                 }
367               }
368
369             }
370
371             if (det == 0) FillRawsData(4,xx,yy);
372             else if (det == 1) FillRawsData(5,xx,yy);
373             
374           }
375
376         pmdddlcont->Delete();
377     }
378
379
380     delete pmdddlcont;
381     pmdddlcont = 0x0;
382     
383     if(nddlp != 0)   FillRawsData(8,nddlp);
384     if(nddlc != 0)   FillRawsData(9,nddlc);
385    
386
387     ResetRawsData(6);
388     ResetRawsData(7);
389     for (Int_t i = 0; i < 10; i++) {
390      
391       Float_t prerC =  nDdl11[i]/nDdl22[i];
392       Float_t cpvrC = nDdl33[i]/nDdl44[i];
393
394       Float_t prera =(nDdl1[i]/nDdl1a[i])/(nDdl2[i]/nDdl2a[i]);
395       Float_t cpvra =(nDdl3[i]/nDdl3a[i])/(nDdl4[i]/nDdl4a[i]);
396       
397      
398       FillRawsData(6,i+2,prera);
399       FillRawsData(6,i+14,cpvra);
400       FillRawsData(7,i+2,prerC);
401       FillRawsData(7,i+14,cpvrC);
402     }
403
404     IncEvCountCycleRaws();
405     IncEvCountTotalRaws();
406  
407     //
408 }
409 //____________________________________________________________________________
410 void AliPMDQADataMakerRec::MakeDigits()
411 {
412   // makes data from Digits
413   
414    Int_t cpvmul = 0, premul = 0;
415   
416   TIter next(fDigitsArray) ; 
417   AliPMDdigit * digit ; 
418   while ( (digit = dynamic_cast<AliPMDdigit *>(next())) )
419     {
420     if(digit->GetDetector() == 0)
421       {
422             FillDigitsData(0, digit->GetADC()) ;
423             premul++;
424       }
425     if(digit->GetDetector() == 1)
426       {
427             FillDigitsData(1, digit->GetADC());
428             cpvmul++;
429       }
430     }  
431   
432   if (premul > 0) FillDigitsData(2,premul);
433   if (cpvmul > 0) FillDigitsData(3,cpvmul);
434   
435   
436 }
437
438 //____________________________________________________________________________
439 void AliPMDQADataMakerRec::MakeDigits(TTree * digitTree)
440 {
441   // makes data from Digit Tree
442   
443   if (fDigitsArray) 
444     fDigitsArray->Clear() ; 
445   else
446     fDigitsArray = new TClonesArray("AliPMDdigit", 1000) ; 
447   
448   TBranch * branch = digitTree->GetBranch("PMDDigit") ;
449   if ( ! branch ) {AliWarning("PMD branch in Digit Tree not found"); return;} 
450   //
451   branch->SetAddress(&fDigitsArray) ;
452   for (Int_t ient = 0; ient < branch->GetEntries(); ient++) {
453     branch->GetEntry(ient) ; 
454     MakeDigits() ; 
455   }
456   //
457   IncEvCountCycleDigits();
458   IncEvCountTotalDigits();
459   //
460 }
461
462 //____________________________________________________________________________
463 void AliPMDQADataMakerRec::MakeRecPoints(TTree * clustersTree)
464 {
465     // makes data from RecPoints
466
467   Int_t multDdl0 = 0, multDdl1 = 0, multDdl2 = 0;
468   Int_t multDdl3 = 0, multDdl4 = 0, multDdl5 = 0;
469   
470   AliPMDrecpoint1 * recpoint; 
471   
472   if (fRecPointsArray) 
473     fRecPointsArray->Clear() ; 
474   else 
475     fRecPointsArray = new TClonesArray("AliPMDrecpoint1", 1000) ; 
476   
477   TBranch * branch = clustersTree->GetBranch("PMDRecpoint") ;
478
479   if ( ! branch )
480     {
481       AliWarning("PMD branch in Recpoints Tree not found") ; 
482     }
483   else
484     {
485       branch->SetAddress(&fRecPointsArray) ;
486       
487       for (Int_t imod = 0; imod < branch->GetEntries(); imod++)
488         {
489           branch->GetEntry(imod) ;
490           
491           TIter next(fRecPointsArray) ; 
492           
493           while ( (recpoint = dynamic_cast<AliPMDrecpoint1 *>(next())) )
494             {
495               //Float_t xpos = recpoint->GetClusX();
496               //Float_t ypos = recpoint->GetClusY();
497               //Int_t smn = recpoint->GetSMNumber();
498               
499               if(recpoint->GetDetector() == 0)
500                 {
501                   if(recpoint->GetSMNumber() >= 0 && recpoint->GetSMNumber() < 6)
502                     {
503                       FillRecPointsData(0,recpoint->GetClusCells());
504                       multDdl0++;
505                     }
506                   if(recpoint->GetSMNumber() >= 6 && recpoint->GetSMNumber() < 12)
507                     {
508                       FillRecPointsData(1,recpoint->GetClusCells());
509                       multDdl1++;
510                     }
511                   if(recpoint->GetSMNumber() >= 12 && recpoint->GetSMNumber() < 18)
512                     {
513                       FillRecPointsData(2,recpoint->GetClusCells());
514                       multDdl2++;
515                     }
516                   if(recpoint->GetSMNumber() >= 18 && recpoint->GetSMNumber() < 24)
517                     {
518                       FillRecPointsData(3,recpoint->GetClusCells());
519                       multDdl3++;
520                     }
521                 }
522               
523               if(recpoint->GetDetector() == 1)
524                 {
525                   if((recpoint->GetSMNumber() >= 0 && recpoint->GetSMNumber() < 6) || 
526                      (recpoint->GetSMNumber() >= 18 && recpoint->GetSMNumber() < 24))
527                     {
528                       FillRecPointsData(4,recpoint->GetClusCells());
529                       multDdl4++;
530                     }
531                   if(recpoint->GetSMNumber() >= 6 && recpoint->GetSMNumber() < 18 )
532                     {
533                       FillRecPointsData(5,recpoint->GetClusCells());
534                       multDdl5++;
535                     }
536                 }
537             } 
538         }
539     }
540   
541   FillRecPointsData(6,multDdl0,multDdl1);
542   FillRecPointsData(7,multDdl2,multDdl3);
543   FillRecPointsData(8,multDdl4,multDdl5);
544   //
545   IncEvCountCycleRecPoints();
546   IncEvCountTotalRecPoints();  
547   //
548 }
549
550 //____________________________________________________________________________
551
552 void AliPMDQADataMakerRec::MakeESDs(AliESDEvent * esd)
553 {
554   // make QA data from ESDs
555
556   Int_t premul = 0, cpvmul = 0;
557
558   for (Int_t icl = 0; icl < esd->GetNumberOfPmdTracks(); icl++)
559     {
560       AliESDPmdTrack *pmdtr = esd->GetPmdTrack(icl);
561       
562       //Int_t   det   = pmdtr->GetDetector(); 
563       //Float_t clsX  = pmdtr->GetClusterX();
564       //Float_t clsY  = pmdtr->GetClusterY();
565       //Float_t clsZ  = pmdtr->GetClusterZ();
566       //Float_t ncell = pmdtr->GetClusterCells();
567       Float_t adc   = pmdtr->GetClusterADC();
568       //Float_t pid   = pmdtr->GetClusterPID();
569
570       if (pmdtr->GetDetector() == 0)
571         {
572           FillESDsData(0,adc);
573           premul++;
574         }
575       if (pmdtr->GetDetector() == 1)
576         {
577           FillESDsData(1,adc) ;
578           cpvmul++;
579         }
580     }
581   
582   FillESDsData(2,cpvmul,premul) ;
583   //FillESDsData(3,cpvmul) ;  
584   //
585   IncEvCountCycleESDs();
586   IncEvCountTotalESDs(); 
587   //
588 }
589
590 //____________________________________________________________________________ 
591
592 void AliPMDQADataMakerRec::StartOfDetectorCycle()
593 {
594   //Detector specific actions at start of cycle
595   
596 }
597 //____________________________________________________________________________ 
598 void AliPMDQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
599 {
600   //Detector specific actions at end of cycle
601   // do the QA checking
602   AliQAChecker::Instance()->Run(AliQAv1::kPMD, task, list) ;  
603 }