]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDClusterFinder.cxx
findisocell method is removed, isolated cells are tagged while doing the clustering...
[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 "AliPMDisocell.h"
44 #include "AliPMDRawStream.h"
45 #include "AliPMDCalibData.h"
46 #include "AliPMDPedestal.h"
47 #include "AliPMDddldata.h"
48 #include "AliPMDHotData.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   fRecoParam(0x0),
67   fTreeD(0),
68   fTreeR(0),
69   fDigits(new TClonesArray("AliPMDdigit", 1000)),
70   fRecpoints(new TClonesArray("AliPMDrecpoint1", 1000)),
71   fRechits(new TClonesArray("AliPMDrechit", 1000)),
72   fNpoint(0),
73   fNhit(0),
74   fDetNo(0),
75   fEcut(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   fRecoParam(0x0),
89   fTreeD(0),
90   fTreeR(0),
91   fDigits(new TClonesArray("AliPMDdigit", 1000)),
92   fRecpoints(new TClonesArray("AliPMDrecpoint1", 1000)),
93   fRechits(new TClonesArray("AliPMDrechit", 1000)),
94   fNpoint(0),
95   fNhit(0),
96   fDetNo(0),
97   fEcut(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   fRecoParam(0x0),
112   fTreeD(0),
113   fTreeR(0),
114   fDigits(NULL),
115   fRecpoints(NULL),
116   fRechits(NULL),
117   fNpoint(0),
118   fNhit(0),
119   fDetNo(0),
120   fEcut(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(Int_t ievt)
153 {
154   // Converts digits to recpoints after running clustering
155   // algorithm on CPV plane and PREshower plane
156   //
157
158   Int_t    det = 0,smn = 0;
159   Int_t    xpos,ypos;
160   Float_t  adc;
161   Int_t    ismn;
162   Int_t    idet;
163   Float_t  clusdata[6];
164
165   AliPMDisocell *pmdiso = 0x0;
166
167   TObjArray *pmdcont = new TObjArray();
168   TObjArray *pmdisocell = new TObjArray();
169
170   AliPMDClustering *pmdclust = new AliPMDClusteringV1();
171
172   // fetch the recoparam object from database
173   fRecoParam = AliPMDReconstructor::GetRecoParam();
174
175   if(fRecoParam == 0x0)
176     {
177        AliFatal("No Reco Param found for PMD!!!");
178     }
179
180
181   fRunLoader->GetEvent(ievt);
182
183
184   fTreeD = fPMDLoader->TreeD();
185   if (fTreeD == 0x0)
186     {
187       AliFatal("AliPMDClusterFinder: Can not get TreeD");
188
189     }
190   AliPMDdigit  *pmddigit;
191   TBranch *branch = fTreeD->GetBranch("PMDDigit");
192   branch->SetAddress(&fDigits);
193
194   ResetRecpoint();
195
196   fTreeR = fPMDLoader->TreeR();
197   if (fTreeR == 0x0)
198     {
199       fPMDLoader->MakeTree("R");
200       fTreeR = fPMDLoader->TreeR();
201     }
202
203   Int_t bufsize = 16000;
204   TBranch * branch1 = fTreeR->Branch("PMDRecpoint", &fRecpoints, bufsize); 
205   TBranch * branch2 = fTreeR->Branch("PMDRechit", &fRechits, bufsize); 
206
207   Int_t nmodules = (Int_t) fTreeD->GetEntries();
208
209   for (Int_t imodule = 0; imodule < nmodules; imodule++)
210     {
211       ResetCellADC();
212       fTreeD->GetEntry(imodule); 
213       Int_t nentries = fDigits->GetLast();
214       for (Int_t ient = 0; ient < nentries+1; ient++)
215         {
216           pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
217           
218           det    = pmddigit->GetDetector();
219           smn    = pmddigit->GetSMNumber();
220           xpos   = pmddigit->GetRow();
221           ypos   = pmddigit->GetColumn();
222           adc    = pmddigit->GetADC();
223           if(xpos < 0 || xpos > 48 || ypos < 0 || ypos > 96)
224             {
225               AliError(Form("*Row %d and Column NUMBER %d NOT Valid *",
226                               xpos, ypos));
227               continue; 
228             }
229
230           // Hot cell - set the cell adc = 0
231           Float_t hotflag = fCalibHot->GetHotChannel(det,smn,xpos,ypos);
232           if (hotflag == 1) adc = 0;
233           // CALIBRATION
234           Float_t gain = fCalibGain->GetGainFact(det,smn,xpos,ypos);
235           // printf("adc = %d gain = %f\n",adc,gain);
236           
237           adc = adc*gain;
238
239           //Int_t trno   = pmddigit->GetTrackNumber();
240           fCellTrack[xpos][ypos] = pmddigit->GetTrackNumber();
241           fCellPid[xpos][ypos]   = pmddigit->GetTrackPid();
242           fCellADC[xpos][ypos]   = (Double_t) adc;
243         }
244
245       idet = det;
246       ismn = smn;
247
248       // Set the minimum noise cut per module before clustering
249
250       Int_t imod = idet*24 + ismn;
251       fEcut = fRecoParam->GetNoiseCut(imod);   // default
252       // fEcut = fRecoParam->GetPbPbParam()->GetNoiseCut(imod);
253       // fEcut = fRecoParam->GetPPParam()->GetNoiseCut(imod);
254       // fEcut = fRecoParam->GetCosmicParam()->GetNoiseCut(imod);
255
256       pmdclust->SetEdepCut(fEcut);
257
258
259       pmdclust->DoClust(idet,ismn,fCellTrack,fCellPid,fCellADC,
260                         pmdisocell,pmdcont);
261
262       Int_t nentries1 = pmdcont->GetEntries();
263
264       AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
265
266       for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
267         {
268           AliPMDcluster *pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
269           idet        = pmdcl->GetDetector();
270           ismn        = pmdcl->GetSMN();
271           clusdata[0] = pmdcl->GetClusX();
272           clusdata[1] = pmdcl->GetClusY();
273           clusdata[2] = pmdcl->GetClusADC();
274           clusdata[3] = pmdcl->GetClusCells();
275           clusdata[4] = pmdcl->GetClusSigmaX();
276           clusdata[5] = pmdcl->GetClusSigmaY();
277
278           AddRecPoint(idet,ismn,clusdata);
279
280           Int_t ncell = (Int_t) clusdata[3];
281           for(Int_t ihit = 0; ihit < ncell; ihit++)
282             {
283               Int_t celldataX = pmdcl->GetClusCellX(ihit);
284               Int_t celldataY = pmdcl->GetClusCellY(ihit);
285               Int_t celldataTr = pmdcl->GetClusCellTrack(ihit);
286               Int_t celldataPid   = pmdcl->GetClusCellPid(ihit);
287               Float_t celldataAdc = pmdcl->GetClusCellAdc(ihit);
288               AddRecHit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
289             }
290           branch2->Fill();
291           ResetRechit();
292         }
293       pmdcont->Delete();
294
295           // Added single isolated cell for offline gain calibration
296           nentries1 = pmdisocell->GetEntries();
297           AliDebug(1,Form("Total number of isolated single cell clusters = %d",nentries1));
298           
299           for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
300             {
301               pmdiso = (AliPMDisocell*)pmdisocell->UncheckedAt(ient1);
302               idet = pmdiso->GetDetector();
303               ismn = pmdiso->GetSmn();
304               clusdata[0] = (Float_t) pmdiso->GetRow();
305               clusdata[1] = (Float_t) pmdiso->GetCol();
306               clusdata[2] = pmdiso->GetADC();
307               clusdata[3] = 1.;
308               clusdata[4] = -99.;
309               clusdata[5] = -99.;
310       
311               AddRecPoint(idet,ismn,clusdata);
312             }
313           pmdisocell->Delete();
314       
315       branch1->Fill();
316       ResetRecpoint();
317
318     } // modules
319
320   ResetCellADC();
321   fPMDLoader = fRunLoader->GetLoader("PMDLoader");  
322   fPMDLoader->WriteRecPoints("OVERWRITE");
323
324   //   delete the pointers
325   delete pmdclust;
326   delete pmdcont;
327   delete pmdisocell;
328     
329 }
330 // ------------------------------------------------------------------------- //
331
332 void AliPMDClusterFinder::Digits2RecPoints(TTree *digitsTree,
333                                            TTree *clustersTree)
334 {
335   // Converts digits to recpoints after running clustering
336   // algorithm on CPV plane and PREshower plane
337   //
338   // This algorithm is called during the reconstruction from digits
339
340   Int_t    det = 0,smn = 0;
341   Int_t    xpos,ypos;
342   Float_t  adc;
343   Int_t    ismn;
344   Int_t    idet;
345   Float_t  clusdata[6];
346
347   AliPMDcluster *pmdcl = 0x0;
348   AliPMDisocell *pmdiso = 0x0;
349
350   TObjArray *pmdcont = new TObjArray();
351   TObjArray *pmdisocell = new TObjArray();
352   AliPMDClustering *pmdclust = new AliPMDClusteringV1();
353
354   // Fetch the reco param object
355
356   fRecoParam = AliPMDReconstructor::GetRecoParam();
357   if(fRecoParam == 0x0)
358     {
359        AliFatal("No Reco Param found for PMD!!!");
360     }
361
362
363   AliPMDdigit  *pmddigit;
364   TBranch *branch = digitsTree->GetBranch("PMDDigit");
365   branch->SetAddress(&fDigits);
366
367   ResetRecpoint();
368
369   Int_t bufsize = 16000;
370   TBranch * branch1 = clustersTree->Branch("PMDRecpoint", &fRecpoints, bufsize); 
371   TBranch * branch2 = clustersTree->Branch("PMDRechit", &fRechits, bufsize); 
372
373   Int_t nmodules = (Int_t) digitsTree->GetEntries();
374
375   for (Int_t imodule = 0; imodule < nmodules; imodule++)
376     {
377
378       Int_t totADCMod = 0;
379       ResetCellADC();
380       digitsTree->GetEntry(imodule); 
381       Int_t nentries = fDigits->GetLast();
382       for (Int_t ient = 0; ient < nentries+1; ient++)
383         {
384           pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
385           
386           det    = pmddigit->GetDetector();
387           smn    = pmddigit->GetSMNumber();
388           xpos   = pmddigit->GetRow();
389           ypos   = pmddigit->GetColumn();
390           adc    = pmddigit->GetADC();
391
392           if(det < 0 || det > 1)
393             {
394               AliError(Form("*CPV/PRE NUMBER WRONG %d *",det));
395               continue; 
396             }
397           if(smn == -1 || smn > 23)
398             {
399               AliError(Form("*MODULE NUMBER WRONG %d *",smn));
400               continue; 
401             }
402
403           if(xpos < 0 || xpos > 47 || ypos < 0 || ypos > 95)
404             {
405               AliError(Form("*Row %d and Column NUMBER %d NOT Valid *",
406                             xpos, ypos));
407               continue; 
408             }
409           
410           // Pedestal Subtraction
411           Int_t   pedmeanrms = fCalibPed->GetPedMeanRms(det,smn,xpos,ypos);
412           Int_t   pedrms1    = (Int_t) pedmeanrms%100;
413           Float_t pedrms     = (Float_t)pedrms1/10.;
414           Float_t pedmean    = (Float_t) (pedmeanrms - pedrms1)/1000.0;
415           //printf("%f %f\n",pedmean, pedrms);
416
417           Float_t adc1 = adc - (pedmean + 3.0*pedrms);
418
419           // Hot cell - set the cell adc = 0
420           Float_t hotflag = fCalibHot->GetHotChannel(det,smn,xpos,ypos);
421           if (hotflag == 1.) adc1 = 0;
422
423           // CALIBRATION
424           Float_t gain = fCalibGain->GetGainFact(det,smn,xpos,ypos);
425           // printf("adc = %d gain = %f\n",adc,gain);
426
427           adc = adc1*gain;
428
429           fCellTrack[xpos][ypos] = pmddigit->GetTrackNumber();
430           fCellPid[xpos][ypos] = pmddigit->GetTrackPid();
431           fCellADC[xpos][ypos] = (Double_t) adc;
432
433           totADCMod += (Int_t) adc;
434
435         }
436
437       idet = det;
438       ismn = smn;
439
440       if (totADCMod <= 0) continue;
441
442       // Set the minimum noise cut per module before clustering
443
444       Int_t imod = idet*24 + ismn;
445       fEcut = fRecoParam->GetNoiseCut(imod);       // default
446       // fEcut = fRecoParam->GetPbPbParam()->GetNoiseCut(imod);
447       // fEcut = fRecoParam->GetPPParam()->GetNoiseCut(imod);
448       // fEcut = fRecoParam->GetCosmicParam()->GetNoiseCut(imod);
449
450
451       pmdclust->SetEdepCut(fEcut);
452
453       pmdclust->DoClust(idet,ismn,fCellTrack,fCellPid,fCellADC,
454                         pmdisocell,pmdcont);
455       
456       Int_t nentries1 = pmdcont->GetEntries();
457
458       AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
459
460       for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
461         {
462           pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
463           idet        = pmdcl->GetDetector();
464           ismn        = pmdcl->GetSMN();
465           clusdata[0] = pmdcl->GetClusX();
466           clusdata[1] = pmdcl->GetClusY();
467           clusdata[2] = pmdcl->GetClusADC();
468           clusdata[3] = pmdcl->GetClusCells();
469           clusdata[4] = pmdcl->GetClusSigmaX();
470           clusdata[5] = pmdcl->GetClusSigmaY();
471
472           AddRecPoint(idet,ismn,clusdata);
473
474           Int_t ncell = (Int_t) clusdata[3];
475           for(Int_t ihit = 0; ihit < ncell; ihit++)
476             {
477               Int_t celldataX = pmdcl->GetClusCellX(ihit);
478               Int_t celldataY = pmdcl->GetClusCellY(ihit);
479               Int_t celldataTr = pmdcl->GetClusCellTrack(ihit);
480               Int_t celldataPid   = pmdcl->GetClusCellPid(ihit);
481               Float_t celldataAdc = pmdcl->GetClusCellAdc(ihit);
482               AddRecHit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
483             }
484           branch2->Fill();
485           ResetRechit();
486         }
487       pmdcont->Delete();
488
489       // Added single isolated cell for offline gain calibration
490       nentries1 = pmdisocell->GetEntries();
491       AliDebug(1,Form("Total number of isolated single cell clusters = %d",nentries1));
492       
493       for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
494         {
495           pmdiso = (AliPMDisocell*)pmdisocell->UncheckedAt(ient1);
496           idet = pmdiso->GetDetector();
497           ismn = pmdiso->GetSmn();
498           clusdata[0] = (Float_t) pmdiso->GetRow();
499           clusdata[1] = (Float_t) pmdiso->GetCol();
500           clusdata[2] = pmdiso->GetADC();
501           clusdata[3] = 1.;
502           clusdata[4] = -99.;
503           clusdata[5] = -99.;
504           
505           AddRecPoint(idet,ismn,clusdata);
506         }
507       pmdisocell->Delete();
508       
509       branch1->Fill();
510       ResetRecpoint();
511
512     } // modules
513
514
515   ResetCellADC();
516
517   //   delete the pointers
518   delete pmdclust;
519   delete pmdcont;
520   delete pmdisocell;
521 }
522 // ------------------------------------------------------------------------- //
523
524 void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
525                                            TTree *clustersTree)
526 {
527   // Converts RAW data to recpoints after running clustering
528   // algorithm on CPV and PREshower plane
529   //
530   // This method is called at the time of reconstruction from RAW data
531
532
533   AliPMDddldata *pmdddl = 0x0;
534   AliPMDcluster *pmdcl  = 0x0;
535   AliPMDisocell *pmdiso = 0x0;
536
537
538   Float_t  clusdata[6];
539   TObjArray pmdddlcont;
540
541   TObjArray *pmdcont = new TObjArray();
542   TObjArray *pmdisocell = new TObjArray();
543
544   AliPMDClustering *pmdclust = new AliPMDClusteringV1();
545
546   // open the ddl file info to know the module
547   TString ddlinfofileName(gSystem->Getenv("ALICE_ROOT"));
548   ddlinfofileName += "/PMD/PMD_ddl_info.dat";
549   
550   ifstream infileddl;
551   infileddl.open(ddlinfofileName.Data(), ios::in); // ascii file
552   if(!infileddl) AliError("Could not read the ddl info file");
553
554   Int_t ddlno;
555   Int_t modno;
556   Int_t modulePerDDL;
557   Int_t moduleddl[6];
558
559   for(Int_t jddl = 0; jddl < 6; jddl++)
560     {
561       if (infileddl.eof()) break;
562       infileddl >> ddlno >> modulePerDDL;
563       moduleddl[jddl] = modulePerDDL;
564
565       if (modulePerDDL == 0) continue;
566       for (Int_t im = 0; im < modulePerDDL; im++)
567         {
568           infileddl >> modno;
569         }
570     }
571
572   infileddl.close();
573
574   // Set the minimum noise cut per module before clustering
575
576   fRecoParam = AliPMDReconstructor::GetRecoParam();
577
578   if(fRecoParam == 0x0)
579     {
580        AliFatal("No Reco Param found for PMD!!!");
581     }
582
583   ResetRecpoint();
584
585   Int_t bufsize = 16000;
586   TBranch *branch1 = clustersTree->Branch("PMDRecpoint", &fRecpoints, bufsize); 
587
588   TBranch * branch2 = clustersTree->Branch("PMDRechit", &fRechits, bufsize); 
589
590   const Int_t kRow = 48;
591   const Int_t kCol = 96;
592
593   Int_t idet = 0;
594   Int_t iSMN = 0;
595
596   Int_t indexDDL = -1;
597   AliPMDRawStream pmdinput(rawReader);
598
599   while ((indexDDL = pmdinput.DdlData(&pmdddlcont)) >=0)
600     {
601       iSMN = moduleddl[indexDDL];
602
603       Int_t ***precpvADC;
604       precpvADC = new int **[iSMN];
605       for (Int_t i=0; i<iSMN; i++) precpvADC[i] = new int *[kRow];
606       for (Int_t i=0; i<iSMN;i++)
607         {
608           for (Int_t j=0; j<kRow; j++) precpvADC[i][j] = new int [kCol];
609         }
610       for (Int_t i = 0; i < iSMN; i++)
611         {
612           for (Int_t j = 0; j < kRow; j++)
613             {
614               for (Int_t k = 0; k < kCol; k++)
615                 {
616                   precpvADC[i][j][k] = 0;
617                 }
618             }
619         }
620       ResetCellADC();
621
622       Int_t indexsmn = 0;
623       Int_t ientries = pmdddlcont.GetEntries();
624       for (Int_t ient = 0; ient < ientries; ient++)
625         {
626           pmdddl = (AliPMDddldata*)pmdddlcont.UncheckedAt(ient);
627           
628           Int_t det = pmdddl->GetDetector();
629           Int_t smn = pmdddl->GetSMN();
630           //Int_t mcm = pmdddl->GetMCM();
631           //Int_t chno = pmdddl->GetChannel();
632           Int_t row = pmdddl->GetRow();
633           Int_t col = pmdddl->GetColumn();
634           Int_t sig = pmdddl->GetSignal();
635
636
637           if(det < 0 || det > 1)
638             {
639               AliError(Form("*CPV/PRE NUMBER WRONG %d *",det));
640               continue; 
641             }
642           if(smn < 0 || smn > 23)
643             {
644               AliError(Form("*MODULE NUMBER WRONG %d *",smn));
645               continue; 
646             }
647           if(row < 0 || row > 47 || col < 0 || col > 95)
648             {
649               AliError(Form("*Row %d and Column NUMBER %d NOT Valid *",
650                             row, col));
651
652               continue; 
653             }
654
655
656           // Pedestal Subtraction
657           Int_t   pedmeanrms = fCalibPed->GetPedMeanRms(det,smn,row,col);
658           Int_t   pedrms1    = (Int_t) pedmeanrms%100;
659           Float_t pedrms     = (Float_t)pedrms1/10.;
660           Float_t pedmean    = (Float_t) (pedmeanrms - pedrms1)/1000.0;
661
662           //printf("%f %f\n",pedmean, pedrms);
663
664           // Float_t sig1 = (Float_t) sig;
665           Float_t sig1 = (Float_t) sig - (pedmean + 3.0*pedrms);
666
667
668           // Hot cell - set the cell adc = 0
669           Float_t hotflag = fCalibHot->GetHotChannel(det,smn,row,col);
670           if (hotflag == 1.) sig1 = 0;
671
672           // CALIBRATION
673           Float_t gain = fCalibGain->GetGainFact(det,smn,row,col);
674           //printf("sig = %d gain = %f\n",sig,gain);
675           sig = (Int_t) (sig1*gain);
676
677           if (indexDDL == 0)
678             {
679               if (det != 0)
680                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
681                               indexDDL, det));
682               if (iSMN == 6)
683                 {
684                   indexsmn = smn;
685                 }
686               else if (iSMN == 12)
687                 {
688                   if (smn < 6)
689                     indexsmn = smn;
690                   else if (smn >= 18 && smn < 24)
691                     indexsmn = smn-12;
692                 }
693             }
694           else if (indexDDL >= 1 && indexDDL < 4)
695             {
696               if (det != 0)
697                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
698                               indexDDL, det));
699               indexsmn = smn - indexDDL * 6;
700             }
701           else if (indexDDL == 4)
702             {
703               if (det != 1)
704                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
705                               indexDDL, det));
706               if (smn < 6)
707                 {
708                   indexsmn = smn;
709                 }
710               else if (smn >= 18 && smn < 24)
711                 {
712                   indexsmn = smn - 12;
713                 }
714             }
715           else if (indexDDL == 5)
716             {
717               if (det != 1)
718                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
719                               indexDDL, det));
720               if (smn >= 6 && smn < 18)
721                 {
722                   indexsmn = smn - 6;
723                 }
724             }         
725
726           precpvADC[indexsmn][row][col] = sig;
727         }
728       
729       pmdddlcont.Delete();
730
731       Int_t totAdcMod = 0;
732
733       Int_t ismn = 0;
734       for (indexsmn = 0; indexsmn < iSMN; indexsmn++)
735         {
736           ResetCellADC();
737           totAdcMod = 0;
738           for (Int_t irow = 0; irow < kRow; irow++)
739             {
740               for (Int_t icol = 0; icol < kCol; icol++)
741                 {
742                   fCellTrack[irow][icol] = -1;
743                   fCellPid[irow][icol]   = -1;
744
745                   fCellADC[irow][icol] = 
746                     (Double_t) precpvADC[indexsmn][irow][icol];
747                   totAdcMod += precpvADC[indexsmn][irow][icol];
748                 } // row
749             }     // col
750           
751           if (indexDDL == 0)
752             {
753               if (iSMN == 6)
754                 {
755                   ismn = indexsmn;
756                 }
757               else if (iSMN == 12)
758                 {
759                   
760                   if (indexsmn < 6)
761                     ismn = indexsmn;
762                   else if (indexsmn >= 6 && indexsmn < 12)
763                     ismn = indexsmn + 12;
764                 }
765               idet = 0;
766             }
767           else if (indexDDL >= 1 && indexDDL < 4)
768             {
769               ismn = indexsmn + indexDDL * 6;
770               idet = 0;
771             }
772           else if (indexDDL == 4)
773             {
774               if (indexsmn < 6)
775                 {
776                   ismn = indexsmn;
777                 }
778               else if (indexsmn >= 6 && indexsmn < 12)
779                 {
780                   ismn = indexsmn + 12;
781                 }
782               idet = 1;
783             }
784           else if (indexDDL == 5)
785             {
786               ismn = indexsmn + 6;
787               idet = 1;
788             }
789
790           if (totAdcMod <= 0) continue;
791
792           Int_t imod = idet*24 + ismn;
793
794           fEcut = fRecoParam->GetNoiseCut(imod);       // default
795           // fEcut = fRecoParam->GetPbPbParam()->GetNoiseCut(imod);
796           // fEcut = fRecoParam->GetPPParam()->GetNoiseCut(imod);
797           // fEcut = fRecoParam->GetCosmicParam()->GetNoiseCut(imod);
798
799
800           pmdclust->SetEdepCut(fEcut);
801
802           pmdclust->DoClust(idet,ismn,fCellTrack,fCellPid,fCellADC,
803                             pmdisocell,pmdcont);
804
805           Int_t nentries1 = pmdcont->GetEntries();
806
807           AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
808
809           for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
810             {
811               pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
812               idet        = pmdcl->GetDetector();
813               ismn        = pmdcl->GetSMN();
814               clusdata[0] = pmdcl->GetClusX();
815               clusdata[1] = pmdcl->GetClusY();
816               clusdata[2] = pmdcl->GetClusADC();
817               clusdata[3] = pmdcl->GetClusCells();
818               clusdata[4] = pmdcl->GetClusSigmaX();
819               clusdata[5] = pmdcl->GetClusSigmaY();
820
821               AddRecPoint(idet,ismn,clusdata);
822
823               Int_t ncell = (Int_t) clusdata[3];
824               for(Int_t ihit = 0; ihit < ncell; ihit++)
825                 {
826                   Int_t celldataX = pmdcl->GetClusCellX(ihit);
827                   Int_t celldataY = pmdcl->GetClusCellY(ihit);
828                   Int_t celldataTr = pmdcl->GetClusCellTrack(ihit);
829                   Int_t celldataPid   = pmdcl->GetClusCellPid(ihit);
830                   Float_t celldataAdc = pmdcl->GetClusCellAdc(ihit);
831                   AddRecHit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
832                 }
833               branch2->Fill();
834               ResetRechit();
835
836             }
837           pmdcont->Delete();
838
839
840           // Added single isolated cell for offline gain calibration
841           nentries1 = pmdisocell->GetEntries();
842           AliDebug(1,Form("Total number of isolated single cell clusters = %d",nentries1));
843           
844           for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
845             {
846               pmdiso = (AliPMDisocell*)pmdisocell->UncheckedAt(ient1);
847               idet = pmdiso->GetDetector();
848               ismn = pmdiso->GetSmn();
849               clusdata[0] = (Float_t) pmdiso->GetRow();
850               clusdata[1] = (Float_t) pmdiso->GetCol();
851               clusdata[2] = pmdiso->GetADC();
852               clusdata[3] = 1.;
853               clusdata[4] = -99.;
854               clusdata[5] = -99.;
855       
856               AddRecPoint(idet,ismn,clusdata);
857             }
858           pmdisocell->Delete();
859           
860           branch1->Fill();
861           ResetRecpoint();
862
863
864         } // smn
865
866       for (Int_t i=0; i<iSMN; i++)
867         {
868           for (Int_t j=0; j<kRow; j++) delete [] precpvADC[i][j];
869         }
870       for (Int_t i=0; i<iSMN; i++) delete [] precpvADC[i];
871       delete [] precpvADC;
872
873     } // DDL Loop
874
875   
876   ResetCellADC();
877   
878   //   delete the pointers
879   delete pmdclust;
880   delete pmdcont;
881   delete pmdisocell;
882
883 }
884 // ------------------------------------------------------------------------- //
885
886 void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
887 {
888   // Converts RAW data to recpoints after running clustering
889   // algorithm on CPV and PREshower plane
890   //
891
892   Float_t  clusdata[6];
893
894   TObjArray pmdddlcont;
895
896   AliPMDcluster *pmdcl  = 0x0;
897   AliPMDisocell *pmdiso  = 0x0;
898
899
900   TObjArray *pmdcont = new TObjArray();
901   TObjArray *pmdisocell = new TObjArray();
902
903   AliPMDClustering *pmdclust = new AliPMDClusteringV1();
904
905   // open the ddl file info to know the module
906   TString ddlinfofileName(gSystem->Getenv("ALICE_ROOT"));
907   ddlinfofileName += "/PMD/PMD_ddl_info.dat";
908   
909   ifstream infileddl;
910   infileddl.open(ddlinfofileName.Data(), ios::in); // ascii file
911   if(!infileddl) AliError("Could not read the ddl info file");
912
913   Int_t ddlno;
914   Int_t modno;
915   Int_t modulePerDDL;
916   Int_t moduleddl[6];
917
918   for(Int_t jddl = 0; jddl < 6; jddl++)
919     {
920       if (infileddl.eof()) break;
921       infileddl >> ddlno >> modulePerDDL;
922       moduleddl[jddl] = modulePerDDL;
923
924       if (modulePerDDL == 0) continue;
925       for (Int_t im = 0; im < modulePerDDL; im++)
926         {
927           infileddl >> modno;
928         }
929     }
930
931   infileddl.close();
932
933   // Set the minimum noise cut per module before clustering
934
935   fRecoParam = AliPMDReconstructor::GetRecoParam();
936
937   if(fRecoParam == 0x0)
938     {
939        AliFatal("No Reco Param found for PMD!!!");
940     }
941
942
943   fRunLoader->GetEvent(ievt);
944
945   ResetRecpoint();
946
947   fTreeR = fPMDLoader->TreeR();
948   if (fTreeR == 0x0)
949     {
950       fPMDLoader->MakeTree("R");
951       fTreeR = fPMDLoader->TreeR();
952     }
953   Int_t bufsize = 16000;
954   TBranch *branch1 = fTreeR->Branch("PMDRecpoint", &fRecpoints, bufsize); 
955   TBranch *branch2 = fTreeR->Branch("PMDRechit", &fRechits, bufsize); 
956
957   const Int_t kRow = 48;
958   const Int_t kCol = 96;
959
960   Int_t idet = 0;
961   Int_t iSMN = 0;
962
963   AliPMDRawStream pmdinput(rawReader);
964   Int_t indexDDL = -1;
965
966   while ((indexDDL = pmdinput.DdlData(&pmdddlcont)) >=0)
967     {
968       
969       iSMN = moduleddl[indexDDL];
970
971       Int_t ***precpvADC;
972       precpvADC = new int **[iSMN];
973       for (Int_t i=0; i<iSMN; i++) precpvADC[i] = new int *[kRow];
974       for (Int_t i=0; i<iSMN;i++)
975         {
976           for (Int_t j=0; j<kRow; j++) precpvADC[i][j] = new int [kCol];
977         }
978       for (Int_t i = 0; i < iSMN; i++)
979         {
980           for (Int_t j = 0; j < kRow; j++)
981             {
982               for (Int_t k = 0; k < kCol; k++)
983                 {
984                   precpvADC[i][j][k] = 0;
985                 }
986             }
987         }
988       ResetCellADC();
989
990     
991       Int_t indexsmn = 0;
992       Int_t ientries = pmdddlcont.GetEntries();
993       for (Int_t ient = 0; ient < ientries; ient++)
994         {
995           AliPMDddldata *pmdddl = (AliPMDddldata*)pmdddlcont.UncheckedAt(ient);
996           
997           Int_t det = pmdddl->GetDetector();
998           Int_t smn = pmdddl->GetSMN();
999           //Int_t mcm = pmdddl->GetMCM();
1000           //Int_t chno = pmdddl->GetChannel();
1001           Int_t row = pmdddl->GetRow();
1002           Int_t col = pmdddl->GetColumn();
1003           Int_t sig = pmdddl->GetSignal();
1004           if(row < 0 || row > 48 || col < 0 || col > 96)
1005             {
1006               AliError(Form("*Row %d and Column NUMBER %d NOT Valid *",
1007                               row, col));
1008               continue; 
1009             }
1010           // Pedestal Subtraction
1011           Int_t   pedmeanrms = fCalibPed->GetPedMeanRms(det,smn,row,col);
1012           Int_t   pedrms1    = (Int_t) pedmeanrms%100;
1013           Float_t pedrms     = (Float_t)pedrms1/10.;
1014           Float_t pedmean    = (Float_t) (pedmeanrms - pedrms1)/1000.0;
1015           //printf("%f %f\n",pedmean, pedrms);
1016           Float_t sig1 = (Float_t) sig - (pedmean + 3.0*pedrms);
1017
1018           // Hot cell - set the cell adc = 0
1019           Float_t hotflag = fCalibHot->GetHotChannel(det,smn,row,col);
1020           if (hotflag == 1) sig1 = 0;
1021
1022           // CALIBRATION
1023           Float_t gain = fCalibGain->GetGainFact(det,smn,row,col);
1024
1025           //printf("sig = %d gain = %f\n",sig,gain);
1026           sig = (Int_t) (sig1*gain);
1027
1028           if (indexDDL == 0)
1029             {
1030               if (det != 0)
1031                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
1032                               indexDDL, det));
1033               if (iSMN == 6)
1034                 {
1035                   indexsmn = smn;
1036                 }
1037               else if (iSMN == 12)
1038                 {
1039                   if (smn < 6)
1040                     indexsmn = smn;
1041                   else if (smn >= 18 && smn < 24)
1042                     indexsmn = smn-12;
1043                 }
1044             }
1045           else if (indexDDL >= 1 && indexDDL < 4)
1046             {
1047               if (det != 0)
1048                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
1049                               indexDDL, det));
1050               indexsmn = smn - indexDDL * 6;
1051             }
1052           else if (indexDDL == 4)
1053             {
1054               if (det != 1)
1055                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
1056                               indexDDL, det));
1057               if (smn < 6)
1058                 {
1059                   indexsmn = smn;
1060                 }
1061               else if (smn >= 18 && smn < 24)
1062                 {
1063                   indexsmn = smn - 12;
1064                 }
1065             }
1066           else if (indexDDL == 5)
1067             {
1068               if (det != 1)
1069                 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
1070                               indexDDL, det));
1071               if (smn >= 6 && smn < 18)
1072                 {
1073                   indexsmn = smn - 6;
1074                 }
1075             }         
1076           
1077           precpvADC[indexsmn][row][col] = sig;
1078
1079         }
1080       
1081       pmdddlcont.Delete();
1082
1083       Int_t ismn = 0;
1084       for (indexsmn = 0; indexsmn < iSMN; indexsmn++)
1085         {
1086           ResetCellADC();
1087           for (Int_t irow = 0; irow < kRow; irow++)
1088             {
1089               for (Int_t icol = 0; icol < kCol; icol++)
1090                 {
1091                   fCellTrack[irow][icol] = -1;
1092                   fCellPid[irow][icol]   = -1;
1093                   fCellADC[irow][icol] = 
1094                     (Double_t) precpvADC[indexsmn][irow][icol];
1095                 } // row
1096             }     // col
1097
1098
1099           if (indexDDL == 0)
1100             {
1101               if (iSMN == 6)
1102                 {
1103                   ismn = indexsmn;
1104                 }
1105               else if (iSMN == 12)
1106                 {
1107                   
1108                   if (indexsmn < 6)
1109                     ismn = indexsmn;
1110                   else if (indexsmn >= 6 && indexsmn < 12)
1111                     ismn = indexsmn + 12;
1112                 }
1113               idet = 0;
1114             }
1115           else if (indexDDL >= 1 && indexDDL < 4)
1116             {
1117               ismn = indexsmn + indexDDL * 6;
1118               idet = 0;
1119             }
1120           else if (indexDDL == 4)
1121             {
1122               if (indexsmn < 6)
1123                 {
1124                   ismn = indexsmn;
1125                 }
1126               else if (indexsmn >= 6 && indexsmn < 12)
1127                 {
1128                   ismn = indexsmn + 12;
1129                 }
1130               idet = 1;
1131             }
1132           else if (indexDDL == 5)
1133             {
1134               ismn = indexsmn + 6;
1135               idet = 1;
1136             }
1137
1138           Int_t imod = idet*24 + ismn;
1139           fEcut = fRecoParam->GetNoiseCut(imod);       // default
1140           // fEcut = fRecoParam->GetPbPbParam()->GetNoiseCut(imod);
1141           // fEcut = fRecoParam->GetPPParam()->GetNoiseCut(imod);
1142           // fEcut = fRecoParam->GetCosmicParam()->GetNoiseCut(imod);
1143
1144           pmdclust->SetEdepCut(fEcut);
1145
1146           pmdclust->DoClust(idet,ismn,fCellTrack,fCellPid,fCellADC,
1147                             pmdisocell,pmdcont);
1148
1149           Int_t nentries1 = pmdcont->GetEntries();
1150
1151           AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
1152
1153           for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
1154             {
1155               pmdcl       = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
1156               idet        = pmdcl->GetDetector();
1157               ismn        = pmdcl->GetSMN();
1158               clusdata[0] = pmdcl->GetClusX();
1159               clusdata[1] = pmdcl->GetClusY();
1160               clusdata[2] = pmdcl->GetClusADC();
1161               clusdata[3] = pmdcl->GetClusCells();
1162               clusdata[4] = pmdcl->GetClusSigmaX();
1163               clusdata[5] = pmdcl->GetClusSigmaY();
1164
1165               AddRecPoint(idet,ismn,clusdata);
1166
1167               Int_t ncell = (Int_t) clusdata[3];
1168               for(Int_t ihit = 0; ihit < ncell; ihit++)
1169                 {
1170                   Int_t celldataX = pmdcl->GetClusCellX(ihit);
1171                   Int_t celldataY = pmdcl->GetClusCellY(ihit);
1172                   Int_t celldataTr = pmdcl->GetClusCellTrack(ihit);
1173                   Int_t celldataPid = pmdcl->GetClusCellPid(ihit);
1174                   Float_t celldataAdc = pmdcl->GetClusCellAdc(ihit);
1175                   AddRecHit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
1176                 }
1177               branch2->Fill();
1178               ResetRechit();
1179
1180             }
1181           pmdcont->Delete();
1182
1183           // Added single isolated cell for offline gain calibration
1184           nentries1 = pmdisocell->GetEntries();
1185           AliDebug(1,Form("Total number of isolated single cell clusters = %d",nentries1));
1186           
1187           for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
1188             {
1189               pmdiso = (AliPMDisocell*)pmdisocell->UncheckedAt(ient1);
1190               idet = pmdiso->GetDetector();
1191               ismn = pmdiso->GetSmn();
1192               clusdata[0] = (Float_t) pmdiso->GetRow();
1193               clusdata[1] = (Float_t) pmdiso->GetCol();
1194               clusdata[2] = pmdiso->GetADC();
1195               clusdata[3] = 1.;
1196               clusdata[4] = -99.;
1197               clusdata[5] = -99.;
1198       
1199               AddRecPoint(idet,ismn,clusdata);
1200             }
1201           pmdisocell->Delete();
1202           
1203           branch1->Fill();
1204           ResetRecpoint();
1205
1206
1207         } // smn
1208
1209       for (Int_t i=0; i<iSMN; i++)
1210         {
1211           for (Int_t j=0; j<kRow; j++) delete [] precpvADC[i][j];
1212         }
1213       for (Int_t i=0; i<iSMN; i++) delete [] precpvADC[i];
1214       delete precpvADC;
1215     } // DDL Loop
1216
1217
1218   ResetCellADC();
1219   
1220   fPMDLoader = fRunLoader->GetLoader("PMDLoader");  
1221   fPMDLoader->WriteRecPoints("OVERWRITE");
1222
1223   //   delete the pointers
1224   delete pmdclust;
1225   delete pmdcont;
1226   delete pmdisocell;
1227 }
1228 // ------------------------------------------------------------------------- //
1229 void AliPMDClusterFinder::SetCellEdepCut(Float_t ecut)
1230 {
1231   fEcut = ecut;
1232 }
1233 // ------------------------------------------------------------------------- //
1234 void AliPMDClusterFinder::AddRecPoint(Int_t idet,Int_t ismn,Float_t *clusdata)
1235 {
1236   // Add Reconstructed points
1237   //
1238   TClonesArray &lrecpoints = *fRecpoints;
1239   AliPMDrecpoint1 *newrecpoint;
1240   newrecpoint = new AliPMDrecpoint1(idet, ismn, clusdata);
1241   new(lrecpoints[fNpoint++]) AliPMDrecpoint1(newrecpoint);
1242   delete newrecpoint;
1243 }
1244 // ------------------------------------------------------------------------- //
1245 void AliPMDClusterFinder::AddRecHit(Int_t celldataX,Int_t celldataY,
1246                                     Int_t celldataTr, Int_t celldataPid,
1247                                     Float_t celldataAdc)
1248 {
1249   // Add associated cell hits to the Reconstructed points
1250   //
1251   TClonesArray &lrechits = *fRechits;
1252   AliPMDrechit *newrechit;
1253   newrechit = new AliPMDrechit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
1254   new(lrechits[fNhit++]) AliPMDrechit(newrechit);
1255   delete newrechit;
1256 }
1257 // ------------------------------------------------------------------------- //
1258 void AliPMDClusterFinder::ResetCellADC()
1259 {
1260   // Reset the individual cell ADC value to zero
1261   //
1262   for(Int_t irow = 0; irow < fgkRow; irow++)
1263     {
1264       for(Int_t icol = 0; icol < fgkCol; icol++)
1265         {
1266           fCellTrack[irow][icol] = -1;
1267           fCellPid[irow][icol]   = -1;
1268           fCellADC[irow][icol]   = 0.;
1269         }
1270     }
1271 }
1272 // ------------------------------------------------------------------------- //
1273 void AliPMDClusterFinder::ResetRecpoint()
1274 {
1275   // Clear the list of reconstructed points
1276   fNpoint = 0;
1277   if (fRecpoints) fRecpoints->Clear();
1278 }
1279 // ------------------------------------------------------------------------- //
1280 void AliPMDClusterFinder::ResetRechit()
1281 {
1282   // Clear the list of reconstructed points
1283   fNhit = 0;
1284   if (fRechits) fRechits->Clear();
1285 }
1286 // ------------------------------------------------------------------------- //
1287 void AliPMDClusterFinder::Load()
1288 {
1289   // Load all the *.root files
1290   //
1291   fPMDLoader->LoadDigits("READ");
1292   fPMDLoader->LoadRecPoints("recreate");
1293 }
1294 // ------------------------------------------------------------------------- //
1295 void AliPMDClusterFinder::LoadClusters()
1296 {
1297   // Load all the *.root files
1298   //
1299   fPMDLoader->LoadRecPoints("recreate");
1300 }
1301 // ------------------------------------------------------------------------- //
1302 void AliPMDClusterFinder::UnLoad()
1303 {
1304   // Unload all the *.root files
1305   //
1306   fPMDLoader->UnloadDigits();
1307   fPMDLoader->UnloadRecPoints();
1308 }
1309 // ------------------------------------------------------------------------- //
1310 void AliPMDClusterFinder::UnLoadClusters()
1311 {
1312   // Unload all the *.root files
1313   //
1314   fPMDLoader->UnloadRecPoints();
1315 }
1316 // ------------------------------------------------------------------------- //
1317 AliPMDCalibData* AliPMDClusterFinder::GetCalibGain() const
1318 {
1319   // The run number will be centralized in AliCDBManager,
1320   // you don't need to set it here!
1321   // Added by ZA
1322   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("PMD/Calib/Gain");
1323   
1324   if(!entry)  AliFatal("Calibration object retrieval failed! ");
1325   
1326   AliPMDCalibData *calibdata=0;
1327   if (entry) calibdata = (AliPMDCalibData*) entry->GetObject();
1328   
1329   if (!calibdata)  AliFatal("No calibration data from calibration database !");
1330   
1331   return calibdata;
1332 }
1333 // ------------------------------------------------------------------------- //
1334 AliPMDPedestal* AliPMDClusterFinder::GetCalibPed() const
1335 {
1336   // The run number will be centralized in AliCDBManager,
1337   // you don't need to set it here!
1338   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("PMD/Calib/Ped");
1339   
1340   if(!entry) AliFatal("Pedestal object retrieval failed!");
1341     
1342   AliPMDPedestal *pedestal = 0;
1343   if (entry) pedestal = (AliPMDPedestal*) entry->GetObject();
1344   
1345   if (!pedestal)  AliFatal("No pedestal data from pedestal database !");
1346   
1347   return pedestal;
1348 }
1349 //--------------------------------------------------------------------//
1350 AliPMDHotData* AliPMDClusterFinder::GetCalibHot() const
1351 {
1352   // The run number will be centralized in AliCDBManager,
1353   // you don't need to set it here!
1354   AliCDBEntry  *entry = AliCDBManager::Instance()->Get("PMD/Calib/Hot");
1355   
1356   if(!entry) AliFatal("HotData object retrieval failed!");
1357   
1358   AliPMDHotData *hot = 0;
1359   if (entry) hot = (AliPMDHotData*) entry->GetObject();
1360   
1361   if (!hot)  AliFatal("No hot data from  database !");
1362   
1363   return hot;
1364 }
1365