]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDClusterFinder.cxx
Coding conventions (F.Carminati)
[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::GetDefaultEventFolderName(),
77                                   "UPDATE");
78   
79   if (!fRunLoader)
80    {
81      Error("Open","Can not open session for file %s.",file);
82    }
83   
84   fRunLoader->LoadgAlice();
85   gAlice = fRunLoader->GetAliRun();
86   
87   if (gAlice)
88     {
89       printf("<AliPMDdigitizer::Open> ");
90       printf("AliRun object found on file.\n");
91     }
92   else
93     {
94       printf("<AliPMDdigitizer::Open> ");
95       printf("Could not find AliRun object.\n");
96     }
97   fPMDLoader = fRunLoader->GetLoader("PMDLoader");
98   if (fPMDLoader == 0x0)
99     {
100       cerr<<"OpengAlice : Can not find PMD or PMDLoader\n";
101     }
102
103   const char *cDR = strstr(option,"DR");
104
105   if (cDR)
106     {
107       fPMDLoader->LoadDigits("READ");
108       fPMDLoader->LoadRecPoints("recreate");
109     }
110 }
111
112 void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
113 {
114   // Converts digits to recpoints after running clustering
115   // algorithm on CPV plane and PREshower plane
116   //
117   Int_t    det = 0,smn = 0;
118   Int_t    xpos,ypos;
119   Float_t  adc;
120   Int_t    isup;
121   Int_t    idet;
122   Float_t  clusdata[7];
123
124   TObjArray *pmdcont = new TObjArray();
125   AliPMDcluster  *pmdcl  = new AliPMDcluster;
126   AliPMDClustering *pmdclust = new AliPMDClustering();
127   pmdclust->SetDebug(fDebug);
128   pmdclust->SetEdepCut(fEcut);
129
130   fRunLoader->GetEvent(ievt);
131   //cout << " ***** Beginning::Digits2RecPoints *****" << endl;
132   fTreeD = fPMDLoader->TreeD();
133   if (fTreeD == 0x0)
134     {
135       cout << " Can not get TreeD" << endl;
136     }
137   AliPMDdigit  *pmddigit;
138   TBranch *branch = fTreeD->GetBranch("PMDDigit");
139   branch->SetAddress(&fDigits);
140
141   ResetRecpoint();
142   fTreeR = fPMDLoader->TreeR();
143   if (fTreeR == 0x0)
144     {
145       fPMDLoader->MakeTree("R");
146       fTreeR = fPMDLoader->TreeR();
147     }
148
149   Int_t bufsize = 16000;
150   fTreeR->Branch("PMDRecpoint", &fRecpoints, bufsize); 
151
152   Int_t nmodules = (Int_t) fTreeD->GetEntries();
153   
154   for (Int_t imodule = 0; imodule < nmodules; imodule++)
155     {
156       ResetCellADC();
157       fTreeD->GetEntry(imodule); 
158       Int_t nentries = fDigits->GetLast();
159       for (Int_t ient = 0; ient < nentries+1; ient++)
160         {
161           pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
162           
163           det    = pmddigit->GetDetector();
164           smn    = pmddigit->GetSMNumber();
165           xpos   = pmddigit->GetRow();
166           ypos   = pmddigit->GetColumn();
167           adc    = pmddigit->GetADC();
168           //Int_t trno   = pmddigit->GetTrackNumber();
169
170           fCellADC[xpos][ypos] = (Double_t) adc;
171         }
172
173       idet = det;
174       isup = smn;
175       pmdclust->DoClust(fCellADC,pmdcont);
176       
177       Int_t nentries1 = pmdcont->GetEntries();
178       cout << " nentries1 = " << nentries1 << endl;
179       for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
180         {
181           clusdata[0] = (Float_t) idet;
182           clusdata[1] = (Float_t) isup;
183               
184           pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
185               
186           clusdata[2] = pmdcl->GetClusX();
187           clusdata[3] = pmdcl->GetClusY();
188           clusdata[4] = pmdcl->GetClusADC();
189           clusdata[5] = pmdcl->GetClusCells();
190           clusdata[6] = pmdcl->GetClusRadius();
191           
192           AddRecPoint(clusdata);
193         }
194       pmdcont->Clear();
195       
196       fTreeR->Fill();
197       ResetRecpoint();
198
199     } // modules
200
201   ResetCellADC();
202   fPMDLoader = fRunLoader->GetLoader("PMDLoader");  
203   fPMDLoader->WriteRecPoints("OVERWRITE");
204
205   //   delete the pointers
206   delete pmdclust;
207   delete pmdcont;
208     
209   //  cout << " ***** End::Digits2RecPoints *****" << endl;
210 }
211
212 void AliPMDClusterFinder::SetCellEdepCut(Float_t ecut)
213 {
214   fEcut = ecut;
215 }
216 void AliPMDClusterFinder::SetDebug(Int_t idebug)
217 {
218   fDebug = idebug;
219 }
220
221 void AliPMDClusterFinder::AddRecPoint(Float_t *clusdata)
222 {
223   // Add Reconstructed points
224   //
225   TClonesArray &lrecpoints = *fRecpoints;
226   AliPMDrecpoint1 *newrecpoint;
227   newrecpoint = new AliPMDrecpoint1(clusdata);
228   new(lrecpoints[fNpoint++]) AliPMDrecpoint1(newrecpoint);
229   delete newrecpoint;
230 }
231 void AliPMDClusterFinder::ResetCellADC()
232 {
233   // Reset the individual cell ADC value to zero
234   //
235   for(Int_t irow = 0; irow < fgkRow; irow++)
236     {
237       for(Int_t icol = 0; icol < fgkCol; icol++)
238         {
239           fCellADC[irow][icol] = 0.;
240         }
241     }
242 }
243
244 void AliPMDClusterFinder::ResetRecpoint()
245 {
246   // Clear the list of reconstructed points
247   fNpoint = 0;
248   if (fRecpoints) fRecpoints->Clear();
249 }
250 void AliPMDClusterFinder::UnLoad(Option_t *option)
251 {
252   // Unload all the *.root files
253   //
254   const char *cR = strstr(option,"R");
255
256   fRunLoader->UnloadgAlice();
257
258   if (cR)
259     {
260       fPMDLoader->UnloadDigits();
261       fPMDLoader->UnloadRecPoints();
262     }
263 }