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