]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MFT/AliMFTClusterFinder.cxx
Modified init of AliMuonForwardTrack
[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(const 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]->Delete();
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   TClonesArray *myDigitList = 0;
109
110   for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
111
112     AliDebug(1, Form("Plane %02d", iPlane));
113
114     myDigitList = (TClonesArray*) pDigitList->At(iPlane);
115
116     AliDebug(1, Form("myDigitList->GetEntries() = %d", myDigitList->GetEntries()));
117
118     Int_t cycleOverDigits = 0;
119     Double_t myCutForAvailableDigits = fCutForAvailableDigits;
120     
121     while (myDigitList->GetEntries()) {
122
123       for (Int_t iDig=0; iDig<myDigitList->GetEntries(); iDig++) {
124
125         AliDebug(1, Form("Check %d: Digit %5d of %5d", cycleOverDigits, iDig, myDigitList->GetEntries()));
126
127         fCurrentDigit = (AliMFTDigit*) myDigitList->At(iDig);
128         isDigAvailableForNewCluster = kTRUE;
129
130         for (Int_t iCluster=0; iCluster<fClustersPerPlane[iPlane]->GetEntries(); iCluster++) {
131           fCurrentCluster = (AliMFTCluster*) fClustersPerPlane[iPlane]->At(iCluster);
132           AliDebug(2, Form("Distance between cluster and digit = %f",fCurrentCluster->GetDistanceFromPixel(fCurrentDigit))); 
133           if (fCurrentCluster->GetDistanceFromPixel(fCurrentDigit) < fCutForAttachingDigits) { 
134             fCurrentCluster->AddPixel(fCurrentDigit); 
135             myDigitList->Remove(fCurrentDigit);
136             myDigitList->Compress();
137             iDig--;
138             isDigAvailableForNewCluster = kFALSE;
139             break; 
140           }
141           if (fCurrentCluster->GetDistanceFromPixel(fCurrentDigit) < myCutForAvailableDigits) isDigAvailableForNewCluster=kFALSE;
142         }
143
144         if (isDigAvailableForNewCluster) {
145           AliMFTCluster *newCluster = new AliMFTCluster();
146           newCluster->AddPixel(fCurrentDigit);
147           myDigitList->Remove(fCurrentDigit);
148           myDigitList->Compress();
149           iDig--;
150           new ((*fClustersPerPlane[iPlane])[fClustersPerPlane[iPlane]->GetEntries()]) AliMFTCluster(*newCluster);
151         }
152
153       }   // end of cycle over the digits
154
155       if (cycleOverDigits) myCutForAvailableDigits -= 0.5;
156       cycleOverDigits++;
157
158     }   // no more digits to check in current plane!
159
160     for (Int_t iCluster=0; iCluster<fClustersPerPlane[iPlane]->GetEntries(); iCluster++) {
161       fCurrentCluster = (AliMFTCluster*) fClustersPerPlane[iPlane]->At(iCluster);
162       fCurrentCluster -> TerminateCluster();
163     }
164
165     AliDebug(1, Form("Found %d clusters in plane %02d", fClustersPerPlane[iPlane]->GetEntries(), iPlane));
166
167     myDigitList -> Delete();
168
169   }  // end of cycle over the planes
170
171 }
172
173 //====================================================================================================================================================
174
175 void AliMFTClusterFinder::MakeClusterBranch(TTree *treeCluster) {
176
177   // Creating the cluster branches, one for each plane (see AliMFTReconstructor::Reconstruct)
178
179   AliDebug(1, "Making Cluster Branch");
180
181   CreateClusters();
182
183   if (treeCluster) {
184     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
185       AliDebug(1, Form("Setting Branch Plane_%02d for treeCluster",iPlane));
186       if (treeCluster->GetBranch(Form("Plane_%02d",iPlane))) continue;
187       AliDebug(1, Form("Branch Plane_%02d does not exist, creating!",iPlane));
188       treeCluster->Branch(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane])); 
189     }
190   }
191
192 }
193
194 //====================================================================================================================================================
195
196 void AliMFTClusterFinder::SetClusterTreeAddress(TTree *treeCluster) {
197
198   // Addressing the cluster branches, one for each plane (see AliMFTReconstructor::Reconstruct)
199
200   if (treeCluster && treeCluster->GetBranch("Plane_00")) {
201     CreateClusters();
202     for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
203       if (treeCluster->GetBranch(Form("Plane_%02d",iPlane))) {
204         treeCluster->SetBranchAddress(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane]));
205       }
206       else AliError(Form("No branch available with name Plane_%02d", iPlane));
207     }
208   }
209
210 }
211
212 //====================================================================================================================================================
213
214 void AliMFTClusterFinder::CreateClusters() {
215
216   // create cluster list
217
218   AliDebug(1, Form("Creating clusters list: nPlanes = %d",fNPlanes));
219
220   if (fClustersPerPlane[0]) return;
221   
222   for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
223     AliDebug(1, Form("plane %02d", iPlane));
224     fClustersPerPlane[iPlane] = new TClonesArray("AliMFTCluster");
225   }
226
227 }
228
229 //====================================================================================================================================================