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