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