]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDClusterFinder.cxx
clustering V2 added
[u/mrichter/AliRoot.git] / PMD / AliPMDClusterFinder.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 //           Date   : August 05 2003                   //
19 //  This reads the file PMD.digits.root(TreeD),        //
20 //  calls the Clustering algorithm and stores the      //
21 //  clustering output in PMD.RecPoints.root(TreeR)     // 
22 //                                                     //
23 //-----------------------------------------------------//
24
25 #include <Riostream.h>
26 #include <TTree.h>
27 #include <TObjArray.h>
28 #include <TClonesArray.h>
29
30 #include "AliLog.h"
31 #include "AliRunLoader.h"
32 #include "AliLoader.h"
33 #include "AliRawReader.h"
34
35 #include "AliPMDdigit.h"
36 #include "AliPMDClusterFinder.h"
37 #include "AliPMDClustering.h"
38 #include "AliPMDClusteringV1.h"
39 #include "AliPMDcluster.h"
40 #include "AliPMDrecpoint1.h"
41 #include "AliPMDrechit.h"
42 #include "AliPMDRawStream.h"
43 #include "AliPMDCalibData.h"
44 #include "AliPMDddldata.h"
45
46 #include "AliCDBManager.h"
47 #include "AliCDBEntry.h"
48
49
50
51 ClassImp(AliPMDClusterFinder)
52
53 AliPMDClusterFinder::AliPMDClusterFinder():
54   fRunLoader(0),
55   fPMDLoader(0),
56   fTreeD(0),
57   fTreeR(0),
58   fDigits(new TClonesArray("AliPMDdigit", 1000)),
59   fRecpoints(new TClonesArray("AliPMDrecpoint1", 1000)),
60   fRechits(new TClonesArray("AliPMDrechit", 1000)),
61   fNpoint(0),
62   fNhit(0),
63   fEcut(0.)
64 {
65 //
66 // Constructor
67 //
68   fCalibData = GetCalibData();
69 }
70 // ------------------------------------------------------------------------- //
71 AliPMDClusterFinder::AliPMDClusterFinder(AliRunLoader* runLoader):
72   fRunLoader(runLoader),
73   fPMDLoader(runLoader->GetLoader("PMDLoader")),
74   fTreeD(0),
75   fTreeR(0),
76   fDigits(new TClonesArray("AliPMDdigit", 1000)),
77   fRecpoints(new TClonesArray("AliPMDrecpoint1", 1000)),
78   fRechits(new TClonesArray("AliPMDrechit", 1000)),
79   fNpoint(0),
80   fNhit(0),
81   fEcut(0.)
82 {
83 //
84 // Constructor
85 //
86   fCalibData = GetCalibData();
87 }
88 // ------------------------------------------------------------------------- //
89 AliPMDClusterFinder::~AliPMDClusterFinder()
90 {
91   // Destructor
92   if (fDigits)
93     {
94       fDigits->Delete();
95       delete fDigits;
96       fDigits=0;
97     }
98   if (fRecpoints)
99     {
100       fRecpoints->Delete();
101       delete fRecpoints;
102       fRecpoints=0;
103     }
104   if (fRechits)
105     {
106       fRechits->Delete();
107       delete fRechits;
108       fRechits=0;
109     }
110 }
111 // ------------------------------------------------------------------------- //
112
113 void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
114 {
115   // Converts digits to recpoints after running clustering
116   // algorithm on CPV plane and PREshower plane
117   //
118   Int_t    det = 0,smn = 0;
119   Int_t    xpos,ypos;
120   Float_t  adc;
121   Int_t    ismn;
122   Int_t    idet;
123   Float_t  clusdata[6];
124
125   TObjArray *pmdcont = new TObjArray();
126   AliPMDClustering *pmdclust = new AliPMDClusteringV1();
127
128   pmdclust->SetEdepCut(fEcut);
129
130   fRunLoader->GetEvent(ievt);
131
132
133   fTreeD = fPMDLoader->TreeD();
134   if (fTreeD == 0x0)
135     {
136       AliFatal("AliPMDClusterFinder: Can not get TreeD");
137
138     }
139   AliPMDdigit  *pmddigit;
140   TBranch *branch = fTreeD->GetBranch("PMDDigit");
141   branch->SetAddress(&fDigits);
142
143   ResetRecpoint();
144
145   fTreeR = fPMDLoader->TreeR();
146   if (fTreeR == 0x0)
147     {
148       fPMDLoader->MakeTree("R");
149       fTreeR = fPMDLoader->TreeR();
150     }
151
152   Int_t bufsize = 16000;
153   TBranch * branch1 = fTreeR->Branch("PMDRecpoint", &fRecpoints, bufsize); 
154   TBranch * branch2 = fTreeR->Branch("PMDRechit", &fRechits, bufsize); 
155
156   Int_t nmodules = (Int_t) fTreeD->GetEntries();
157
158   for (Int_t imodule = 0; imodule < nmodules; imodule++)
159     {
160       ResetCellADC();
161       fTreeD->GetEntry(imodule); 
162       Int_t nentries = fDigits->GetLast();
163       for (Int_t ient = 0; ient < nentries+1; ient++)
164         {
165           pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
166           
167           det    = pmddigit->GetDetector();
168           smn    = pmddigit->GetSMNumber();
169           xpos   = pmddigit->GetRow();
170           ypos   = pmddigit->GetColumn();
171           adc    = pmddigit->GetADC();
172           
173           // CALIBRATION
174           Float_t gain = fCalibData->GetGainFact(det,smn,xpos,ypos);
175           
176          // printf("adc = %d gain = %f\n",adc,gain);
177           
178           adc = adc*gain;
179
180
181           //Int_t trno   = pmddigit->GetTrackNumber();
182           fCellADC[xpos][ypos] = (Double_t) adc;
183         }
184
185       idet = det;
186       ismn = smn;
187       pmdclust->DoClust(idet,ismn,fCellADC,pmdcont);
188       
189       Int_t nentries1 = pmdcont->GetEntries();
190
191       AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
192
193       for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
194         {
195           AliPMDcluster *pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
196           idet        = pmdcl->GetDetector();
197           ismn        = pmdcl->GetSMN();
198           clusdata[0] = pmdcl->GetClusX();
199           clusdata[1] = pmdcl->GetClusY();
200           clusdata[2] = pmdcl->GetClusADC();
201           clusdata[3] = pmdcl->GetClusCells();
202           clusdata[4] = pmdcl->GetClusSigmaX();
203           clusdata[5] = pmdcl->GetClusSigmaY();
204
205           AddRecPoint(idet,ismn,clusdata);
206
207           Int_t ncell = (Int_t) clusdata[3];
208           for(Int_t ihit = 0; ihit < ncell; ihit++)
209             {
210               Int_t celldataX = pmdcl->GetClusCellX(ihit);
211               Int_t celldataY = pmdcl->GetClusCellY(ihit);
212               AddRecHit(celldataX, celldataY);
213             }
214           branch2->Fill();
215           ResetRechit();
216         }
217       pmdcont->Clear();
218       
219       branch1->Fill();
220       ResetRecpoint();
221
222     } // modules
223
224   ResetCellADC();
225   fPMDLoader = fRunLoader->GetLoader("PMDLoader");  
226   fPMDLoader->WriteRecPoints("OVERWRITE");
227
228   //   delete the pointers
229   delete pmdclust;
230   delete pmdcont;
231     
232 }
233 // ------------------------------------------------------------------------- //
234
235 void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
236                                            TTree *clustersTree)
237 {
238   // Converts RAW data to recpoints after running clustering
239   // algorithm on CPV and PREshower plane
240   //
241
242   Float_t  clusdata[6];
243   TObjArray pmdddlcont;
244
245   TObjArray *pmdcont = new TObjArray();
246   AliPMDClustering *pmdclust = new AliPMDClusteringV1();
247
248   pmdclust->SetEdepCut(fEcut);
249
250   ResetRecpoint();
251
252   Int_t bufsize = 16000;
253   TBranch *branch1 = clustersTree->Branch("PMDRecpoint", &fRecpoints, bufsize); 
254
255   TBranch * branch2 = clustersTree->Branch("PMDRechit", &fRechits, bufsize); 
256
257   const Int_t kDDL = 6;
258   const Int_t kRow = 48;
259   const Int_t kCol = 96;
260
261   Int_t idet = 0;
262   Int_t iSMN = 0;
263   
264   for (Int_t indexDDL = 0; indexDDL < kDDL; indexDDL++)
265     {
266       if (indexDDL < 4)
267         {
268           iSMN = 6;
269         }
270       else if (indexDDL >= 4)
271         {
272           iSMN = 12;
273         }
274       Int_t ***precpvADC;
275       precpvADC = new int **[iSMN];
276       for (Int_t i=0; i<iSMN; i++) precpvADC[i] = new int *[kRow];
277       for (Int_t i=0; i<iSMN;i++)
278         {
279           for (Int_t j=0; j<kRow; j++) precpvADC[i][j] = new int [kCol];
280         }
281       for (Int_t i = 0; i < iSMN; i++)
282         {
283           for (Int_t j = 0; j < kRow; j++)
284             {
285               for (Int_t k = 0; k < kCol; k++)
286                 {
287                   precpvADC[i][j][k] = 0;
288                 }
289             }
290         }
291       ResetCellADC();
292       rawReader->Reset();
293       AliPMDRawStream pmdinput(rawReader);
294       rawReader->Select("PMD", indexDDL, indexDDL);
295       pmdinput.DdlData(&pmdddlcont);
296       Int_t indexsmn = 0;
297       Int_t ientries = pmdddlcont.GetEntries();
298       for (Int_t ient = 0; ient < ientries; ient++)
299         {
300           AliPMDddldata *pmdddl = (AliPMDddldata*)pmdddlcont.UncheckedAt(ient);
301           
302           Int_t det = pmdddl->GetDetector();
303           Int_t smn = pmdddl->GetSMN();
304           //Int_t mcm = pmdddl->GetMCM();
305           //Int_t chno = pmdddl->GetChannel();
306           Int_t row = pmdddl->GetRow();
307           Int_t col = pmdddl->GetColumn();
308           Int_t sig = pmdddl->GetSignal();
309
310           Float_t sig1 = (Float_t) sig;
311           // CALIBRATION
312           Float_t gain = fCalibData->GetGainFact(det,smn,row,col);
313           //printf("sig = %d gain = %f\n",sig,gain);
314           sig = (Int_t) (sig1*gain);
315
316           if (indexDDL < 4)
317             {
318               if (det != 0)
319                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
320                               indexDDL, det));
321               indexsmn = smn - indexDDL * 6;
322             }
323           else if (indexDDL == 4)
324             {
325               if (det != 1)
326                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
327                               indexDDL, det));
328               if (smn < 6)
329                 {
330                   indexsmn = smn;
331                 }
332               else if (smn >= 18 && smn < 24)
333                 {
334                   indexsmn = smn - 12;
335                 }
336             }
337           else if (indexDDL == 5)
338             {
339               if (det != 1)
340                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
341                               indexDDL, det));
342               if (smn >= 6 && smn < 12)
343                 {
344                   indexsmn = smn - 6;
345                 }
346               else if (smn >= 12 && smn < 18)
347                 {
348                   indexsmn = smn - 6;
349                 }
350             }         
351           precpvADC[indexsmn][row][col] = sig;
352         }
353       
354       pmdddlcont.Clear();
355
356       Int_t ismn = 0;
357       for (Int_t indexsmn = 0; indexsmn < iSMN; indexsmn++)
358         {
359           ResetCellADC();
360           for (Int_t irow = 0; irow < kRow; irow++)
361             {
362               for (Int_t icol = 0; icol < kCol; icol++)
363                 {
364                   fCellADC[irow][icol] = 
365                     (Double_t) precpvADC[indexsmn][irow][icol];
366                 } // row
367             }     // col
368           
369           if (indexDDL < 4)
370             {
371               ismn = indexsmn + indexDDL * 6;
372               idet = 0;
373             }
374           else if (indexDDL == 4)
375             {
376               if (indexsmn < 6)
377                 {
378                   ismn = indexsmn;
379                 }
380               else if (indexsmn >= 6 && indexsmn < 12)
381                 {
382                   ismn = indexsmn + 12;
383                 }
384               idet = 1;
385             }
386           else if (indexDDL == 5)
387             {
388               if (indexsmn < 6)
389                 {
390                   ismn = indexsmn + 6;
391                 }
392               else if (indexsmn >= 6 && indexsmn < 12)
393                 {
394                   ismn = indexsmn + 6;
395                 }
396               idet = 1;
397             }
398
399           pmdclust->DoClust(idet,ismn,fCellADC,pmdcont);
400           Int_t nentries1 = pmdcont->GetEntries();
401
402           AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
403
404           for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
405             {
406               AliPMDcluster *pmdcl = 
407                 (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
408               idet        = pmdcl->GetDetector();
409               ismn        = pmdcl->GetSMN();
410               clusdata[0] = pmdcl->GetClusX();
411               clusdata[1] = pmdcl->GetClusY();
412               clusdata[2] = pmdcl->GetClusADC();
413               clusdata[3] = pmdcl->GetClusCells();
414               clusdata[4] = pmdcl->GetClusSigmaX();
415               clusdata[5] = pmdcl->GetClusSigmaY();
416
417               AddRecPoint(idet,ismn,clusdata);
418
419               Int_t ncell = (Int_t) clusdata[3];
420               for(Int_t ihit = 0; ihit < ncell; ihit++)
421                 {
422                   Int_t celldataX = pmdcl->GetClusCellX(ihit);
423                   Int_t celldataY = pmdcl->GetClusCellY(ihit);
424                   AddRecHit(celldataX, celldataY);
425                 }
426               branch2->Fill();
427               ResetRechit();
428
429             }
430           pmdcont->Clear();
431           
432           branch1->Fill();
433           ResetRecpoint();
434
435
436         } // smn
437
438       for (Int_t i=0; i<iSMN; i++)
439         {
440           for (Int_t j=0; j<kRow; j++) delete [] precpvADC[i][j];
441         }
442       for (Int_t i=0; i<iSMN; i++) delete [] precpvADC[i];
443       delete precpvADC;
444     } // DDL Loop
445   
446   ResetCellADC();
447   
448   //   delete the pointers
449   delete pmdclust;
450   delete pmdcont;
451
452 }
453 // ------------------------------------------------------------------------- //
454
455 void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
456 {
457   // Converts RAW data to recpoints after running clustering
458   // algorithm on CPV and PREshower plane
459   //
460
461   Float_t  clusdata[6];
462   TObjArray pmdddlcont;
463   TObjArray *pmdcont = new TObjArray();
464
465   AliPMDClustering *pmdclust = new AliPMDClusteringV1();
466
467   pmdclust->SetEdepCut(fEcut);
468
469   fRunLoader->GetEvent(ievt);
470
471   ResetRecpoint();
472
473   fTreeR = fPMDLoader->TreeR();
474   if (fTreeR == 0x0)
475     {
476       fPMDLoader->MakeTree("R");
477       fTreeR = fPMDLoader->TreeR();
478     }
479   Int_t bufsize = 16000;
480   TBranch *branch1 = fTreeR->Branch("PMDRecpoint", &fRecpoints, bufsize); 
481   TBranch *branch2 = fTreeR->Branch("PMDRechit", &fRechits, bufsize); 
482
483   const Int_t kDDL = 6;
484   const Int_t kRow = 48;
485   const Int_t kCol = 96;
486
487   Int_t idet = 0;
488   Int_t iSMN = 0;
489   
490   for (Int_t indexDDL = 0; indexDDL < kDDL; indexDDL++)
491     {
492       if (indexDDL < 4)
493         {
494           iSMN = 6;
495         }
496       else if (indexDDL >= 4)
497         {
498           iSMN = 12;
499         }
500       Int_t ***precpvADC;
501       precpvADC = new int **[iSMN];
502       for (Int_t i=0; i<iSMN; i++) precpvADC[i] = new int *[kRow];
503       for (Int_t i=0; i<iSMN;i++)
504         {
505           for (Int_t j=0; j<kRow; j++) precpvADC[i][j] = new int [kCol];
506         }
507       for (Int_t i = 0; i < iSMN; i++)
508         {
509           for (Int_t j = 0; j < kRow; j++)
510             {
511               for (Int_t k = 0; k < kCol; k++)
512                 {
513                   precpvADC[i][j][k] = 0;
514                 }
515             }
516         }
517       ResetCellADC();
518       rawReader->Reset();
519       AliPMDRawStream pmdinput(rawReader);
520       rawReader->Select("PMD", indexDDL, indexDDL);
521       
522       pmdinput.DdlData(&pmdddlcont);
523     
524       Int_t indexsmn = 0;
525       Int_t ientries = pmdddlcont.GetEntries();
526       for (Int_t ient = 0; ient < ientries; ient++)
527         {
528           AliPMDddldata *pmdddl = (AliPMDddldata*)pmdddlcont.UncheckedAt(ient);
529           
530           Int_t det = pmdddl->GetDetector();
531           Int_t smn = pmdddl->GetSMN();
532           //Int_t mcm = pmdddl->GetMCM();
533           //Int_t chno = pmdddl->GetChannel();
534           Int_t row = pmdddl->GetRow();
535           Int_t col = pmdddl->GetColumn();
536           Int_t sig = pmdddl->GetSignal();
537
538           Float_t sig1 = (Float_t) sig;
539           // CALIBRATION
540           Float_t gain = fCalibData->GetGainFact(det,smn,row,col);
541           //printf("sig = %d gain = %f\n",sig,gain);
542           sig = (Int_t) (sig1*gain);
543
544
545           if (indexDDL < 4)
546             {
547               if (det != 0)
548                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
549                               indexDDL, det));
550               indexsmn = smn - indexDDL * 6;
551             }
552           else if (indexDDL == 4)
553             {
554               if (det != 1)
555                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
556                               indexDDL, det));
557               if (smn < 6)
558                 {
559                   indexsmn = smn;
560                 }
561               else if (smn >= 18 && smn < 24)
562                 {
563                   indexsmn = smn - 12;
564                 }
565             }
566           else if (indexDDL == 5)
567             {
568               if (det != 1)
569                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
570                               indexDDL, det));
571               if (smn >= 6 && smn < 12)
572                 {
573                   indexsmn = smn - 6;
574                 }
575               else if (smn >= 12 && smn < 18)
576                 {
577                   indexsmn = smn - 6;
578                 }
579             }         
580           precpvADC[indexsmn][row][col] = sig;
581
582         }
583       
584       pmdddlcont.Clear();
585
586       Int_t ismn = 0;
587       for (Int_t indexsmn = 0; indexsmn < iSMN; indexsmn++)
588         {
589           ResetCellADC();
590           for (Int_t irow = 0; irow < kRow; irow++)
591             {
592               for (Int_t icol = 0; icol < kCol; icol++)
593                 {
594                   fCellADC[irow][icol] = 
595                     (Double_t) precpvADC[indexsmn][irow][icol];
596                 } // row
597             }     // col
598
599           
600           if (indexDDL < 4)
601             {
602               ismn = indexsmn + indexDDL * 6;
603               idet = 0;
604             }
605           else if (indexDDL == 4)
606             {
607               if (indexsmn < 6)
608                 {
609                   ismn = indexsmn;
610                 }
611               else if (indexsmn >= 6 && indexsmn < 12)
612                 {
613                   ismn = indexsmn + 12;
614                 }
615               idet = 1;
616             }
617           else if (indexDDL == 5)
618             {
619               if (indexsmn < 6)
620                 {
621                   ismn = indexsmn + 6;
622                 }
623               else if (indexsmn >= 6 && indexsmn < 12)
624                 {
625                   ismn = indexsmn + 6;
626                 }
627               idet = 1;
628             }
629
630           pmdclust->DoClust(idet,ismn,fCellADC,pmdcont);
631           Int_t nentries1 = pmdcont->GetEntries();
632
633           AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
634
635           for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
636             {
637               AliPMDcluster *pmdcl = 
638                 (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
639               idet        = pmdcl->GetDetector();
640               ismn        = pmdcl->GetSMN();
641               clusdata[0] = pmdcl->GetClusX();
642               clusdata[1] = pmdcl->GetClusY();
643               clusdata[2] = pmdcl->GetClusADC();
644               clusdata[3] = pmdcl->GetClusCells();
645               clusdata[4] = pmdcl->GetClusSigmaX();
646               clusdata[5] = pmdcl->GetClusSigmaY();
647
648               AddRecPoint(idet,ismn,clusdata);
649
650               Int_t ncell = (Int_t) clusdata[3];
651               for(Int_t ihit = 0; ihit < ncell; ihit++)
652                 {
653                   Int_t celldataX = pmdcl->GetClusCellX(ihit);
654                   Int_t celldataY = pmdcl->GetClusCellY(ihit);
655                   AddRecHit(celldataX, celldataY);
656                 }
657               branch2->Fill();
658               ResetRechit();
659
660             }
661           pmdcont->Clear();
662           
663           branch1->Fill();
664           ResetRecpoint();
665
666
667         } // smn
668
669       for (Int_t i=0; i<iSMN; i++)
670         {
671           for (Int_t j=0; j<kRow; j++) delete [] precpvADC[i][j];
672         }
673       for (Int_t i=0; i<iSMN; i++) delete [] precpvADC[i];
674       delete precpvADC;
675     } // DDL Loop
676   
677   ResetCellADC();
678   
679   fPMDLoader = fRunLoader->GetLoader("PMDLoader");  
680   fPMDLoader->WriteRecPoints("OVERWRITE");
681
682   //   delete the pointers
683   delete pmdclust;
684   delete pmdcont;
685
686 }
687 // ------------------------------------------------------------------------- //
688 void AliPMDClusterFinder::SetCellEdepCut(Float_t ecut)
689 {
690   fEcut = ecut;
691 }
692 // ------------------------------------------------------------------------- //
693 void AliPMDClusterFinder::AddRecPoint(Int_t idet,Int_t ismn,Float_t *clusdata)
694 {
695   // Add Reconstructed points
696   //
697   TClonesArray &lrecpoints = *fRecpoints;
698   AliPMDrecpoint1 *newrecpoint;
699   newrecpoint = new AliPMDrecpoint1(idet, ismn, clusdata);
700   new(lrecpoints[fNpoint++]) AliPMDrecpoint1(newrecpoint);
701   delete newrecpoint;
702 }
703 // ------------------------------------------------------------------------- //
704 void AliPMDClusterFinder::AddRecHit(Int_t celldataX,Int_t celldataY)
705 {
706   // Add associated cell hits to the Reconstructed points
707   //
708   TClonesArray &lrechits = *fRechits;
709   AliPMDrechit *newrechit;
710   newrechit = new AliPMDrechit(celldataX, celldataY);
711   new(lrechits[fNhit++]) AliPMDrechit(newrechit);
712   delete newrechit;
713 }
714 // ------------------------------------------------------------------------- //
715 void AliPMDClusterFinder::ResetCellADC()
716 {
717   // Reset the individual cell ADC value to zero
718   //
719   for(Int_t irow = 0; irow < fgkRow; irow++)
720     {
721       for(Int_t icol = 0; icol < fgkCol; icol++)
722         {
723           fCellADC[irow][icol] = 0.;
724         }
725     }
726 }
727 // ------------------------------------------------------------------------- //
728
729 void AliPMDClusterFinder::ResetRecpoint()
730 {
731   // Clear the list of reconstructed points
732   fNpoint = 0;
733   if (fRecpoints) fRecpoints->Clear();
734 }
735 // ------------------------------------------------------------------------- //
736 void AliPMDClusterFinder::ResetRechit()
737 {
738   // Clear the list of reconstructed points
739   fNhit = 0;
740   if (fRechits) fRechits->Clear();
741 }
742 // ------------------------------------------------------------------------- //
743 void AliPMDClusterFinder::Load()
744 {
745   // Load all the *.root files
746   //
747   fPMDLoader->LoadDigits("READ");
748   fPMDLoader->LoadRecPoints("recreate");
749 }
750 // ------------------------------------------------------------------------- //
751 void AliPMDClusterFinder::LoadClusters()
752 {
753   // Load all the *.root files
754   //
755   fPMDLoader->LoadRecPoints("recreate");
756 }
757 // ------------------------------------------------------------------------- //
758 void AliPMDClusterFinder::UnLoad()
759 {
760   // Unload all the *.root files
761   //
762   fPMDLoader->UnloadDigits();
763   fPMDLoader->UnloadRecPoints();
764 }
765 // ------------------------------------------------------------------------- //
766 void AliPMDClusterFinder::UnLoadClusters()
767 {
768   // Unload all the *.root files
769   //
770   fPMDLoader->UnloadRecPoints();
771 }
772 // ------------------------------------------------------------------------- //
773
774 AliPMDCalibData* AliPMDClusterFinder::GetCalibData() const
775 {
776   // The run number will be centralized in AliCDBManager,
777   // you don't need to set it here!
778   // Added by ZA
779   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("PMD/Calib/Data");
780   
781   if(!entry){
782     AliWarning("Calibration object retrieval failed! Dummy calibration will be used.");
783     
784     // this just remembers the actual default storage. No problem if it is null.
785     AliCDBStorage *origStorage = AliCDBManager::Instance()->GetDefaultStorage();
786     AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
787     
788     entry = AliCDBManager::Instance()->Get("PMD/Calib/Data");
789     
790     // now reset the original default storage to AliCDBManager...
791     AliCDBManager::Instance()->SetDefaultStorage(origStorage);  
792   }
793   
794   AliPMDCalibData *calibdata=0;
795   if (entry) calibdata = (AliPMDCalibData*) entry->GetObject();
796   
797   if (!calibdata)  AliError("No calibration data from calibration database !");
798   
799   return calibdata;
800 }