]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDClusterFinder.cxx
cellno is replaced by row and column
[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 <TMath.h>
27 #include <TBRIK.h>
28 #include <TNode.h>
29 #include <TTree.h>
30 #include <TGeometry.h>
31 #include <TObjArray.h>
32 #include <TClonesArray.h>
33 #include <TFile.h>
34 #include <TNtuple.h>
35 #include <TParticle.h>
36
37 #include "AliRun.h"
38 #include "AliPMD.h"
39 #include "AliDetector.h"
40 #include "AliRunLoader.h"
41 #include "AliLoader.h"
42 #include "AliHeader.h"
43
44 #include "AliPMDdigit.h"
45 #include "AliPMDClusterFinder.h"
46 #include "AliPMDClustering.h"
47 #include "AliPMDcluster.h"
48 #include "AliPMDrecpoint1.h"
49
50
51 ClassImp(AliPMDClusterFinder)
52
53 AliPMDClusterFinder::AliPMDClusterFinder()
54 {
55 //
56 // Default Constructor
57 //
58   if (!fRecpoints) fRecpoints = new TClonesArray("AliPMDrecpoint1", 1000);  
59   fNpoint = 0;
60
61   fDebug = 0;
62   fEcut  = 0.;
63
64 }
65 AliPMDClusterFinder::~AliPMDClusterFinder()
66 {
67   // Destructor
68   delete fRecpoints;
69 }
70
71 void AliPMDClusterFinder::OpengAliceFile(Char_t *file, Option_t *option)
72 {
73   // Loads galice.root file and corresponding header, kinematics
74   // hits and sdigits or digits depending on the option
75   //
76   fRunLoader = AliRunLoader::Open(file,AliConfig::fgkDefaultEventFolderName,
77                                   "UPDATE");
78   
79   if (!fRunLoader)
80    {
81      Error("Open","Can not open session for file %s.",file);
82    }
83   
84   fRunLoader->LoadgAlice();
85   fRunLoader->LoadHeader();
86   gAlice = fRunLoader->GetAliRun();
87   
88   if (gAlice)
89     {
90       printf("<AliPMDdigitizer::Open> ");
91       printf("AliRun object found on file.\n");
92     }
93   else
94     {
95       printf("<AliPMDdigitizer::Open> ");
96       printf("Could not find AliRun object.\n");
97     }
98   fPMD  = (AliPMD*)gAlice->GetDetector("PMD");
99   fPMDLoader = fRunLoader->GetLoader("PMDLoader");
100   if (fPMDLoader == 0x0)
101     {
102       cerr<<"OpengAlice : Can not find PMD or PMDLoader\n";
103     }
104
105   const char *cDR = strstr(option,"DR");
106
107   if (cDR)
108     {
109       fPMDLoader->LoadDigits("READ");
110       fPMDLoader->LoadRecPoints("recreate");
111     }
112 }
113
114 void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
115 {
116   // Converts digits to recpoints after running clustering
117   // algorithm on CPV plane and PREshower plane
118   //
119   Int_t    det = 0,smn = 0;
120   Int_t    xpos,ypos;
121   Float_t  adc;
122   Int_t    isup;
123   Int_t    idet;
124   Float_t  clusdata[7];
125
126   TObjArray *pmdcont = new TObjArray();
127   AliPMDcluster  *pmdcl  = new AliPMDcluster;
128   AliPMDClustering *pmdclust = new AliPMDClustering();
129   pmdclust->SetDebug(fDebug);
130   pmdclust->SetEdepCut(fEcut);
131
132   fRunLoader->GetEvent(ievt);
133   //cout << " ***** Beginning::Digits2RecPoints *****" << endl;
134   fTreeD = fPMDLoader->TreeD();
135   if (fTreeD == 0x0)
136     {
137       cout << " Can not get TreeD" << endl;
138     }
139   AliPMDdigit  *pmddigit;
140   TBranch *branch = fTreeD->GetBranch("PMDDigit");
141   branch->SetAddress(&fDigits);
142
143   ResetRecpoint();
144   fTreeR = fPMDLoader->TreeR();
145   if (fTreeR == 0x0)
146     {
147       fPMDLoader->MakeTree("R");
148       fTreeR = fPMDLoader->TreeR();
149     }
150
151   Int_t bufsize = 16000;
152   fTreeR->Branch("PMDRecpoint", &fRecpoints, bufsize); 
153
154   Int_t nmodules = (Int_t) fTreeD->GetEntries();
155   
156   for (Int_t imodule = 0; imodule < nmodules; imodule++)
157     {
158       ResetCellADC();
159       fTreeD->GetEntry(imodule); 
160       Int_t nentries = fDigits->GetLast();
161       for (Int_t ient = 0; ient < nentries+1; ient++)
162         {
163           pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
164           
165           det    = pmddigit->GetDetector();
166           smn    = pmddigit->GetSMNumber();
167           xpos   = pmddigit->GetRow();
168           ypos   = pmddigit->GetColumn();
169           adc    = pmddigit->GetADC();
170           //Int_t trno   = pmddigit->GetTrackNumber();
171
172           fCellADC[xpos][ypos] = (Double_t) adc;
173         }
174
175       idet = det;
176       isup = smn;
177       pmdclust->DoClust(fCellADC,pmdcont);
178       
179       Int_t nentries1 = pmdcont->GetEntries();
180       cout << " nentries1 = " << nentries1 << endl;
181       for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
182         {
183           clusdata[0] = (Float_t) idet;
184           clusdata[1] = (Float_t) isup;
185               
186           pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
187               
188           clusdata[2] = pmdcl->GetClusX();
189           clusdata[3] = pmdcl->GetClusY();
190           clusdata[4] = pmdcl->GetClusADC();
191           clusdata[5] = pmdcl->GetClusCells();
192           clusdata[6] = pmdcl->GetClusRadius();
193           
194           AddRecPoint(clusdata);
195         }
196       pmdcont->Clear();
197       
198       fTreeR->Fill();
199       ResetRecpoint();
200
201     } // modules
202
203   ResetCellADC();
204   
205   fPMDLoader->WriteRecPoints("OVERWRITE");
206
207   //   delete the pointers
208   delete pmdclust;
209   delete pmdcont;
210     
211   //  cout << " ***** End::Digits2RecPoints *****" << endl;
212 }
213
214 void AliPMDClusterFinder::SetCellEdepCut(Float_t ecut)
215 {
216   fEcut = ecut;
217 }
218 void AliPMDClusterFinder::SetDebug(Int_t idebug)
219 {
220   fDebug = idebug;
221 }
222
223 void AliPMDClusterFinder::AddRecPoint(Float_t *clusdata)
224 {
225   // Add Reconstructed points
226   //
227   TClonesArray &lrecpoints = *fRecpoints;
228   AliPMDrecpoint1 *newrecpoint;
229   newrecpoint = new AliPMDrecpoint1(clusdata);
230   new(lrecpoints[fNpoint++]) AliPMDrecpoint1(newrecpoint);
231   delete newrecpoint;
232 }
233 void AliPMDClusterFinder::ResetCellADC()
234 {
235   // Reset the individual cell ADC value to zero
236   //
237   for(Int_t irow = 0; irow < fgkRow; irow++)
238     {
239       for(Int_t icol = 0; icol < fgkCol; icol++)
240         {
241           fCellADC[irow][icol] = 0.;
242         }
243     }
244 }
245
246 void AliPMDClusterFinder::ResetRecpoint()
247 {
248   // Clear the list of reconstructed points
249   fNpoint = 0;
250   if (fRecpoints) fRecpoints->Clear();
251 }
252 void AliPMDClusterFinder::UnLoad(Option_t *option)
253 {
254   // Unload all the *.root files
255   //
256   const char *cR = strstr(option,"R");
257
258   fRunLoader->UnloadgAlice();
259   fRunLoader->UnloadHeader();
260   fRunLoader->UnloadKinematics();
261
262   if (cR)
263     {
264       fPMDLoader->UnloadDigits();
265     }
266 }