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