]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MFT/AliMFTClusterFinder.cxx
Compilation warnings
[u/mrichter/AliRoot.git] / MFT / AliMFTClusterFinder.cxx
CommitLineData
820b4d9e 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// Class for finding and building the clusters of the ALICE Muon Forward Tracker
19//
20// Contact author: antonio.uras@cern.ch
21//
22//====================================================================================================================================================
23
24#include "AliLog.h"
25#include "TObjArray.h"
26#include "TClonesArray.h"
27#include "AliMFTDigit.h"
28#include "AliMFTCluster.h"
29#include "AliMFTSegmentation.h"
30#include "TTree.h"
31#include "TMath.h"
32#include "AliMFTClusterFinder.h"
33
34ClassImp(AliMFTClusterFinder)
35
36//====================================================================================================================================================
37
38AliMFTClusterFinder::AliMFTClusterFinder() :
39 TObject(),
40 fDigitsInCluster(0),
41 fCurrentDig(0),
42 fSegmentation(0),
43 fNPlanes(0)
44{
45
46 // Default constructor
47
48 for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) fClustersPerPlane[iPlane] = NULL;
49 fDigitsInCluster = new TClonesArray("AliMFTDigit", fNMaxDigitsPerCluster);
50
51}
52
53//====================================================================================================================================================
54
55AliMFTClusterFinder::~AliMFTClusterFinder() {
56
57 AliDebug(1, "Deleting AliMFTClusterFinder...");
58
59 for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) delete fClustersPerPlane[iPlane];
60
61 AliDebug(1, "... done!");
62
63}
64
65//====================================================================================================================================================
66
58af8e48 67void AliMFTClusterFinder::Init(const Char_t *nameGeomFile) {
820b4d9e 68
69 fSegmentation = new AliMFTSegmentation(nameGeomFile);
70 fNPlanes = fSegmentation -> GetNPlanes();
71
72}
73
74//====================================================================================================================================================
75
76void AliMFTClusterFinder::StartEvent() {
77
78 // Cleaning up and preparation for the clustering procedure
79
80 AliDebug(1, "Starting Event...");
81
82 for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
83 fClustersPerPlane[iPlane]->Clear();
84 }
85
86 AliDebug(1, "... done!");
87
88}
89
90//====================================================================================================================================================
91
92void AliMFTClusterFinder::DigitsToClusters(const TObjArray *pDigitList) {
93
94 // where the clusterization is performed
95
96 AliInfo("Starting Clusterization for MFT");
97 AliDebug(1, Form("nPlanes = %d", fNPlanes));
98
99 StartEvent();
100
101 for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
102
103 AliDebug(1, Form("Plane %02d", iPlane));
104
105 TClonesArray *myDigitList = (TClonesArray*) pDigitList->At(iPlane);
106 // myDigitList->Sort();
107
108 AliDebug(1, Form("myDigitList->GetEntries() = %d", myDigitList->GetEntries()));
109
110 while (myDigitList->GetEntries()) {
111
112 fDigitsInCluster->Clear();
113
114 Bool_t clusterUpdated=kTRUE;
115
116 //---------------------------------------------------------------------------------------------------------
117
118 while (clusterUpdated) { // repeat the loop on the digits until no new compatible digit is found
119
120 clusterUpdated = kFALSE;
121
122 for (Int_t iDig=0; iDig<myDigitList->GetEntries(); iDig++) {
123 fCurrentDig = (AliMFTDigit*) myDigitList->At(iDig);
124 if (fDigitsInCluster->GetEntries()<fNMaxDigitsPerCluster) {
125 if (fDigitsInCluster->GetEntries()==0) {
126 new ((*fDigitsInCluster)[fDigitsInCluster->GetEntries()]) AliMFTDigit(*fCurrentDig);
127 myDigitList->Remove(fCurrentDig);
128 clusterUpdated = kTRUE;
129 }
130 else if (IsCurrentDigitCompatible()) {
131 new ((*fDigitsInCluster)[fDigitsInCluster->GetEntries()]) AliMFTDigit(*fCurrentDig);
132 myDigitList->Remove(fCurrentDig);
133 clusterUpdated = kTRUE;
134 }
135 }
136 }
137
138 if (clusterUpdated) myDigitList->Compress();
139
140 }
141
142 //---------------------------------------------------------------------------------------------------------
143
144 AliDebug(1, "Building new cluster");
145
146 BuildNewCluster(iPlane); // here the new cluster is built
147
148 }
149
150 }
151
152}
153
154//====================================================================================================================================================
155
156Bool_t AliMFTClusterFinder::IsCurrentDigitCompatible() {
157
158 // where it is decided if the current digit (fCurrentDig) is compatible with the current digits array (fDigitsInCluster)
159
160 for (Int_t iDigit=0; iDigit<fDigitsInCluster->GetEntries(); iDigit++) {
161 AliMFTDigit *tmpDig = (AliMFTDigit*) fDigitsInCluster->At(iDigit);
162 Int_t distX = TMath::Abs(tmpDig->GetPixelX() - fCurrentDig->GetPixelX());
163 Int_t distY = TMath::Abs(tmpDig->GetPixelY() - fCurrentDig->GetPixelY());
164 if (distX<=1 && distY<=1) return kTRUE;
165 }
166
167 return kFALSE;
168
169}
170
171//====================================================================================================================================================
172
173void AliMFTClusterFinder::BuildNewCluster(Int_t plane) {
174
175 // where a new cluster is built, starting from the array of compatible digits (fDigitsInCluster)
176
177 AliDebug(1, Form("Starting cluster building from %d digits", fDigitsInCluster->GetEntries()));
178
179 AliMFTCluster *newCluster = new AliMFTCluster();
180
181 Double_t xCenters[fNMaxDigitsPerCluster] = {0};
182 Double_t yCenters[fNMaxDigitsPerCluster] = {0};
183 Double_t nElectrons = 0.;
184
185 for (Int_t iDigit=0; iDigit<fDigitsInCluster->GetEntries(); iDigit++) {
186 AliMFTDigit *tmpDig = (AliMFTDigit*) fDigitsInCluster->At(iDigit);
187 xCenters[iDigit] = tmpDig->GetPixelCenterX();
188 yCenters[iDigit] = tmpDig->GetPixelCenterY();
189 nElectrons += tmpDig->GetNElectrons();
190 for (Int_t iTrack=0; iTrack<tmpDig->GetNMCTracks(); iTrack++) newCluster->AddMCLabel(tmpDig->GetMCLabel(iTrack));
191 }
192
193 newCluster -> SetX(TMath::Mean(fDigitsInCluster->GetEntries(), xCenters));
194 newCluster -> SetY(TMath::Mean(fDigitsInCluster->GetEntries(), yCenters));
195 newCluster -> SetZ(((AliMFTDigit*) fDigitsInCluster->At(0))->GetPixelCenterZ());
196
197 Double_t minErrX = ((AliMFTDigit*) fDigitsInCluster->At(0))->GetPixelWidthX() / TMath::Sqrt(12.);
198 Double_t minErrY = ((AliMFTDigit*) fDigitsInCluster->At(0))->GetPixelWidthY() / TMath::Sqrt(12.);
199 Double_t minErrZ = ((AliMFTDigit*) fDigitsInCluster->At(0))->GetPixelWidthZ() / TMath::Sqrt(12.);
200 newCluster -> SetErrX( TMath::Max(TMath::RMS(fDigitsInCluster->GetEntries(), xCenters), minErrX) );
201 newCluster -> SetErrY( TMath::Max(TMath::RMS(fDigitsInCluster->GetEntries(), yCenters), minErrY) );
202 newCluster -> SetErrZ( minErrZ );
203
204 newCluster -> SetNElectrons(nElectrons);
205 newCluster -> SetPlane(plane);
206
207 newCluster -> SetSize(fDigitsInCluster->GetEntries());
208
209 AliDebug(1, Form("Adding cluster #%02d to tree (%f, %f, %f)",
210 fClustersPerPlane[plane]->GetEntries(), newCluster->GetX(), newCluster->GetY(), newCluster->GetZ()));
211
212 new ((*fClustersPerPlane[plane])[fClustersPerPlane[plane]->GetEntries()]) AliMFTCluster(*newCluster);
213
214}
215
216//====================================================================================================================================================
217
218void AliMFTClusterFinder::MakeClusterBranch(TTree *treeCluster) {
219
220 // Creating the cluster branches, one for each plane (see AliMFTReconstructor::Reconstruct)
221
222 AliDebug(1, "Making Cluster Branch");
223
224 CreateClusters();
225
226 if (treeCluster) {
227 for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
228 AliDebug(1, Form("Setting Branch Plane_%02d for treeCluster",iPlane));
229 TBranch *branch = treeCluster->GetBranch(Form("Plane_%02d",iPlane));
230 if (branch) continue;
231 AliDebug(1, Form("Branch Plane_%02d does not exist, creating!",iPlane));
232 branch = treeCluster->Branch(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane]));
233 }
234 }
235
236}
237
238//====================================================================================================================================================
239
240void AliMFTClusterFinder::SetClusterTreeAddress(TTree *treeCluster) {
241
242 // Addressing the cluster branches, one for each plane (see AliMFTReconstructor::Reconstruct)
243
244 if (treeCluster && treeCluster->GetBranch("Plane_00")) {
245 CreateClusters();
246 for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
247 if (treeCluster->GetBranch(Form("Plane_%02d",iPlane))) {
248 treeCluster->SetBranchAddress(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane]));
249 }
250 else AliError(Form("No branch available with name Plane_%02d", iPlane));
251 }
252 }
253
254}
255
256//====================================================================================================================================================
257
258void AliMFTClusterFinder::CreateClusters() {
259
260 // create cluster list
261
262 AliDebug(1, Form("Creating clusters list: nPlanes = %d",fNPlanes));
263
264 if (fClustersPerPlane[0]) return;
265
266 for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
267 AliDebug(1, Form("plane %02d", iPlane));
268 fClustersPerPlane[iPlane] = new TClonesArray("AliMFTCluster");
269 }
270
271}
272
273//====================================================================================================================================================