]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDClusterFinder.cxx
removed pmdcontainer
[u/mrichter/AliRoot.git] / PMD / AliPMDClusterFinder.cxx
1 //-----------------------------------------------------//
2 //                                                     //
3 //           Date   : August 05 2003                   //
4 //  This reads the file PMD.digits.root(TreeD),        //
5 //  calls the Clustering algorithm and stores the      //
6 //  clustering output in PMD.RecPoints.root(TreeR)     // 
7 //                                                     //
8 //-----------------------------------------------------//
9
10 #include <Riostream.h>
11 #include <TBRIK.h>
12 #include <TNode.h>
13 #include <TTree.h>
14 #include <TGeometry.h>
15 #include <TObjArray.h>
16 #include <TClonesArray.h>
17 #include <TFile.h>
18 #include <TNtuple.h>
19 #include <TParticle.h>
20
21 #include "AliRun.h"
22 #include "AliPMD.h"
23 #include "AliDetector.h"
24 #include "AliRunLoader.h"
25 #include "AliLoader.h"
26 #include "AliHeader.h"
27
28 #include "AliPMDdigit.h"
29 #include "AliPMDClusterFinder.h"
30 #include "AliPMDClustering.h"
31 #include "AliPMDcluster.h"
32 #include "AliPMDrecpoint.h"
33
34
35 ClassImp(AliPMDClusterFinder)
36 //
37 // Constructor
38 //
39 AliPMDClusterFinder::AliPMDClusterFinder()
40 {
41   if (!fRecpoints) fRecpoints = new TClonesArray("AliPMDrecpoint", 1000);  
42   fNpoint = 0;
43
44   for (Int_t i = 0; i < fTotSM; i++)
45     {
46       for (Int_t j = 0; j < fNCell; j++)
47         {
48           for (Int_t k = 0; k < fNCell; k++)
49             {
50               fCPV[i][j][k] = 0.; 
51               fPMD[i][j][k] = 0.; 
52             }
53         }
54     }
55
56 }
57 AliPMDClusterFinder::~AliPMDClusterFinder()
58 {
59   delete fRecpoints;
60 }
61 //
62 // Member functions
63 //
64 void AliPMDClusterFinder::OpengAliceFile(Char_t *file, Option_t *option)
65 {
66
67   fRunLoader = AliRunLoader::Open(file,AliConfig::fgkDefaultEventFolderName,
68                                   "UPDATE");
69   
70   if (!fRunLoader)
71    {
72      Error("Open","Can not open session for file %s.",file);
73    }
74   
75   fRunLoader->LoadgAlice();
76   fRunLoader->LoadHeader();
77   fRunLoader->LoadKinematics();
78
79   gAlice = fRunLoader->GetAliRun();
80   
81   if (gAlice)
82     {
83       printf("<AliPMDdigitizer::Open> ");
84       printf("AliRun object found on file.\n");
85     }
86   else
87     {
88       printf("<AliPMDdigitizer::Open> ");
89       printf("Could not find AliRun object.\n");
90     }
91   PMD  = (AliPMD*)gAlice->GetDetector("PMD");
92   pmdloader = fRunLoader->GetLoader("PMDLoader");
93   if (pmdloader == 0x0)
94     {
95       cerr<<"OpengAlice : Can not find PMD or PMDLoader\n";
96     }
97
98   const char *cDR = strstr(option,"DR");
99
100   if (cDR)
101     {
102       pmdloader->LoadDigits("READ");
103       pmdloader->LoadRecPoints("recreate");
104     }
105 }
106
107 void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
108 {
109   Int_t    det,smn;
110   Int_t    cellno;
111   Int_t    xpos,ypos;
112   Float_t  adc;
113   Int_t    isup, ix, iy;
114   Int_t    idet;
115   Double_t d[72][72];
116   Float_t  clusdata[7];
117   Int_t    fMessage = 1;
118
119   fRunLoader->GetEvent(ievt);
120   //cout << " ***** Beginning::Digits2RecPoints *****" << endl;
121   treeD = pmdloader->TreeD();
122   if (treeD == 0x0)
123     {
124       cout << " Can not get TreeD" << endl;
125     }
126   AliPMDdigit  *pmddigit;
127   TBranch *branch = treeD->GetBranch("PMDDigit");
128   branch->SetAddress(&fDigits);
129
130   ResetRecpoint();
131   treeR = pmdloader->TreeR();
132   if (treeR == 0x0)
133     {
134       pmdloader->MakeTree("R");
135       treeR = pmdloader->TreeR();
136     }
137
138   Int_t bufsize = 16000;
139   treeR->Branch("PMDRecpoint", &fRecpoints, bufsize); 
140
141   Int_t nmodules = (Int_t) treeD->GetEntries();
142   
143   for (Int_t imodule = 0; imodule < nmodules; imodule++)
144     {
145       treeD->GetEntry(imodule); 
146       Int_t nentries = fDigits->GetLast();
147       for (Int_t ient = 0; ient < nentries+1; ient++)
148         {
149           pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
150           
151           det    = pmddigit->GetDetector();
152           smn    = pmddigit->GetSMNumber();
153           cellno = pmddigit->GetCellNumber();
154           adc    = pmddigit->GetADC();
155
156           ypos = cellno/fNCell;
157           xpos = cellno - ypos*fNCell;
158
159           if (det == 0)
160             {
161               fPMD[smn][xpos][ypos] = adc;
162             }
163           else if(det == 1)
164             {
165               fCPV[smn][xpos][ypos] = adc;
166             }
167         }
168     } // modules
169
170   //
171   // Clustering started
172   //
173
174   TObjArray *pmdcont = new TObjArray();
175
176   AliPMDcluster  *pmdcl  = new AliPMDcluster;
177
178   AliPMDClustering *pmdclust = new AliPMDClustering();
179   pmdclust->SetMessage(fMessage);
180   
181   for (idet = 0; idet < 2; idet++)
182     {
183       for (isup = 0; isup < fTotSM; isup++)
184         {
185           for (ix = 0; ix < fNCell; ix++)
186             {
187               for (iy = 0; iy < fNCell; iy++)
188                 {
189                   if (idet == 0)
190                     {
191                       d[ix][iy] = (Double_t) fPMD[isup][ix][iy];
192                     }
193                   else if (idet == 1)
194                     {
195                       d[ix][iy] = (Double_t) fCPV[isup][ix][iy];
196                     }
197                     
198                 }
199             }
200           pmdclust->DoClust(idet,isup,d,pmdcont);
201           
202           Int_t nentries = pmdcont->GetEntries();
203           cout << " nentries = " << nentries << endl;
204           for (Int_t ient = 0; ient < nentries; ient++)
205             {
206               clusdata[0] = (Float_t) idet;
207               clusdata[1] = (Float_t) isup;
208               
209               pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient);
210               
211               clusdata[2] = pmdcl->GetClusX();
212               clusdata[3] = pmdcl->GetClusY();
213               clusdata[4] = pmdcl->GetClusADC();
214               clusdata[5] = pmdcl->GetClusCells();
215               clusdata[6] = pmdcl->GetClusRadius();
216               
217               AddRecPoint(clusdata);
218             }
219           pmdcont->Clear();
220           
221           treeR->Fill();
222           ResetRecpoint();
223         }  // SuperModule
224     }  // Detector
225   
226   ResetCellADC();
227   
228   pmdloader->WriteRecPoints("OVERWRITE");
229
230   //   delete the pointers
231   delete pmdclust;
232   delete pmdcont;
233     
234   //  cout << " ***** End::Digits2RecPoints *****" << endl;
235 }
236
237
238 void AliPMDClusterFinder::AddRecPoint(Float_t *clusdata)
239 {
240   TClonesArray &lrecpoints = *fRecpoints;
241   AliPMDrecpoint *newrecpoint;
242   newrecpoint = new AliPMDrecpoint(clusdata);
243   new(lrecpoints[fNpoint++]) AliPMDrecpoint(newrecpoint);
244   delete newrecpoint;
245 }
246 void AliPMDClusterFinder::ResetCellADC()
247 {
248   for (Int_t i = 0; i < fTotSM; i++)
249     {
250       for (Int_t j = 0; j < fNCell; j++)
251         {
252           for (Int_t k = 0; k < fNCell; k++)
253             {
254               fCPV[i][j][k] = 0.; 
255               fPMD[i][j][k] = 0.; 
256             }
257         }
258     }
259 }
260
261 void AliPMDClusterFinder::ResetRecpoint()
262 {
263   fNpoint = 0;
264   if (fRecpoints) fRecpoints->Clear();
265 }
266 void AliPMDClusterFinder::UnLoad(Option_t *option)
267 {
268   const char *cR = strstr(option,"R");
269
270   fRunLoader->UnloadgAlice();
271   fRunLoader->UnloadHeader();
272   fRunLoader->UnloadKinematics();
273
274   if (cR)
275     {
276       pmdloader->UnloadDigits();
277     }
278 }