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