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