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