]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliITSMultReconstructor.cxx
Added script to draw calibrated raw data
[u/mrichter/AliRoot.git] / ITS / AliITSMultReconstructor.cxx
CommitLineData
7ca4655f 1/**************************************************************************
eefb3acc 2 * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
7ca4655f 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/* $Id$ */
17
7284b2b2 18//_________________________________________________________________________
ac903f1b 19//
7284b2b2 20// Implementation of the ITS-SPD trackleter class
ac903f1b 21//
fa9ed8e9 22// It retrieves clusters in the pixels (theta and phi) and finds tracklets.
23// These can be used to extract charged particle multiplicity from the ITS.
ac903f1b 24//
fa9ed8e9 25// A tracklet consists of two ITS clusters, one in the first pixel layer and
26// one in the second. The clusters are associated if the differences in
27// Phi (azimuth) and Theta (polar angle) are within fiducial windows.
28// In case of multiple candidates the candidate with minimum
29// distance is selected.
968e8539 30//
fa9ed8e9 31// Two methods return the number of tracklets and the number of unassociated
7284b2b2 32// clusters (i.e. not used in any tracklet) in the first SPD layer
33// (GetNTracklets and GetNSingleClusters)
34//
35// The cuts on phi and theta depend on the interacting system (p-p or Pb-Pb)
36// and can be set via AliITSRecoParam class
37// (SetPhiWindow and SetThetaWindow)
ac903f1b 38//
7284b2b2 39// Origin: Tiziano Virgili
40//
41// Current support and development:
42// Domenico Elia, Maria Nicassio (INFN Bari)
43// Domenico.Elia@ba.infn.it, Maria.Nicassio@ba.infn.it
44//
45// Most recent updates:
46// - multiple association forbidden (fOnlyOneTrackletPerC2 = kTRUE)
f606f16a 47// - phi definition changed to ALICE convention (0,2*TMath::pi())
48// - cluster coordinates taken with GetGlobalXYZ()
9b373e9a 49// - fGeometry removed
50// - number of fired chips on the two layers
fa9ed8e9 51// - option to cut duplicates in the overlaps
7b116aa1 52// - options and fiducial cuts via AliITSRecoParam
fa9ed8e9 53// - move from DeltaZeta to DeltaTheta cut
54// - update to the new algorithm by Mariella and Jan Fiete
55// - store also DeltaTheta in the ESD
56// - less new and delete calls when creating the needed arrays
7284b2b2 57//_________________________________________________________________________
ac903f1b 58
7ca4655f 59#include <TClonesArray.h>
60#include <TH1F.h>
61#include <TH2F.h>
62#include <TTree.h>
7284b2b2 63#include "TArrayI.h"
ac903f1b 64
7ca4655f 65#include "AliITSMultReconstructor.h"
7b116aa1 66#include "AliITSReconstructor.h"
9b373e9a 67#include "AliITSsegmentationSPD.h"
b51872de 68#include "AliITSRecPoint.h"
ac903f1b 69#include "AliITSgeom.h"
70#include "AliLog.h"
fa9ed8e9 71#include "TGeoGlobalMagField.h"
72#include "AliMagF.h"
ac903f1b 73
74//____________________________________________________________________
0762f3a8 75ClassImp(AliITSMultReconstructor)
ac903f1b 76
3ef75756 77
ac903f1b 78//____________________________________________________________________
7537d03c 79AliITSMultReconstructor::AliITSMultReconstructor():
58e8dc31 80TObject(),
7537d03c 81fClustersLay1(0),
82fClustersLay2(0),
7b116aa1 83fDetectorIndexClustersLay1(0),
84fDetectorIndexClustersLay2(0),
85fOverlapFlagClustersLay1(0),
86fOverlapFlagClustersLay2(0),
7537d03c 87fTracklets(0),
968e8539 88fSClusters(0),
7537d03c 89fNClustersLay1(0),
90fNClustersLay2(0),
91fNTracklets(0),
968e8539 92fNSingleCluster(0),
7537d03c 93fPhiWindow(0),
7284b2b2 94fThetaWindow(0),
fa9ed8e9 95fPhiShift(0),
7b116aa1 96fRemoveClustersFromOverlaps(0),
97fPhiOverlapCut(0),
98fZetaOverlapCut(0),
7537d03c 99fHistOn(0),
100fhClustersDPhiAcc(0),
101fhClustersDThetaAcc(0),
7537d03c 102fhClustersDPhiAll(0),
103fhClustersDThetaAll(0),
7537d03c 104fhDPhiVsDThetaAll(0),
105fhDPhiVsDThetaAcc(0),
7537d03c 106fhetaTracklets(0),
107fhphiTracklets(0),
108fhetaClustersLay1(0),
109fhphiClustersLay1(0){
9b373e9a 110
111 fNFiredChips[0] = 0;
112 fNFiredChips[1] = 0;
113
3ef75756 114 // Method to reconstruct the charged particles multiplicity with the
115 // SPD (tracklets).
ac903f1b 116
ac903f1b 117
118 SetHistOn();
ac903f1b 119
7b116aa1 120 if(AliITSReconstructor::GetRecoParam()) {
7b116aa1 121 SetPhiWindow(AliITSReconstructor::GetRecoParam()->GetTrackleterPhiWindow());
7284b2b2 122 SetThetaWindow(AliITSReconstructor::GetRecoParam()->GetTrackleterThetaWindow());
fa9ed8e9 123 SetPhiShift(AliITSReconstructor::GetRecoParam()->GetTrackleterPhiShift());
7b116aa1 124 SetRemoveClustersFromOverlaps(AliITSReconstructor::GetRecoParam()->GetTrackleterRemoveClustersFromOverlaps());
125 SetPhiOverlapCut(AliITSReconstructor::GetRecoParam()->GetTrackleterPhiOverlapCut());
126 SetZetaOverlapCut(AliITSReconstructor::GetRecoParam()->GetTrackleterZetaOverlapCut());
127 } else {
7b116aa1 128 SetPhiWindow();
7284b2b2 129 SetThetaWindow();
fa9ed8e9 130 SetPhiShift();
7b116aa1 131 SetRemoveClustersFromOverlaps();
132 SetPhiOverlapCut();
133 SetZetaOverlapCut();
134 }
135
136
fa9ed8e9 137 fClustersLay1 = 0;
138 fClustersLay2 = 0;
139 fDetectorIndexClustersLay1 = 0;
140 fDetectorIndexClustersLay2 = 0;
141 fOverlapFlagClustersLay1 = 0;
142 fOverlapFlagClustersLay2 = 0;
143 fTracklets = 0;
144 fSClusters = 0;
ac903f1b 145
146 // definition of histograms
fa9ed8e9 147 Bool_t oldStatus = TH1::AddDirectoryStatus();
148 TH1::AddDirectory(kFALSE);
149
7284b2b2 150 fhClustersDPhiAcc = new TH1F("dphiacc", "dphi", 100,-0.1,0.1);
ddced3c8 151 fhClustersDThetaAcc = new TH1F("dthetaacc","dtheta",100,-0.1,0.1);
ddced3c8 152
7284b2b2 153 fhDPhiVsDThetaAcc = new TH2F("dphiVsDthetaAcc","",100,-0.1,0.1,100,-0.1,0.1);
ac903f1b 154
02a95988 155 fhClustersDPhiAll = new TH1F("dphiall", "dphi", 100,0.0,0.5);
7284b2b2 156 fhClustersDThetaAll = new TH1F("dthetaall","dtheta",100,0.0,0.5);
ddced3c8 157
7284b2b2 158 fhDPhiVsDThetaAll = new TH2F("dphiVsDthetaAll","",100,0.,0.5,100,0.,0.5);
ddced3c8 159
160 fhetaTracklets = new TH1F("etaTracklets", "eta", 100,-2.,2.);
f606f16a 161 fhphiTracklets = new TH1F("phiTracklets", "phi", 100, 0., 2*TMath::Pi());
ddced3c8 162 fhetaClustersLay1 = new TH1F("etaClustersLay1", "etaCl1", 100,-2.,2.);
f606f16a 163 fhphiClustersLay1 = new TH1F("phiClustersLay1", "phiCl1", 100, 0., 2*TMath::Pi());
fa9ed8e9 164
165 TH1::AddDirectory(oldStatus);
ac903f1b 166}
ddced3c8 167
3ef75756 168//______________________________________________________________________
7537d03c 169AliITSMultReconstructor::AliITSMultReconstructor(const AliITSMultReconstructor &mr) : TObject(mr),
7537d03c 170fClustersLay1(mr.fClustersLay1),
171fClustersLay2(mr.fClustersLay2),
7b116aa1 172fDetectorIndexClustersLay1(mr.fDetectorIndexClustersLay1),
173fDetectorIndexClustersLay2(mr.fDetectorIndexClustersLay2),
174fOverlapFlagClustersLay1(mr.fOverlapFlagClustersLay1),
175fOverlapFlagClustersLay2(mr.fOverlapFlagClustersLay2),
7537d03c 176fTracklets(mr.fTracklets),
968e8539 177fSClusters(mr.fSClusters),
7537d03c 178fNClustersLay1(mr.fNClustersLay1),
179fNClustersLay2(mr.fNClustersLay2),
180fNTracklets(mr.fNTracklets),
968e8539 181fNSingleCluster(mr.fNSingleCluster),
7537d03c 182fPhiWindow(mr.fPhiWindow),
7284b2b2 183fThetaWindow(mr.fThetaWindow),
fa9ed8e9 184fPhiShift(mr.fPhiShift),
7b116aa1 185fRemoveClustersFromOverlaps(mr.fRemoveClustersFromOverlaps),
186fPhiOverlapCut(mr.fPhiOverlapCut),
187fZetaOverlapCut(mr.fZetaOverlapCut),
7537d03c 188fHistOn(mr.fHistOn),
189fhClustersDPhiAcc(mr.fhClustersDPhiAcc),
190fhClustersDThetaAcc(mr.fhClustersDThetaAcc),
7537d03c 191fhClustersDPhiAll(mr.fhClustersDPhiAll),
192fhClustersDThetaAll(mr.fhClustersDThetaAll),
7537d03c 193fhDPhiVsDThetaAll(mr.fhDPhiVsDThetaAll),
194fhDPhiVsDThetaAcc(mr.fhDPhiVsDThetaAcc),
7537d03c 195fhetaTracklets(mr.fhetaTracklets),
196fhphiTracklets(mr.fhphiTracklets),
197fhetaClustersLay1(mr.fhetaClustersLay1),
198fhphiClustersLay1(mr.fhphiClustersLay1) {
3ef75756 199 // Copy constructor
7537d03c 200
3ef75756 201}
202
203//______________________________________________________________________
7537d03c 204AliITSMultReconstructor& AliITSMultReconstructor::operator=(const AliITSMultReconstructor& mr){
3ef75756 205 // Assignment operator
7537d03c 206 this->~AliITSMultReconstructor();
207 new(this) AliITSMultReconstructor(mr);
3ef75756 208 return *this;
209}
210
211//______________________________________________________________________
212AliITSMultReconstructor::~AliITSMultReconstructor(){
213 // Destructor
1ba5b31c 214
215 // delete histograms
216 delete fhClustersDPhiAcc;
217 delete fhClustersDThetaAcc;
1ba5b31c 218 delete fhClustersDPhiAll;
219 delete fhClustersDThetaAll;
1ba5b31c 220 delete fhDPhiVsDThetaAll;
221 delete fhDPhiVsDThetaAcc;
1ba5b31c 222 delete fhetaTracklets;
223 delete fhphiTracklets;
224 delete fhetaClustersLay1;
225 delete fhphiClustersLay1;
226
227 // delete arrays
fa9ed8e9 228 for(Int_t i=0; i<fNClustersLay1; i++)
1ba5b31c 229 delete [] fClustersLay1[i];
fa9ed8e9 230
231 for(Int_t i=0; i<fNClustersLay2; i++)
1ba5b31c 232 delete [] fClustersLay2[i];
fa9ed8e9 233
234 for(Int_t i=0; i<fNTracklets; i++)
1ba5b31c 235 delete [] fTracklets[i];
fa9ed8e9 236
237 for(Int_t i=0; i<fNSingleCluster; i++)
968e8539 238 delete [] fSClusters[i];
fa9ed8e9 239
1ba5b31c 240 delete [] fClustersLay1;
241 delete [] fClustersLay2;
7b116aa1 242 delete [] fDetectorIndexClustersLay1;
243 delete [] fDetectorIndexClustersLay2;
244 delete [] fOverlapFlagClustersLay1;
245 delete [] fOverlapFlagClustersLay2;
1ba5b31c 246 delete [] fTracklets;
968e8539 247 delete [] fSClusters;
ddced3c8 248}
ac903f1b 249
250//____________________________________________________________________
7284b2b2 251void AliITSMultReconstructor::Reconstruct(TTree* clusterTree, Float_t* vtx, Float_t* /* vtxRes*/) {
ac903f1b 252 //
253 // - calls LoadClusterArray that finds the position of the clusters
254 // (in global coord)
255 // - convert the cluster coordinates to theta, phi (seen from the
7284b2b2 256 // interaction vertex).
ac903f1b 257 // - makes an array of tracklets
258 //
259 // After this method has been called, the clusters of the two layers
260 // and the tracklets can be retrieved by calling the Get'er methods.
261
ac903f1b 262 // reset counters
263 fNClustersLay1 = 0;
264 fNClustersLay2 = 0;
265 fNTracklets = 0;
7284b2b2 266 fNSingleCluster = 0;
267
ac903f1b 268 // loading the clusters
269 LoadClusterArrays(clusterTree);
3ef75756 270
7284b2b2 271 const Double_t pi = TMath::Pi();
fa9ed8e9 272
273 // dPhi shift is field dependent
274 // get average magnetic field
275 Float_t bz = 0;
276 AliMagF* field = 0;
277 if (TGeoGlobalMagField::Instance())
278 field = dynamic_cast<AliMagF*>(TGeoGlobalMagField::Instance()->GetField());
279 if (!field)
280 {
281 AliError("Could not retrieve magnetic field. Assuming no field. Delta Phi shift will be deactivated in AliITSMultReconstructor.")
282 }
283 else
284 bz = TMath::Abs(field->SolenoidField());
285
286 const Double_t dPhiShift = fPhiShift / 5 * bz;
287 AliDebug(1, Form("Using phi shift of %f", dPhiShift));
288
289 const Double_t dPhiWindow2 = fPhiWindow * fPhiWindow;
290 const Double_t dThetaWindow2 = fThetaWindow * fThetaWindow;
291
7284b2b2 292 Int_t* partners = new Int_t[fNClustersLay2];
293 Float_t* minDists = new Float_t[fNClustersLay2];
294 Int_t* associatedLay1 = new Int_t[fNClustersLay1];
295 TArrayI** blacklist = new TArrayI*[fNClustersLay1];
296
297 for (Int_t i=0; i<fNClustersLay2; i++) {
298 partners[i] = -1;
299 minDists[i] = 2;
300 }
301 for (Int_t i=0; i<fNClustersLay1; i++)
302 associatedLay1[i] = 0;
303 for (Int_t i=0; i<fNClustersLay1; i++)
304 blacklist[i] = 0;
305
ac903f1b 306 // find the tracklets
307 AliDebug(1,"Looking for tracklets... ");
fa9ed8e9 308
ac903f1b 309 //###########################################################
310 // Loop on layer 1 : finding theta, phi and z
311 for (Int_t iC1=0; iC1<fNClustersLay1; iC1++) {
312 Float_t x = fClustersLay1[iC1][0] - vtx[0];
313 Float_t y = fClustersLay1[iC1][1] - vtx[1];
314 Float_t z = fClustersLay1[iC1][2] - vtx[2];
ddced3c8 315
fa9ed8e9 316 Float_t r = TMath::Sqrt(x*x + y*y + z*z);
ac903f1b 317
eefb3acc 318 fClustersLay1[iC1][0] = TMath::ACos(z/r); // Store Theta
319 fClustersLay1[iC1][1] = TMath::Pi() + TMath::ATan2(-y,-x); // Store Phi
fa9ed8e9 320
ddced3c8 321 if (fHistOn) {
322 Float_t eta=fClustersLay1[iC1][0];
323 eta= TMath::Tan(eta/2.);
324 eta=-TMath::Log(eta);
325 fhetaClustersLay1->Fill(eta);
de4c520e 326 fhphiClustersLay1->Fill(fClustersLay1[iC1][1]);
ddced3c8 327 }
96c2c35d 328 }
ac903f1b 329
330 // Loop on layer 2 : finding theta, phi and r
331 for (Int_t iC2=0; iC2<fNClustersLay2; iC2++) {
332 Float_t x = fClustersLay2[iC2][0] - vtx[0];
333 Float_t y = fClustersLay2[iC2][1] - vtx[1];
334 Float_t z = fClustersLay2[iC2][2] - vtx[2];
ddced3c8 335
fa9ed8e9 336 Float_t r = TMath::Sqrt(x*x + y*y + z*z);
ac903f1b 337
eefb3acc 338 fClustersLay2[iC2][0] = TMath::ACos(z/r); // Store Theta
339 fClustersLay2[iC2][1] = TMath::Pi() + TMath::ATan2(-y,-x); // Store Phi
ac903f1b 340 }
341
342 //###########################################################
7284b2b2 343 Int_t found = 1;
344 while (found > 0) {
7284b2b2 345 found = 0;
346
347 // Step1: find all tracklets allowing double assocation
348 // Loop on layer 1
349 for (Int_t iC1=0; iC1<fNClustersLay1; iC1++) {
350
351 // already used or in the overlap ?
352 if (associatedLay1[iC1] != 0 || fOverlapFlagClustersLay1[iC1]) continue;
ac903f1b 353
7284b2b2 354 found++;
355
356 // reset of variables for multiple candidates
357 Int_t iC2WithBestDist = -1; // reset
358 Double_t minDist = 2; // reset
7b116aa1 359
7284b2b2 360 // Loop on layer 2
361 for (Int_t iC2=0; iC2<fNClustersLay2; iC2++) {
362
363 // in the overlap ?
364 if (fOverlapFlagClustersLay2[iC2]) continue;
365
366 if (blacklist[iC1]) {
367 Bool_t blacklisted = kFALSE;
368 for (Int_t i=0; i<blacklist[iC1]->GetSize(); i++) {
369 if (blacklist[iC1]->At(i) == iC2) {
370 blacklisted = kTRUE;
371 break;
372 }
373 }
374 if (blacklisted) continue;
375 }
376
ac903f1b 377 // find the difference in angles
7284b2b2 378 Double_t dTheta = TMath::Abs(fClustersLay2[iC2][0] - fClustersLay1[iC1][0]);
379 Double_t dPhi = TMath::Abs(fClustersLay2[iC2][1] - fClustersLay1[iC1][1]);
02a95988 380 // take into account boundary condition
7284b2b2 381 if (dPhi>pi) dPhi=2.*pi-dPhi;
fa9ed8e9 382
ddced3c8 383 if (fHistOn) {
7284b2b2 384 fhClustersDPhiAll->Fill(dPhi);
ddced3c8 385 fhClustersDThetaAll->Fill(dTheta);
ac903f1b 386 fhDPhiVsDThetaAll->Fill(dTheta, dPhi);
387 }
fa9ed8e9 388
389 dPhi -= dPhiShift;
390
7284b2b2 391 // make "elliptical" cut in Phi and Theta!
fa9ed8e9 392 Float_t d = dPhi*dPhi/dPhiWindow2 + dTheta*dTheta/dThetaWindow2;
7284b2b2 393
394 // look for the minimum distance: the minimum is in iC2WithBestDist
fa9ed8e9 395 if (d<1 && d<minDist) {
7284b2b2 396 minDist=d;
ddced3c8 397 iC2WithBestDist = iC2;
ac903f1b 398 }
7284b2b2 399 } // end of loop over clusters in layer 2
ac903f1b 400
7284b2b2 401 if (minDist<1) { // This means that a cluster in layer 2 was found that matches with iC1
402
403 if (minDists[iC2WithBestDist] > minDist) {
404 Int_t oldPartner = partners[iC2WithBestDist];
405 partners[iC2WithBestDist] = iC1;
406 minDists[iC2WithBestDist] = minDist;
407
408 // mark as assigned
409 associatedLay1[iC1] = 1;
410
411 if (oldPartner != -1) {
fa9ed8e9 412 // redo partner search for cluster in L0 (oldPartner), putting this one (iC2WithBestDist) on its blacklist
7284b2b2 413 if (blacklist[oldPartner] == 0) {
414 blacklist[oldPartner] = new TArrayI(1);
415 } else blacklist[oldPartner]->Set(blacklist[oldPartner]->GetSize()+1);
416
417 blacklist[oldPartner]->AddAt(iC2WithBestDist, blacklist[oldPartner]->GetSize()-1);
418
419 // mark as free
fa9ed8e9 420 associatedLay1[oldPartner] = 0;
7284b2b2 421 }
422 } else {
423 // try again to find a cluster without considering iC2WithBestDist
424 if (blacklist[iC1] == 0) {
425 blacklist[iC1] = new TArrayI(1);
426 }
fa9ed8e9 427 else
428 blacklist[iC1]->Set(blacklist[iC1]->GetSize()+1);
7284b2b2 429
430 blacklist[iC1]->AddAt(iC2WithBestDist, blacklist[iC1]->GetSize()-1);
431 }
de4c520e 432
7284b2b2 433 } else // cluster has no partner; remove
434 associatedLay1[iC1] = 2;
435 } // end of loop over clusters in layer 1
436 }
437
438 // Step2: store tracklets; remove used clusters
439 for (Int_t iC2=0; iC2<fNClustersLay2; iC2++) {
de4c520e 440
7284b2b2 441 if (partners[iC2] == -1) continue;
7b116aa1 442
7284b2b2 443 if (fOverlapFlagClustersLay1[partners[iC2]] || fOverlapFlagClustersLay2[iC2]) continue;
444 if (fRemoveClustersFromOverlaps) FlagClustersInOverlapRegions (partners[iC2],iC2);
445
fa9ed8e9 446 fTracklets[fNTracklets] = new Float_t[6];
447
7284b2b2 448 // use the theta from the clusters in the first layer
449 fTracklets[fNTracklets][0] = fClustersLay1[partners[iC2]][0];
450 // use the phi from the clusters in the first layer
451 fTracklets[fNTracklets][1] = fClustersLay1[partners[iC2]][1];
452 // store the difference between phi1 and phi2
453 fTracklets[fNTracklets][2] = fClustersLay1[partners[iC2]][1] - fClustersLay2[iC2][1];
454
455 // define dphi in the range [0,pi] with proper sign (track charge correlated)
456 if (fTracklets[fNTracklets][2] > TMath::Pi())
457 fTracklets[fNTracklets][2] = fTracklets[fNTracklets][2]-2.*TMath::Pi();
458 if (fTracklets[fNTracklets][2] < -TMath::Pi())
459 fTracklets[fNTracklets][2] = fTracklets[fNTracklets][2]+2.*TMath::Pi();
7b116aa1 460
7284b2b2 461 // store the difference between theta1 and theta2
462 fTracklets[fNTracklets][3] = fClustersLay1[partners[iC2]][0] - fClustersLay2[iC2][0];
463
464 if (fHistOn) {
465 fhClustersDPhiAcc->Fill(fTracklets[fNTracklets][2]);
466 fhClustersDThetaAcc->Fill(fTracklets[fNTracklets][3]);
467 fhDPhiVsDThetaAcc->Fill(fTracklets[fNTracklets][3],fTracklets[fNTracklets][2]);
ac903f1b 468 }
3ef75756 469
7284b2b2 470 // find label
471 // if equal label in both clusters found this label is assigned
472 // if no equal label can be found the first labels of the L1 AND L2 cluster are assigned
473 Int_t label1 = 0;
474 Int_t label2 = 0;
475 while (label2 < 3) {
476 if ((Int_t) fClustersLay1[partners[iC2]][3+label1] != -2 && (Int_t) fClustersLay1[partners[iC2]][3+label1] == (Int_t) fClustersLay2[iC2][3+label2])
477 break;
478 label1++;
479 if (label1 == 3) {
480 label1 = 0;
481 label2++;
482 }
483 }
484 if (label2 < 3) {
485 AliDebug(AliLog::kDebug, Form("Found label %d == %d for tracklet candidate %d\n", (Int_t) fClustersLay1[partners[iC2]][3+label1], (Int_t) fClustersLay2[iC2][3+label2], fNTracklets));
486 fTracklets[fNTracklets][4] = fClustersLay1[partners[iC2]][3+label1];
487 fTracklets[fNTracklets][5] = fClustersLay2[iC2][3+label2];
488 } else {
489 AliDebug(AliLog::kDebug, Form("Did not find label %d %d %d %d %d %d for tracklet candidate %d\n", (Int_t) fClustersLay1[partners[iC2]][3], (Int_t) fClustersLay1[partners[iC2]][4], (Int_t) fClustersLay1[partners[iC2]][5], (Int_t) fClustersLay2[iC2][3], (Int_t) fClustersLay2[iC2][4], (Int_t) fClustersLay2[iC2][5], fNTracklets));
490 fTracklets[fNTracklets][4] = fClustersLay1[partners[iC2]][3];
491 fTracklets[fNTracklets][5] = fClustersLay2[iC2][3];
492 }
493
494 if (fHistOn) {
495 Float_t eta=fTracklets[fNTracklets][0];
496 eta= TMath::Tan(eta/2.);
497 eta=-TMath::Log(eta);
498 fhetaTracklets->Fill(eta);
499 fhphiTracklets->Fill(fTracklets[fNTracklets][1]);
500 }
3ef75756 501
7284b2b2 502 AliDebug(1,Form(" Adding tracklet candidate %d ", fNTracklets));
503 AliDebug(1,Form(" Cl. %d of Layer 1 and %d of Layer 2", partners[iC2], iC2));
504 fNTracklets++;
3ef75756 505
7284b2b2 506 associatedLay1[partners[iC2]] = 1;
507 }
508
509 // Delete the following else if you do not want to save Clusters!
510 // store the cluster
511 for (Int_t iC1=0; iC1<fNClustersLay1; iC1++) {
512 if (associatedLay1[iC1]==2||associatedLay1[iC1]==0) {
fa9ed8e9 513 fSClusters[fNSingleCluster] = new Float_t[2];
968e8539 514 fSClusters[fNSingleCluster][0] = fClustersLay1[iC1][0];
515 fSClusters[fNSingleCluster][1] = fClustersLay1[iC1][1];
de4c520e 516 AliDebug(1,Form(" Adding a single cluster %d (cluster %d of layer 1)",
7284b2b2 517 fNSingleCluster, iC1));
968e8539 518 fNSingleCluster++;
3ef75756 519 }
7284b2b2 520 }
521
522 delete[] partners;
523 delete[] minDists;
524
525 for (Int_t i=0; i<fNClustersLay1; i++)
526 if (blacklist[i])
527 delete blacklist[i];
528 delete[] blacklist;
529
ac903f1b 530 AliDebug(1,Form("%d tracklets found", fNTracklets));
531}
532
533//____________________________________________________________________
534void
535AliITSMultReconstructor::LoadClusterArrays(TTree* itsClusterTree) {
536 // This method
537 // - gets the clusters from the cluster tree
538 // - convert them into global coordinates
539 // - store them in the internal arrays
9b373e9a 540 // - count the number of cluster-fired chips
ac903f1b 541
9b373e9a 542 AliDebug(1,"Loading clusters and cluster-fired chips ...");
ac903f1b 543
544 fNClustersLay1 = 0;
545 fNClustersLay2 = 0;
9b373e9a 546 fNFiredChips[0] = 0;
547 fNFiredChips[1] = 0;
ac903f1b 548
9b373e9a 549 AliITSsegmentationSPD *seg = new AliITSsegmentationSPD();
550
b51872de 551 TClonesArray* itsClusters = new TClonesArray("AliITSRecPoint");
552 TBranch* itsClusterBranch=itsClusterTree->GetBranch("ITSRecPoints");
ddced3c8 553
ac903f1b 554 itsClusterBranch->SetAddress(&itsClusters);
ddced3c8 555
ac903f1b 556 Int_t nItsSubs = (Int_t)itsClusterTree->GetEntries();
f606f16a 557 Float_t cluGlo[3]={0.,0.,0.};
ddced3c8 558
fa9ed8e9 559
560 // count clusters
561 // loop over the its subdetectors
562 for (Int_t iIts=0; iIts < nItsSubs; iIts++) {
563 if (!itsClusterTree->GetEvent(iIts))
564 continue;
565
566 Int_t nClusters = itsClusters->GetEntriesFast();
567 // loop over clusters
568 while(nClusters--) {
569 AliITSRecPoint* cluster = (AliITSRecPoint*)itsClusters->UncheckedAt(nClusters);
570
571 Int_t layer = cluster->GetLayer();
572 if (layer == 0)
573 fNClustersLay1++;
574 else if (layer == 1)
575 fNClustersLay2++;
576 }
577 }
578
579 // create arrays
580 fClustersLay1 = new Float_t*[fNClustersLay1];
581 fDetectorIndexClustersLay1 = new Int_t[fNClustersLay1];
582 fOverlapFlagClustersLay1 = new Bool_t[fNClustersLay1];
583
584 fClustersLay2 = new Float_t*[fNClustersLay2];
585 fDetectorIndexClustersLay2 = new Int_t[fNClustersLay2];
586 fOverlapFlagClustersLay2 = new Bool_t[fNClustersLay2];
587
588 // no double association allowed
589 fTracklets = new Float_t*[TMath::Min(fNClustersLay1, fNClustersLay2)];
590 fSClusters = new Float_t*[fNClustersLay1];
591
592 for (Int_t i=0; i<fNClustersLay1; i++) {
593 fClustersLay1[i] = new Float_t[6];
594 fOverlapFlagClustersLay1[i] = kFALSE;
595 fSClusters[i] = 0;
596 }
597
598 for (Int_t i=0; i<fNClustersLay2; i++) {
599 fClustersLay2[i] = new Float_t[6];
600 fOverlapFlagClustersLay2[i] = kFALSE;
601 }
602
603 for (Int_t i=0; i<TMath::Min(fNClustersLay1, fNClustersLay2); i++)
604 fTracklets[i] = 0;
605
606 // fill clusters
ac903f1b 607 // loop over the its subdetectors
fa9ed8e9 608 fNClustersLay1 = 0; // reset to 0
609 fNClustersLay2 = 0;
ac903f1b 610 for (Int_t iIts=0; iIts < nItsSubs; iIts++) {
611
612 if (!itsClusterTree->GetEvent(iIts))
613 continue;
614
615 Int_t nClusters = itsClusters->GetEntriesFast();
9b373e9a 616
617 // number of clusters in each chip of the current module
618 Int_t nClustersInChip[5] = {0,0,0,0,0};
619 Int_t layer = 0;
ac903f1b 620
ac903f1b 621 // loop over clusters
622 while(nClusters--) {
de4c520e 623 AliITSRecPoint* cluster = (AliITSRecPoint*)itsClusters->UncheckedAt(nClusters);
ac903f1b 624
9b373e9a 625 layer = cluster->GetLayer();
626 if (layer>1) continue;
ac903f1b 627
f606f16a 628 cluster->GetGlobalXYZ(cluGlo);
629 Float_t x = cluGlo[0];
630 Float_t y = cluGlo[1];
631 Float_t z = cluGlo[2];
9b373e9a 632
633 // find the chip for the current cluster
634 Float_t locz = cluster->GetDetLocalZ();
635 Int_t iChip = seg->GetChipFromLocal(0,locz);
636 nClustersInChip[iChip]++;
ac903f1b 637
9b373e9a 638 if (layer==0) {
ac903f1b 639 fClustersLay1[fNClustersLay1][0] = x;
640 fClustersLay1[fNClustersLay1][1] = y;
641 fClustersLay1[fNClustersLay1][2] = z;
7b116aa1 642
643 fDetectorIndexClustersLay1[fNClustersLay1]=iIts;
644
de4c520e 645 for (Int_t i=0; i<3; i++)
646 fClustersLay1[fNClustersLay1][3+i] = cluster->GetLabel(i);
ac903f1b 647 fNClustersLay1++;
648 }
9b373e9a 649 if (layer==1) {
ac903f1b 650 fClustersLay2[fNClustersLay2][0] = x;
651 fClustersLay2[fNClustersLay2][1] = y;
652 fClustersLay2[fNClustersLay2][2] = z;
7b116aa1 653
654 fDetectorIndexClustersLay2[fNClustersLay2]=iIts;
655
de4c520e 656 for (Int_t i=0; i<3; i++)
657 fClustersLay2[fNClustersLay2][3+i] = cluster->GetLabel(i);
ac903f1b 658 fNClustersLay2++;
659 }
660
661 }// end of cluster loop
9b373e9a 662
663 // get number of fired chips in the current module
664 if(layer<2)
665 for(Int_t ifChip=0; ifChip<5; ifChip++) {
666 if(nClustersInChip[ifChip] >= 1) fNFiredChips[layer]++;
667 }
668
ac903f1b 669 } // end of its "subdetector" loop
9b373e9a 670
cedf398d 671 if (itsClusters) {
672 itsClusters->Delete();
673 delete itsClusters;
9b373e9a 674 delete seg;
cedf398d 675 itsClusters = 0;
676 }
ac903f1b 677 AliDebug(1,Form("(clusters in layer 1 : %d, layer 2: %d)",fNClustersLay1,fNClustersLay2));
9b373e9a 678 AliDebug(1,Form("(cluster-fired chips in layer 1 : %d, layer 2: %d)",fNFiredChips[0],fNFiredChips[1]));
679}
680//____________________________________________________________________
681void
682AliITSMultReconstructor::LoadClusterFiredChips(TTree* itsClusterTree) {
683 // This method
684 // - gets the clusters from the cluster tree
685 // - counts the number of (cluster)fired chips
686
687 AliDebug(1,"Loading cluster-fired chips ...");
688
689 fNFiredChips[0] = 0;
690 fNFiredChips[1] = 0;
691
692 AliITSsegmentationSPD *seg = new AliITSsegmentationSPD();
693
694 TClonesArray* itsClusters = new TClonesArray("AliITSRecPoint");
695 TBranch* itsClusterBranch=itsClusterTree->GetBranch("ITSRecPoints");
696
697 itsClusterBranch->SetAddress(&itsClusters);
698
699 Int_t nItsSubs = (Int_t)itsClusterTree->GetEntries();
700
701 // loop over the its subdetectors
702 for (Int_t iIts=0; iIts < nItsSubs; iIts++) {
703
704 if (!itsClusterTree->GetEvent(iIts))
705 continue;
706
707 Int_t nClusters = itsClusters->GetEntriesFast();
708
709 // number of clusters in each chip of the current module
710 Int_t nClustersInChip[5] = {0,0,0,0,0};
711 Int_t layer = 0;
712
713 // loop over clusters
714 while(nClusters--) {
715 AliITSRecPoint* cluster = (AliITSRecPoint*)itsClusters->UncheckedAt(nClusters);
716
717 layer = cluster->GetLayer();
718 if (layer>1) continue;
719
720 // find the chip for the current cluster
721 Float_t locz = cluster->GetDetLocalZ();
722 Int_t iChip = seg->GetChipFromLocal(0,locz);
723 nClustersInChip[iChip]++;
724
725 }// end of cluster loop
726
727 // get number of fired chips in the current module
728 if(layer<2)
729 for(Int_t ifChip=0; ifChip<5; ifChip++) {
730 if(nClustersInChip[ifChip] >= 1) fNFiredChips[layer]++;
731 }
732
733 } // end of its "subdetector" loop
734
735 if (itsClusters) {
736 itsClusters->Delete();
737 delete itsClusters;
738 delete seg;
739 itsClusters = 0;
740 }
741 AliDebug(1,Form("(cluster-fired chips in layer 1 : %d, layer 2: %d)",fNFiredChips[0],fNFiredChips[1]));
ac903f1b 742}
743//____________________________________________________________________
744void
745AliITSMultReconstructor::SaveHists() {
3ef75756 746 // This method save the histograms on the output file
747 // (only if fHistOn is TRUE).
ac903f1b 748
749 if (!fHistOn)
750 return;
751
ddced3c8 752 fhClustersDPhiAll->Write();
753 fhClustersDThetaAll->Write();
ac903f1b 754 fhDPhiVsDThetaAll->Write();
ddced3c8 755
756 fhClustersDPhiAcc->Write();
757 fhClustersDThetaAcc->Write();
ac903f1b 758 fhDPhiVsDThetaAcc->Write();
ddced3c8 759
760 fhetaTracklets->Write();
761 fhphiTracklets->Write();
762 fhetaClustersLay1->Write();
763 fhphiClustersLay1->Write();
ac903f1b 764}
7b116aa1 765
766//____________________________________________________________________
767void
768AliITSMultReconstructor::FlagClustersInOverlapRegions (Int_t iC1, Int_t iC2WithBestDist) {
769
770 Float_t distClSameMod=0.;
771 Float_t distClSameModMin=0.;
772 Int_t iClOverlap =0;
773 Float_t meanRadiusLay1 = 3.99335; // average radius inner layer
774 Float_t meanRadiusLay2 = 7.37935; // average radius outer layer;
775
776 Float_t zproj1=0.;
777 Float_t zproj2=0.;
778 Float_t deZproj=0.;
779
780 // Loop on inner layer clusters
781 for (Int_t iiC1=0; iiC1<fNClustersLay1; iiC1++) {
782 if (!fOverlapFlagClustersLay1[iiC1]) {
783 // only for adjacent modules
784 if ((TMath::Abs(fDetectorIndexClustersLay1[iC1]-fDetectorIndexClustersLay1[iiC1])==4)||
785 (TMath::Abs(fDetectorIndexClustersLay1[iC1]-fDetectorIndexClustersLay1[iiC1])==76)) {
786 Float_t dePhi=TMath::Abs(fClustersLay1[iiC1][1]-fClustersLay1[iC1][1]);
787 if (dePhi>TMath::Pi()) dePhi=2.*TMath::Pi()-dePhi;
788
789 zproj1=meanRadiusLay1/TMath::Tan(fClustersLay1[iC1][0]);
790 zproj2=meanRadiusLay1/TMath::Tan(fClustersLay1[iiC1][0]);
791
792 deZproj=TMath::Abs(zproj1-zproj2);
793
794 distClSameMod = TMath::Sqrt(TMath::Power(deZproj/fZetaOverlapCut,2)+TMath::Power(dePhi/fPhiOverlapCut,2));
795 if (distClSameMod<=1.) fOverlapFlagClustersLay1[iiC1]=kTRUE;
796
797// if (distClSameMod<=1.) {
798// if (distClSameModMin==0. || distClSameMod<distClSameModMin) {
799// distClSameModMin=distClSameMod;
800// iClOverlap=iiC1;
801// }
802// }
803
804
805 } // end adjacent modules
806 }
807 } // end Loop on inner layer clusters
808
809// if (distClSameModMin!=0.) fOverlapFlagClustersLay1[iClOverlap]=kTRUE;
810
811 distClSameMod=0.;
812 distClSameModMin=0.;
813 iClOverlap =0;
814 // Loop on outer layer clusters
815 for (Int_t iiC2=0; iiC2<fNClustersLay2; iiC2++) {
816 if (!fOverlapFlagClustersLay2[iiC2]) {
817 // only for adjacent modules
818 if ((TMath::Abs(fDetectorIndexClustersLay2[iC2WithBestDist]-fDetectorIndexClustersLay2[iiC2])==4) ||
819 (TMath::Abs(fDetectorIndexClustersLay2[iC2WithBestDist]-fDetectorIndexClustersLay2[iiC2])==156)) {
820 Float_t dePhi=TMath::Abs(fClustersLay2[iiC2][1]-fClustersLay2[iC2WithBestDist][1]);
821 if (dePhi>TMath::Pi()) dePhi=2.*TMath::Pi()-dePhi;
822
823 zproj1=meanRadiusLay2/TMath::Tan(fClustersLay2[iC2WithBestDist][0]);
824 zproj2=meanRadiusLay2/TMath::Tan(fClustersLay2[iiC2][0]);
825
826 deZproj=TMath::Abs(zproj1-zproj2);
827 distClSameMod = TMath::Sqrt(TMath::Power(deZproj/fZetaOverlapCut,2)+TMath::Power(dePhi/fPhiOverlapCut,2));
828 if (distClSameMod<=1.) fOverlapFlagClustersLay2[iiC2]=kTRUE;
829
830// if (distClSameMod<=1.) {
831// if (distClSameModMin==0. || distClSameMod<distClSameModMin) {
832// distClSameModMin=distClSameMod;
833// iClOverlap=iiC2;
834// }
835// }
836
837 } // end adjacent modules
838 }
839 } // end Loop on outer layer clusters
840
841// if (distClSameModMin!=0.) fOverlapFlagClustersLay2[iClOverlap]=kTRUE;
842
843}