]>
Commit | Line | Data |
---|---|---|
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 "AliMFTConstants.h" | |
33 | #include "AliMFTClusterFinder.h" | |
34 | #include "TRandom.h" | |
35 | ||
36 | const Double_t AliMFTClusterFinder::fCutForAvailableDigits = AliMFTConstants::fCutForAvailableDigits; | |
37 | const Double_t AliMFTClusterFinder::fCutForAttachingDigits = AliMFTConstants::fCutForAttachingDigits; | |
38 | const Double_t AliMFTClusterFinder::fMisalignmentMagnitude = AliMFTConstants::fMisalignmentMagnitude; | |
39 | ||
40 | ClassImp(AliMFTClusterFinder) | |
41 | ||
42 | //==================================================================================================================================================== | |
43 | ||
44 | AliMFTClusterFinder::AliMFTClusterFinder() : | |
45 | TObject(), | |
46 | fDigitsInCluster(0), | |
47 | fCurrentDigit(0), | |
48 | fCurrentCluster(0), | |
49 | fSegmentation(0), | |
50 | fNPlanes(0), | |
51 | fApplyMisalignment(kFALSE), | |
52 | fStopWatch(0) | |
53 | { | |
54 | ||
55 | // Default constructor | |
56 | ||
57 | for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) fClustersPerPlane[iPlane] = NULL; | |
58 | fDigitsInCluster = new TClonesArray("AliMFTDigit", fNMaxDigitsPerCluster); | |
59 | fDigitsInCluster -> SetOwner(kTRUE); | |
60 | fStopWatch = new TStopwatch(); | |
61 | fStopWatch -> Reset(); | |
62 | ||
63 | } | |
64 | ||
65 | //==================================================================================================================================================== | |
66 | ||
67 | AliMFTClusterFinder::~AliMFTClusterFinder() { | |
68 | ||
69 | AliDebug(1, "Deleting AliMFTClusterFinder..."); | |
70 | ||
71 | for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) { | |
72 | if (fClustersPerPlane[iPlane]) fClustersPerPlane[iPlane]->Delete(); delete fClustersPerPlane[iPlane]; fClustersPerPlane[iPlane] = 0x0; | |
73 | } | |
74 | if (fDigitsInCluster) fDigitsInCluster->Delete(); delete fDigitsInCluster; fDigitsInCluster = NULL; | |
75 | delete fSegmentation; fSegmentation=NULL; | |
76 | ||
77 | delete fStopWatch; | |
78 | ||
79 | AliDebug(1, "... done!"); | |
80 | ||
81 | } | |
82 | ||
83 | //==================================================================================================================================================== | |
84 | ||
85 | void AliMFTClusterFinder::Clear(const Option_t* /*opt*/) { | |
86 | ||
87 | AliDebug(1, "Clearing AliMFTClusterFinder..."); | |
88 | ||
89 | for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) { | |
90 | if(fClustersPerPlane[iPlane]) fClustersPerPlane[iPlane]->Clear("C"); | |
91 | } | |
92 | if(fDigitsInCluster) fDigitsInCluster->Clear("C"); | |
93 | fSegmentation->Clear("C"); | |
94 | ||
95 | AliDebug(1, "... done!"); | |
96 | ||
97 | } | |
98 | ||
99 | //==================================================================================================================================================== | |
100 | ||
101 | void AliMFTClusterFinder::Init(const Char_t *nameGeomFile) { | |
102 | ||
103 | fSegmentation = new AliMFTSegmentation(nameGeomFile); | |
104 | fNPlanes = fSegmentation -> GetNPlanes(); | |
105 | ||
106 | } | |
107 | ||
108 | //==================================================================================================================================================== | |
109 | ||
110 | void AliMFTClusterFinder::StartEvent() { | |
111 | ||
112 | // Cleaning up and preparation for the clustering procedure | |
113 | ||
114 | AliDebug(1, "Starting Event..."); | |
115 | ||
116 | for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) { | |
117 | fClustersPerPlane[iPlane]->Delete(); | |
118 | } | |
119 | ||
120 | AliDebug(1, "... done!"); | |
121 | ||
122 | } | |
123 | ||
124 | //==================================================================================================================================================== | |
125 | ||
126 | void AliMFTClusterFinder::DigitsToClusters(const TObjArray *pDigitList) { | |
127 | ||
128 | // where the clusterization is performed | |
129 | ||
130 | AliInfo("Starting Clusterization for MFT"); | |
131 | AliDebug(1, Form("nPlanes = %d", fNPlanes)); | |
132 | ||
133 | if (!fStopWatch) fStopWatch = new TStopwatch(); | |
134 | fStopWatch -> Reset(); | |
135 | ||
136 | StartEvent(); | |
137 | Bool_t isDigAvailableForNewCluster = kTRUE; | |
138 | ||
139 | TClonesArray *myDigitList = 0; | |
140 | ||
141 | for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) { | |
142 | ||
143 | AliDebug(1, Form("Plane %02d", iPlane)); | |
144 | ||
145 | Int_t nDetElem = fSegmentation->GetPlane(iPlane)->GetNActiveElements(); | |
146 | TClonesArray *clustersPerDetElem[fNMaxDetElemPerPlane] = {0}; | |
147 | for (Int_t iDetElem=0; iDetElem<nDetElem; iDetElem++) clustersPerDetElem[iDetElem] = new TClonesArray("AliMFTCluster"); | |
148 | ||
149 | myDigitList = (TClonesArray*) pDigitList->At(iPlane); | |
150 | ||
151 | AliDebug(1, Form("myDigitList->GetEntries() = %d", myDigitList->GetEntries())); | |
152 | ||
153 | Int_t cycleOverDigits = 0; | |
154 | Double_t myCutForAvailableDigits = 0; | |
155 | ||
156 | Int_t currentDetElem = -1; | |
157 | Int_t currentDetElemLocal = -1; | |
158 | Bool_t areThereSkippedDigits = kFALSE; | |
159 | ||
160 | fStopWatch -> Start(); | |
161 | ||
162 | while (myDigitList->GetEntries()) { | |
163 | ||
164 | for (Int_t iDig=0; iDig<myDigitList->GetEntries(); iDig++) { | |
165 | ||
166 | fCurrentDigit = (AliMFTDigit*) myDigitList->At(iDig); | |
167 | ||
168 | if (!iDig) { | |
169 | if (fCurrentDigit->GetDetElemID() != currentDetElem) { | |
170 | // first iteration over the digits of a specific detection element | |
171 | currentDetElem = fCurrentDigit->GetDetElemID(); | |
172 | currentDetElemLocal = fSegmentation->GetDetElemLocalID(currentDetElem); | |
173 | cycleOverDigits = 0; | |
174 | myCutForAvailableDigits = fCutForAvailableDigits; | |
175 | } | |
176 | else if (fCurrentDigit->GetDetElemID()==currentDetElem && areThereSkippedDigits) { | |
177 | // second (or further) iteration over the digits of a specific detection element | |
178 | cycleOverDigits++; | |
179 | myCutForAvailableDigits -= 0.5; | |
180 | } | |
181 | areThereSkippedDigits = kFALSE; | |
182 | } | |
183 | else { | |
184 | areThereSkippedDigits = kTRUE; | |
185 | if (fCurrentDigit->GetDetElemID() != currentDetElem) break; | |
186 | } | |
187 | ||
188 | isDigAvailableForNewCluster = kTRUE; | |
189 | ||
190 | for (Int_t iCluster=0; iCluster<clustersPerDetElem[currentDetElemLocal]->GetEntries(); iCluster++) { | |
191 | fCurrentCluster = (AliMFTCluster*) clustersPerDetElem[currentDetElemLocal]->At(iCluster); | |
192 | if (fCurrentCluster->GetDistanceFromPixel(fCurrentDigit) < fCutForAttachingDigits) { | |
193 | fCurrentCluster->AddPixel(fCurrentDigit); | |
194 | myDigitList->Remove(fCurrentDigit); | |
195 | myDigitList->Compress(); | |
196 | iDig--; | |
197 | isDigAvailableForNewCluster = kFALSE; | |
198 | break; | |
199 | } | |
200 | if (fCurrentCluster->GetDistanceFromPixel(fCurrentDigit) < myCutForAvailableDigits) isDigAvailableForNewCluster=kFALSE; | |
201 | } | |
202 | ||
203 | if (isDigAvailableForNewCluster) { | |
204 | AliMFTCluster *newCluster = new AliMFTCluster(); | |
205 | newCluster->AddPixel(fCurrentDigit); | |
206 | myDigitList->Remove(fCurrentDigit); | |
207 | myDigitList->Compress(); | |
208 | iDig--; | |
209 | new ((*clustersPerDetElem[currentDetElemLocal])[clustersPerDetElem[currentDetElemLocal]->GetEntries()]) AliMFTCluster(*newCluster); | |
210 | delete newCluster; | |
211 | } | |
212 | ||
213 | } // end of cycle over the digits | |
214 | ||
215 | } // no more digits to check in current plane! | |
216 | ||
217 | fStopWatch -> Print("m"); | |
218 | ||
219 | printf("Plane %d: clusters found in %f seconds\n",iPlane,fStopWatch->CpuTime()); | |
220 | fStopWatch->Start(); | |
221 | ||
222 | // Now we merge the cluster lists coming from each detection element, to build the cluster list of the plane | |
223 | ||
224 | Double_t misalignmentPhi = 2.*TMath::Pi()*gRandom->Rndm(); | |
225 | Double_t misalignmentX = fMisalignmentMagnitude*TMath::Cos(misalignmentPhi); | |
226 | Double_t misalignmentY = fMisalignmentMagnitude*TMath::Sin(misalignmentPhi); | |
227 | ||
228 | AliMFTCluster *newCluster = NULL; | |
229 | for (Int_t iDetElem=0; iDetElem<nDetElem; iDetElem++) { | |
230 | for (Int_t iCluster=0; iCluster<clustersPerDetElem[iDetElem]->GetEntries(); iCluster++) { | |
231 | newCluster = (AliMFTCluster*) (clustersPerDetElem[iDetElem]->At(iCluster)); | |
232 | newCluster -> TerminateCluster(); | |
233 | newCluster -> SetClusterEditable(kTRUE); | |
234 | if (TMath::Abs(newCluster->GetZ())<TMath::Abs(fSegmentation->GetPlane(iPlane)->GetZCenter())) newCluster->SetClusterFront(kTRUE); | |
235 | else newCluster->SetClusterFront(kFALSE); | |
236 | if (fApplyMisalignment) { | |
237 | newCluster -> SetX(newCluster->GetX()+misalignmentX); | |
238 | newCluster -> SetY(newCluster->GetY()+misalignmentY); | |
239 | } | |
240 | newCluster -> SetClusterEditable(kFALSE); | |
241 | ||
242 | new ((*fClustersPerPlane[iPlane])[fClustersPerPlane[iPlane]->GetEntries()]) AliMFTCluster(*newCluster); | |
243 | } | |
244 | } | |
245 | ||
246 | printf("%d Clusters in plane %02d merged in %f seconds\n", fClustersPerPlane[iPlane]->GetEntries(), iPlane, fStopWatch->CpuTime()); | |
247 | ||
248 | for (Int_t iDetElem=0; iDetElem<nDetElem; iDetElem++) { | |
249 | clustersPerDetElem[iDetElem] -> Delete(); | |
250 | delete clustersPerDetElem[iDetElem]; | |
251 | } | |
252 | ||
253 | myDigitList -> Delete(); | |
254 | ||
255 | } // end of cycle over the planes | |
256 | ||
257 | } | |
258 | ||
259 | //==================================================================================================================================================== | |
260 | ||
261 | void AliMFTClusterFinder::MakeClusterBranch(TTree *treeCluster) { | |
262 | ||
263 | // Creating the cluster branches, one for each plane (see AliMFTReconstructor::Reconstruct) | |
264 | ||
265 | AliDebug(1, "Making Cluster Branch"); | |
266 | ||
267 | CreateClusters(); | |
268 | ||
269 | if (treeCluster) { | |
270 | for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) { | |
271 | AliDebug(1, Form("Setting Branch Plane_%02d for treeCluster",iPlane)); | |
272 | if (treeCluster->GetBranch(Form("Plane_%02d",iPlane))) continue; | |
273 | AliDebug(1, Form("Branch Plane_%02d does not exist, creating!",iPlane)); | |
274 | treeCluster->Branch(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane])); | |
275 | } | |
276 | } | |
277 | ||
278 | } | |
279 | ||
280 | //==================================================================================================================================================== | |
281 | ||
282 | void AliMFTClusterFinder::SetClusterTreeAddress(TTree *treeCluster) { | |
283 | ||
284 | // Addressing the cluster branches, one for each plane (see AliMFTReconstructor::Reconstruct) | |
285 | ||
286 | if (treeCluster && treeCluster->GetBranch("Plane_00")) { | |
287 | CreateClusters(); | |
288 | for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) { | |
289 | if (treeCluster->GetBranch(Form("Plane_%02d",iPlane))) { | |
290 | treeCluster->SetBranchAddress(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane])); | |
291 | } | |
292 | else AliError(Form("No branch available with name Plane_%02d", iPlane)); | |
293 | } | |
294 | } | |
295 | ||
296 | } | |
297 | ||
298 | //==================================================================================================================================================== | |
299 | ||
300 | void AliMFTClusterFinder::CreateClusters() { | |
301 | ||
302 | // create cluster list | |
303 | ||
304 | AliDebug(1, Form("Creating clusters list: nPlanes = %d",fNPlanes)); | |
305 | ||
306 | if (fClustersPerPlane[0]) return; | |
307 | ||
308 | for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) { | |
309 | AliDebug(1, Form("plane %02d", iPlane)); | |
310 | fClustersPerPlane[iPlane] = new TClonesArray("AliMFTCluster"); | |
311 | fClustersPerPlane[iPlane] -> SetOwner(kTRUE); | |
312 | ||
313 | } | |
314 | ||
315 | } | |
316 | ||
317 | //==================================================================================================================================================== |