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