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