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