1 /***************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //-----------------------------------------------------//
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) //
23 //-----------------------------------------------------//
25 #include <Riostream.h>
30 #include <TGeometry.h>
31 #include <TObjArray.h>
32 #include <TClonesArray.h>
35 #include <TParticle.h>
39 #include "AliDetector.h"
40 #include "AliRunLoader.h"
41 #include "AliLoader.h"
42 #include "AliHeader.h"
44 #include "AliPMDdigit.h"
45 #include "AliPMDClusterFinder.h"
46 #include "AliPMDClustering.h"
47 #include "AliPMDcluster.h"
48 #include "AliPMDrecpoint1.h"
51 ClassImp(AliPMDClusterFinder)
53 AliPMDClusterFinder::AliPMDClusterFinder()
56 // Default Constructor
58 if (!fRecpoints) fRecpoints = new TClonesArray("AliPMDrecpoint1", 1000);
65 AliPMDClusterFinder::~AliPMDClusterFinder()
71 void AliPMDClusterFinder::OpengAliceFile(Char_t *file, Option_t *option)
73 // Loads galice.root file and corresponding header, kinematics
74 // hits and sdigits or digits depending on the option
76 fRunLoader = AliRunLoader::Open(file,AliConfig::fgkDefaultEventFolderName,
81 Error("Open","Can not open session for file %s.",file);
84 fRunLoader->LoadgAlice();
85 gAlice = fRunLoader->GetAliRun();
89 printf("<AliPMDdigitizer::Open> ");
90 printf("AliRun object found on file.\n");
94 printf("<AliPMDdigitizer::Open> ");
95 printf("Could not find AliRun object.\n");
97 fPMDLoader = fRunLoader->GetLoader("PMDLoader");
98 if (fPMDLoader == 0x0)
100 cerr<<"OpengAlice : Can not find PMD or PMDLoader\n";
103 const char *cDR = strstr(option,"DR");
107 fPMDLoader->LoadDigits("READ");
108 fPMDLoader->LoadRecPoints("recreate");
112 void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
114 // Converts digits to recpoints after running clustering
115 // algorithm on CPV plane and PREshower plane
117 Int_t det = 0,smn = 0;
124 TObjArray *pmdcont = new TObjArray();
125 AliPMDcluster *pmdcl = new AliPMDcluster;
126 AliPMDClustering *pmdclust = new AliPMDClustering();
127 pmdclust->SetDebug(fDebug);
128 pmdclust->SetEdepCut(fEcut);
130 fRunLoader->GetEvent(ievt);
131 //cout << " ***** Beginning::Digits2RecPoints *****" << endl;
132 fTreeD = fPMDLoader->TreeD();
135 cout << " Can not get TreeD" << endl;
137 AliPMDdigit *pmddigit;
138 TBranch *branch = fTreeD->GetBranch("PMDDigit");
139 branch->SetAddress(&fDigits);
142 fTreeR = fPMDLoader->TreeR();
145 fPMDLoader->MakeTree("R");
146 fTreeR = fPMDLoader->TreeR();
149 Int_t bufsize = 16000;
150 fTreeR->Branch("PMDRecpoint", &fRecpoints, bufsize);
152 Int_t nmodules = (Int_t) fTreeD->GetEntries();
154 for (Int_t imodule = 0; imodule < nmodules; imodule++)
157 fTreeD->GetEntry(imodule);
158 Int_t nentries = fDigits->GetLast();
159 for (Int_t ient = 0; ient < nentries+1; ient++)
161 pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
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();
170 fCellADC[xpos][ypos] = (Double_t) adc;
175 pmdclust->DoClust(fCellADC,pmdcont);
177 Int_t nentries1 = pmdcont->GetEntries();
178 cout << " nentries1 = " << nentries1 << endl;
179 for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
181 clusdata[0] = (Float_t) idet;
182 clusdata[1] = (Float_t) isup;
184 pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
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();
192 AddRecPoint(clusdata);
202 fPMDLoader = fRunLoader->GetLoader("PMDLoader");
203 fPMDLoader->WriteRecPoints("OVERWRITE");
205 // delete the pointers
209 // cout << " ***** End::Digits2RecPoints *****" << endl;
212 void AliPMDClusterFinder::SetCellEdepCut(Float_t ecut)
216 void AliPMDClusterFinder::SetDebug(Int_t idebug)
221 void AliPMDClusterFinder::AddRecPoint(Float_t *clusdata)
223 // Add Reconstructed points
225 TClonesArray &lrecpoints = *fRecpoints;
226 AliPMDrecpoint1 *newrecpoint;
227 newrecpoint = new AliPMDrecpoint1(clusdata);
228 new(lrecpoints[fNpoint++]) AliPMDrecpoint1(newrecpoint);
231 void AliPMDClusterFinder::ResetCellADC()
233 // Reset the individual cell ADC value to zero
235 for(Int_t irow = 0; irow < fgkRow; irow++)
237 for(Int_t icol = 0; icol < fgkCol; icol++)
239 fCellADC[irow][icol] = 0.;
244 void AliPMDClusterFinder::ResetRecpoint()
246 // Clear the list of reconstructed points
248 if (fRecpoints) fRecpoints->Clear();
250 void AliPMDClusterFinder::UnLoad(Option_t *option)
252 // Unload all the *.root files
254 const char *cR = strstr(option,"R");
256 fRunLoader->UnloadgAlice();
260 fPMDLoader->UnloadDigits();
261 fPMDLoader->UnloadRecPoints();