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