]>
Commit | Line | Data |
---|---|---|
4c904d93 | 1 | #ifndef ALIANALYSISTASKLINKTOMC_H |
2 | #define ALIANALYSISTASKLINKTOMC_H | |
3 | /* This file is property of and copyright by the ALICE HLT Project * | |
4 | * ALICE Experiment at CERN, All rights reserved. * | |
5 | * See cxx source for full Copyright notice */ | |
6 | ||
27de2dfb | 7 | /* $Id$ */ |
8 | ||
4c904d93 | 9 | /// |
10 | /// @file AliAnalysisTaskLinkToMC.h | |
11 | /// @author Artur Szostak <artursz@iafrica.com> | |
12 | /// @date 27 Oct 2008 | |
13 | /// @brief Definition of the AliAnalysisTaskLinkToMC task to connect reconstructed tracks with MC tracks. | |
14 | /// | |
15 | ||
16 | #include "AliAnalysisTaskSE.h" | |
17 | ||
18 | class TList; | |
19 | class TH2D; | |
20 | class TMap; | |
21 | class AliESDMuonTrack; | |
22 | class AliMCParticle; | |
23 | ||
24 | /** | |
25 | * \class AliAnalysisTaskLinkToMC | |
26 | * Connects reconstructed tracks found in the ESD to their corresponding Monte | |
27 | * Carlo tracks via the proximity of their hit coordinates. | |
28 | */ | |
29 | class AliAnalysisTaskLinkToMC : public AliAnalysisTaskSE | |
30 | { | |
31 | public: | |
32 | AliAnalysisTaskLinkToMC(); | |
33 | AliAnalysisTaskLinkToMC(const char* name); | |
34 | ||
35 | virtual void UserCreateOutputObjects(); | |
36 | virtual void UserExec(Option_t* option); | |
37 | virtual void Terminate(Option_t* option); | |
38 | ||
39 | /// Returns the absolute maximum distance in centimetres between the | |
40 | /// cluster and track reference that is allowed in the X direction. | |
41 | Double_t HardCutLimitX() const { return fHardCutLimitX; } | |
42 | ||
43 | /// Returns the absolute maximum distance in centimetres between the | |
44 | /// cluster and track reference that is allowed in the Y direction. | |
45 | Double_t HardCutLimitY() const { return fHardCutLimitY; } | |
46 | ||
47 | /// Returns the absolute maximum distance in centimetres between the | |
48 | /// cluster and track reference that is allowed in the Z direction. | |
49 | Double_t HardCutLimitZ() const { return fHardCutLimitZ; } | |
50 | ||
51 | void HardCutLimitX(Double_t value); | |
52 | void HardCutLimitY(Double_t value); | |
53 | void HardCutLimitZ(Double_t value); | |
54 | ||
55 | Double_t SigmaCut() const; | |
56 | ||
57 | /// Sets the maximum sigma value to allow for each X and Y dimension | |
58 | /// between clusters and track references. | |
59 | void SigmaCut(Double_t value) { fSigmaCut2 = value*value; } | |
60 | ||
61 | /// Returns the minimum number of clusters that must match between the | |
62 | /// ESD track clusters and track references. | |
63 | Int_t MinClusters() const { return fMinClusters; } | |
64 | ||
65 | /// Sets the minimum number of clusters that must match between the ESD | |
66 | /// track clusters and track references. | |
67 | void MinClusters(Int_t value) { fMinClusters = value; } | |
68 | ||
69 | /// Returns the minimum number of clusters that must match between the | |
70 | /// ESD track clusters and track references in stations 1 and 2. | |
71 | Int_t MinClustersInSt12() const { return fMinClustersInSt12; } | |
72 | ||
73 | /// Sets the minimum number of clusters that must match between the ESD | |
74 | /// track clusters and track references in stations 1 and 2. | |
75 | void MinClustersInSt12(Int_t value) { fMinClustersInSt12 = value; } | |
76 | ||
77 | /// Returns the minimum number of clusters that must match between the | |
78 | /// ESD track clusters and track references in stations 4 and 5. | |
79 | Int_t MinClustersInSt45() const { return fMinClustersInSt45; } | |
80 | ||
81 | /// Sets the minimum number of clusters that must match between the ESD | |
82 | /// track clusters and track references in stations 4 and 5. | |
83 | void MinClustersInSt45(Int_t value) { fMinClustersInSt45 = value; } | |
84 | ||
85 | /// Returns the flag indicating if histograms should be generated for | |
86 | /// checking the quality of track matching. | |
87 | Bool_t GenerateHistograms() const { return fMakeHists; } | |
88 | ||
89 | /// Set the flag indicating if histograms should be generated for checking | |
90 | /// the quality of track matching. | |
91 | void GenerateHistograms(Bool_t value) { fMakeHists = value; } | |
92 | ||
93 | /// Returns the flag indicating if the Z coordinate is to be used for | |
94 | /// track matching or not. | |
95 | Bool_t UseZCoordinate() const { return fUseZCoordinate; } | |
96 | ||
97 | /// Set the flag indicating if the Z coordinate should be used for track | |
98 | /// matching. | |
99 | void UseZCoordinate(Bool_t value) { fUseZCoordinate = value; } | |
100 | ||
101 | Bool_t ChamberMustMatch(Int_t chamber) const; | |
102 | void ChamberMustMatch(Int_t chamber, Bool_t value); | |
103 | Bool_t StationMustMatch(Int_t station) const; | |
104 | void StationMustMatch(Int_t station, Bool_t value); | |
105 | ||
106 | private: | |
107 | ||
108 | /// Do not allow copying of this object. Method is not implemented. | |
109 | AliAnalysisTaskLinkToMC(const AliAnalysisTaskLinkToMC& obj); | |
110 | /// Do not allow copying of this object. Method is not implemented. | |
111 | AliAnalysisTaskLinkToMC& operator = (const AliAnalysisTaskLinkToMC& obj); | |
112 | ||
113 | void FindLinksToMC(TMap& links); | |
114 | void CalculateTrackCompatibility( | |
115 | AliESDMuonTrack* esdTrack, AliMCParticle* mcTrack, | |
116 | bool& tracksCompatible, Double_t& chi2, Int_t& noOfClusters | |
117 | ); | |
118 | ||
119 | bool IsFindable(AliMCParticle* track) const; | |
120 | ||
121 | void FillHistsFromMC(); | |
122 | void FillHistsFromLinks(TMap& links); | |
123 | void CreateAODTracks(TMap& links); | |
124 | ||
125 | TList* fHistos; //! List of all histograms produced by this analysis task. | |
126 | TH2D* fFindableHist; //! Histogram of reconstructible tracks (filled with MC information). | |
127 | TH2D* fFoundHistMC; //! Histogram of reconstructed tracks (filled with MC information). | |
128 | TH2D* fFoundHist; //! Histogram of reconstructed tracks (filled with reconstructed information). | |
129 | TH2D* fFakeHist; //! Histogram of fake tracks (filled with reconstructed information). | |
130 | ||
131 | Double_t fHardCutLimitX; //! The absolute maximum difference between a cluster and reference track point in X direction in centimetres. | |
132 | Double_t fHardCutLimitY; //! The absolute maximum difference between a cluster and reference track point in Y direction in centimetres. | |
133 | Double_t fHardCutLimitZ; //! The absolute maximum difference between a cluster and reference track point in Z direction in centimetres. | |
134 | Double_t fVarX; //! The variance to use in the X direction if error value in the cluster was not set. | |
135 | Double_t fVarY; //! The variance to use in the Y direction if error value in the cluster was not set. | |
136 | Double_t fSigmaCut2; //! The cut on the chi2 X/Y direction component to make for each cluster fit. | |
137 | Int_t fMinClusters; //! Minimum number of clusters that must match between an ESD and MC track for it to be considered compatible. | |
138 | Int_t fMinClustersInSt12; //! The minimum number of clusters that must match in stations 1 and 2. | |
139 | Int_t fMinClustersInSt45; //! The minimum number of clusters that must match in stations 4 and 5. | |
140 | Bool_t fMakeHists; //! Flag indicating if histogram filling and generation should be performed. | |
141 | Bool_t fUseZCoordinate; //! Flag indicating if the Z coordinate should be used for comparrison. | |
142 | Bool_t fChamberMustMatch[14]; //! A flag for each chamber which must have a compatible hit matched. | |
143 | Bool_t fStationMustMatch[7]; //! A flag for each station which must have at least one compatible hit matched. | |
144 | ||
805e576c | 145 | ClassDef(AliAnalysisTaskLinkToMC, 0) // Analysis task to connect reconstructed tracks to their Monte Carlo counterparts. |
4c904d93 | 146 | }; |
147 | ||
148 | #endif // ALIANALYSISTASKLINKTOMC_H | |
149 |