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