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