1 //-----------------------------------------------------//
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) //
8 //-----------------------------------------------------//
10 #include <Riostream.h>
14 #include <TGeometry.h>
15 #include <TObjArray.h>
16 #include <TClonesArray.h>
19 #include <TParticle.h>
23 #include "AliDetector.h"
24 #include "AliRunLoader.h"
25 #include "AliLoader.h"
26 #include "AliHeader.h"
28 #include "AliPMDdigit.h"
29 #include "AliPMDClusterFinder.h"
30 #include "AliPMDClustering.h"
31 #include "AliPMDContainer.h"
32 #include "AliPMDcluster.h"
33 #include "AliPMDrecpoint.h"
36 ClassImp(AliPMDClusterFinder)
40 AliPMDClusterFinder::AliPMDClusterFinder()
42 if (!fRecpoints) fRecpoints = new TClonesArray("AliPMDrecpoint", 1000);
45 for (Int_t i = 0; i < fTotSM; i++)
47 for (Int_t j = 0; j < fNCell; j++)
49 for (Int_t k = 0; k < fNCell; k++)
58 AliPMDClusterFinder::~AliPMDClusterFinder()
65 void AliPMDClusterFinder::OpengAliceFile(Char_t *file, Option_t *option)
68 fRunLoader = AliRunLoader::Open(file,AliConfig::fgkDefaultEventFolderName,
73 Error("Open","Can not open session for file %s.",file);
76 fRunLoader->LoadgAlice();
77 fRunLoader->LoadHeader();
78 fRunLoader->LoadKinematics();
80 gAlice = fRunLoader->GetAliRun();
84 printf("<AliPMDdigitizer::Open> ");
85 printf("AliRun object found on file.\n");
89 printf("<AliPMDdigitizer::Open> ");
90 printf("Could not find AliRun object.\n");
92 PMD = (AliPMD*)gAlice->GetDetector("PMD");
93 pmdloader = fRunLoader->GetLoader("PMDLoader");
96 cerr<<"OpengAlice : Can not find PMD or PMDLoader\n";
99 const char *cDR = strstr(option,"DR");
103 pmdloader->LoadDigits("READ");
104 pmdloader->LoadRecPoints("recreate");
108 void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
120 fRunLoader->GetEvent(ievt);
121 //cout << " ***** Beginning::Digits2RecPoints *****" << endl;
122 treeD = pmdloader->TreeD();
125 cout << " Can not get TreeD" << endl;
127 AliPMDdigit *pmddigit;
128 TBranch *branch = treeD->GetBranch("PMDDigit");
129 branch->SetAddress(&fDigits);
132 treeR = pmdloader->TreeR();
135 pmdloader->MakeTree("R");
136 treeR = pmdloader->TreeR();
139 Int_t bufsize = 16000;
140 treeR->Branch("PMDRecpoint", &fRecpoints, bufsize);
142 Int_t nmodules = (Int_t) treeD->GetEntries();
144 for (Int_t imodule = 0; imodule < nmodules; imodule++)
146 treeD->GetEntry(imodule);
147 Int_t nentries = fDigits->GetLast();
148 for (Int_t ient = 0; ient < nentries+1; ient++)
150 pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
152 det = pmddigit->GetDetector();
153 smn = pmddigit->GetSMNumber();
154 cellno = pmddigit->GetCellNumber();
155 adc = pmddigit->GetADC();
157 ypos = cellno/fNCell;
158 xpos = cellno - ypos*fNCell;
162 fPMD[smn][xpos][ypos] = adc;
166 fCPV[smn][xpos][ypos] = adc;
172 // Clustering started
175 // AliPMDContainer *pmdcont = new AliPMDContainer;
177 TObjArray *pmdcont = new TObjArray();
179 AliPMDcluster *pmdcl = new AliPMDcluster;
181 AliPMDClustering *pmdclust = new AliPMDClustering();
182 pmdclust->SetMessage(fMessage);
184 for (idet = 0; idet < 2; idet++)
186 for (isup = 0; isup < fTotSM; isup++)
188 for (ix = 0; ix < fNCell; ix++)
190 for (iy = 0; iy < fNCell; iy++)
194 d[ix][iy] = (Double_t) fPMD[isup][ix][iy];
198 d[ix][iy] = (Double_t) fCPV[isup][ix][iy];
203 pmdclust->DoClust(idet,isup,d,pmdcont);
205 Int_t nentries = pmdcont->GetEntries();
206 cout << " nentries = " << nentries << endl;
207 for (Int_t ient = 0; ient < nentries; ient++)
209 clusdata[0] = (Float_t) idet;
210 clusdata[1] = (Float_t) isup;
212 pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient);
214 clusdata[2] = pmdcl->GetClusX();
215 clusdata[3] = pmdcl->GetClusY();
216 clusdata[4] = pmdcl->GetClusADC();
217 clusdata[5] = pmdcl->GetClusCells();
218 clusdata[6] = pmdcl->GetClusRadius();
220 AddRecPoint(clusdata);
231 pmdloader->WriteRecPoints("OVERWRITE");
233 // delete the pointers
237 // cout << " ***** End::Digits2RecPoints *****" << endl;
241 void AliPMDClusterFinder::AddRecPoint(Float_t *clusdata)
243 TClonesArray &lrecpoints = *fRecpoints;
244 AliPMDrecpoint *newrecpoint;
245 newrecpoint = new AliPMDrecpoint(clusdata);
246 new(lrecpoints[fNpoint++]) AliPMDrecpoint(newrecpoint);
249 void AliPMDClusterFinder::ResetCellADC()
251 for (Int_t i = 0; i < fTotSM; i++)
253 for (Int_t j = 0; j < fNCell; j++)
255 for (Int_t k = 0; k < fNCell; k++)
264 void AliPMDClusterFinder::ResetRecpoint()
267 if (fRecpoints) fRecpoints->Clear();
269 void AliPMDClusterFinder::UnLoad(Option_t *option)
271 const char *cR = strstr(option,"R");
273 fRunLoader->UnloadgAlice();
274 fRunLoader->UnloadHeader();
275 fRunLoader->UnloadKinematics();
279 pmdloader->UnloadDigits();