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