]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MFT/AliMuonForwardTrackFinder.h
K0s analysis update (Matt Steinpreis)
[u/mrichter/AliRoot.git] / MFT / AliMuonForwardTrackFinder.h
CommitLineData
820b4d9e 1#ifndef AliMuonForwardTrackFinder_H
2#define AliMuonForwardTrackFinder_H
3
4// ROOT includes
5#include "TObject.h"
6#include "TClonesArray.h"
7#include "TH2D.h"
8#include "TH1D.h"
9#include "TFile.h"
10#include "TGeoManager.h"
11#include "TMatrixD.h"
12#include "TParticle.h"
13#include "TMath.h"
14#include "TGraph.h"
15#include "TEllipse.h"
16#include "TCanvas.h"
17#include "TString.h"
18#include "TLatex.h"
19#include "TMarker.h"
20#include "TNtuple.h"
21#include "TRandom.h"
22#include "TIterator.h"
23
24// STEER includes
25#include "AliLog.h"
26#include "AliRun.h"
27#include "AliRunLoader.h"
28#include "AliLoader.h"
29#include "AliHeader.h"
30#include "AliMC.h"
31#include "AliStack.h"
32#include "AliMagF.h"
33#include "AliTracker.h"
34#include "AliGRPObject.h"
35#include "AliRunInfo.h"
36
37// MUON includes
38#include "AliMUONConstants.h"
39#include "AliMUONTrack.h"
40#include "AliMUONRecoCheck.h"
41#include "AliMUONTrackParam.h"
42#include "AliMUONTrackExtrap.h"
43#include "AliMUONVTrackStore.h"
44#include "AliMUONVCluster.h"
45
46// MFT includes
47#include "AliMuonForwardTrack.h"
48#include "AliMFTCluster.h"
49#include "AliMFT.h"
50#include "AliMFTSegmentation.h"
d4643a10 51#include "AliMFTConstants.h"
820b4d9e 52
53//====================================================================================================================================================
54//
55// Class for the creation of the "global muon tracks" built from the clusters in the
56// muon spectrometer and the clusters of the Muon Forward Tracker. QA histograms are also created
57//
58// Contact author: antonio.uras@cern.ch
59//
60//====================================================================================================================================================
61
62class AliMuonForwardTrackFinder : public TObject {
63
64public:
65
66 enum matchingOption {kRealMatching, kIdealMatching};
67
68 AliMuonForwardTrackFinder();
69
d4643a10 70 virtual ~AliMuonForwardTrackFinder();
820b4d9e 71
72 enum {kAllClusters, kClustersGoodChi2, kClusterOfTrack, kClusterCorrectMC};
73
74 void Init(Int_t nRun, Char_t *readDir, Char_t *outDir, Int_t nEventsToAnalyze = -1);
75 Bool_t LoadNextEvent();
76 Int_t LoadNextTrack();
77 Int_t GetNDF(Int_t nClusters);
78 void PDGNameConverter(const Char_t *nameIn, Char_t *nameOut);
79 void DrawPlanes();
80 void Terminate();
81 void WriteHistos();
82
83 void SetRun(Int_t nRun) { fRun = nRun; }
84 void SetNEventsToAnalyze(Int_t nEventsToAnalyze) { fNEventsToAnalyze = nEventsToAnalyze; }
85 void SetSigmaSpectrometerCut(Double_t sigmaSpectrometerCut) { fSigmaSpectrometerCut = sigmaSpectrometerCut; }
86 void SetSigmaClusterCut(Double_t sigmaClusterCut) { fSigmaClusterCut = sigmaClusterCut; }
87 void SetChi2GlobalCut(Double_t chi2GlobalCut) { fChi2GlobalCut = chi2GlobalCut; }
88 void SetReadDir(Char_t *dirName) { fReadDir = dirName; }
89 void SetOutDir(Char_t *dirName) { fOutDir = dirName; }
90 void SetDraw(Bool_t drawOption);
91 void SetRAbsorberCut(Double_t rAbsorberCut) { fRAbsorberCut = rAbsorberCut; }
92 void SetLowPtCut(Double_t lowPtCut) { fLowPtCut = lowPtCut; }
93 void SetNFinalCandidatesCut(Int_t nFinalCandidatesCut) { fNFinalCandidatesCut = nFinalCandidatesCut; }
b03d16a3 94 void SetVertexError(Double_t xErr, Double_t yErr, Double_t zErr) { fVertexErrorX=xErr; fVertexErrorY=yErr; fVertexErrorZ=zErr; }
820b4d9e 95
96 Int_t GetRun() { return fRun; }
97 Int_t GetNEvents() { return fNEventsToAnalyze; }
98 Double_t GetSigmaSpectrometerCut() { return fSigmaSpectrometerCut; }
99 Double_t GetSigmaClusterCut() { return fSigmaClusterCut; }
100 Double_t GetChi2GlobalCut() { return fChi2GlobalCut; }
101 Double_t GetRAbsorberCut() { return fRAbsorberCut; }
102 Double_t GetLowPtCut() { return fLowPtCut; }
820b4d9e 103 Int_t GetNPlanesMFT() { return fNPlanesMFT; }
104 Int_t GetNFinalCandidatesCut() { return fNFinalCandidatesCut; }
105 Int_t GetCurrentEvent() { return fEv; }
106
107 Int_t GetNRealTracksAnalyzed() const { return fCountRealTracksAnalyzed; }
108 Int_t GetNRealTracksAnalyzedOfEvent() const { return fCountRealTracksAnalyzedOfEvent; }
109 Int_t GetNRealTracksWithRefMC() const { return fCountRealTracksWithRefMC; }
110 Int_t GetNRealTracksWithRefMC_andTrigger() const { return fCountRealTracksWithRefMC_andTrigger; }
111 Int_t GetNRealTracksWithRefMC_andTrigger_andGoodPt() const { return fCountRealTracksWithRefMC_andTrigger_andGoodPt; }
112 Int_t GetNRealTracksWithRefMC_andTrigger_andGoodPt_andGoodTheta() const { return fCountRealTracksWithRefMC_andTrigger_andGoodPt_andGoodTheta; }
113
114 void SetNPlanesMFT(Int_t nPlanesMFT) { fNPlanesMFT = nPlanesMFT; }
115 void SeparateFrontBackClusters();
116 void SetNMaxMissingMFTClusters(Int_t nMaxMissingMFTClusters) { fNMaxMissingMFTClusters = nMaxMissingMFTClusters; }
d4643a10 117 void SetMandatoryPlane(Int_t iPlane) { if (0<=iPlane && iPlane<AliMFTConstants::fNMaxPlanes) fIsPlaneMandatory[iPlane] = kTRUE; }
820b4d9e 118
a2b7dc2a 119 Int_t FindClusterInPlane(Int_t planeId);
820b4d9e 120 void AttachGoodClusterInPlane(Int_t planeId);
d4643a10 121 void FillPlanesWithTrackHistory();
820b4d9e 122 Double_t TryOneCluster(const AliMUONTrackParam &trackParam, AliMFTCluster *cluster);
123 void BookHistos();
124 void SetTitleHistos();
125 void BookPlanes();
126 void ResetPlanes();
127 void PrintParticleHistory();
128
129 Bool_t IsCorrectMatch(AliMFTCluster *cluster);
130 void CheckCurrentMuonTrackable();
df4bbbc3 131 Bool_t IsMother(const Char_t *nameMother);
820b4d9e 132
133 void SetMatchingMode(Int_t matchingMode) { fMatchingMode = matchingMode; }
b5ab1ac4 134 void SetMinResearchRadiusAtPlane(Int_t plane, Double_t radius) { if (plane>=0 && plane<fNMaxPlanes) fMinResearchRadiusAtPlane[plane] = radius; }
820b4d9e 135
136 void FillOutputTree();
137 void WriteOutputTree();
138
139 Bool_t InitGRP();
140 Bool_t SetRunNumber();
141
d4643a10 142 void SetMaxNTracksToBeAnalyzed(Int_t nTracks) { fMaxNTracksToBeAnalyzed = nTracks; }
7e3dd1af 143 void SetBransonCorrection(Bool_t correction) { fBransonCorrection = correction; }
d4643a10 144
820b4d9e 145private:
146
147 AliMuonForwardTrackFinder(const AliMuonForwardTrackFinder& obj);
148 AliMuonForwardTrackFinder& operator=(const AliMuonForwardTrackFinder& other);
149
150protected:
151
d4643a10 152 static const Int_t fNMaxPlanes = AliMFTConstants::fNMaxPlanes; // max number of MFT planes
153 static const Double_t fRadLengthSi;
a2b7dc2a 154 static const Int_t fMaxNCandidates = 1000;
d4643a10 155
820b4d9e 156 Int_t fRun;
157 Int_t fNEventsToAnalyze; // events to analyze
158 Double_t fSigmaClusterCut; // to select the clusters in the MFT planes which are compatible with the extrapolated muon track
a2b7dc2a 159 Double_t fScaleSigmaClusterCut; // to tune the cut on the compatible clusters in case of too many candidates
160 Bool_t fGlobalTrackingDiverged; // to keep memory of a possible divergence in the global tracking finding
820b4d9e 161 Double_t fChi2GlobalCut; // cut on the final chi2 of the global muon track
162 Double_t fSigmaSpectrometerCut; // for the selection of the tracks in the muon spectrometer
b03d16a3 163 Double_t fVertexErrorX; // uncertainty on the x position of the muon's origin
164 Double_t fVertexErrorY; // uncertainty on the y position of the muon's origin
165 Double_t fVertexErrorZ; // uncertainty on the z position of the muon's origin
820b4d9e 166 Int_t fNFinalCandidatesCut; // cut on the number of the final candidates for the muon track
167 TString fReadDir;
168 TString fOutDir;
169 Bool_t fDrawOption;
170
171 Double_t fDistanceFromGoodClusterAndTrackAtLastPlane;
172 Double_t fDistanceFromBestClusterAndTrackAtLastPlane;
173
d4643a10 174 TClonesArray *fMFTClusterArray[fNMaxPlanes]; //! array of clusters for the planes of the MFT
175 TClonesArray *fMFTClusterArrayFront[fNMaxPlanes]; //! array of front clusters for the planes of the MFT
176 TClonesArray *fMFTClusterArrayBack[fNMaxPlanes]; //! array of back clusters for the planes of the MFT
820b4d9e 177
178 Double_t fRAbsorberCut; // in cm, corresponds to the radial position of a 3 degrees track at the end of the absorber (-503 cm)
179 Double_t fLowPtCut; // in GeV/c, the lower limit for the pt of a track in the muon spectrometer
180 Int_t fNPlanesMFT; // number of planes of the Vertex Telescope (Muon Internal Tracker) -> This should be taken from the new version of AliVZERO2
181 Int_t fNPlanesMFTAnalyzed;
182 Int_t fNMaxMissingMFTClusters; // max. number of MFT clusters which can be missed in the global fit procedure
d4643a10 183 Bool_t fIsPlaneMandatory[fNMaxPlanes]; // specifies which MFT planes cannot be missed in the global fit procedure
820b4d9e 184
185 Int_t fEv; // current event being analyzed
186 Int_t fLabelMC; // MC label of the muon track reconstructed in the spectrometer
187
188 Bool_t fIsClusterCompatible[10]; // here the clusters in the Muon Spectrometer are concerned
189
d4643a10 190 Double_t fZPlane[fNMaxPlanes]; // z-position of the MFT planes (center of the support)
191 Double_t fRPlaneMax[fNMaxPlanes]; // max radius of the MFT planes (the support)
192 Double_t fRPlaneMin[fNMaxPlanes]; // min radius of the MFT planes (the support)
193
d4643a10 194 TH1D *fHistRadiusEndOfAbsorber, *fHistNTracksAfterExtrapolation[fNMaxPlanes]; //
195 TH1D *fHistNGoodClustersForFinalTracks, *fHistResearchRadius[fNMaxPlanes]; //
196 TH1D *fHistDistanceGoodClusterFromTrackMinusDistanceBestClusterFromTrackAtLastPlane; //
197 TH1D *fHistDistanceGoodClusterFromTrackAtLastPlane; //
198 TH1D *fHistChi2Cluster_GoodCluster[fNMaxPlanes], *fHistChi2Cluster_BadCluster[fNMaxPlanes]; //
199 TH1D *fHistGlobalChi2AtPlaneFor_GOOD_CandidatesOfTrackableMuons[fNMaxPlanes]; //
200 TH1D *fHistGlobalChi2AtPlaneFor_BAD_CandidatesOfTrackableMuons[fNMaxPlanes]; //
201
202 TNtuple *fNtuFinalCandidates;
203 TNtuple *fNtuFinalBestCandidates;
204
205 TGraph *fGrMFTPlane[4][20]; //!
206 TEllipse *fCircleExt[fNMaxPlanes], *fCircleInt[fNMaxPlanes]; //!
207 TCanvas *fCanvas; //!
208
209 TLatex *fTxtMuonHistory, *fTxtTrackGoodClusters, *fTxtTrackChi2[fNMaxPlanes]; //!
210 TLatex *fTxtTrackFinalChi2, *fTxtTrackMomentum, *fTxtFinalCandidates, *fTxtDummy; //!
211 TLatex *fTxtAllClust, *fTxtClustGoodChi2, *fTxtClustMC, *fTxtClustOfTrack; //!
212 TMarker *fMrkAllClust, *fMrkClustGoodChi2, *fMrkClustMC, *fMrkClustOfTrack; //!
213
214 Int_t fCountRealTracksAnalyzed;
215 Int_t fMaxNTracksToBeAnalyzed;
820b4d9e 216 Int_t fCountRealTracksWithRefMC;
217 Int_t fCountRealTracksWithRefMC_andTrigger;
218 Int_t fCountRealTracksWithRefMC_andTrigger_andGoodPt;
219 Int_t fCountRealTracksWithRefMC_andTrigger_andGoodPt_andGoodTheta;
220 Int_t fCountRealTracksAnalyzedOfEvent;
221 Int_t fCountRealTracksAnalyzedWithFinalCandidates;
222
d4643a10 223 TFile *fFileCluster; //!
224 TFile *fFileESD; //!
225 TFile *fFile_gAlice; //!
820b4d9e 226
d4643a10 227 AliRunLoader *fRunLoader; //!
228 AliLoader *fMFTLoader; //!
229 AliMUONRecoCheck *fMuonRecoCheck; //!
820b4d9e 230
d4643a10 231 TTree *fMFTClusterTree; //!
820b4d9e 232
d4643a10 233 AliMUONTrack *fMuonTrackReco; //! muon track being analyzed
234 AliMuonForwardTrack *fCurrentTrack; //! muon extrapolated track being tested
235 AliMuonForwardTrack *fFinalBestCandidate; //! best final candidate (if any)
820b4d9e 236 Bool_t fIsCurrentMuonTrackable;
d4643a10 237 Bool_t fIsGoodClusterInPlane[fNMaxPlanes];
820b4d9e 238
d4643a10 239 TClonesArray *fCandidateTracks; //! array of track we are going to build (starting from fMuonTrackReco)
820b4d9e 240
d4643a10 241 AliMUONVTrackStore *fTrackStore; //! list of reconstructed MUON tracks
242 AliMUONVTrackStore *fTrackRefStore; //! list of reconstructible MUON tracks
820b4d9e 243
244 TIterator *fNextTrack; //! Iterator for reading the MUON tracks
245
d4643a10 246 AliStack *fStack; //!
820b4d9e 247
d4643a10 248 AliMFT *fMFT; //!
249 AliMFTSegmentation *fSegmentation; //!
820b4d9e 250
d4643a10 251 TFile *fOutputTreeFile, *fOutputQAFile; //
252 TTree *fOutputEventTree; //!
820b4d9e 253
d4643a10 254 TClonesArray *fMuonForwardTracks; //! array of AliMuonForwardTrack
820b4d9e 255
256 Int_t fMatchingMode;
b5ab1ac4 257 Double_t fMinResearchRadiusAtPlane[fNMaxPlanes];
820b4d9e 258
d4643a10 259 AliGRPObject *fGRPData; //! Data from the GRP/GRP/Data CDB folder
260 AliRunInfo *fRunInfo; //!
7e3dd1af 261
262 Bool_t fBransonCorrection; // if TRUE, Branson Correction is applied when extrapolating the MUON tracks to the vertex region
820b4d9e 263
264 ClassDef(AliMuonForwardTrackFinder, 1);
265
266};
267
7e3dd1af 268//====================================================================================================================================================
820b4d9e 269
270#endif
271
272