]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskJetMatching.h
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / UserTasks / AliAnalysisTaskJetMatching.h
1 #ifndef AliAnalysisTaskJetMatching_H
2 #define AliAnalysisTaskJetMatching_H
3
4 //#define DEBUGTASK
5
6 #include <AliAnalysisTaskEmcalJet.h>
7 #include <AliEmcalJet.h>
8 #include <AliVTrack.h>
9 #include <TClonesArray.h>
10 #include <TMath.h>
11 #include <TRandom3.h>
12
13 class TF1;
14 class TH1F;
15 class TH2F;
16 class TH3F;
17 class TProfile;
18 class AliRhoParameter;
19 class AliLocalRhoParameter;
20 class TClonesArray;
21
22 class AliAnalysisTaskJetMatching : public AliAnalysisTaskEmcalJet
23 {
24     public:
25         // enumerators
26         enum matchingSceme      {kGeoEtaPhi, kGeoR, kDiJet};
27         enum sourceBKG          {kNoSourceBKG, kSourceRho, kSourceLocalRho};
28         enum targetBKG          {kNoTargetBKG, kTargetRho, kTargetLocalRho};
29         // constructors, destructor
30                                 AliAnalysisTaskJetMatching();
31                                 AliAnalysisTaskJetMatching(const char *name);
32         virtual                 ~AliAnalysisTaskJetMatching();
33         // setting up the task and technical aspects
34         void                    ExecOnce();
35         virtual Bool_t          Notify();
36         virtual void            UserCreateOutputObjects();
37         TH1F*                   BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Bool_t append = kTRUE);
38         TH2F*                   BookTH2F(const char* name, const char* x, const char* y, Int_t binsx, Double_t minx, Double_t maxx, Int_t binsy, Double_t miny, Double_t maxy, Bool_t append = kTRUE);
39         TH3F*                   BookTH3F(const char* name, const char* x, const char* y, const char* z, Int_t binsx, Double_t minx, Double_t maxx, Int_t binsy, Double_t miny, Double_t maxy, Int_t binsz, Double_t minz, Double_t maxz, Bool_t append = kTRUE);
40
41         virtual Bool_t          Run();
42         /* inline */    Double_t PhaseShift(Double_t x) const {  
43             while (x>=TMath::TwoPi())x-=TMath::TwoPi();
44             while (x<0.)x+=TMath::TwoPi();
45             return x; }
46         /* inline */    Double_t PhaseShift(Double_t x, Double_t n) const {
47             x = PhaseShift(x);
48             if(TMath::Nint(n)==2) while (x>TMath::Pi()) x-=TMath::Pi();
49             if(TMath::Nint(n)==3) {
50                 if(x>2.*TMath::TwoPi()/n) x = TMath::TwoPi() - x;
51                 if(x>TMath::TwoPi()/n) x = TMath::TwoPi()-(x+TMath::TwoPi()/n);
52             }
53             return x; }
54         /* inline */    Bool_t CompareTracks(AliVParticle* a, AliVParticle* b) const {
55             return (TMath::AreEqualAbs(a->Pt(), b->Pt(), 1e-2) && TMath::AreEqualAbs(a->Eta(), b->Eta(), 1e-2) && TMath::AreEqualAbs(a->Phi(), b->Phi(), 1e-2)) ? kTRUE : kFALSE; }
56         /* inline */    Double_t GetR(AliVParticle* a, AliVParticle* b) const {
57                if(!(a&&b)) return 999;
58                Double_t phiA(a->Phi()), phiB(b->Phi()), etaA(a->Eta()), etaB(b->Eta());
59                if(TMath::Abs(phiA-phiB) > TMath::Abs(phiA-phiB + TMath::TwoPi())) phiA+=TMath::TwoPi();
60                if(TMath::Abs(phiA-phiB) > TMath::Abs(phiA-phiB - TMath::TwoPi())) phiA-=TMath::TwoPi();
61                return TMath::Sqrt(TMath::Abs((etaA-etaB)*(etaA-etaB)+(phiA-phiB)*(phiA-phiB))); }
62
63         // setters - setup how to run
64         void                    SetMatchingScheme(matchingSceme m)              {fMatchingScheme = m;}
65         void                    SetMatchConstituents(Bool_t m)                  {fMatchConstituents = m;}
66         void                    SetMinFracRecoveredConstituents(Float_t f)      {fMinFracRecoveredConstituents = f;}
67         void                    SetMinFracRecoveredConstituentPt(Float_t f)     {fMinFracRecoveredConstituentPt = f;}
68         void                    SetUseEmcalBaseJetCuts(Bool_t b)                {fUseEmcalBaseJetCuts = b;}
69         void                    SetGetBijection(Bool_t b)                       {fGetBijection = b;}
70         void                    SetSourceBKG(sourceBKG b)                       {fSourceBKG = b;}
71         void                    SetTargetBKG(targetBKG b)                       {fTargetBKG = b;}
72         void                    SetSourceJetsName(const char* n)                {fSourceJetsName = n;}
73         void                    SetTargetJetsName(const char* n)                {fTargetJetsName = n; }
74         void                    SetMatchedJetsName(const char* n)               {fMatchedJetsName = n;}
75         void                    SetSourceLocalRhoName(const char* n)            {fSourceRhoName = n;}
76         void                    SetTargetLocalRhoName(const char* n)            {fTargetRhoName = n;}
77         void                    SetSourceRadius(Float_t r)                      {fSourceRadius = r;}
78         void                    SetTargetRadius(Float_t r)                      {fTargetRadius = r;}
79         void                    SetMatchEta(Float_t f)                          {fMatchEta = f;}
80         void                    SetMatchPhi(Float_t f)                          {fMatchPhi = f;}
81         void                    SetMatchR(Float_t f)                            {fMatchR = f;}
82         void                    SetDoDetectorResponse(Bool_t d)                 {fDoDetectorResponse = d;}
83         // methods
84         void                    DoGeometricMatchingEtaPhi();
85         void                    DoGeometricMatchingR();
86         void                    DoDiJetMatching();
87         void                    DoConstituentMatching();
88         void                    GetBijection();
89         void                    PostMatchedJets();
90         // fill histogrmas
91         void                    FillAnalysisSummaryHistogram() const;
92         void                    FillMatchedJetHistograms();
93         // setters - analysis details
94         /* inline */    Bool_t PassesCuts(const AliVTrack* track) const {
95             return (!track || TMath::Abs(track->Eta()) > 0.9 || track->Phi() < 0 || track->Phi() > TMath::TwoPi()) ? kFALSE : kTRUE; }
96         /* inline */    Bool_t PassesCuts(AliEmcalJet* jet, Int_t c)  { if(!fUseEmcalBaseJetCuts) return (jet) ? kTRUE : kFALSE; 
97             else return AcceptJet(jet, c);}
98         /* inline */    void ClearMatchedJetsCache() {
99             for (Int_t i(0); i < fNoMatchedJets; i++) {
100                 fMatchedJetContainer[i][0] = 0x0; fMatchedJetContainer[i][1] = 0x0; }
101             fNoMatchedJets = 0;
102         }
103         AliEmcalJet*            GetLeadingJet(TClonesArray* source, Int_t &leadingJetIndex, Double_t etaMin = -.7, Double_t etaMax = .7);
104         AliEmcalJet*            GetSubLeadingJet(TClonesArray* source, Int_t leadingJetIndex, Int_t &subLeadingJetIndex);
105         void                    PrintInfo() const;
106         virtual void            Terminate(Option_t* option);
107
108     private: 
109         TClonesArray*           fSourceJets;            //! array with source jets
110         TString                 fSourceJetsName;        // name of array with source jets
111         TClonesArray*           fTargetJets;            //! array with target jets
112         TString                 fTargetJetsName;        // name of array with target jets
113         TClonesArray*           fMatchedJets;           //! final list of matched jets which is added to event
114         TString                 fMatchedJetsName;       // name of list of matched jets
115         AliLocalRhoParameter*   fSourceRho;             //! source rho
116         TString                 fSourceRhoName;         // source rho  name
117         AliLocalRhoParameter*   fTargetRho;             //! target rho
118         TString                 fTargetRhoName;         // target rho name
119         Bool_t                  fUseScaledRho;          // use scaled rho
120         Float_t                 fSourceRadius;          // source radius 
121         Float_t                 fTargetRadius;          // target radius
122         matchingSceme           fMatchingScheme;        // select your favorite matching algorithm
123         Bool_t                  fUseEmcalBaseJetCuts;   // use the emcal jet base class for jet cuts
124         sourceBKG               fSourceBKG;             // subtracted background of source jets
125         targetBKG               fTargetBKG;             // subtracted background of target jets
126         // transient object pointers
127         TList*                  fOutputList;            //! output list
128         TH1F*                   fHistUnsortedCorrelation;       //! dphi correlation of unsorted jets
129         TH1F*                   fHistMatchedCorrelation;        //! dphi correlation of matched jets
130         TH1F*                   fHistSourceJetPt;       //! pt of source jets
131         TH1F*                   fHistMatchedSourceJetPt;//! pt of matched source jets
132         TH1F*                   fHistTargetJetPt;       //! pt of target jets
133         TH1F*                   fHistMatchedJetPt;      //! pt of matched jets
134         TH2F*                   fHistSourceMatchedJetPt;//! pt or source vs matched jets
135         TH2F*                   fHistDetectorResponseProb;      //! det prob
136         TH1F*                   fHistNoConstSourceJet;  //! number of constituents of source jet
137         TH1F*                   fHistNoConstTargetJet;  //! number of constituents of target jet
138         TH1F*                   fHistNoConstMatchJet;   //! number of constituents of matched jets
139         TProfile*               fProfFracPtMatched;     //! sum pt fraction for matched tracks / jet
140         TProfile*               fProfFracPtJets;        //! sum pt fraction for matched jets
141         TProfile*               fProfFracNoMatched;     //! no of constituents fraction found / jet
142         TProfile*               fProfFracNoJets;        //! no of consstituents fraction jet / jet
143         TH3F*                   fHistDiJet;             //! matched dijet eta, phi
144         TH3F*                   fHistDiJetLeadingJet;   //! leading jet (for dijet) eta, phi
145         TH2F*                   fHistDiJetDPhi;         //! dijet dphi
146         TH2F*                   fHistDiJetDPt;          //! dijet dpt
147         TH1F*                   fHistAnalysisSummary;   //! flags
148         TProfile*               fProfQAMatched;         //! QA spreads of matched jets
149         TProfile*               fProfQA;                //! QA spreads of source and target jets
150         Int_t                   fNoMatchedJets;         //! number of matched jets
151         AliEmcalJet*            fMatchedJetContainer[200][2];   //! pointers to matched jets
152         // geometric matching parameters
153         Float_t                 fMatchEta;              // max eta distance between centers of matched jets
154         Float_t                 fMatchPhi;              // max phi distance between centers of matched jets
155         Float_t                 fMatchR;                // max distance between matched jets
156         Bool_t                  fDoDetectorResponse;    // get detector response matrix (not normalized)
157         Bool_t                  fMatchConstituents;     // match constituents
158         Float_t                 fMinFracRecoveredConstituents;  // minimium fraction of constituents that needs to be found
159         Float_t                 fMinFracRecoveredConstituentPt; // minimium fraction of constituent pt that needs to be recovered
160         Bool_t                  fGetBijection;          // get bijection of source and matched jets
161         // pythia ntrials xsec bookkeeping
162         TH1F*                   fh1Trials;              //! sum of trails (in Notify)
163         TH1F*                   fh1AvgTrials;           //! average sum of trials (in Run)
164         TProfile*               fh1Xsec;                //! pythia cross section and trials
165         Float_t                 fAvgTrials;             // average number of trails per event
166
167         AliAnalysisTaskJetMatching(const AliAnalysisTaskJetMatching&);                  // not implemented
168         AliAnalysisTaskJetMatching& operator=(const AliAnalysisTaskJetMatching&);       // not implemented
169
170         ClassDef(AliAnalysisTaskJetMatching, 6);
171 };
172
173 #endif