ROOT logo
#ifndef ALIANALYSISTASKLOCALRHO_H
#define ALIANALYSISTASKLOCALRHO_H

// $Id$

#include <AliAnalysisTaskEmcalJet.h>
#include <AliEmcalJet.h>
#include <AliVEvent.h>
#include <AliVTrack.h>
#include <AliVCluster.h>
#include <TClonesArray.h>
#include <TMath.h>
#include <TRandom3.h>
#include <AliLog.h>
#include <AliJetContainer.h>

class TF1;
class THF1;
class THF2;
class TProfile;
class AliLocalRhoParameter;

class AliAnalysisTaskLocalRho : public AliAnalysisTaskEmcalJet {
 public:
  // enumerators
  enum fitModulationType  { kNoFit, kV2, kV3, kCombined, kFourierSeries, kIntegratedFlow, kQC2, kQC4 }; // fit type
  enum detectorType       { kTPC, kVZEROA, kVZEROC, kVZEROComb};  // detector that was used
  enum qcRecovery         { kFixedRho, kNegativeVn, kTryFit };    // how to deal with negative cn value for qcn value
  enum runModeType        { kLocal, kGrid };                      // run mode type
  // constructors, destructor
  AliAnalysisTaskLocalRho();
  AliAnalysisTaskLocalRho(const char *name, runModeType type);
  virtual                 ~AliAnalysisTaskLocalRho();
  // setting up the task and technical aspects
  void                    ExecOnce();
  Bool_t                  InitializeAnalysis();
  virtual void            UserCreateOutputObjects();
  TH1F*                   BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Int_t c = -1, Bool_t append = kTRUE);
  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, Int_t c = -1, Bool_t append = kTRUE);
  virtual Bool_t          Run();
  /* inline */   Double_t PhaseShift(Double_t x) const {  
    while (x>=TMath::TwoPi())x-=TMath::TwoPi();
    while (x<0.)x+=TMath::TwoPi();
    return x; }
  /* inline */   Double_t PhaseShift(Double_t x, Double_t n) const {
    x = PhaseShift(x);
    if(TMath::Nint(n)==2) while (x>TMath::Pi()) x-=TMath::Pi();
    if(TMath::Nint(n)==3) {
      if(x>2.*TMath::TwoPi()/n) x = TMath::TwoPi() - x;
      if(x>TMath::TwoPi()/n) x = TMath::TwoPi()-(x+TMath::TwoPi()/n);
    }
    return x; }
  /* inline */   Double_t ChiSquarePDF(Int_t ndf, Double_t x) const {
    Double_t n(ndf/2.), denom(TMath::Power(2, n)*TMath::Gamma(n));
    if (denom!=0)  return ((1./denom)*TMath::Power(x, n-1)*TMath::Exp(-x/2.)); 
    return -999; }
  // the cdf of the chisquare distribution is the normalized lower incomplete gamma function
  /* inline */    Double_t ChiSquareCDF(Int_t ndf, Double_t x) const { return TMath::Gamma(ndf/2., x/2.); }
  // setters - setup how to run
  void                    SetDebugMode(Int_t d)                           {fDebug = d;}
  void                    SetCentralityClasses(TArrayI* c)                {fCentralityClasses = c;}
  void                    SetAttachToEvent(Bool_t a)                      {fAttachToEvent = a;}
  void                    SetUseScaledRho(Bool_t s)                       {fUseScaledRho = s;}
  void                    SetFillHistograms(Bool_t b)                     {fFillHistograms = b;}
  // setters - analysis details
  void                    SetNoEventWeightsForQC(Bool_t e)                {fNoEventWeightsForQC = e;}
  void                    SetIntegratedFlow(TH1F* i, TH1F* j)             {fUserSuppliedV2 = i; fUserSuppliedV3 = j; }
  void                    SetOnTheFlyResCorrection(TH1F* r2, TH1F* r3)    {fUserSuppliedR2 = r2; fUserSuppliedR3 = r3; }
  void                    SetModulationFit(TF1* fit);
  void                    SetModulationFitMinMaxP(Float_t m, Float_t n)   {fMinPvalue = m; fMaxPvalue = n; }
  void                    SetModulationFitType(fitModulationType type)    {fFitModulationType = type; }
  void                    SetQCnRecoveryType(qcRecovery type)             {fQCRecovery = type; }
  void                    SetModulationFitOptions(TString opt)            {fFitModulationOptions = opt; }
  void                    SetReferenceDetector(detectorType type)         {fDetectorType = type; }
  void                    SetUsePtWeight(Bool_t w)                        {fUsePtWeight = w; }
  void                    SetUsePtWeightErrorPropagation(Bool_t w)        {fUsePtWeightErrorPropagation = w;}
  void                    SetRunModeType(runModeType type)                {fRunModeType = type; }
  void                    SetForceAbsVnHarmonics(Bool_t f)                {fAbsVnHarmonics = f; }
  void                    SetExcludeLeadingJetsFromFit(Float_t n)         {fExcludeLeadingJetsFromFit = n; }
  void                    SetRebinSwapHistoOnTheFly(Bool_t r)             {fRebinSwapHistoOnTheFly = r; }
  void                    SetSaveThisPercentageOfFits(Float_t p)          {fPercentageOfFits = p; }
  void                    SetUseV0EventPlaneFromHeader(Bool_t h)          {fUseV0EventPlaneFromHeader = h;}
  void                    SetSoftTrackMinMaxPt(Float_t min, Float_t max)  {fSoftTrackMinPt = min; fSoftTrackMaxPt = max;}
  // getters
  TString                 GetLocalRhoName() const                         {return fLocalRhoName; }
  // numerical evaluations
  void                    CalculateEventPlaneVZERO(Double_t vzero[2][2]) const;
  void                    CalculateEventPlaneTPC(Double_t* tpc);
  void                    CalculateEventPlaneCombinedVZERO(Double_t* comb) const;
  Double_t                CalculateQC2(Int_t harm);
  Double_t                CalculateQC4(Int_t harm);
  // helper calculations for the q-cumulant analysis, also used by AliAnalyisTaskJetFlow
  void                    QCnQnk(Int_t n, Int_t k, Double_t &reQ, Double_t &imQ);
  Double_t                QCnS(Int_t i, Int_t j);
  Double_t                QCnM();
  Double_t                QCnM11();
  Double_t                QCnM1111();
  Bool_t                  QCnRecovery(Double_t psi2, Double_t psi3);
  // analysis details
  Bool_t                  CorrectRho(Double_t psi2, Double_t psi3);
  void                    FillEventPlaneHistograms(Double_t psi2, Double_t psi3) const;
  void                    FillAnalysisSummaryHistogram() const;
  // track selection
  /* inline */     Bool_t PassesCuts(AliVTrack* track) const { return AcceptTrack(track, 0);}
  /* inline */     Bool_t PassesCuts(AliEmcalJet* jet) { return AcceptJet(jet, 0);}
  virtual void            Terminate(Option_t* option);

 private: 
  Int_t                   fDebug;                 // debug level (0 none, 1 fcn calls, 2 verbose)
  Bool_t                  fInitialized;           //! is the analysis initialized?
  Bool_t                  fAttachToEvent;         // attach local rho to the event
  Bool_t                  fFillHistograms;        // fill qa histograms
  Bool_t                  fNoEventWeightsForQC;   // don't store event weights for qc analysis
  Bool_t                  fUseScaledRho;          // use scaled rho
  TArrayI*                fCentralityClasses;     // centrality classes (maximum 10) used for QA
  TH1F*                   fUserSuppliedV2;        // histo with integrated v2
  TH1F*                   fUserSuppliedV3;        // histo with integrated v3
  TH1F*                   fUserSuppliedR2;        // correct the extracted v2 with this r
  TH1F*                   fUserSuppliedR3;        // correct the extracted v3 with this r
  Int_t                   fNAcceptedTracks;       //! number of accepted tracks
  Int_t                   fNAcceptedTracksQCn;    //! accepted tracks for QCn
  Int_t                   fInCentralitySelection; //! centrality bin, only for QA plots
  fitModulationType       fFitModulationType;     // fit modulation type
  qcRecovery              fQCRecovery;            // recovery type for e-by-e qc method
  Bool_t                  fUsePtWeight;           // use dptdphi instead of dndphi
  Bool_t                  fUsePtWeightErrorPropagation; // recalculate the bin error on the dpt dphi histogram
  detectorType            fDetectorType;          // type of detector used for modulation fit
  TString                 fFitModulationOptions;  // fit options for modulation fit
  runModeType             fRunModeType;           // run mode type 
  TF1*                    fFitModulation;         // modulation fit for rho
  Float_t                 fMinPvalue;             // minimum value of p
  Float_t                 fMaxPvalue;             // maximum value of p
  // additional jet cuts (most are inherited)
  Float_t                 fLocalJetMinEta;        // local eta cut for jets
  Float_t                 fLocalJetMaxEta;        // local eta cut for jets
  Float_t                 fLocalJetMinPhi;        // local phi cut for jets
  Float_t                 fLocalJetMaxPhi;        // local phi cut for jets
  Float_t                 fSoftTrackMinPt;        // min pt for soft tracks
  Float_t                 fSoftTrackMaxPt;        // max pt for soft tracks
  // general qa histograms
  TH1F*                   fHistPvalueCDF;         //! cdf value of chisquare p
  TH2F*                   fHistRhoStatusCent;     //! status of rho vs centrality
  // general settings
  Bool_t                  fAbsVnHarmonics;        // force postive local rho
  Float_t                 fExcludeLeadingJetsFromFit;    // exclude n leading jets from fit
  Bool_t                  fRebinSwapHistoOnTheFly;       // rebin swap histo on the fly
  Float_t                 fPercentageOfFits;      // save this percentage of fits
  Bool_t                  fUseV0EventPlaneFromHeader;    // use the vzero event plane from the header
  // transient object pointers
  TList*                  fOutputList;            //! output list
  TList*                  fOutputListGood;        //! output list for local analysis
  TList*                  fOutputListBad;         //! output list for local analysis
  TH1F*                   fHistSwap;              //! swap histogram
  TH1F*                   fHistAnalysisSummary;   //! flags
  TProfile*               fProfV2;                //! extracted v2
  TProfile*               fProfV2Cumulant;        //! v2 cumulant
  TProfile*               fProfV3;                //! extracted v3
  TProfile*               fProfV3Cumulant;        //! v3 cumulant
  TH1F*                   fHistPsi2[10];          //! psi 2
  TH1F*                   fHistPsi3[10];          //! psi 3

  AliAnalysisTaskLocalRho(const AliAnalysisTaskLocalRho&);                  // not implemented
  AliAnalysisTaskLocalRho& operator=(const AliAnalysisTaskLocalRho&);       // not implemented

  ClassDef(AliAnalysisTaskLocalRho, 5);
};
#endif
 AliAnalysisTaskLocalRho.h:1
 AliAnalysisTaskLocalRho.h:2
 AliAnalysisTaskLocalRho.h:3
 AliAnalysisTaskLocalRho.h:4
 AliAnalysisTaskLocalRho.h:5
 AliAnalysisTaskLocalRho.h:6
 AliAnalysisTaskLocalRho.h:7
 AliAnalysisTaskLocalRho.h:8
 AliAnalysisTaskLocalRho.h:9
 AliAnalysisTaskLocalRho.h:10
 AliAnalysisTaskLocalRho.h:11
 AliAnalysisTaskLocalRho.h:12
 AliAnalysisTaskLocalRho.h:13
 AliAnalysisTaskLocalRho.h:14
 AliAnalysisTaskLocalRho.h:15
 AliAnalysisTaskLocalRho.h:16
 AliAnalysisTaskLocalRho.h:17
 AliAnalysisTaskLocalRho.h:18
 AliAnalysisTaskLocalRho.h:19
 AliAnalysisTaskLocalRho.h:20
 AliAnalysisTaskLocalRho.h:21
 AliAnalysisTaskLocalRho.h:22
 AliAnalysisTaskLocalRho.h:23
 AliAnalysisTaskLocalRho.h:24
 AliAnalysisTaskLocalRho.h:25
 AliAnalysisTaskLocalRho.h:26
 AliAnalysisTaskLocalRho.h:27
 AliAnalysisTaskLocalRho.h:28
 AliAnalysisTaskLocalRho.h:29
 AliAnalysisTaskLocalRho.h:30
 AliAnalysisTaskLocalRho.h:31
 AliAnalysisTaskLocalRho.h:32
 AliAnalysisTaskLocalRho.h:33
 AliAnalysisTaskLocalRho.h:34
 AliAnalysisTaskLocalRho.h:35
 AliAnalysisTaskLocalRho.h:36
 AliAnalysisTaskLocalRho.h:37
 AliAnalysisTaskLocalRho.h:38
 AliAnalysisTaskLocalRho.h:39
 AliAnalysisTaskLocalRho.h:40
 AliAnalysisTaskLocalRho.h:41
 AliAnalysisTaskLocalRho.h:42
 AliAnalysisTaskLocalRho.h:43
 AliAnalysisTaskLocalRho.h:44
 AliAnalysisTaskLocalRho.h:45
 AliAnalysisTaskLocalRho.h:46
 AliAnalysisTaskLocalRho.h:47
 AliAnalysisTaskLocalRho.h:48
 AliAnalysisTaskLocalRho.h:49
 AliAnalysisTaskLocalRho.h:50
 AliAnalysisTaskLocalRho.h:51
 AliAnalysisTaskLocalRho.h:52
 AliAnalysisTaskLocalRho.h:53
 AliAnalysisTaskLocalRho.h:54
 AliAnalysisTaskLocalRho.h:55
 AliAnalysisTaskLocalRho.h:56
 AliAnalysisTaskLocalRho.h:57
 AliAnalysisTaskLocalRho.h:58
 AliAnalysisTaskLocalRho.h:59
 AliAnalysisTaskLocalRho.h:60
 AliAnalysisTaskLocalRho.h:61
 AliAnalysisTaskLocalRho.h:62
 AliAnalysisTaskLocalRho.h:63
 AliAnalysisTaskLocalRho.h:64
 AliAnalysisTaskLocalRho.h:65
 AliAnalysisTaskLocalRho.h:66
 AliAnalysisTaskLocalRho.h:67
 AliAnalysisTaskLocalRho.h:68
 AliAnalysisTaskLocalRho.h:69
 AliAnalysisTaskLocalRho.h:70
 AliAnalysisTaskLocalRho.h:71
 AliAnalysisTaskLocalRho.h:72
 AliAnalysisTaskLocalRho.h:73
 AliAnalysisTaskLocalRho.h:74
 AliAnalysisTaskLocalRho.h:75
 AliAnalysisTaskLocalRho.h:76
 AliAnalysisTaskLocalRho.h:77
 AliAnalysisTaskLocalRho.h:78
 AliAnalysisTaskLocalRho.h:79
 AliAnalysisTaskLocalRho.h:80
 AliAnalysisTaskLocalRho.h:81
 AliAnalysisTaskLocalRho.h:82
 AliAnalysisTaskLocalRho.h:83
 AliAnalysisTaskLocalRho.h:84
 AliAnalysisTaskLocalRho.h:85
 AliAnalysisTaskLocalRho.h:86
 AliAnalysisTaskLocalRho.h:87
 AliAnalysisTaskLocalRho.h:88
 AliAnalysisTaskLocalRho.h:89
 AliAnalysisTaskLocalRho.h:90
 AliAnalysisTaskLocalRho.h:91
 AliAnalysisTaskLocalRho.h:92
 AliAnalysisTaskLocalRho.h:93
 AliAnalysisTaskLocalRho.h:94
 AliAnalysisTaskLocalRho.h:95
 AliAnalysisTaskLocalRho.h:96
 AliAnalysisTaskLocalRho.h:97
 AliAnalysisTaskLocalRho.h:98
 AliAnalysisTaskLocalRho.h:99
 AliAnalysisTaskLocalRho.h:100
 AliAnalysisTaskLocalRho.h:101
 AliAnalysisTaskLocalRho.h:102
 AliAnalysisTaskLocalRho.h:103
 AliAnalysisTaskLocalRho.h:104
 AliAnalysisTaskLocalRho.h:105
 AliAnalysisTaskLocalRho.h:106
 AliAnalysisTaskLocalRho.h:107
 AliAnalysisTaskLocalRho.h:108
 AliAnalysisTaskLocalRho.h:109
 AliAnalysisTaskLocalRho.h:110
 AliAnalysisTaskLocalRho.h:111
 AliAnalysisTaskLocalRho.h:112
 AliAnalysisTaskLocalRho.h:113
 AliAnalysisTaskLocalRho.h:114
 AliAnalysisTaskLocalRho.h:115
 AliAnalysisTaskLocalRho.h:116
 AliAnalysisTaskLocalRho.h:117
 AliAnalysisTaskLocalRho.h:118
 AliAnalysisTaskLocalRho.h:119
 AliAnalysisTaskLocalRho.h:120
 AliAnalysisTaskLocalRho.h:121
 AliAnalysisTaskLocalRho.h:122
 AliAnalysisTaskLocalRho.h:123
 AliAnalysisTaskLocalRho.h:124
 AliAnalysisTaskLocalRho.h:125
 AliAnalysisTaskLocalRho.h:126
 AliAnalysisTaskLocalRho.h:127
 AliAnalysisTaskLocalRho.h:128
 AliAnalysisTaskLocalRho.h:129
 AliAnalysisTaskLocalRho.h:130
 AliAnalysisTaskLocalRho.h:131
 AliAnalysisTaskLocalRho.h:132
 AliAnalysisTaskLocalRho.h:133
 AliAnalysisTaskLocalRho.h:134
 AliAnalysisTaskLocalRho.h:135
 AliAnalysisTaskLocalRho.h:136
 AliAnalysisTaskLocalRho.h:137
 AliAnalysisTaskLocalRho.h:138
 AliAnalysisTaskLocalRho.h:139
 AliAnalysisTaskLocalRho.h:140
 AliAnalysisTaskLocalRho.h:141
 AliAnalysisTaskLocalRho.h:142
 AliAnalysisTaskLocalRho.h:143
 AliAnalysisTaskLocalRho.h:144
 AliAnalysisTaskLocalRho.h:145
 AliAnalysisTaskLocalRho.h:146
 AliAnalysisTaskLocalRho.h:147
 AliAnalysisTaskLocalRho.h:148
 AliAnalysisTaskLocalRho.h:149
 AliAnalysisTaskLocalRho.h:150
 AliAnalysisTaskLocalRho.h:151
 AliAnalysisTaskLocalRho.h:152
 AliAnalysisTaskLocalRho.h:153
 AliAnalysisTaskLocalRho.h:154
 AliAnalysisTaskLocalRho.h:155
 AliAnalysisTaskLocalRho.h:156
 AliAnalysisTaskLocalRho.h:157
 AliAnalysisTaskLocalRho.h:158
 AliAnalysisTaskLocalRho.h:159
 AliAnalysisTaskLocalRho.h:160
 AliAnalysisTaskLocalRho.h:161
 AliAnalysisTaskLocalRho.h:162
 AliAnalysisTaskLocalRho.h:163
 AliAnalysisTaskLocalRho.h:164
 AliAnalysisTaskLocalRho.h:165
 AliAnalysisTaskLocalRho.h:166
 AliAnalysisTaskLocalRho.h:167
 AliAnalysisTaskLocalRho.h:168