]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MFT/AliMFTClusterFinder.cxx
Updated version of the MFT code (Antonio)
[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 {
51
52   // Default constructor
53   
54   for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) fClustersPerPlane[iPlane] = NULL;
55   fDigitsInCluster = new TClonesArray("AliMFTDigit", fNMaxDigitsPerCluster);
56
57 }
58
59 //====================================================================================================================================================
60
61 AliMFTClusterFinder::~AliMFTClusterFinder() {
62
63   AliDebug(1, "Deleting AliMFTClusterFinder...");
64
65   for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) delete fClustersPerPlane[iPlane];
66
67   AliDebug(1, "... done!");
68
69 }
70
71 //====================================================================================================================================================
72
73 void AliMFTClusterFinder::Init(Char_t *nameGeomFile) {
74
75   fSegmentation = new AliMFTSegmentation(nameGeomFile);
76   fNPlanes = fSegmentation -> GetNPlanes();
77   
78 }
79
80 //====================================================================================================================================================
81
82 void AliMFTClusterFinder::StartEvent() {
83
84   // Cleaning up and preparation for the clustering procedure
85
86   AliDebug(1, "Starting Event...");
87   
88   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
89     fClustersPerPlane[iPlane]->Clear();
90   }
91
92   AliDebug(1, "... done!");
93
94 }
95
96 //====================================================================================================================================================
97
98 void AliMFTClusterFinder::DigitsToClusters(const TObjArray *pDigitList) {
99
100   // where the clusterization is performed
101
102   AliInfo("Starting Clusterization for MFT");
103   AliDebug(1, Form("nPlanes = %d", fNPlanes));
104
105   StartEvent(); 
106   Bool_t isDigAvailableForNewCluster = kTRUE;
107   
108   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
109
110     AliDebug(1, Form("Plane %02d", iPlane));
111
112     TClonesArray *myDigitList = (TClonesArray*) pDigitList->At(iPlane);
113
114     AliDebug(1, Form("myDigitList->GetEntries() = %d", myDigitList->GetEntries()));
115
116     Int_t cycleOverDigits = 0;
117     Double_t myCutForAvailableDigits = fCutForAvailableDigits;
118     
119     while (myDigitList->GetEntries()) {
120
121       for (Int_t iDig=0; iDig<myDigitList->GetEntries(); iDig++) {
122
123         AliDebug(1, Form("Check %d: Digit %5d of %5d", cycleOverDigits, iDig, myDigitList->GetEntries()));
124
125         fCurrentDigit = (AliMFTDigit*) myDigitList->At(iDig);
126         isDigAvailableForNewCluster = kTRUE;
127
128         for (Int_t iCluster=0; iCluster<fClustersPerPlane[iPlane]->GetEntries(); iCluster++) {
129           fCurrentCluster = (AliMFTCluster*) fClustersPerPlane[iPlane]->At(iCluster);
130           AliDebug(2, Form("Distance between cluster and digit = %f",fCurrentCluster->GetDistanceFromPixel(fCurrentDigit))); 
131           if (fCurrentCluster->GetDistanceFromPixel(fCurrentDigit) < fCutForAttachingDigits) { 
132             fCurrentCluster->AddPixel(fCurrentDigit); 
133             myDigitList->Remove(fCurrentDigit);
134             myDigitList->Compress();
135             iDig--;
136             isDigAvailableForNewCluster = kFALSE;
137             break; 
138           }
139           if (fCurrentCluster->GetDistanceFromPixel(fCurrentDigit) < myCutForAvailableDigits) isDigAvailableForNewCluster=kFALSE;
140         }
141
142         if (isDigAvailableForNewCluster) {
143           AliMFTCluster *newCluster = new AliMFTCluster();
144           newCluster->AddPixel(fCurrentDigit);
145           myDigitList->Remove(fCurrentDigit);
146           myDigitList->Compress();
147           iDig--;
148           new ((*fClustersPerPlane[iPlane])[fClustersPerPlane[iPlane]->GetEntries()]) AliMFTCluster(*newCluster);
149         }
150
151       }   // end of cycle over the digits
152
153       if (cycleOverDigits) myCutForAvailableDigits -= 0.5;
154       cycleOverDigits++;
155
156     }   // no more digits to check in current plane!
157
158     for (Int_t iCluster=0; iCluster<fClustersPerPlane[iPlane]->GetEntries(); iCluster++) {
159       fCurrentCluster = (AliMFTCluster*) fClustersPerPlane[iPlane]->At(iCluster);
160       fCurrentCluster -> TerminateCluster();
161     }
162
163     AliDebug(1, Form("Found %d clusters in plane %02d", fClustersPerPlane[iPlane]->GetEntries(), iPlane));
164
165   }  // end of cycle over the planes
166
167 }
168
169 //====================================================================================================================================================
170
171 void AliMFTClusterFinder::MakeClusterBranch(TTree *treeCluster) {
172
173   // Creating the cluster branches, one for each plane (see AliMFTReconstructor::Reconstruct)
174
175   AliDebug(1, "Making Cluster Branch");
176
177   CreateClusters();
178
179   if (treeCluster) {
180     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
181       AliDebug(1, Form("Setting Branch Plane_%02d for treeCluster",iPlane));
182       TBranch *branch = treeCluster->GetBranch(Form("Plane_%02d",iPlane));
183       if (branch) continue;
184       AliDebug(1, Form("Branch Plane_%02d does not exist, creating!",iPlane));
185       branch = treeCluster->Branch(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane])); 
186     }
187   }
188
189 }
190
191 //====================================================================================================================================================
192
193 void AliMFTClusterFinder::SetClusterTreeAddress(TTree *treeCluster) {
194
195   // Addressing the cluster branches, one for each plane (see AliMFTReconstructor::Reconstruct)
196
197   if (treeCluster && treeCluster->GetBranch("Plane_00")) {
198     CreateClusters();
199     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
200       if (treeCluster->GetBranch(Form("Plane_%02d",iPlane))) {
201         treeCluster->SetBranchAddress(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane]));
202       }
203       else AliError(Form("No branch available with name Plane_%02d", iPlane));
204     }
205   }
206
207 }
208
209 //====================================================================================================================================================
210
211 void AliMFTClusterFinder::CreateClusters() {
212
213   // create cluster list
214
215   AliDebug(1, Form("Creating clusters list: nPlanes = %d",fNPlanes));
216
217   if (fClustersPerPlane[0]) return;
218   
219   for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
220     AliDebug(1, Form("plane %02d", iPlane));
221     fClustersPerPlane[iPlane] = new TClonesArray("AliMFTCluster");
222   }
223
224 }
225
226 //====================================================================================================================================================