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 //====================================================================================================================================================
20 // Class for the creation of the "global muon tracks" built from the tracks reconstructed in the
21 // muon spectrometer and the clusters of the Muon Forward Tracker
23 // Contact author: antonio.uras@cern.ch
25 //====================================================================================================================================================
29 #include "AliGeomManager.h"
30 #include "AliESDEvent.h"
31 #include "AliESDMuonTrack.h"
32 #include "AliESDMuonGlobalTrack.h"
33 #include "AliMFTTrackerMU.h"
36 #include "AliRunLoader.h"
37 #include "AliHeader.h"
38 #include "AliGenEventHeader.h"
41 #include "AliMUONTrackExtrap.h"
42 #include "AliMUONTrack.h"
43 #include "AliMUONESDInterface.h"
44 #include "AliMuonForwardTrack.h"
45 #include "AliMUONConstants.h"
48 ClassImp(AliMFTTrackerMU)
50 const Double_t AliMFTTrackerMU::fRadLengthSi = AliMFTConstants::fRadLengthSi;
52 //====================================================================================================================================================
54 AliMFTTrackerMU::AliMFTTrackerMU() :
60 fNPlanesMFTAnalyzed(0),
62 fScaleSigmaClusterCut(1.),
63 fNMaxMissingMFTClusters(0),
64 fGlobalTrackingDiverged(kFALSE),
68 fFinalBestCandidate(0),
72 fXExtrapVertexError(0),
73 fYExtrapVertexError(0),
74 fBransonCorrection(kFALSE)
77 //--------------------------------------------------------------------
78 // This is the AliMFTTrackerMU constructor
79 //--------------------------------------------------------------------
81 fMFT = (AliMFT*) gAlice->GetDetector("MFT");
82 fSegmentation = fMFT->GetSegmentation();
83 SetNPlanesMFT(fSegmentation->GetNPlanes());
84 AliMUONTrackExtrap::SetField(); // set the magnetic field for track extrapolations
86 AliInfo(Form("fMFT = %p, fSegmentation = %p", fMFT, fSegmentation));
88 for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
89 fMFTClusterArray[iPlane] = new TClonesArray("AliMFTCluster");
90 fMFTClusterArrayFront[iPlane] = new TClonesArray("AliMFTCluster");
91 fMFTClusterArrayBack[iPlane] = new TClonesArray("AliMFTCluster");
92 fMFTClusterArray[iPlane] -> SetOwner(kTRUE);
93 fMFTClusterArrayFront[iPlane] -> SetOwner(kTRUE);
94 fMFTClusterArrayBack[iPlane] -> SetOwner(kTRUE);
95 fMinResearchRadiusAtPlane[iPlane] = 0.;
98 fCandidateTracks = new TClonesArray("AliMuonForwardTrack",50000);
102 //====================================================================================================================================================
104 AliMFTTrackerMU::~AliMFTTrackerMU() {
108 for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
109 fMFTClusterArray[iPlane] -> Delete();
110 delete fMFTClusterArray[iPlane];
111 delete fMFTClusterArrayFront[iPlane];
112 delete fMFTClusterArrayBack[iPlane];
115 delete fCandidateTracks;
119 //====================================================================================================================================================
121 Int_t AliMFTTrackerMU::LoadClusters(TTree *cTree) {
123 //--------------------------------------------------------------------
124 // This function loads the MFT clusters
125 //--------------------------------------------------------------------
127 for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
128 AliDebug(1, Form("Setting Address for Branch Plane_%02d", iPlane));
129 cTree->SetBranchAddress(Form("Plane_%02d",iPlane), &fMFTClusterArray[iPlane]);
132 if (!cTree->GetEvent()) return kFALSE;
133 for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
134 AliInfo(Form("plane %02d: nClusters = %d\n", iPlane, fMFTClusterArray[iPlane]->GetEntries()));
137 AddClustersFromUnderlyingEvent();
138 AddClustersFromPileUpEvents();
140 SeparateFrontBackClusters();
146 //====================================================================================================================================================
148 void AliMFTTrackerMU::UnloadClusters() {
150 //--------------------------------------------------------------------
151 // This function unloads MFT clusters
152 //--------------------------------------------------------------------
154 for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
155 fMFTClusterArray[iPlane] -> Clear("C");
156 fMFTClusterArrayFront[iPlane] -> Clear("C");
157 fMFTClusterArrayBack[iPlane] -> Clear("C");
162 //====================================================================================================================================================
164 Int_t AliMFTTrackerMU::Clusters2Tracks(AliESDEvent *event) {
166 //--------------------------------------------------------------------
167 // This functions reconstructs the Muon Forward Tracks
168 // The clusters must be already loaded !
169 //--------------------------------------------------------------------
171 // Tree of AliMuonForwardTrack objects. Created outside the ESD framework for cross-check purposes
173 TFile *outputFileMuonGlobalTracks = new TFile("MuonGlobalTracks.root", "update");
174 TTree *outputTreeMuonGlobalTracks = new TTree("AliMuonForwardTracks", "Tree of AliMuonForwardTracks");
175 TClonesArray *muonForwardTracks = new TClonesArray("AliMuonForwardTrack");
176 outputTreeMuonGlobalTracks -> Branch("tracks", &muonForwardTracks);
178 //--------------------------------------------------------------------
185 fXExtrapVertexError = AliMFTConstants::fXVertexTolerance;
186 fYExtrapVertexError = AliMFTConstants::fYVertexTolerance;
188 // Imposing a fixed primary vertex, in order to keep memory of its position. Taking the primary vertex would imply the risk
189 // of loosing memory of its position when passing from ESD to AOD, due to possible refitting
190 const AliESDVertex* esdVert = fESD->GetVertex();
191 if (esdVert->GetNContributors() > 0 || !strcmp(esdVert->GetTitle(),"vertexer: smearMC")) {
192 fXExtrapVertex = esdVert->GetX();
193 fYExtrapVertex = esdVert->GetY();
194 fZExtrapVertex = esdVert->GetZ();
195 fXExtrapVertexError = TMath::Max(AliMFTConstants::fXVertexTolerance, esdVert->GetXRes());
196 fYExtrapVertexError = TMath::Max(AliMFTConstants::fYVertexTolerance, esdVert->GetYRes());
197 AliInfo(Form("Found ESD vertex from %d contributors (%f +/- %f, %f +/- %f, %f)",
198 esdVert->GetNContributors(),fXExtrapVertex,fXExtrapVertexError,fYExtrapVertex,fYExtrapVertexError,fZExtrapVertex));
200 else GetVertexFromMC();
202 //----------- Read ESD MUON tracks -------------------
204 Int_t nTracksMUON = event->GetNumberOfMuonTracks();
206 AliInfo(Form("Number of ESD MUON tracks: %d\n", nTracksMUON));
209 while (iTrack<nTracksMUON) {
211 fNPlanesMFTAnalyzed = 0;
213 AliInfo("****************************************************************************************");
214 AliInfo(Form("*************************** MUON TRACK %3d/%d ***************************************", iTrack, nTracksMUON));
215 AliInfo("****************************************************************************************");
217 fCandidateTracks -> Delete();
219 fNPlanesMFTAnalyzed = 0;
221 const AliESDMuonTrack *esdTrack = event->GetMuonTrack(iTrack);
222 if (fMUONTrack) delete fMUONTrack;
223 fMUONTrack = new AliMUONTrack();
225 AliMUONESDInterface::ESDToMUON(*esdTrack, *fMUONTrack, kFALSE);
227 if (!fMUONTrack->GetTrackParamAtCluster()->First()) {
228 AliInfo("Skipping track, no parameters available!!!");
233 // the track we are going to build, starting from fMUONTrack and adding the MFT clusters
234 AliMuonForwardTrack *track = new ((*fCandidateTracks)[0]) AliMuonForwardTrack();
235 track -> SetMUONTrack(new AliMUONTrack(*fMUONTrack));
236 track -> SetMCLabel(fMUONTrack->GetMCLabel());
237 track -> SetMatchTrigger(fMUONTrack->GetMatchTrigger());
239 // track parameters linearly extrapolated from the first tracking station to the end of the absorber
240 AliMUONTrackParam trackParamEndOfAbsorber(*((AliMUONTrackParam*)(fMUONTrack->GetTrackParamAtCluster()->First())));
241 AliMUONTrackExtrap::ExtrapToZCov(&trackParamEndOfAbsorber, AliMUONConstants::AbsZEnd()); // absorber extends from -90 to -503 cm
242 Double_t xEndOfAbsorber = trackParamEndOfAbsorber.GetNonBendingCoor();
243 Double_t yEndOfAbsorber = trackParamEndOfAbsorber.GetBendingCoor();
244 Double_t rAbsorber = TMath::Sqrt(xEndOfAbsorber*xEndOfAbsorber + yEndOfAbsorber*yEndOfAbsorber);
245 track -> SetRAtAbsorberEnd(rAbsorber);
247 //------------------------- NOW THE CYCLE OVER THE MFT PLANES STARTS ---------------------------------------
249 for (Int_t iPlane=fNPlanesMFT-1; iPlane>=0; iPlane--) { /* *** do not reverse the order of this cycle!!!
250 *** this reflects the fact that the extrapolation is performed
251 *** starting from the last MFT plane back to the origin */
253 // --------- updating the array of candidates according to the clusters available in the i-th plane ---------
255 fNPlanesMFTAnalyzed++;
257 Int_t nCandidates = fCandidateTracks->GetEntriesFast();
258 for (Int_t iCandidate=0; iCandidate<nCandidates; iCandidate++) {
259 fCurrentTrack = (AliMuonForwardTrack*) fCandidateTracks->UncheckedAt(iCandidate);
261 // if the old track is compatible with the new cluster, the track is updated and inserted as new track in the array
262 // (several new tracks can be created for one old track)
263 if (FindClusterInPlane(iPlane) == kDiverged) {
264 fGlobalTrackingDiverged = kTRUE;
268 if ((fNPlanesMFTAnalyzed-fCurrentTrack->GetNMFTClusters())>fNMaxMissingMFTClusters || fIsPlaneMandatory[iPlane]) {
269 fCandidateTracks->Remove(fCurrentTrack); // the old track is removed after the check;
272 if (fGlobalTrackingDiverged) {
273 if (fScaleSigmaClusterCut>0) fScaleSigmaClusterCut -= 0.1;
277 fCandidateTracks->Compress();
281 // -------------------------- END OF THE CYCLE OVER THE MFT PLANES --------------------------------------------
283 fGlobalTrackingDiverged = kFALSE;
284 fScaleSigmaClusterCut = 1.0;
286 AliDebug(1, "Finished cycle over planes");
290 // If we have several final tracks, we must find the best candidate:
292 Int_t nFinalTracks = fCandidateTracks->GetEntriesFast();
293 AliInfo(Form("nFinalTracks = %d", nFinalTracks));
295 Int_t nGoodClustersBestCandidate = 0;
296 Int_t idBestCandidate = 0;
297 Double_t bestChi2 = -1.; // variable defining the best candidate
299 for (Int_t iFinalCandidate=0; iFinalCandidate<nFinalTracks; iFinalCandidate++) {
301 AliMuonForwardTrack *finalTrack = (AliMuonForwardTrack*) fCandidateTracks->UncheckedAt(iFinalCandidate);
302 Int_t nMFTClusters = finalTrack->GetNMFTClusters();
305 for (Int_t iCluster=0; iCluster<nMFTClusters; iCluster++) {
306 AliMFTCluster *localCluster = finalTrack->GetMFTCluster(iCluster);
307 chi2 += localCluster->GetLocalChi2();
309 chi2 /= nMFTClusters;
311 // now comparing the tracks in order to find the best one
313 if (chi2<bestChi2 || bestChi2<0) {
315 idBestCandidate = iFinalCandidate;
322 AliMuonForwardTrack *newTrack = (AliMuonForwardTrack*) fCandidateTracks->UncheckedAt(idBestCandidate);
323 newTrack -> SetNWrongClustersMC(newTrack->GetNMFTClusters() - nGoodClustersBestCandidate);
325 new ((*muonForwardTracks)[muonForwardTracks->GetEntries()]) AliMuonForwardTrack(*newTrack);
327 //----------------------- Save the information to the AliESDMuonGlobalTrack object
329 newTrack -> EvalKinem(AliMFTConstants::fZEvalKinem);
331 AliDebug(2,"Creating a new Muon Global Track");
332 AliESDMuonGlobalTrack *myESDTrack = event->NewMuonGlobalTrack();
333 myESDTrack -> SetPxPyPz(newTrack->Px(), newTrack->Py(), newTrack->Pz());
335 myESDTrack -> SetLabel(newTrack->GetMCLabel());
336 myESDTrack -> SetChi2OverNdf(newTrack->GetChi2OverNdf());
337 myESDTrack -> SetCharge(newTrack->GetCharge());
338 myESDTrack -> SetMatchTrigger(newTrack->GetMatchTrigger());
339 myESDTrack -> SetNMFTClusters(newTrack->GetNMFTClusters());
340 myESDTrack -> SetNWrongMFTClustersMC(newTrack->GetNWrongClustersMC());
341 myESDTrack -> SetFirstTrackingPoint(newTrack->GetMFTCluster(0)->GetX(), newTrack->GetMFTCluster(0)->GetY(), newTrack->GetMFTCluster(0)->GetZ());
342 myESDTrack -> SetXYAtVertex(newTrack->GetOffsetX(0., AliMFTConstants::fZEvalKinem), newTrack->GetOffsetY(0., AliMFTConstants::fZEvalKinem));
343 myESDTrack -> SetRAtAbsorberEnd(newTrack->GetRAtAbsorberEnd());
344 myESDTrack -> SetCovariances(newTrack->GetTrackParamAtMFTCluster(0)->GetCovariances());
345 myESDTrack -> SetChi2MatchTrigger(esdTrack->GetChi2MatchTrigger());
346 myESDTrack -> SetMuonClusterMap(esdTrack->GetMuonClusterMap());
347 myESDTrack -> SetHitsPatternInTrigCh(esdTrack->GetHitsPatternInTrigCh());
348 myESDTrack -> SetHitsPatternInTrigChTrk(esdTrack->GetHitsPatternInTrigChTrk());
349 myESDTrack -> Connected(esdTrack->IsConnected());
351 ULong_t mftClusterPattern = 0;
352 for (Int_t iCluster=0; iCluster<newTrack->GetNMFTClusters(); iCluster++) {
353 AliMFTCluster *localCluster = newTrack->GetMFTCluster(iCluster);
354 mftClusterPattern |= 1 << localCluster->GetPlane();
355 mftClusterPattern |= IsCorrectMatch(localCluster, newTrack->GetMCLabel()) << (fNMaxPlanes + localCluster->GetPlane());
357 myESDTrack -> SetMFTClusterPattern(mftClusterPattern);
359 //---------------------------------------------------------------------------------
363 fFinalBestCandidate = NULL;
367 outputTreeMuonGlobalTracks -> Fill();
370 while (outputFileMuonGlobalTracks->cd(Form("Event%d",myEventID))) myEventID++;
371 outputFileMuonGlobalTracks -> mkdir(Form("Event%d",myEventID));
372 outputFileMuonGlobalTracks -> cd(Form("Event%d",myEventID));
373 outputTreeMuonGlobalTracks -> Write();
374 outputFileMuonGlobalTracks -> Close();
376 muonForwardTracks -> Delete();
377 delete muonForwardTracks;
383 //=========================================================================================================================================
385 void AliMFTTrackerMU::SeparateFrontBackClusters() {
387 for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
388 fMFTClusterArrayFront[iPlane]->Delete();
389 fMFTClusterArrayBack[iPlane] ->Delete();
390 for (Int_t iCluster=0; iCluster<fMFTClusterArray[iPlane]->GetEntries(); iCluster++) {
391 AliMFTCluster *cluster = (AliMFTCluster*) fMFTClusterArray[iPlane]->At(iCluster);
392 if (TMath::Abs(cluster->GetZ())<TMath::Abs(fSegmentation->GetPlane(iPlane)->GetZCenter())) {
393 new ((*fMFTClusterArrayFront[iPlane])[fMFTClusterArrayFront[iPlane]->GetEntries()]) AliMFTCluster(*cluster);
396 new ((*fMFTClusterArrayBack[iPlane])[fMFTClusterArrayBack[iPlane]->GetEntries()]) AliMFTCluster(*cluster);
403 //==========================================================================================================================================
405 Int_t AliMFTTrackerMU::FindClusterInPlane(Int_t planeId) {
407 // !!!!!!!!! coordinates and errors on the interaction vertex should be taken from the event itself (ITS) if available
409 // propagate track to plane #planeId (both to front and back active sensors)
410 // look for compatible clusters
411 // update TrackParam at found cluster (if any) using Kalman Filter
413 AliMUONTrackParam currentParamFront, currentParamBack, currentParamForResearchFront, currentParamForResearchBack;
415 if (planeId == fNPlanesMFT-1) { // last plane of the telecope
416 currentParamFront = (*((AliMUONTrackParam*)(fMUONTrack->GetTrackParamAtCluster()->First())));
417 currentParamBack = (*((AliMUONTrackParam*)(fMUONTrack->GetTrackParamAtCluster()->First())));
418 currentParamForResearchFront = currentParamFront;
419 currentParamForResearchBack = currentParamBack;
420 if (fBransonCorrection) {
421 AliMUONTrackExtrap::ExtrapToVertex(¤tParamFront, fXExtrapVertex, fYExtrapVertex, fZExtrapVertex, fXExtrapVertexError, fYExtrapVertexError);
422 AliMUONTrackExtrap::ExtrapToVertex(¤tParamBack, fXExtrapVertex, fYExtrapVertex, fZExtrapVertex, fXExtrapVertexError, fYExtrapVertexError);
425 AliMUONTrackExtrap::ExtrapToVertexWithoutBranson(¤tParamFront, fZExtrapVertex);
426 AliMUONTrackExtrap::ExtrapToVertexWithoutBranson(¤tParamBack, fZExtrapVertex);
428 AliMUONTrackExtrap::ExtrapToVertex(¤tParamForResearchFront, fXExtrapVertex, fYExtrapVertex, fZExtrapVertex, fXExtrapVertexError, fYExtrapVertexError);
429 AliMUONTrackExtrap::ExtrapToVertex(¤tParamForResearchBack, fXExtrapVertex, fYExtrapVertex, fZExtrapVertex, fXExtrapVertexError, fYExtrapVertexError);
431 else { // MFT planes others than the last one: mult. scattering correction because of the upstream MFT planes is performed
432 currentParamFront = (*((AliMUONTrackParam*)(fCurrentTrack->GetTrackParamAtCluster()->First())));
433 currentParamBack = (*((AliMUONTrackParam*)(fCurrentTrack->GetTrackParamAtCluster()->First())));
434 currentParamForResearchFront = currentParamFront;
435 currentParamForResearchBack = currentParamBack;
436 AliMUONTrackExtrap::AddMCSEffect(¤tParamFront, (fSegmentation->GetPlane(planeId+1)->GetEquivalentSilicon()+
437 fSegmentation->GetPlane(planeId)->GetEquivalentSiliconBeforeFront())/fRadLengthSi,-1.);
438 AliMUONTrackExtrap::AddMCSEffect(¤tParamForResearchFront,(fSegmentation->GetPlane(planeId+1)->GetEquivalentSilicon()+
439 fSegmentation->GetPlane(planeId)->GetEquivalentSiliconBeforeFront())/fRadLengthSi,-1.);
440 AliMUONTrackExtrap::AddMCSEffect(¤tParamBack, (fSegmentation->GetPlane(planeId+1)->GetEquivalentSilicon()+
441 fSegmentation->GetPlane(planeId)->GetEquivalentSiliconBeforeBack())/fRadLengthSi,-1.);
442 AliMUONTrackExtrap::AddMCSEffect(¤tParamForResearchBack, (fSegmentation->GetPlane(planeId+1)->GetEquivalentSilicon()+
443 fSegmentation->GetPlane(planeId)->GetEquivalentSiliconBeforeBack())/fRadLengthSi,-1.);
445 // for all planes: extrapolation to the Z of the plane
446 AliMUONTrackExtrap::ExtrapToZCov(¤tParamFront, -1.*fSegmentation->GetPlane(planeId)->GetZCenterActiveFront());
447 AliMUONTrackExtrap::ExtrapToZCov(¤tParamForResearchFront, -1.*fSegmentation->GetPlane(planeId)->GetZCenterActiveFront());
448 AliMUONTrackExtrap::ExtrapToZCov(¤tParamBack, -1.*fSegmentation->GetPlane(planeId)->GetZCenterActiveBack());
449 AliMUONTrackExtrap::ExtrapToZCov(¤tParamForResearchBack, -1.*fSegmentation->GetPlane(planeId)->GetZCenterActiveBack());
451 //---------------------------------------------------------------------------------------
453 TMatrixD covFront(5,5); covFront = currentParamForResearchFront.GetCovariances();
454 TMatrixD covBack(5,5); covBack = currentParamForResearchBack.GetCovariances();
456 Double_t squaredError_X_Front = covFront(0,0);
457 Double_t squaredError_Y_Front = covFront(2,2);
458 Double_t squaredError_X_Back = covBack(0,0);
459 Double_t squaredError_Y_Back = covBack(2,2);
461 Double_t corrFact = 1.0;
463 Double_t researchRadiusFront = TMath::Sqrt(squaredError_X_Front + squaredError_Y_Front);
464 Double_t researchRadiusBack = TMath::Sqrt(squaredError_X_Back + squaredError_Y_Back);
465 if (0.5*(researchRadiusFront+researchRadiusBack)<fMinResearchRadiusAtPlane[planeId]) {
466 corrFact = fMinResearchRadiusAtPlane[planeId]/(0.5*(researchRadiusFront+researchRadiusBack));
469 //---------------------------------------------------------------------------------------
471 Double_t chi2cut = 2.*fScaleSigmaClusterCut*fScaleSigmaClusterCut*fSigmaClusterCut*fSigmaClusterCut; // depends on the number of variables (here, 2)
473 // Analyizing the clusters: FRONT ACTIVE ELEMENTS
475 Int_t nClustersFront = fMFTClusterArrayFront[planeId]->GetEntries();
476 AliDebug(2, Form("There are %3d clusters in plane %02d FRONT\n", nClustersFront, planeId));
478 for (Int_t iCluster=0; iCluster<nClustersFront; iCluster++) {
480 Bool_t isGoodChi2 = kFALSE;
482 AliMFTCluster *cluster = (AliMFTCluster*) fMFTClusterArrayFront[planeId]->At(iCluster);
483 Double_t chi2 = (1./(corrFact*corrFact)) * TryOneCluster(currentParamForResearchFront, cluster); // describes the compatibility between the track and the cluster
484 if (chi2<chi2cut) isGoodChi2 = kTRUE;
487 AliDebug(3, Form("accepting cluster: chi2=%f (cut = %f)\n", chi2, chi2cut));
488 AliMuonForwardTrack *newTrack = new ((*fCandidateTracks)[fCandidateTracks->GetEntriesFast()]) AliMuonForwardTrack(*fCurrentTrack);
489 if (fCandidateTracks->GetEntriesFast() > fMaxNCandidates) return kDiverged;
490 newTrack->AddTrackParamAtMFTCluster(currentParamFront, *cluster); // creating new track param and attaching the cluster
491 AliDebug(2, Form("After plane %02d: newTrack->GetNMFTClusters() = %d (fCurrentTrack->GetNMFTClusters() = %d)",
492 planeId, newTrack->GetNMFTClusters(), fCurrentTrack->GetNMFTClusters()));
493 newTrack->SetPlaneExists(planeId);
495 else AliDebug(3, Form("discarding cluster: chi2=%f (cut = %f)\n", chi2, chi2cut));
499 // Analyizing the clusters: BACK ACTIVE ELEMENTS
501 Int_t nClustersBack = fMFTClusterArrayBack[planeId]->GetEntries();
502 AliDebug(2, Form("There are %3d clusters in plane %02d BACK\n", nClustersBack, planeId));
504 for (Int_t iCluster=0; iCluster<nClustersBack; iCluster++) {
506 Bool_t isGoodChi2 = kFALSE;
508 AliMFTCluster *cluster = (AliMFTCluster*) fMFTClusterArrayBack[planeId]->At(iCluster);
509 Double_t chi2 = (1./(corrFact*corrFact)) * TryOneCluster(currentParamForResearchBack, cluster); // describes the compatibility between the track and the cluster
510 if (chi2<chi2cut) isGoodChi2 = kTRUE;
513 AliDebug(3,Form("accepting cluster: chi2=%f (cut = %f)\n", chi2, chi2cut));
514 AliMuonForwardTrack *newTrack = new ((*fCandidateTracks)[fCandidateTracks->GetEntriesFast()]) AliMuonForwardTrack(*fCurrentTrack);
515 if (fCandidateTracks->GetEntriesFast() > fMaxNCandidates) return kDiverged;
516 newTrack->AddTrackParamAtMFTCluster(currentParamBack, *cluster); // creating new track param and attaching the cluster
517 AliDebug(2, Form("After plane %02d: newTrack->GetNMFTClusters() = %d (fCurrentTrack->GetNMFTClusters() = %d)",
518 planeId, newTrack->GetNMFTClusters(), fCurrentTrack->GetNMFTClusters()));
519 newTrack->SetPlaneExists(planeId);
521 else AliDebug(3,Form("discarding cluster: chi2=%f (cut = %f)\n", chi2, chi2cut));
525 //---------------------------------------------------------------------------------------------
531 //==========================================================================================================================================
533 Double_t AliMFTTrackerMU::TryOneCluster(const AliMUONTrackParam &trackParam, AliMFTCluster *cluster) {
535 // Test the compatibility between the track and the cluster (using trackParam's covariance matrix):
536 // return the corresponding Chi2
537 // assume the track parameters are given at the Z of the cluster
539 // Set differences between trackParam and cluster in the bending and non bending directions
540 Double_t dX = cluster->GetX() - trackParam.GetNonBendingCoor();
541 Double_t dY = cluster->GetY() - trackParam.GetBendingCoor();
542 AliDebug(3,Form("dX = %f, dY = %f\n", dX, dY));
544 // Calculate errors and covariances
545 const TMatrixD& kParamCov = trackParam.GetCovariances();
546 Double_t sigmaX2 = kParamCov(0,0) + cluster->GetErrX2();
547 Double_t sigmaY2 = kParamCov(2,2) + cluster->GetErrY2();
548 AliDebug(3, Form("dX2 = %f, dY2 = %f\n", sigmaX2, sigmaY2));
549 Double_t covXY = kParamCov(0,2);
550 Double_t det = sigmaX2 * sigmaY2 - covXY * covXY;
553 if (det==0.) return 1.e10;
554 return (dX*dX*sigmaY2 + dY*dY*sigmaX2 - 2.*dX*dY*covXY) / det;
558 //=========================================================================================================================================
560 Bool_t AliMFTTrackerMU::IsCorrectMatch(AliMFTCluster *cluster, Int_t labelMC) {
562 Bool_t result = kFALSE;
564 // check if the cluster belongs to the correct MC track
566 for (Int_t iTrack=0; iTrack<cluster->GetNMCTracks(); iTrack++) {
567 if (cluster->GetMCLabel(iTrack)==labelMC) {
577 //======================================================================================================================================
579 void AliMFTTrackerMU::GetVertexFromMC() {
581 AliRunLoader *runLoader = AliRunLoader::Open("galice.root");
583 AliError("no run loader found in file galice.root");
587 runLoader->CdGAFile();
588 runLoader->LoadgAlice();
589 runLoader->LoadHeader();
590 runLoader->GetEvent(gAlice->GetEvNumber());
593 runLoader->GetHeader()->GenEventHeader()->PrimaryVertex(vtx);
594 AliInfo(Form("Primary vertex from MC found in (%f, %f, %f)\n",vtx[0], vtx[1], vtx[2]));
596 fXExtrapVertex = gRandom->Gaus(vtx[0], AliMFTConstants::fPrimaryVertexResX);
597 fYExtrapVertex = gRandom->Gaus(vtx[1], AliMFTConstants::fPrimaryVertexResY);
598 fZExtrapVertex = gRandom->Gaus(vtx[2], AliMFTConstants::fPrimaryVertexResZ);
599 fXExtrapVertexError = AliMFTConstants::fXVertexTolerance;
600 fYExtrapVertexError = AliMFTConstants::fYVertexTolerance;
601 AliInfo(Form("Set ESD vertex from MC (%f +/- %f, %f +/- %f, %f)",
602 fXExtrapVertex,fXExtrapVertexError,fYExtrapVertex,fYExtrapVertexError,fZExtrapVertex));
606 //======================================================================================================================================
608 void AliMFTTrackerMU::AddClustersFromUnderlyingEvent() {
610 AliInfo("Adding clusters from underlying event");
614 TGrid::Connect("alien://");
616 TFile* fileWithClustersToAdd = TFile::Open(fMFT->GetFileNameForUnderlyingEvent());
617 if (!fileWithClustersToAdd) return;
618 if (!(fileWithClustersToAdd->IsOpen())) return;
619 if (!(fileWithClustersToAdd->cd(Form("Event%d",fMFT->GetUnderlyingEventID())))) return;
621 TClonesArray *recPointsPerPlaneToAdd[AliMFTConstants::fNMaxPlanes] = {0};
624 for (Int_t iPlane=0; iPlane<AliMFTConstants::fNMaxPlanes; iPlane++) recPointsPerPlaneToAdd[iPlane] = new TClonesArray("AliMFTCluster");
626 treeIn = (TTree*) gDirectory->Get("TreeR");
628 for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
629 if (!(treeIn->GetBranch(Form("Plane_%02d",iPlane)))) {
630 for (Int_t jPlane=0; jPlane<AliMFTConstants::fNMaxPlanes; jPlane++) delete recPointsPerPlaneToAdd[jPlane];
633 else treeIn->SetBranchAddress(Form("Plane_%02d",iPlane), &(recPointsPerPlaneToAdd[iPlane]));
636 treeIn -> GetEntry(0);
638 for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
639 printf("plane %d -> before = %d ",iPlane,fMFTClusterArray[iPlane]->GetEntries());
640 Int_t nClusters = recPointsPerPlaneToAdd[iPlane]->GetEntries();
641 for (Int_t iCluster=0; iCluster<nClusters; iCluster++) {
642 AliMFTCluster *newCluster = (AliMFTCluster*) recPointsPerPlaneToAdd[iPlane]->At(iCluster);
643 for (Int_t iTrack=0; iTrack<newCluster->GetNMCTracks(); iTrack++) newCluster->SetMCLabel(iTrack, newCluster->GetMCLabel(iTrack)+AliMFTConstants::fLabelOffsetMC);
644 new ((*fMFTClusterArray[iPlane])[fMFTClusterArray[iPlane]->GetEntries()]) AliMFTCluster(*newCluster);
646 printf("after = %d\n",fMFTClusterArray[iPlane]->GetEntries());
649 for (Int_t jPlane=0; jPlane<AliMFTConstants::fNMaxPlanes; jPlane++) delete recPointsPerPlaneToAdd[jPlane];
653 //======================================================================================================================================
655 void AliMFTTrackerMU::AddClustersFromPileUpEvents() {
657 AliInfo("Adding clusters from pile-up event(s)");
661 TGrid::Connect("alien://");
663 TFile* fileWithClustersToAdd = TFile::Open(fMFT->GetFileNameForPileUpEvents());
664 if (!fileWithClustersToAdd) return;
665 if (!(fileWithClustersToAdd->IsOpen())) return;
667 TClonesArray *recPointsPerPlaneToAdd[AliMFTConstants::fNMaxPlanes] = {0};
670 for (Int_t iPileUp=0; iPileUp<AliMFTConstants::fNMaxPileUpEvents; iPileUp++) {
672 if (!(fileWithClustersToAdd->cd(Form("Event%d",fMFT->GetPileUpEventID(iPileUp))))) continue;
674 for (Int_t iPlane=0; iPlane<AliMFTConstants::fNMaxPlanes; iPlane++) recPointsPerPlaneToAdd[iPlane] = new TClonesArray("AliMFTCluster");
676 treeIn = (TTree*) gDirectory->Get("TreeR");
678 for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
679 if (!(treeIn->GetBranch(Form("Plane_%02d",iPlane)))) {
680 for (Int_t jPlane=0; jPlane<AliMFTConstants::fNMaxPlanes; jPlane++) delete recPointsPerPlaneToAdd[jPlane];
683 else treeIn->SetBranchAddress(Form("Plane_%02d",iPlane), &(recPointsPerPlaneToAdd[iPlane]));
686 treeIn -> GetEntry(0);
688 for (Int_t iPlane=0; iPlane<fNPlanesMFT; iPlane++) {
689 printf("plane %d -> before = %d ",iPlane,fMFTClusterArray[iPlane]->GetEntries());
690 Int_t nClusters = recPointsPerPlaneToAdd[iPlane]->GetEntries();
691 for (Int_t iCluster=0; iCluster<nClusters; iCluster++) {
692 AliMFTCluster *newCluster = (AliMFTCluster*) recPointsPerPlaneToAdd[iPlane]->At(iCluster);
693 for (Int_t iTrack=0; iTrack<newCluster->GetNMCTracks(); iTrack++) newCluster->SetMCLabel(iTrack, newCluster->GetMCLabel(iTrack)+AliMFTConstants::fLabelOffsetMC);
694 new ((*fMFTClusterArray[iPlane])[fMFTClusterArray[iPlane]->GetEntries()]) AliMFTCluster(*newCluster);
696 printf("after = %d\n",fMFTClusterArray[iPlane]->GetEntries());
699 for (Int_t jPlane=0; jPlane<AliMFTConstants::fNMaxPlanes; jPlane++) delete recPointsPerPlaneToAdd[jPlane];
705 //======================================================================================================================================