1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //====================================================================================================================================================
18 // Class for finding and building the clusters of the ALICE Muon Forward Tracker
20 // Contact author: antonio.uras@cern.ch
22 //====================================================================================================================================================
25 #include "TObjArray.h"
26 #include "TClonesArray.h"
27 #include "AliMFTDigit.h"
28 #include "AliMFTCluster.h"
29 #include "AliMFTSegmentation.h"
32 #include "AliMFTConstants.h"
33 #include "AliMFTClusterFinder.h"
35 const Double_t AliMFTClusterFinder::fCutForAvailableDigits = AliMFTConstants::fCutForAvailableDigits;
36 const Double_t AliMFTClusterFinder::fCutForAttachingDigits = AliMFTConstants::fCutForAttachingDigits;
39 ClassImp(AliMFTClusterFinder)
41 //====================================================================================================================================================
43 AliMFTClusterFinder::AliMFTClusterFinder() :
52 // Default constructor
54 for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) fClustersPerPlane[iPlane] = NULL;
55 fDigitsInCluster = new TClonesArray("AliMFTDigit", fNMaxDigitsPerCluster);
59 //====================================================================================================================================================
61 AliMFTClusterFinder::~AliMFTClusterFinder() {
63 AliDebug(1, "Deleting AliMFTClusterFinder...");
65 for (Int_t iPlane=0; iPlane<fNMaxPlanes; iPlane++) delete fClustersPerPlane[iPlane];
67 AliDebug(1, "... done!");
71 //====================================================================================================================================================
73 void AliMFTClusterFinder::Init(Char_t *nameGeomFile) {
75 fSegmentation = new AliMFTSegmentation(nameGeomFile);
76 fNPlanes = fSegmentation -> GetNPlanes();
80 //====================================================================================================================================================
82 void AliMFTClusterFinder::StartEvent() {
84 // Cleaning up and preparation for the clustering procedure
86 AliDebug(1, "Starting Event...");
88 for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
89 fClustersPerPlane[iPlane]->Clear();
92 AliDebug(1, "... done!");
96 //====================================================================================================================================================
98 void AliMFTClusterFinder::DigitsToClusters(const TObjArray *pDigitList) {
100 // where the clusterization is performed
102 AliInfo("Starting Clusterization for MFT");
103 AliDebug(1, Form("nPlanes = %d", fNPlanes));
106 Bool_t isDigAvailableForNewCluster = kTRUE;
108 for (Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
110 AliDebug(1, Form("Plane %02d", iPlane));
112 TClonesArray *myDigitList = (TClonesArray*) pDigitList->At(iPlane);
114 AliDebug(1, Form("myDigitList->GetEntries() = %d", myDigitList->GetEntries()));
116 Int_t cycleOverDigits = 0;
117 Double_t myCutForAvailableDigits = fCutForAvailableDigits;
119 while (myDigitList->GetEntries()) {
121 for (Int_t iDig=0; iDig<myDigitList->GetEntries(); iDig++) {
123 AliDebug(1, Form("Check %d: Digit %5d of %5d", cycleOverDigits, iDig, myDigitList->GetEntries()));
125 fCurrentDigit = (AliMFTDigit*) myDigitList->At(iDig);
126 isDigAvailableForNewCluster = kTRUE;
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();
136 isDigAvailableForNewCluster = kFALSE;
139 if (fCurrentCluster->GetDistanceFromPixel(fCurrentDigit) < myCutForAvailableDigits) isDigAvailableForNewCluster=kFALSE;
142 if (isDigAvailableForNewCluster) {
143 AliMFTCluster *newCluster = new AliMFTCluster();
144 newCluster->AddPixel(fCurrentDigit);
145 myDigitList->Remove(fCurrentDigit);
146 myDigitList->Compress();
148 new ((*fClustersPerPlane[iPlane])[fClustersPerPlane[iPlane]->GetEntries()]) AliMFTCluster(*newCluster);
151 } // end of cycle over the digits
153 if (cycleOverDigits) myCutForAvailableDigits -= 0.5;
156 } // no more digits to check in current plane!
158 for (Int_t iCluster=0; iCluster<fClustersPerPlane[iPlane]->GetEntries(); iCluster++) {
159 fCurrentCluster = (AliMFTCluster*) fClustersPerPlane[iPlane]->At(iCluster);
160 fCurrentCluster -> TerminateCluster();
163 AliDebug(1, Form("Found %d clusters in plane %02d", fClustersPerPlane[iPlane]->GetEntries(), iPlane));
165 } // end of cycle over the planes
169 //====================================================================================================================================================
171 void AliMFTClusterFinder::MakeClusterBranch(TTree *treeCluster) {
173 // Creating the cluster branches, one for each plane (see AliMFTReconstructor::Reconstruct)
175 AliDebug(1, "Making Cluster Branch");
180 for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
181 AliDebug(1, Form("Setting Branch Plane_%02d for treeCluster",iPlane));
182 if (treeCluster->GetBranch(Form("Plane_%02d",iPlane))) continue;
183 AliDebug(1, Form("Branch Plane_%02d does not exist, creating!",iPlane));
184 treeCluster->Branch(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane]));
190 //====================================================================================================================================================
192 void AliMFTClusterFinder::SetClusterTreeAddress(TTree *treeCluster) {
194 // Addressing the cluster branches, one for each plane (see AliMFTReconstructor::Reconstruct)
196 if (treeCluster && treeCluster->GetBranch("Plane_00")) {
198 for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
199 if (treeCluster->GetBranch(Form("Plane_%02d",iPlane))) {
200 treeCluster->SetBranchAddress(Form("Plane_%02d",iPlane), &(fClustersPerPlane[iPlane]));
202 else AliError(Form("No branch available with name Plane_%02d", iPlane));
208 //====================================================================================================================================================
210 void AliMFTClusterFinder::CreateClusters() {
212 // create cluster list
214 AliDebug(1, Form("Creating clusters list: nPlanes = %d",fNPlanes));
216 if (fClustersPerPlane[0]) return;
218 for(Int_t iPlane=0; iPlane<fNPlanes; iPlane++) {
219 AliDebug(1, Form("plane %02d", iPlane));
220 fClustersPerPlane[iPlane] = new TClonesArray("AliMFTCluster");
225 //====================================================================================================================================================