]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/correlationHF/AliHFCorrelationFDsubtraction.h
Adding description
[u/mrichter/AliRoot.git] / PWGHF / correlationHF / AliHFCorrelationFDsubtraction.h
1 #ifndef ALIHFCORRELATIONFDSUBTRACTION_H
2 #define ALIHFCORRELATIONFDSUBTRACTION_H
3 /* Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id: $ */
7
8 /////////////////////////////////////////////////////////////
9 // Class for subtracting D-hadron correlations from secondary D from B meson decay from D-hadron correlations of prompt+secondary D mesons
10 //       n.b. requires the evaluation of the fraction of prompt D mesons using the same receipes used for D meson cross-section and RAA in D2H pag
11 //       (-> fprompt obtained from FONLL predictions of prompt and secondary D, prompt and secondary D meson efficiencies and a range of hypotheses for
12 //       RAA(DfromB)/RAA(promptD) (if needed)
13 //
14 // Author: A. Rossi, andrea.rossi@cern.ch
15 /////////////////////////////////////////////////////////////
16
17
18 #include <TH1D.h>
19 #include <TGraphAsymmErrors.h>
20 #include <TCanvas.h>
21
22 class AliHFCorrelationFDsubtraction : public TNamed {
23   
24  public:
25
26   AliHFCorrelationFDsubtraction();
27   AliHFCorrelationFDsubtraction(const char* name, const char* title);
28   ~AliHFCorrelationFDsubtraction();
29
30   void SetUncorrectedHistogram(TH1D *h){fhDataUncorrected=(TH1D*)h->Clone();}
31   Bool_t AddTemplateHisto(TH1D *h);
32   void SetDptRange(Double_t ptmin,Double_t ptmax){fptmin=ptmin;fptmax=ptmax;}
33   Bool_t Init();
34   void SetNTemplates(Int_t ntempl){fnTemplates=ntempl;}
35   void SetFpromptGraphFc(TGraphAsymmErrors *gr){fgrConservativeFc=(TGraphAsymmErrors*)gr->Clone();}
36   void SetFpromptGraphNb(TGraphAsymmErrors *gr){fgrConservativeNb=(TGraphAsymmErrors*)gr->Clone();}
37   void SetMethod(Int_t method){fmethod=method;}
38   void SubtractFeedDown(TH1D *hFDTempl);
39   TGraphAsymmErrors* GetUncGraphFromHistos(TH1D *hRef,TH1D *hMin,TH1D *hMax);
40   void CalculateEnvelope();
41   TGraphAsymmErrors* GetGraphEnvelope(){return fgrEnvelope;}
42   TGraphAsymmErrors* GetGraphEnvelopeRatio(){return fgrEnvelopeRatio;}
43   TH1D* GetHistoEnvelopeRatioMin(){return fhEnvelopeMinRatio;}
44   TH1D* GetHistoEnvelopeRatioMax(){return fhEnvelopeMaxRatio;}
45   TH1D* GetHistoRelSystUncMin();
46   TH1D* GetHistoRelSystUncMax();
47   TH1D* GetHistoEnvelopeMin(){return fhEnvelopeMin;}
48   TH1D* GetHistoEnvelopeMax(){return fhEnvelopeMax;}
49   TH1D* GetCentralSubtractedPlot(){return fhDataCorrected[0];}
50   TH1D* GetTemplate(Int_t i){       
51     if (i>=fLastTemplAdded){Printf("Get Template: Error"); return 0;}
52     else return fhTemplates[i];
53   }
54   
55  private:
56
57   Double_t fptmin;                                      //  min pt (D meson pt range)
58   Double_t fptmax;                                      //  max pt (D meson pt range)
59   TH1D* fhDataUncorrected;                              // Input correlation histogram
60   Int_t fmethod;                                        // method: 1= fc, 2= Nb, 3= both
61   TGraphAsymmErrors* fgrConservativeFc;                 // fc graph, fc method
62   TGraphAsymmErrors* fgrConservativeNb;                 // fc graph, Nb method
63   Int_t fnTemplates;                                   // maximum number of templates that will be included
64   Int_t fLastTemplAdded;                                // counter of templates included
65   TH1D **fhTemplates;                                  //  Array with template histo, size fnTemplates
66   TH1D **fhDataCorrected;                               // Array with corrected histo, size 3*fnTemplates
67   TH1D **fhRatioSameTemplate;                            // Array with ratio of histo based on the same template
68   TH1D **fhRatioAllToCentralTempl;                                // array with ratio of all histos with respect to the central hypo
69   TCanvas **fcCompare;                                            // array of canvases with comparison (template by template)
70   TCanvas **fcRatio;                                          // array of canvases with ratio (templ by templ)
71   Int_t fCountTempl;                                               // internal counter
72   TCanvas *fcAllRatio;                                         // canva with all ratios with respect to central hypo
73   TH1D* fhEnvelopeMax;                                         // envelope with max variation
74   TH1D* fhEnvelopeMin;                                          // envelope with min variation
75   TH1D* fhEnvelopeMaxRatio;                                     // envelope with max relative variation
76   TH1D* fhEnvelopeMinRatio;                                      // envelope with min relative variation
77   TGraphAsymmErrors* fgrEnvelope;                                // graph with envelope
78   TGraphAsymmErrors* fgrEnvelopeRatio;                           // graph with envelope of ratios
79   
80   ClassDef(AliHFCorrelationFDsubtraction,1);
81   
82 };
83 #endif