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(AliRunLoader* runLoader):
54 fRunLoader(runLoader),
55 fPMDLoader(runLoader->GetLoader("PMDLoader")),
58 fDigits(new TClonesArray("AliPMDdigit", 1000)),
59 fRecpoints(new TClonesArray("AliPMDrecpoint1", 1000)),
68 // ------------------------------------------------------------------------- //
69 AliPMDClusterFinder::~AliPMDClusterFinder()
85 // ------------------------------------------------------------------------- //
87 void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
89 // Converts digits to recpoints after running clustering
90 // algorithm on CPV plane and PREshower plane
92 Int_t det = 0,smn = 0;
99 TObjArray *pmdcont = new TObjArray();
100 AliPMDClustering *pmdclust = new AliPMDClustering();
101 pmdclust->SetDebug(fDebug);
102 pmdclust->SetEdepCut(fEcut);
104 fRunLoader->GetEvent(ievt);
105 //cout << " ***** Beginning::Digits2RecPoints *****" << endl;
106 fTreeD = fPMDLoader->TreeD();
109 cout << " Can not get TreeD" << endl;
111 AliPMDdigit *pmddigit;
112 TBranch *branch = fTreeD->GetBranch("PMDDigit");
113 branch->SetAddress(&fDigits);
116 fTreeR = fPMDLoader->TreeR();
119 fPMDLoader->MakeTree("R");
120 fTreeR = fPMDLoader->TreeR();
123 Int_t bufsize = 16000;
124 fTreeR->Branch("PMDRecpoint", &fRecpoints, bufsize);
126 Int_t nmodules = (Int_t) fTreeD->GetEntries();
128 for (Int_t imodule = 0; imodule < nmodules; imodule++)
131 fTreeD->GetEntry(imodule);
132 Int_t nentries = fDigits->GetLast();
133 for (Int_t ient = 0; ient < nentries+1; ient++)
135 pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
137 det = pmddigit->GetDetector();
138 smn = pmddigit->GetSMNumber();
139 xpos = pmddigit->GetRow();
140 ypos = pmddigit->GetColumn();
141 adc = pmddigit->GetADC();
142 //Int_t trno = pmddigit->GetTrackNumber();
144 fCellADC[xpos][ypos] = (Double_t) adc;
149 pmdclust->DoClust(idet,ismn,fCellADC,pmdcont);
151 Int_t nentries1 = pmdcont->GetEntries();
152 // cout << " nentries1 = " << nentries1 << endl;
153 for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
155 AliPMDcluster *pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
156 idet = pmdcl->GetDetector();
157 ismn = pmdcl->GetSMN();
158 clusdata[0] = pmdcl->GetClusX();
159 clusdata[1] = pmdcl->GetClusY();
160 clusdata[2] = pmdcl->GetClusADC();
161 clusdata[3] = pmdcl->GetClusCells();
162 clusdata[4] = pmdcl->GetClusRadius();
164 AddRecPoint(idet,ismn,clusdata);
174 fPMDLoader = fRunLoader->GetLoader("PMDLoader");
175 fPMDLoader->WriteRecPoints("OVERWRITE");
177 // delete the pointers
181 // cout << " ***** End::Digits2RecPoints *****" << endl;
183 // ------------------------------------------------------------------------- //
184 void AliPMDClusterFinder::SetCellEdepCut(Float_t ecut)
188 // ------------------------------------------------------------------------- //
189 void AliPMDClusterFinder::SetDebug(Int_t idebug)
193 // ------------------------------------------------------------------------- //
194 void AliPMDClusterFinder::AddRecPoint(Int_t idet,Int_t ismn,Float_t *clusdata)
196 // Add Reconstructed points
198 TClonesArray &lrecpoints = *fRecpoints;
199 AliPMDrecpoint1 *newrecpoint;
200 newrecpoint = new AliPMDrecpoint1(idet, ismn, clusdata);
201 new(lrecpoints[fNpoint++]) AliPMDrecpoint1(newrecpoint);
204 // ------------------------------------------------------------------------- //
205 void AliPMDClusterFinder::ResetCellADC()
207 // Reset the individual cell ADC value to zero
209 for(Int_t irow = 0; irow < fgkRow; irow++)
211 for(Int_t icol = 0; icol < fgkCol; icol++)
213 fCellADC[irow][icol] = 0.;
217 // ------------------------------------------------------------------------- //
219 void AliPMDClusterFinder::ResetRecpoint()
221 // Clear the list of reconstructed points
223 if (fRecpoints) fRecpoints->Clear();
225 // ------------------------------------------------------------------------- //
226 void AliPMDClusterFinder::Load()
228 // Load all the *.root files
230 fPMDLoader->LoadDigits("READ");
231 fPMDLoader->LoadRecPoints("recreate");
233 // ------------------------------------------------------------------------- //
234 void AliPMDClusterFinder::UnLoad()
236 // Unload all the *.root files
238 fPMDLoader->UnloadDigits();
239 fPMDLoader->UnloadRecPoints();
241 // ------------------------------------------------------------------------- //