]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PMD/AliPMDClusterFinder.cxx
Remove minor warning
[u/mrichter/AliRoot.git] / PMD / AliPMDClusterFinder.cxx
CommitLineData
ed228cbc 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
01709453 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 <TBRIK.h>
27#include <TNode.h>
28#include <TTree.h>
29#include <TGeometry.h>
30#include <TObjArray.h>
31#include <TClonesArray.h>
32#include <TFile.h>
33#include <TNtuple.h>
34#include <TParticle.h>
35
36#include "AliRun.h"
37#include "AliPMD.h"
38#include "AliDetector.h"
39#include "AliRunLoader.h"
40#include "AliLoader.h"
41#include "AliHeader.h"
42
43#include "AliPMDdigit.h"
44#include "AliPMDClusterFinder.h"
45#include "AliPMDClustering.h"
01709453 46#include "AliPMDcluster.h"
96377d57 47#include "AliPMDrecpoint1.h"
01709453 48
49
50ClassImp(AliPMDClusterFinder)
51//
52// Constructor
53//
54AliPMDClusterFinder::AliPMDClusterFinder()
55{
ed228cbc 56 if (!fRecpoints) fRecpoints = new TClonesArray("AliPMDrecpoint1", 1000);
01709453 57 fNpoint = 0;
58
ed228cbc 59 fDebug = 0;
60 fEcut = 0.;
01709453 61
62}
63AliPMDClusterFinder::~AliPMDClusterFinder()
64{
65 delete fRecpoints;
66}
67//
68// Member functions
69//
70void AliPMDClusterFinder::OpengAliceFile(Char_t *file, Option_t *option)
71{
72
73 fRunLoader = AliRunLoader::Open(file,AliConfig::fgkDefaultEventFolderName,
74 "UPDATE");
75
76 if (!fRunLoader)
77 {
78 Error("Open","Can not open session for file %s.",file);
79 }
80
81 fRunLoader->LoadgAlice();
82 fRunLoader->LoadHeader();
83 fRunLoader->LoadKinematics();
84
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 PMD = (AliPMD*)gAlice->GetDetector("PMD");
98 pmdloader = fRunLoader->GetLoader("PMDLoader");
99 if (pmdloader == 0x0)
100 {
101 cerr<<"OpengAlice : Can not find PMD or PMDLoader\n";
102 }
103
104 const char *cDR = strstr(option,"DR");
105
106 if (cDR)
107 {
108 pmdloader->LoadDigits("READ");
109 pmdloader->LoadRecPoints("recreate");
110 }
111}
112
113void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
114{
ed228cbc 115 Int_t det = 0,smn = 0;
01709453 116 Int_t cellno;
117 Int_t xpos,ypos;
118 Float_t adc;
ed228cbc 119 Int_t isup;
01709453 120 Int_t idet;
01709453 121 Float_t clusdata[7];
ed228cbc 122
123 TObjArray *pmdcont = new TObjArray();
124 AliPMDcluster *pmdcl = new AliPMDcluster;
125 AliPMDClustering *pmdclust = new AliPMDClustering();
126 pmdclust->SetDebug(fDebug);
127 pmdclust->SetEdepCut(fEcut);
01709453 128
129 fRunLoader->GetEvent(ievt);
130 //cout << " ***** Beginning::Digits2RecPoints *****" << endl;
131 treeD = pmdloader->TreeD();
132 if (treeD == 0x0)
133 {
134 cout << " Can not get TreeD" << endl;
135 }
136 AliPMDdigit *pmddigit;
137 TBranch *branch = treeD->GetBranch("PMDDigit");
138 branch->SetAddress(&fDigits);
139
140 ResetRecpoint();
141 treeR = pmdloader->TreeR();
142 if (treeR == 0x0)
143 {
144 pmdloader->MakeTree("R");
145 treeR = pmdloader->TreeR();
146 }
147
148 Int_t bufsize = 16000;
149 treeR->Branch("PMDRecpoint", &fRecpoints, bufsize);
150
151 Int_t nmodules = (Int_t) treeD->GetEntries();
152
153 for (Int_t imodule = 0; imodule < nmodules; imodule++)
154 {
ed228cbc 155 ResetCellADC();
01709453 156 treeD->GetEntry(imodule);
157 Int_t nentries = fDigits->GetLast();
158 for (Int_t ient = 0; ient < nentries+1; ient++)
159 {
160 pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
161
162 det = pmddigit->GetDetector();
163 smn = pmddigit->GetSMNumber();
164 cellno = pmddigit->GetCellNumber();
165 adc = pmddigit->GetADC();
ed228cbc 166 //Int_t trno = pmddigit->GetTrackNumber();
01709453 167
ed228cbc 168 xpos = cellno/fCol;
169 ypos = cellno - xpos*fCol;
170 fCellADC[xpos][ypos] = (Double_t) adc;
01709453 171 }
01709453 172
ed228cbc 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++)
01709453 180 {
ed228cbc 181 clusdata[0] = (Float_t) idet;
182 clusdata[1] = (Float_t) isup;
01709453 183
ed228cbc 184 pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
01709453 185
ed228cbc 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();
01709453 191
ed228cbc 192 AddRecPoint(clusdata);
193 }
194 pmdcont->Clear();
195
196 treeR->Fill();
197 ResetRecpoint();
198
199 } // modules
200
01709453 201 ResetCellADC();
202
203 pmdloader->WriteRecPoints("OVERWRITE");
204
205 // delete the pointers
206 delete pmdclust;
207 delete pmdcont;
208
209 // cout << " ***** End::Digits2RecPoints *****" << endl;
210}
211
ed228cbc 212void AliPMDClusterFinder::SetCellEdepCut(Float_t ecut)
213{
214 fEcut = ecut;
215}
216void AliPMDClusterFinder::SetDebug(Int_t idebug)
217{
218 fDebug = idebug;
219}
01709453 220
221void AliPMDClusterFinder::AddRecPoint(Float_t *clusdata)
222{
223 TClonesArray &lrecpoints = *fRecpoints;
ed228cbc 224 AliPMDrecpoint1 *newrecpoint;
225 newrecpoint = new AliPMDrecpoint1(clusdata);
226 new(lrecpoints[fNpoint++]) AliPMDrecpoint1(newrecpoint);
01709453 227 delete newrecpoint;
228}
229void AliPMDClusterFinder::ResetCellADC()
230{
ed228cbc 231 for(Int_t irow = 0; irow < fRow; irow++)
01709453 232 {
ed228cbc 233 for(Int_t icol = 0; icol < fCol; icol++)
01709453 234 {
ed228cbc 235 fCellADC[irow][icol] = 0.;
01709453 236 }
237 }
238}
239
240void AliPMDClusterFinder::ResetRecpoint()
241{
242 fNpoint = 0;
243 if (fRecpoints) fRecpoints->Clear();
244}
245void AliPMDClusterFinder::UnLoad(Option_t *option)
246{
247 const char *cR = strstr(option,"R");
248
249 fRunLoader->UnloadgAlice();
250 fRunLoader->UnloadHeader();
251 fRunLoader->UnloadKinematics();
252
253 if (cR)
254 {
255 pmdloader->UnloadDigits();
256 }
257}