AliPhysics  b0c77bb (b0c77bb)
HFPtSpectrumRaa.C
Go to the documentation of this file.
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include "TFile.h"
3 #include "TH1.h"
4 #include "TH1D.h"
5 #include "TH2.h"
6 #include "TH2D.h"
7 #include "TH3.h"
8 #include "TH3D.h"
9 #include "TNtuple.h"
10 #include "TGraphAsymmErrors.h"
11 #include "TMath.h"
12 #include "TCanvas.h"
13 #include "TLegend.h"
14 #include "TROOT.h"
15 #include "TStyle.h"
16 #include "TLine.h"
17 #include "TLatex.h"
18 
19 #include "AliHFSystErr.h"
20 #include <Riostream.h>
21 #endif
22 
23 /* $Id$ */
25 //
26 // Macro to compute the Raa, taking as inputs the output of the corrected yields
27 // and the pp reference
28 //
29 // R_AB = [ ( dsigma/dpt )_AB / sigma_AB ] / <TAB> * ( dsigma/dpt )_pp
30 //
31 //
32 // Parameters:
33 // 1. ppfile = reference pp in the same pt binning
34 // 2. ABfile = corrected AB yields
35 // 3. outfile = output file name
36 // 4. decay = decay as in HFSystErr class
37 // 5. sigmaABCINT1B = cross section for normalization (**USE THE SAME AS ON 2.**)
38 // 6. fdMethod = feed-down subtraction method (kNb, kfc)
39 // 7. cc = centrality class
40 // 8. Energy = colliding energy (k276,k55)
41 // 9. MinHypo = minimum energy loss hypothesis (Default 1./3.)
42 // 10. MaxHypo = maximum energy loss hypothesis (Default 3.0)
43 // 11. MaxRb = maximum Raa(b) hypothesis (Default 6.0, won't do anything)
44 // 12. RbHypo : flag to decide whether the Eloss hypothesis is Rb or Rb/Rc
45 // 13. CentralHypo = central energy loss hypothesis, DEFAULT TO 1.0
46 // 14. isRaavsEP = flag to compute the Raa IN/OUT of plane, divides the reference by 2.0
47 // 15. ScaledAndExtrapRef: flag to tag scaled+reference pp-scaled-data
48 //
49 // Complains to : Zaida Conesa del Valle
50 //
52 
55 enum energy{ k276, k5dot023, k55 };
60 
61 
62 Bool_t printout = false;
67 
68 //____________________________________________________________
69 Bool_t PbPbDataSyst(AliHFSystErr *syst, Double_t pt, Int_t cc, Double_t &dataSystUp, Double_t &dataSystDown);
70 
71 //____________________________________________________________
73  // total^2 = data^2 + fd^2
74  Double_t data2 = total*total - fd*fd ;
75  return TMath::Sqrt( data2 );
76 }
77 
78 //____________________________________________________________
80 {
81  Int_t istart =0;
82  Int_t npoints = gr->GetN();
83  for(Int_t i=0; i<=npoints; i++){
84  Double_t x=0.,y=0.;
85  gr->GetPoint(i,x,y);
86  if ( TMath::Abs ( x - pt ) < 0.4 ) {
87  istart = i;
88  break;
89  }
90  }
91  return istart;
92 }
93 
94 
95 //
96 //
97 // R_AB = [ ( dsigma/dpt )_AB / sigma_AB ] / <TAB> * ( dsigma/dpt )_pp
98 //
99 //
100 //____________________________________________________________
101 void HFPtSpectrumRaa(const char *ppfile="HFPtSpectrum_D0Kpi_method2_rebinnedth_230311_newsigma.root",
102  const char *ABfile="HFPtSpectrum_D0Kpi_PbPbcuts_method2_rebinnedth_230311_newsigma.root",
103  const char *outfile="HFPtSpectrumRaa.root",
104  Int_t decay=1,
105  Double_t sigmaABCINT1B=54.e9,
106  Int_t fdMethod = kNb, Int_t cc=kpp, Int_t Energy=k276,
107  Double_t MinHypo=1./3., Double_t MaxHypo=3.0, Double_t MaxRb=6.0,
108  Bool_t isRbHypo=false, Double_t CentralHypo = 1.0,
109  Int_t ccestimator = kV0M,
110  Bool_t isUseTaaForRaa=true, const char *shadRbcFile="", Int_t nSigmaShad=3.0,
111  Int_t isRaavsEP=kPhiIntegrated, Bool_t isScaledAndExtrapRef=kFALSE,
112  Int_t rapiditySlice=kdefault, Int_t analysisSpeciality=kTopological)
113 {
114 
115  gROOT->Macro("$ALICE_PHYSICS/PWGHF/vertexingHF/macros/LoadLibraries.C");
116 
117  //
118  // Defining the TAB values for the given centrality class
119  //
120  Double_t Tab = 1., TabSyst = 0., A=207.2, B=207.2;
121  if ( Energy!=k276 && Energy!=k5dot023) {
122  printf("\n The Tab values for this cms energy have not yet been implemented, please do it ! \n");
123  return;
124  }
125  if ( cc == kpp ){
126  Tab = 1.;
127  }
128  // Values from Alberica's twiki:
129  // https://twiki.cern.ch/twiki/bin/viewauth/ALICE/CentStudies
130  if( (ccestimator == kV0M) && (Energy==k276) ) {
131  if ( cc == k07half ) {
132  Tab = 24.81; TabSyst = 0.8037;
133  } else if ( cc == k010 ) {
134  Tab = 23.48; TabSyst = 0.97;
135  } else if ( cc == k1020 ) {
136  Tab = 14.4318; TabSyst = 0.5733;
137  } else if ( cc == k020 ) {
138  Tab = 18.93; TabSyst = 0.74;
139  } else if ( cc == k1030 ) {
140  Tab = 11.4; TabSyst = 0.36;
141  } else if ( cc == k2040 ) {
142  Tab = 6.86; TabSyst = 0.28;
143  } else if ( cc == k2030 ) {
144  Tab = 8.73769; TabSyst = 0.370219;
145  } else if ( cc == k3040 ) {
146  Tab = 5.02755; TabSyst = 0.22099;
147  } else if ( cc == k4050 ) {
148  Tab = 2.68327; TabSyst = 0.137073;
149  } else if ( cc == k3050 ) {
150  Tab = 3.87011; TabSyst = 0.183847;
151  } else if ( cc == k4060 ) {
152  Tab = 2.00; TabSyst= 0.11;
153  } else if ( cc == k4080 ) {
154  Tab = 1.20451; TabSyst = 0.071843;
155  } else if ( cc == k5060 ) {
156  Tab = 1.32884; TabSyst = 0.0929536;
157  } else if ( cc == k6080 ) {
158  Tab = 0.419; TabSyst = 0.033;
159  } else if ( cc == k5080 ) {
160  Tab = 0.719; TabSyst = 0.054;
161  } else if ( cc == k80100 ){
162  Tab = 0.0690; TabSyst = 0.0062;
163  }
164  }
165  if( (ccestimator == kV0M) && (Energy==k5dot023) ) {
166  if ( cc == k010 ) {
167  Tab = 23.07; TabSyst = 0.44;
168  } else if ( cc == k3050 ) {
169  Tab = 3.897; TabSyst = 0.11;
170  } else if ( cc == k6080 ) {
171  Tab = 0.4173; TabSyst = 0.014;
172  }
173  }
174 
175 
176  // pPb Glauber (A. Toia)
177  // https://twiki.cern.ch/twiki/bin/viewauth/ALICE/PACentStudies#Glauber_Calculations_with_sigma
178  if( cc == kpPb0100 ){
179  Tab = 0.098334; TabSyst = 0.0070679;
180  A=207.2; B=1.;
181  }
182  else if( ccestimator == kV0A ){
183  if ( cc == kpPb020 ) {
184  Tab = 0.183; TabSyst = 0.006245;
185  } else if ( cc == kpPb2040 ) {
186  Tab = 0.134; TabSyst = 0.004899;
187  } else if ( cc == kpPb4060 ) {
188  Tab = 0.092; TabSyst = 0.004796;
189  } else if ( cc == kpPb60100 ) {
190  Tab = 0.041; TabSyst = 0.008832;
191  }
192  }
193  else if( ccestimator == kZNA ){
194  if ( cc == kpPb010 ) {
195  Tab = 0.1812; TabSyst = 0.01413;
196  } else if ( cc == kpPb020 ) {
197  Tab = 0.1742; TabSyst = 0.00992;
198  } else if ( cc == kpPb1020 ) {
199  Tab = 0.1672; TabSyst = 0.005852;
200  } else if ( cc == kpPb2040 ) {
201  Tab = 0.1453; TabSyst = 0.002615;
202  } else if ( cc == kpPb4060 ) {
203  Tab = 0.1074; TabSyst = 0.004726;
204  } else if ( cc == kpPb60100 ) {
205  Tab = 0.0486; TabSyst = 0.002430;
206  }
207  }
208  else if( ccestimator == kCL1 ){
209  if ( cc == kpPb020 ) {
210  Tab = 0.19; TabSyst = 0.007;
211  } else if ( cc == kpPb2040 ) {
212  Tab = 0.136; TabSyst = 0.005;
213  } else if ( cc == kpPb4060 ) {
214  Tab = 0.088; TabSyst = 0.005;
215  } else if ( cc == kpPb60100 ) {
216  Tab = 0.0369; TabSyst = 0.0085;
217  }
218  }
219 
220  //
221  // Reading the pp file
222  //
223  TFile * ppf = new TFile(ppfile,"read");
224  TH1D * hSigmaPP;
225  TGraphAsymmErrors * gSigmaPPSyst;
226  TGraphAsymmErrors * gSigmaPPSystData = (TGraphAsymmErrors*)ppf->Get("gScaledDataSystData");
227  TGraphAsymmErrors * gSigmaPPSystTheory = (TGraphAsymmErrors*)ppf->Get("gScaledDataSystExtrap");
228  TGraphAsymmErrors * gSigmaPPSystFeedDown = (TGraphAsymmErrors*)ppf->Get("gScaledDataSystFeedDown");
229  TH1I * hCombinedReferenceFlag=0x0;
230  TGraphAsymmErrors * gReferenceFdSyst=0x0;
231  if(isScaledAndExtrapRef){
232  hCombinedReferenceFlag = (TH1I*)ppf->Get("hCombinedReferenceFlag");
233  hSigmaPP = (TH1D*)ppf->Get("hReference");
234  gSigmaPPSyst = (TGraphAsymmErrors*)ppf->Get("gReferenceSyst");
235  gReferenceFdSyst = (TGraphAsymmErrors*)ppf->Get("gReferenceFdSyst");
236  } else {
237  hSigmaPP = (TH1D*)ppf->Get("fhScaledData");
238  gSigmaPPSyst = (TGraphAsymmErrors*)ppf->Get("gScaledData");
239  }
240  Double_t scalePPRefToMatchRapidityBin = 1.0;
241 
242 
243  // Call the systematics uncertainty class for a given decay
244  AliHFSystErr *systematicsPP = 0x0;
245  if(ppf->Get("AliHFSystErr")){
246  systematicsPP = (AliHFSystErr*)ppf->Get("AliHFSystErr");
247  printf(" --> AliHFSystErr object for pp reference read from file. Object title = %s\n",systematicsPP->GetTitle());
248  }else{
249  systematicsPP = new AliHFSystErr();
250  if(analysisSpeciality==kLowPt){
251  systematicsPP->SetIsLowPtAnalysis(true);
252  }
253  if(analysisSpeciality==kPP7TeVPass4){
254  systematicsPP->SetIsPass4Analysis(true);
255  }
256  systematicsPP->Init(decay);
257  printf(" --> AliHFSystErr object for pp reference created based on macro arguments. Object title = %s\n",systematicsPP->GetTitle());
258  }
259  //
260  // Reading the AB collisions file
261  //
262  TFile * ABf = new TFile(ABfile,"read");
263  TH1D *hSigmaAB = (TH1D*)ABf->Get("histoSigmaCorr");
264  // TH2D *hSigmaABRcb = (TH2D*)ABf->Get("histoSigmaCorrRcb");
265  // TGraphAsymmErrors * gSigmaABSyst = (TGraphAsymmErrors*)ABf->Get("gSigmaCorr");
266  TGraphAsymmErrors * gSigmaABSystFeedDown = (TGraphAsymmErrors*)ABf->Get("gSigmaCorrConservative");
267  TNtuple * nSigmaAB = (TNtuple*)ABf->Get("fnSigma");
268  //
269  TH1D *hMassAB = (TH1D*)ABf->Get("hRECpt");
270  TH1D *hDirectEffptAB = (TH1D*)ABf->Get("hDirectEffpt");
271  TH1D *histofcAB = (TH1D*)ABf->Get("histofc");
272  //
273  TH1D* fhStatUncEffcSigmaAB = (TH1D*)ABf->Get("fhStatUncEffcSigma");
274  TH1D* fhStatUncEffbSigmaAB = (TH1D*)ABf->Get("fhStatUncEffbSigma");
275  TH1D* fhStatUncEffcFDAB = (TH1D*)ABf->Get("fhStatUncEffcFD");
276  TH1D* fhStatUncEffbFDAB = (TH1D*)ABf->Get("fhStatUncEffbFD");
277  //
278  TH1D* fhStatUncEffcSigmaAB_Raa = (TH1D*)fhStatUncEffcSigmaAB->Clone("fhStatUncEffcSigmaAB_Raa");
279  TH1D* fhStatUncEffbSigmaAB_Raa = (TH1D*)fhStatUncEffbSigmaAB->Clone("fhStatUncEffbSigmaAB_Raa");
280  TH1D* fhStatUncEffcFDAB_Raa = (TH1D*)fhStatUncEffcFDAB->Clone("fhStatUncEffcFDAB_Raa");
281  TH1D* fhStatUncEffbFDAB_Raa = (TH1D*)fhStatUncEffbFDAB->Clone("fhStatUncEffbFDAB_Raa");
282  fhStatUncEffcSigmaAB_Raa->Reset();
283  fhStatUncEffbSigmaAB_Raa->Reset();
284  fhStatUncEffcFDAB_Raa->Reset();
285  fhStatUncEffbFDAB_Raa->Reset();
286  fhStatUncEffcSigmaAB_Raa->SetName("fhStatUncEffcSigmaAB_Raa");
287  fhStatUncEffbSigmaAB_Raa->SetName("fhStatUncEffbSigmaAB_Raa");
288  fhStatUncEffcFDAB_Raa->SetName("fhStatUncEffcFDAB_Raa");
289  fhStatUncEffbFDAB_Raa->SetName("fhStatUncEffvFDAB_Raa");
290 
291 
292  //
293  // Call the systematics uncertainty class for a given decay
294  AliHFSystErr *systematicsAB = 0x0;
295  if(ABf->Get("AliHFSystErr")){
296  systematicsAB=(AliHFSystErr*)ABf->Get("AliHFSystErr");
297  printf(" --> AliHFSystErr object for A-A read from HFPtSpectrum file. Object title = %s\n",systematicsAB->GetTitle());
298  }else{
299  systematicsAB = new AliHFSystErr();
300  systematicsAB->SetCollisionType(1);
301  systematicsAB->SetRunNumber(2016);// check this
302  if(Energy==k276){
303  systematicsAB->SetRunNumber(11);
304  if ( cc == k07half ) systematicsAB->SetCentrality("07half");
305  else if ( cc == k010 ) systematicsAB->SetCentrality("010");
306  else if ( cc == k1020 ) systematicsAB->SetCentrality("1020");
307  else if ( cc == k020 ) systematicsAB->SetCentrality("020");
308  else if ( cc == k2040 || cc == k2030 || cc == k3040 ) {
309  systematicsAB->SetCentrality("2040");
310  systematicsAB->SetIsPbPb2010EnergyScan(true);
311  }
312  else if ( cc == k4060 || cc == k4050 || cc == k5060 ) systematicsAB->SetCentrality("4060");
313  else if ( cc == k6080 || cc == k5080 ) systematicsAB->SetCentrality("6080");
314  else if ( cc == k4080 ) systematicsAB->SetCentrality("4080");
315  else if ( cc == k3050 ) {
316  if (isRaavsEP == kPhiIntegrated) systematicsAB->SetCentrality("4080");
317  else if (isRaavsEP == kInPlane) systematicsAB->SetCentrality("3050InPlane");
318  else if (isRaavsEP == kOutOfPlane) systematicsAB->SetCentrality("3050OutOfPlane");
319  }
320  } else if (Energy==k5dot023){
321  systematicsAB->SetRunNumber(15);
322  if ( cc == k010 ){
323  systematicsAB->SetCentrality("010");
324  } else if ( cc == k1030 ) {
325  systematicsAB->SetCentrality("3050"); //no systematics available for 10--30
326  } else if ( cc == k3050 ) {
327  systematicsAB->SetCentrality("3050");
328  } else if ( cc == k6080 ) {
329  systematicsAB->SetCentrality("6080");
330  }
331  }
332  //
333  else if ( cc == kpPb0100 || cc == kpPb010 || cc == kpPb020 || cc == kpPb1020 || cc == kpPb2040 || cc == kpPb4060 || cc == kpPb60100 ) {
334  systematicsAB->SetCollisionType(2);
335  // Rapidity slices
336  if(rapiditySlice!=kdefault){
337  systematicsAB->SetIspPb2011RapidityScan(true);
338  TString rapidity="";
339  switch(rapiditySlice) {
340  case k08to04: rapidity="0804"; scalePPRefToMatchRapidityBin=(0.093+0.280)/1.0; break;
341  case k07to04: rapidity="0804"; scalePPRefToMatchRapidityBin=0.280/1.0; break;
342  case k04to01: rapidity="0401"; scalePPRefToMatchRapidityBin=0.284/1.0; break;
343  case k01to01: rapidity="0101"; scalePPRefToMatchRapidityBin=0.191/1.0; break;
344  case k01to04: rapidity="0104"; scalePPRefToMatchRapidityBin=0.288/1.0; break;
345  case k04to07: rapidity="0408"; scalePPRefToMatchRapidityBin=0.288/1.0; break;
346  case k04to08: rapidity="0408"; scalePPRefToMatchRapidityBin=(0.288+0.096)/1.0; break;
347  case k01to05: rapidity="0401"; scalePPRefToMatchRapidityBin=0.4; break;
348  }
349  systematicsAB->SetRapidity(rapidity);
350  }
351  // Centrality slices
352  if(ccestimator==kV0A) {
353  if(cc == kpPb020) systematicsAB->SetCentrality("020V0A");
354  else if(cc == kpPb2040) systematicsAB->SetCentrality("2040V0A");
355  else if(cc == kpPb4060) systematicsAB->SetCentrality("4060V0A");
356  else if(cc == kpPb60100) systematicsAB->SetCentrality("60100V0A");
357  } else if (ccestimator==kZNA) {
358  if(cc == kpPb010) {systematicsAB->SetRunNumber(16); systematicsAB->SetCentrality("010ZNA");}
359  else if(cc == kpPb020) systematicsAB->SetCentrality("020ZNA");
360  else if(cc == kpPb1020) {systematicsAB->SetRunNumber(16); systematicsAB->SetCentrality("1020ZNA");}
361  else if(cc == kpPb2040) {systematicsAB->SetRunNumber(16); systematicsAB->SetCentrality("2040ZNA");}
362  else if(cc == kpPb4060) {systematicsAB->SetRunNumber(16); systematicsAB->SetCentrality("4060ZNA");}
363  else if(cc == kpPb60100) {systematicsAB->SetRunNumber(16); systematicsAB->SetCentrality("60100ZNA");}
364  } else if (ccestimator==kCL1) {
365  if(cc == kpPb020) systematicsAB->SetCentrality("020CL1");
366  else if(cc == kpPb2040) systematicsAB->SetCentrality("2040CL1");
367  else if(cc == kpPb4060) systematicsAB->SetCentrality("4060CL1");
368  else if(cc == kpPb60100) systematicsAB->SetCentrality("60100CL1");
369  }else {
370  if(!(cc == kpPb0100)) {
371  cout <<" Error on the pPb options"<<endl;
372  return;
373  }
374  }
375  }
376  else {
377  cout << " Systematics not yet implemented " << endl;
378  return;
379  }
380  if(analysisSpeciality==kLowPt){
381  systematicsAB->SetIsLowPtAnalysis(true);
382  }
383  else if(analysisSpeciality==kBDT){
384  systematicsAB->SetIsBDTAnalysis(true);
385  }
386  //
387  systematicsAB->Init(decay);
388  printf(" --> AliHFSystErr object for A-A created based on macro arguments. Object title = %s\n",systematicsAB->GetTitle());
389  }
390  //
391  Int_t entries = nSigmaAB->GetEntries();
392  Float_t pt=0., signal=0., Rb=0., Rcb=0., fcAB=0., yieldAB=0., sigmaAB=0., statUncSigmaAB=0., sigmaABMin=0.,sigmaABMax=0.;
393  nSigmaAB->SetBranchAddress("pt",&pt);
394  nSigmaAB->SetBranchAddress("Signal",&signal);
395  if (fdMethod==kNb) nSigmaAB->SetBranchAddress("Rb",&Rb);
396  else if (fdMethod==kfc) nSigmaAB->SetBranchAddress("Rcb",&Rcb);
397  nSigmaAB->SetBranchAddress("fc",&fcAB);
398  nSigmaAB->SetBranchAddress("Yield",&yieldAB);
399  nSigmaAB->SetBranchAddress("Sigma",&sigmaAB);
400  nSigmaAB->SetBranchAddress("SigmaStatUnc",&statUncSigmaAB);
401  nSigmaAB->SetBranchAddress("SigmaMax",&sigmaABMax);
402  nSigmaAB->SetBranchAddress("SigmaMin",&sigmaABMin);
403 
404 
405  // define the binning
406  Int_t nbins = hSigmaAB->GetNbinsX();
407  Double_t binwidth = hSigmaAB->GetBinWidth(1);
408  Double_t *limits = new Double_t[nbins+1];
409  Double_t *binwidths = new Double_t[nbins];
410  Double_t xlow=0.;
411  for (Int_t i=1; i<=nbins; i++) {
412  binwidth = hSigmaAB->GetBinWidth(i);
413  xlow = hSigmaAB->GetBinLowEdge(i);
414  limits[i-1] = xlow;
415  binwidths[i-1] = binwidth;
416  }
417  limits[nbins] = xlow + binwidth;
418 
419 
420  //
421  // Read the shadowing file if given as input
422  //
423  Double_t centralRbcShad[nbins+1], minRbcShad[nbins+1], maxRbcShad[nbins+1];
424  for(Int_t i=0; i<=nbins; i++) { centralRbcShad[i]=1.0; minRbcShad[i]=6.0; maxRbcShad[i]=0.0; }
425  Bool_t isShadHypothesis = false;
426  if( strcmp(shadRbcFile,"")!=0 ) {
427  isShadHypothesis = true;
428  cout<<endl<<">> Beware, using the shadowing prediction file with an "<<nSigmaShad<<"*sigma <<"<<endl<<endl;
429  TFile *fshad = new TFile(shadRbcFile,"read");
430  if(!fshad){ cout <<" >> Shadowing file not properly opened!!!"<<endl<<endl; return;}
431  // TH1D *hRbcShadCentral = (TH1D*)fshad->Get("hDfromBoverPromptD_Shadowing_central");
432  // TH1D *hRbcShadMin = (TH1D*)fshad->Get("hDfromBoverPromptD_Shadowing_upper");
433  // TH1D *hRbcShadMax = (TH1D*)fshad->Get("hDfromBoverPromptD_Shadowing_lower");
434  TH1D *hRbcShadCentral = (TH1D*)fshad->Get("hDfromBoverDfromc_L0");
435  TH1D *hRbcShadMin = (TH1D*)fshad->Get("hDfromBoverDfromc_L0");
436  TH1D *hRbcShadMax = (TH1D*)fshad->Get("hDfromBoverDfromc_L1");
437  if(!hRbcShadCentral || !hRbcShadMin || !hRbcShadMax) {
438  cout<< endl <<">> Shadowing input histograms are not ok !! "<<endl<<endl;
439  return;
440  }
441  // nSigmaShad
442  // nSigmaShad
443  for(Int_t i=1; i<=nbins; i++) {
444  Double_t xpt = hSigmaAB->GetBinCenter(i);
445  if(xpt>24) xpt = 20;
446  centralRbcShad[i] = hRbcShadCentral->GetBinContent( hRbcShadCentral->FindBin(xpt) );
447  Double_t minValue0 = hRbcShadMin->GetBinContent( hRbcShadMin->FindBin(xpt) );
448  Double_t maxValue0 = hRbcShadMax->GetBinContent( hRbcShadMax->FindBin(xpt) );
449  Double_t arrayEl[3] = {minValue0,maxValue0, centralRbcShad[i]};
450  Double_t minValue = TMath::MinElement(3,arrayEl);
451  Double_t maxValue = TMath::MaxElement(3,arrayEl);
452  cout<<">> Shadowing pt="<<xpt<<" central="<<centralRbcShad[i]<<" min="<<minValue<<" max="<<maxValue<<endl;
453  if(minValue>centralRbcShad[i]){ minValue = centralRbcShad[i]; }
454  if(maxValue<centralRbcShad[i]){ maxValue = centralRbcShad[i]; }
455  minRbcShad[i] = centralRbcShad[i] - nSigmaShad*(centralRbcShad[i] - minValue);
456  maxRbcShad[i] = centralRbcShad[i] + nSigmaShad*(maxValue - centralRbcShad[i]);
457  cout<<">> Shadowing hypothesis pt="<<xpt<<" central="<<centralRbcShad[i]<<" min="<<minRbcShad[i]<<" max="<<maxRbcShad[i]<<endl;
458  }
459  }
460 
461  //
462  // define the bins correspondence bw histos/files/graphs
463  //
464  //
465  TH2D * hRABvsRcb = new TH2D("hRABvsRcb"," R_{AB}(c) vs Rcb Eloss hypothesis; p_{T} [GeV/c] ; R_{AB}(c) ; Rcb Eloss hypothesis ",nbins,limits,800,0.,4.);
466  TH2D * hRABvsRb = new TH2D("hRABvsRb"," R_{AB}(c) vs Rb Eloss hypothesis; p_{T} [GeV/c] ; R_{AB}(c) ; Rb Eloss hypothesis ",nbins,limits,800,0.,4.);
467  // TH2D * hRABBeautyvsRCharm = new TH2D("hRABBeautyvsRCharm"," R_{AB}(c) vs Rb Eloss hypothesis; p_{T} [GeV/c] ; R_{AB}(b) ; R_{AB}(c) ",nbins,limits,800,0.,4.);
468  Int_t nbinsHypo=800;//200;
469  Double_t *limitsHypo = new Double_t[nbinsHypo+1];
470  for(Int_t i=1; i<=nbinsHypo+1; i++) limitsHypo[i-1]= i*4./800.;
471  TH3D * hRABCharmVsRBeautyVsPt = new TH3D("hRABCharmVsRBeautyVsPt"," R_{AB}(c) vs Rb vs p_{T} Eloss hypothesis; p_{T} [GeV/c] ; R_{AB}(b) ; R_{AB}(c) ",nbins,limits,nbinsHypo,limitsHypo,nbinsHypo,limitsHypo);
472  TH2D *hRCharmVsRBeauty[nbins+1];
473  for(Int_t i=0; i<=nbins; i++) hRCharmVsRBeauty[i] = new TH2D(Form("hRCharmVsRBeauty_%i",i),Form("RAB(c) vs RAB(b) for pt bin %i ; R_{AB}(b) ; R_{AB}(c)",i),nbinsHypo,limitsHypo,nbinsHypo,limitsHypo);
474  TH2D *hRCharmVsElossHypo[nbins+1];
475  for(Int_t i=0; i<=nbins; i++) hRCharmVsElossHypo[i] = new TH2D(Form("hRCharmVsElossHypo_%i",i),Form("RAB(c) vs ElossHypo for pt bin %i ; Eloss Hypothesis (c/b) ; R_{AB}(c)",i),nbinsHypo,limitsHypo,nbinsHypo,limitsHypo);
476  //
477  TH1D *hRABEloss00= new TH1D("hRABEloss00","hRABEloss00",nbins,limits);
478  TH1D *hRABEloss05= new TH1D("hRABEloss05","hRABEloss05",nbins,limits);
479  TH1D *hRABEloss10= new TH1D("hRABEloss10","hRABEloss10",nbins,limits);
480  TH1D *hRABEloss15= new TH1D("hRABEloss15","hRABEloss15",nbins,limits);
481  TH1D *hRABEloss20= new TH1D("hRABEloss20","hRABEloss20",nbins,limits);
482  //
483  TH2D * hRABvsRbFDlow = new TH2D("hRABvsRbFDlow"," R_{AB}(c) vs Rb Eloss hypothesis (FD low); p_{T} [GeV/c] ; Rb Eloss hypothesis ; R_{AB}(c) ",nbins,limits,800,0.,4.);
484  TH2D * hRABvsRbFDhigh = new TH2D("hRABvsRbFDhigh"," R_{AB}(c) vs Rb Eloss hypothesis (FD high); p_{T} [GeV/c] ; Rb Eloss hypothesis ; R_{AB}(c) ",nbins,limits,800,0.,4.);
485  //
486  TH1D * hRABvsRbFDhigh_proj = new TH1D("hRABvsRbFDhigh_proj","hRABvsRbFDhigh_proj",nbins,limits);
487  TH1D * hRABvsRbFDlow_proj = new TH1D("hRABvsRbFDlow_proj","hRABvsRbFDlow_proj",nbins,limits);
488  //
489  TNtuple *ntupleRAB=0x0 ;
490  if (fdMethod==kNb) {
491  ntupleRAB = new TNtuple("ntupleRAB","ntupleRAB (Nb)","pt:TAB:sigmaPP:sigmaAB:invyieldAB:invyieldABFDHigh:invyieldABFDLow:RABCharm:RABCharmFDHigh:RABCharmFDLow:RABBeauty:fc",100000);
492  } else if (fdMethod==kfc) {
493  ntupleRAB = new TNtuple("ntupleRAB","ntupleRAB (fc)","pt:TAB:sigmaPP:sigmaAB:invyieldAB:invyieldABFDHigh:invyieldABFDLow:Rcb:RABCharm:RABCharmFDHigh:RABCharmFDLow:RABBeauty:RABBeautyFDHigh:RABBeautyFDLow:fc",100000);
494  }
495  if(!ntupleRAB) printf("ERROR: Wrong method option");
496 
497  TH1D * hYieldABvsPt = new TH1D("hYieldABvsPt"," Yield_{AB}(c) vs p_{T} (no Eloss hypothesis); p_{T} [GeV/c] ; Yield_{charm} ",nbins,limits);
498  TH1D * hRABvsPt = new TH1D("hRABvsPt"," R_{AB}(c) vs p_{T} (no Eloss hypothesis); p_{T} [GeV/c] ; R_{charm} ",nbins,limits);
499  TH1D * hRABvsPt_DataSystematics = new TH1D("hRABvsPt_DataSystematics"," Systematics of R_{AB} (c) vs p_{T} (no Eloss hypothesis); p_{T} [GeV/c] ; R_{charm} ",nbins,limits);
500  TGraphAsymmErrors *gRAB_ElossHypothesis = new TGraphAsymmErrors(nbins+1);
501  gRAB_ElossHypothesis->SetNameTitle("gRAB_ElossHypothesis","RAB Eloss systematics");
502  TGraphAsymmErrors *gRAB_FeedDownSystematics = new TGraphAsymmErrors(nbins+1);
503  gRAB_FeedDownSystematics->SetNameTitle("gRAB_FeedDownSystematics","RAB Feed-Down systematics");
504  TGraphAsymmErrors *gRAB_fcFeedDownOnly = new TGraphAsymmErrors(nbins+1);
505  gRAB_fcFeedDownOnly->SetNameTitle("gRAB_fcFeedDownOnly","RAB fc Feed-Down Only");
506  TGraphAsymmErrors *gRAB_FeedDownSystematicsElossHypothesis = new TGraphAsymmErrors(nbins+1);
507  gRAB_FeedDownSystematicsElossHypothesis->SetNameTitle("gRAB_FeedDownSystematicsElossHypothesis","RAB Feed-Down systematics considering Eloss hypothesis");
508  TGraphAsymmErrors *gRAB_DataSystematics = new TGraphAsymmErrors(nbins+1);
509  gRAB_DataSystematics->SetNameTitle("gRAB_DataSystematics","RAB Measurement (no FD, no Eloss) systematics");
510  TGraphAsymmErrors *gRAB_DataSystematicsPP = new TGraphAsymmErrors(nbins+1);
511  gRAB_DataSystematicsPP->SetNameTitle("gRAB_DataSystematicsPP","RAB Measurement PP meas. systematics (data+scaling)");
512  TGraphAsymmErrors *gRAB_DataSystematicsAB = new TGraphAsymmErrors(nbins+1);
513  gRAB_DataSystematicsAB->SetNameTitle("gRAB_DataSystematicsAB","RAB Measurement AB (no FD, no Eloss, no PP data) systematics");
514  TGraphAsymmErrors *gRAB_GlobalSystematics = new TGraphAsymmErrors(nbins+1);
515  gRAB_GlobalSystematics->SetNameTitle("gRAB_GlobalSystematics","RAB Measurement global (data, FD, Eloss) systematics");
516  Double_t ElossMax[nbins+1], ElossMin[nbins+1];
517  for(Int_t i=0; i<=nbins; i++) { ElossMax[i]=0.; ElossMin[i]=6.; }
518  Double_t fcElossMax[nbins+1], fcElossMin[nbins+1];
519  for(Int_t i=0; i<=nbins; i++) { fcElossMax[i]=0.; fcElossMin[i]=6.; }
520  Double_t FDElossMax[nbins+1], FDElossMin[nbins+1];
521  for(Int_t i=0; i<=nbins; i++) { FDElossMax[i]=0.; FDElossMin[i]=6.; }
522 
523  TGraphAsymmErrors *gRAB_Norm = new TGraphAsymmErrors(1);
524  gRAB_Norm->SetNameTitle("gRAB_Norm","RAB Normalization systematics (pp norm + Tab)");
525  Double_t normUnc = TMath::Sqrt ( NormPPUnc*NormPPUnc + (TabSyst/Tab)*(TabSyst/Tab) );
526  if(!isUseTaaForRaa) normUnc = TMath::Sqrt ( NormPPUnc*NormPPUnc + NormABUnc*NormABUnc );
527  gRAB_Norm->SetPoint(1,0.5,1.);
528  gRAB_Norm->SetPointError(1,0.25,0.25,normUnc,normUnc);
529 
530  //
531  // R_AB = ( dN/dpt )_AB / <Ncoll_AB> * ( dN/dpt )_pp ; <Ncoll> = <Tab> * sigma_NN^inel
532  // R_AB = [ ( dsigma/dpt )_AB / sigma_AB ] / <TAB> * ( dsigma/dpt )_pp
533  //
534  Int_t istartPPfd=0, istartPPsyst=0, istartABfd=0, istartPPextr=0;
535  Double_t yPPh=0., yPPl=0., yABh=0., yABl=0.;
536  Double_t RaaCharm =0., RaaBeauty=0.;
537  Double_t RaaCharmFDhigh = 0., RaaCharmFDlow = 0.;
538  Double_t RaaBeautyFDhigh = 0., RaaBeautyFDlow = 0.;
539  Double_t systUp=0., systLow=0., systPPUp=0., systPPLow=0., systABUp=0., systABLow=0.;
540  //
541  //
542  // Search the central value of the energy loss hypothesis Rb = Rc (bin)
543  //
544  Double_t ElossCentral[nbins+1];
545  for(Int_t i=0; i<=nbins; i++) { ElossCentral[i]=0.; }
546  //
547  for(Int_t ientry=0; ientry<=entries; ientry++){
548 
549  nSigmaAB->GetEntry(ientry);
550  // cout << " pt="<< pt<<" sigma-AB="<<sigmaAB<<endl;
551  if ( !(sigmaAB>0.) ) continue;
552  //if(decay==2 && pt<2.) continue;
553 
554  // Compute RAB and the statistical uncertainty
555  Int_t hppbin = hSigmaPP->FindBin( pt );
556  Int_t hABbin = hSigmaAB->FindBin( pt );
557  Double_t sigmapp = hSigmaPP->GetBinContent( hppbin );
558  sigmapp *= scalePPRefToMatchRapidityBin; // scale to the proper rapidity bin width
559  // cout << " pt="<< pt<<", sigma-pp="<< sigmapp<<endl;
560  if (isRaavsEP>0.) sigmapp = 0.5*sigmapp;
561  if ( !(sigmapp>0.) ) continue;
562 
563  RaaCharm = ( sigmaAB / sigmaABCINT1B ) / ((Tab*1e3) * sigmapp *1e-12 ) ;
564  if(!isUseTaaForRaa) {
565  RaaCharm = ( sigmaAB ) / ( (A*B) * sigmapp ) ;
566  }
567 
568  if (fdMethod==kNb) {
569  RaaBeauty = Rb ;
570  }
571  else if (fdMethod==kfc) {
572  RaaBeauty = ( RaaCharm / Rcb ) ;
573  }
574 
575  Double_t ElossHypo = 0.;
576  if (fdMethod==kfc) { ElossHypo = 1. / Rcb; }
577  else { ElossHypo = 1. / (RaaCharm / RaaBeauty) ; }
578  if(isRbHypo) ElossHypo = RaaBeauty;
579 
580  // If using shadowing hypothesis, change the central hypothesis too
581  if(isShadHypothesis) CentralHypo = centralRbcShad[hABbin];
582 
583  // cout <<" pt "<< pt << " Raa charm " << RaaCharm << " Raa beauty " << RaaBeauty << " eloss hypo "<< ElossHypo<<endl;
584  //
585  // Find the bin for the central Eloss hypo
586  //
587  if( TMath::Abs( ElossHypo - CentralHypo ) < 0.075 ){
588  Double_t DeltaIni = TMath::Abs( ElossCentral[ hABbin ] - CentralHypo );
589  Double_t DeltaV = TMath::Abs( ElossHypo - CentralHypo );
590  // cout << " pt " << pt << " ECentral " << ElossCentral[ hABbin ] << " Ehypo "<< ElossHypo ;
591  if ( DeltaV < DeltaIni ) ElossCentral[ hABbin ] = ElossHypo;
592  // cout << " final ECentral " << ElossCentral[ hABbin ] << endl;
593  }
594  }
595  //
596  // Calculation of the Raa and its uncertainties
597  //
598  for(Int_t ientry=0; ientry<entries; ientry++){
599 
600  nSigmaAB->GetEntry(ientry);
601  if ( !(sigmaAB>0.) ) continue;
602  // if ( pt<2 || pt>16) continue;
603 
604 
605  // Compute RAB and the statistical uncertainty
606  Int_t hppbin = hSigmaPP->FindBin( pt );
607  Double_t sigmapp = hSigmaPP->GetBinContent( hppbin );
608  if (isRaavsEP>0.) sigmapp = 0.5*sigmapp;
609  sigmapp *= scalePPRefToMatchRapidityBin; // scale to the proper rapidity bin width
610  if ( !(sigmapp>0.) ) continue;
611 
612  RaaCharm = ( sigmaAB / sigmaABCINT1B ) / ((Tab*1e3) * sigmapp *1e-12 );
613  if(!isUseTaaForRaa) {
614  RaaCharm = ( sigmaAB ) / ( (A*B) * sigmapp ) ;
615  }
616 
617  // Flag to know if it is an scaled or extrapolated point of the pp reference
618  Bool_t isExtrapolatedBin = kFALSE;
619  if(isScaledAndExtrapRef && hCombinedReferenceFlag) isExtrapolatedBin = hCombinedReferenceFlag->GetBinContent( hppbin );
620  istartPPsyst = -1;
621  istartPPsyst = FindGraphBin(gSigmaPPSyst,pt);
622 
623  //
624  // FONLL Feed-Down systematics
625  //
626  istartPPfd = -1;
627  if(!isExtrapolatedBin) istartPPfd = FindGraphBin(gSigmaPPSystFeedDown,pt);
628  istartABfd = -1;
629  istartABfd = FindGraphBin(gSigmaABSystFeedDown,pt);
630 
631  // cout << " Starting bin for pp is "<< istartPPfd <<", for AA is "<<istartABfd << endl;
632  if(isExtrapolatedBin){
633  if(gReferenceFdSyst){
634  Int_t ibinfd = FindGraphBin(gReferenceFdSyst,pt);
635  yPPh = gReferenceFdSyst->GetErrorYhigh(ibinfd);
636  yPPl = gReferenceFdSyst->GetErrorYlow(ibinfd);
637  }
638  } else {
639  yPPh = gSigmaPPSystFeedDown->GetErrorYhigh(istartPPfd);
640  yPPl = gSigmaPPSystFeedDown->GetErrorYlow(istartPPfd);
641  }
642  if (isRaavsEP>0.) {
643  yPPh = yPPh*0.5;
644  yPPl = yPPl*0.5;
645  }
646  yPPh *= scalePPRefToMatchRapidityBin; // scale to the proper rapidity bin width
647  yPPl *= scalePPRefToMatchRapidityBin; // scale to the proper rapidity bin width
648 
649  yABh = gSigmaABSystFeedDown->GetErrorYhigh(istartABfd);
650  yABl = gSigmaABSystFeedDown->GetErrorYlow(istartABfd);
651 
652 
653  RaaCharmFDhigh = ( sigmaABMax / sigmaABCINT1B ) / ((Tab*1e3) * (sigmapp+yPPh) *1e-12 ) ;
654  RaaCharmFDlow = ( sigmaABMin / sigmaABCINT1B ) / ((Tab*1e3) * (sigmapp-yPPl) *1e-12 ) ;
655  if(printout && TMath::Abs(ptprintout-pt)<0.1 ) cout << endl<<" pt "<< pt << " Raa " << RaaCharm <<" high "<< RaaCharmFDhigh << " low "<< RaaCharmFDlow<<endl;
656  if(!isUseTaaForRaa) {
657  RaaCharmFDhigh = ( sigmaABMax ) / ( (A*B)* (sigmapp+yPPh) ) ;
658  RaaCharmFDlow = ( sigmaABMin ) / ( (A*B)* (sigmapp-yPPl) ) ;
659  }
660 
661 
662  if (fdMethod==kNb) {
663  RaaBeauty = Rb ;
664  RaaBeautyFDlow = Rb ;
665  RaaBeautyFDhigh = Rb ;
666  ntupleRAB->Fill( pt, Tab*1e3, sigmapp*1e-12, sigmaAB*1e-12, sigmaAB/sigmaABCINT1B,
667  sigmaABMax / sigmaABCINT1B, sigmaABMin / sigmaABCINT1B,
668  RaaCharm, RaaCharmFDhigh, RaaCharmFDlow, RaaBeauty, fcAB );
669  }
670  else if (fdMethod==kfc) {
671  RaaBeauty = ( RaaCharm / Rcb ) ;
672  RaaBeautyFDlow = ( RaaCharmFDlow / Rcb ) ;
673  RaaBeautyFDhigh = ( RaaCharmFDhigh / Rcb ) ;
674  hRABvsRcb->Fill( pt, RaaCharm, RaaBeauty );
675  ntupleRAB->Fill( pt, Tab*1e3, sigmapp*1e-12, sigmaAB*1e-12, sigmaAB/sigmaABCINT1B,
676  sigmaABMax / sigmaABCINT1B, sigmaABMin / sigmaABCINT1B,
677  Rcb, RaaCharm, RaaCharmFDhigh, RaaCharmFDlow, RaaBeauty, RaaBeautyFDhigh, RaaBeautyFDlow, fcAB );
678  }
679  hRABvsRb->Fill( pt, RaaCharm, RaaBeauty );
680  hRABvsRbFDlow->Fill( pt, RaaCharmFDlow, RaaBeautyFDlow );
681  hRABvsRbFDhigh->Fill( pt, RaaCharmFDhigh, RaaBeautyFDhigh );
682  if(printout && TMath::Abs(ptprintout-pt)<0.1) cout << " pt "<< pt << " Rb " << RaaBeauty <<" high "<< RaaBeautyFDhigh << " low "<< RaaBeautyFDlow <<endl;
683 
684  hRABCharmVsRBeautyVsPt->Fill( pt, RaaBeauty, RaaCharm );
685  Int_t ptbin = hRABvsPt->FindBin( pt );
686  hRCharmVsRBeauty[ptbin]->Fill( RaaBeauty, RaaCharm );
687  hRCharmVsRBeauty[ptbin]->Fill( RaaBeautyFDlow, RaaCharmFDlow );
688  hRCharmVsRBeauty[ptbin]->Fill( RaaBeautyFDhigh, RaaCharmFDhigh );
689 
690 
691  if (fdMethod==kfc) {
692  if( TMath::Abs(Rcb-0.015)<0.009 ) hRABEloss00->Fill(pt,RaaCharm);
693  if( TMath::Abs(Rcb-0.5)<0.009 ) hRABEloss05->Fill(pt,RaaCharm);
694  if( TMath::Abs(Rcb-1.0)<0.009 ) {
695  hRABEloss10->Fill(pt,RaaCharm);
696  hRABvsRbFDhigh_proj->Fill(pt,RaaCharmFDhigh);
697  hRABvsRbFDlow_proj->Fill(pt,RaaCharmFDlow);
698  }
699  if( TMath::Abs(Rcb-1.5)<0.009 ) hRABEloss15->Fill(pt,RaaCharm);
700  if( TMath::Abs(Rcb-2.0)<0.009 ) hRABEloss20->Fill(pt,RaaCharm);
701  }
702  else if (fdMethod==kNb) {
703  if( TMath::Abs(RaaBeauty-0.015)<0.009 ) hRABEloss00->Fill(pt,RaaCharm);
704  if( TMath::Abs(RaaBeauty-0.5)<0.009 ) hRABEloss05->Fill(pt,RaaCharm);
705  if( TMath::Abs(RaaBeauty-1.0)<0.009 ) {
706  hRABEloss10->Fill(pt,RaaCharm);
707  hRABvsRbFDhigh_proj->Fill(pt,RaaCharmFDhigh);
708  hRABvsRbFDlow_proj->Fill(pt,RaaCharmFDlow);
709  }
710  if( TMath::Abs(RaaBeauty-1.5)<0.009 ) hRABEloss15->Fill(pt,RaaCharm);
711  if( TMath::Abs(RaaBeauty-2.0)<0.009 ) hRABEloss20->Fill(pt,RaaCharm);
712  }
713 
714 
715  Int_t hABbin = hMassAB->FindBin( pt );
716  if(isShadHypothesis) CentralHypo = centralRbcShad[hABbin];
717 
718  if(printout && TMath::Abs(ptprintout-pt)<0.1)
719  if ( fdMethod==kNb && TMath::Abs(Rb -CentralHypo)< 0.05) {
720  cout << " pt "<< pt <<", at bin "<<hABbin<<endl;
721  cout<<" entries "<<entries<<", i="<<ientry<<", pt="<<pt<<", Rb="<<Rb<<", Tab="<<Tab<<", sigmaAB="<<sigmaAB<<", sigmapp="<<sigmapp<<", Raacharm="<<RaaCharm<<", RaaBeauty="<<RaaBeauty<<endl;
722  cout << " AB basis: mass "<< hMassAB->GetBinContent(hABbin)<<", eff "<< hDirectEffptAB->GetBinContent(hABbin)<<endl;
723  cout<<" FD low, err low AB "<< (sigmaAB-sigmaABMin)<<" err low PP "<< yPPl<<" Raacharm="<<RaaCharmFDlow<<", RaaBeauty="<<RaaBeautyFDlow<<endl;
724  cout<<" FD high, err high AB "<< (sigmaABMax-sigmaAB)<<" err high PP "<< yPPh<<" Raacharm="<<RaaCharmFDhigh<<", RaaBeauty="<<RaaBeautyFDhigh<<endl;
725  }
726  if(printout && TMath::Abs(ptprintout-pt)<0.1)
727  if ( fdMethod==kfc) if(TMath::Abs(Rcb -CentralHypo)< 0.05 ){
728  cout << " pt "<< pt <<", at bin "<<hABbin<<endl;
729  cout<<" entries "<<entries<<", i="<<ientry<<", pt="<<pt<<", Rcb="<<Rcb<<", Tab="<<Tab<<", sigmaAB="<<sigmaAB<<", sigmapp="<<sigmapp<<", Raacharm="<<RaaCharm<<", RaaBeauty="<<RaaBeauty<<endl;
730  cout << " AB basis: mass "<< hMassAB->GetBinContent(hABbin)<<", eff "<< hDirectEffptAB->GetBinContent(hABbin)<<", fc "<<histofcAB->GetBinContent(hABbin)<< endl;
731  cout<<" FD low, err low AB "<< (sigmaAB-sigmaABMin)<<" err low PP "<< yPPl<<" Raacharm="<<RaaCharmFDlow<<", RaaBeauty="<<RaaBeautyFDlow<<endl;
732  cout<<" FD high, err high AB "<< (sigmaABMax-sigmaAB)<<" err high PP "<< yPPh<<" Raacharm="<<RaaCharmFDhigh<<", RaaBeauty="<<RaaBeautyFDhigh<<endl;
733  }
734 
735 
736  //
737  // Fill in the global properties ?
738  //
739  Double_t ElossHypo = 0.;
740  if (fdMethod==kfc) { ElossHypo = 1./ Rcb; }
741  else { ElossHypo = 1. / (RaaCharm / RaaBeauty); }
742  if(isRbHypo) ElossHypo = RaaBeauty;
743  hRCharmVsElossHypo[ptbin]->Fill( ElossHypo, RaaCharm );
744 
745  // If using shadowing hypothesis, change the limit hypothesis too
746  if(isShadHypothesis) {
747  MinHypo = minRbcShad[ hABbin ];
748  MaxHypo = maxRbcShad[ hABbin ];
749  }
750 
751  // cout <<" pt "<< pt << " Raa charm " << RaaCharm << " Raa beauty " << RaaBeauty << " eloss hypo "<< ElossHypo
752  if(ientry==0) cout<<" pt"<< pt<< " ElossCentral "<< ElossCentral[hABbin] << " min-hypo "<<MinHypo << " max-hypo "<<MaxHypo<<endl;
753 
754  //
755  // Fill in histos charm (null Eloss hypothesis)
756  //
757  Double_t minFdSyst = 0., maxFdSyst = 0.;
758  if ( ElossHypo == ElossCentral[ hABbin ] ) {
759 
760  //
761  // Data stat uncertainty
762  //
763  Double_t sigmappStat = hSigmaPP->GetBinError( hppbin );
764  if (isRaavsEP>0.) sigmappStat = sigmappStat*0.5;
765  sigmappStat *= scalePPRefToMatchRapidityBin; // scale to the proper rapidity bin width
766  Int_t hRABbin = hRABvsPt->FindBin( pt );
767  Double_t stat = RaaCharm * TMath::Sqrt( (statUncSigmaAB/sigmaAB)*(statUncSigmaAB/sigmaAB) +
768  (sigmappStat/sigmapp)*(sigmappStat/sigmapp) ) ;
769  if ( RaaCharm==0 ) stat =0.;
770  if ( RaaCharm>0 ) {
771  hRABvsPt->SetBinContent( hRABbin, RaaCharm );
772  hRABvsPt->SetBinError( hRABbin, stat );
773  hYieldABvsPt->SetBinContent( hRABbin, sigmaAB/sigmaABCINT1B );
774  hYieldABvsPt->SetBinError( hRABbin, statUncSigmaAB/sigmaABCINT1B );
775 
776  cout << "pt="<< pt<< " Raa " << RaaCharm << " stat unc. "<< stat <<
777  " sigma-pp "<< sigmapp <<" sigma-AB "<< sigmaAB<<endl;
778  if(printout && TMath::Abs(ptprintout-pt)<0.1) {
779  cout << " Raa " << RaaCharm << " stat unc. "<< stat << " is "<< stat/RaaCharm * 100. <<
780  "%, stat-pp "<< sigmappStat/sigmapp*100. <<"% stat-AB "<< statUncSigmaAB/sigmaAB*100.<<"%"<<endl;
781  }
782 
783  Double_t errstatEff = fhStatUncEffcSigmaAB->GetBinError( hRABbin );
784  fhStatUncEffcSigmaAB_Raa->SetBinError( hRABbin, errstatEff*RaaCharm );
785  errstatEff = fhStatUncEffbSigmaAB->GetBinError( hRABbin );
786  fhStatUncEffbSigmaAB_Raa->SetBinError( hRABbin, errstatEff*RaaCharm );
787  errstatEff = fhStatUncEffcFDAB->GetBinError( hRABbin );
788  fhStatUncEffcFDAB_Raa->SetBinError( hRABbin, errstatEff*RaaCharm );
789  errstatEff = fhStatUncEffbFDAB->GetBinError( hRABbin );
790  fhStatUncEffbFDAB_Raa->SetBinError( hRABbin, errstatEff*RaaCharm );
791  }
792 
793 
794  //
795  //
796  // Data systematics (sigma syst-but FD + extrap) syst
797  //
798  //
799  // Data syst: a) Syst in p-p
800  //
801  Double_t ptwidth = hSigmaAB->GetBinWidth(hABbin) / 2. ;
802  istartPPextr = -1;
803  if(!isExtrapolatedBin) istartPPextr = FindGraphBin(gSigmaPPSystTheory,pt);
804 
805  Double_t dataPPUp=0., dataPPLow=0.;
806  if(isExtrapolatedBin) {
807  dataPPUp = gSigmaPPSyst->GetErrorYhigh(istartPPsyst);
808  dataPPLow = gSigmaPPSyst->GetErrorYlow(istartPPsyst);
809  systPPUp = dataPPUp;
810  systPPLow = dataPPLow;
811  } else {
812  dataPPUp = ExtractFDSyst( gSigmaPPSystData->GetErrorYhigh(istartPPextr), gSigmaPPSystFeedDown->GetErrorYhigh(istartPPfd) );
813  dataPPLow = ExtractFDSyst( gSigmaPPSystData->GetErrorYlow(istartPPextr), gSigmaPPSystFeedDown->GetErrorYlow(istartPPfd) );
814  systPPUp = TMath::Sqrt( dataPPUp*dataPPUp + gSigmaPPSystTheory->GetErrorYhigh(istartPPextr)*gSigmaPPSystTheory->GetErrorYhigh(istartPPextr) );
815  systPPLow = TMath::Sqrt( dataPPLow*dataPPLow + gSigmaPPSystTheory->GetErrorYlow(istartPPextr)*gSigmaPPSystTheory->GetErrorYlow(istartPPextr) );
816  }
817  if (isRaavsEP>0.) {
818  dataPPUp = dataPPUp*0.5;
819  dataPPLow = dataPPLow*0.5;
820  if(isExtrapolatedBin) {
821  systPPUp = dataPPUp;
822  systPPLow = dataPPLow;
823  } else {
824  systPPUp = TMath::Sqrt( dataPPUp*dataPPUp + 0.5*gSigmaPPSystTheory->GetErrorYhigh(istartPPextr)*0.5*gSigmaPPSystTheory->GetErrorYhigh(istartPPextr) );
825  systPPLow = TMath::Sqrt( dataPPLow*dataPPLow + 0.5*gSigmaPPSystTheory->GetErrorYlow(istartPPextr)*0.5*gSigmaPPSystTheory->GetErrorYlow(istartPPextr) );
826  }
827  }
828  systPPUp *= scalePPRefToMatchRapidityBin; // scale to the proper rapidity bin width
829  systPPLow *= scalePPRefToMatchRapidityBin; // scale to the proper rapidity bin width
830 
831 
832  if(printout && TMath::Abs(ptprintout-pt)<0.1) {
833  cout << " pt : "<< pt<<" Syst-pp-data "<< dataPPUp/sigmapp << "%, ";
834  if(!isExtrapolatedBin){
835  if (isRaavsEP>0.) cout <<" extr unc + "<< 0.5*gSigmaPPSystTheory->GetErrorYhigh(istartPPextr)/sigmapp <<" - "<< 0.5*gSigmaPPSystTheory->GetErrorYlow(istartPPextr)/sigmapp <<" %";
836  else cout <<" extr unc + "<< (gSigmaPPSystTheory->GetErrorYhigh(istartPPextr)*scalePPRefToMatchRapidityBin)/sigmapp <<" - "<< (gSigmaPPSystTheory->GetErrorYlow(istartPPextr)*scalePPRefToMatchRapidityBin)/sigmapp <<" %";
837  }
838  cout << endl;
839  }
840 
841  //
842  // Data syst: b) Syst in PbPb
843  //
844  Double_t dataSystUp=0., dataSystDown=0.;
845  Bool_t PbPbDataSystOk = PbPbDataSyst(systematicsAB,pt,cc,dataSystUp,dataSystDown);
846  if (!PbPbDataSystOk) { cout <<" There is some issue with the PbPb data systematics, please check and rerun"<<endl; return; }
847  systABUp = sigmaAB * TMath::Sqrt( dataSystUp*dataSystUp +
848  (hDirectEffptAB->GetBinError(hABbin)/hDirectEffptAB->GetBinContent(hABbin))*(hDirectEffptAB->GetBinError(hABbin)/hDirectEffptAB->GetBinContent(hABbin)) );
849 
850  systABLow = sigmaAB * TMath::Sqrt( dataSystDown*dataSystDown +
851  (hDirectEffptAB->GetBinError(hABbin)/hDirectEffptAB->GetBinContent(hABbin))*(hDirectEffptAB->GetBinError(hABbin)/hDirectEffptAB->GetBinContent(hABbin)) );
852  //
853  // Data syst : c) combine pp & PbPb
854  //
855  systLow = sigmapp>0. ?
856  RaaCharm * TMath::Sqrt( (systABLow/sigmaAB)*(systABLow/sigmaAB) + (systPPUp/sigmapp)*(systPPUp/sigmapp) )
857  : 0.;
858 
859  systUp = sigmapp>0. ?
860  RaaCharm * TMath::Sqrt( (systABUp/sigmaAB)*(systABUp/sigmaAB) + (systPPLow/sigmapp)*(systPPLow/sigmapp) )
861  : 0.;
862  if ( RaaCharm==0 ) { systPPUp =0.; systPPLow =0.; }
863 
864  // if(printout)
865  cout << " Syst-pp-up "<< systPPUp/sigmapp <<"%, syst-pp-low "<< systPPLow/sigmapp <<"%, syst-AB-up "<<systABUp/sigmaAB<<"%, syst-AB-low "<<systABLow/sigmaAB<<"%, tot-syst-up "<<systUp/RaaCharm<<"%, tot-syst-low "<<systLow/RaaCharm<<"%"<<endl;
866 
867  if ( RaaCharm>0 ) {
868  hRABvsPt_DataSystematics->SetBinContent( hRABbin, RaaCharm );
869  hRABvsPt_DataSystematics->SetBinError( hRABbin, systUp );
870  gRAB_DataSystematics->SetPoint( hABbin, pt, RaaCharm ); // i, x, y
871  gRAB_DataSystematics->SetPointError( hABbin, ptwidth, ptwidth, systLow, systUp );
872  gRAB_DataSystematics->SetPointEXlow(hABbin, 0.4); gRAB_DataSystematics->SetPointEXhigh(hABbin,0.4);
873  gRAB_DataSystematicsPP->SetPoint( hABbin, pt, RaaCharm ); // i, x, y
874  gRAB_DataSystematicsPP->SetPointError( hABbin, ptwidth, ptwidth, RaaCharm *(systPPUp/sigmapp), RaaCharm *systPPLow/sigmapp );
875  gRAB_DataSystematicsAB->SetPoint( hABbin, pt, RaaCharm ); // i, x, y
876  gRAB_DataSystematicsAB->SetPointError( hABbin, ptwidth, ptwidth, RaaCharm *systABLow/sigmaAB, RaaCharm *systABUp/sigmaAB );
877  }
878 
879  //
880  // Feed-down Systematics
881  //
882  Double_t FDL=0., FDH=0.;
883  if ( RaaCharmFDhigh > RaaCharmFDlow ){
884  FDH = RaaCharmFDhigh; FDL = RaaCharmFDlow;
885  } else {
886  FDL = RaaCharmFDhigh; FDH = RaaCharmFDlow;
887  }
888 
889  if(printout && TMath::Abs(ptprintout-pt)<0.1) cout<<" Raa "<<RaaCharm<<", Raa-fd-low "<<RaaCharmFDlow <<", Raa-fd-high "<<RaaCharmFDhigh <<endl;
890  maxFdSyst = TMath::Abs(FDH - RaaCharm);
891  minFdSyst = TMath::Abs(RaaCharm - FDL);
892  if ( RaaCharm>0 ) {
893  gRAB_FeedDownSystematics->SetPoint( hABbin, pt, RaaCharm ); // i, x, y
894  gRAB_FeedDownSystematics->SetPointError( hABbin, 0.3, 0.3, minFdSyst, maxFdSyst ); // i, x, y
895  gRAB_fcFeedDownOnly->SetPoint( hABbin, pt,fcAB );
896  gRAB_fcFeedDownOnly->SetPointError(hABbin, 0.3, 0.3, fcAB-(sigmaABMin/sigmaAB*fcAB), (sigmaABMax/sigmaAB*fcAB)-fcAB );
897  }
898 
899  // if(printout) {
900  cout<<" FD syst +"<< maxFdSyst/RaaCharm <<" - "<<minFdSyst/RaaCharm<<endl;
901  cout<<" fc = "<<fcAB<<", ("<< sigmaABMax/sigmaAB * fcAB <<","<< sigmaABMin/sigmaAB * fcAB <<")"<<endl;
902  // }
903 
904  //
905  // Filling part of the Eloss scenarii information
906  //
907  if(RaaCharm>0 ) {
908  gRAB_ElossHypothesis->SetPoint( hABbin, pt, RaaCharm ); // i, x, y
909  gRAB_ElossHypothesis->SetPointEXlow( hABbin, ptwidth);
910  gRAB_ElossHypothesis->SetPointEXhigh( hABbin, ptwidth);
911  gRAB_FeedDownSystematicsElossHypothesis->SetPoint( hABbin, pt, RaaCharm ); // i, x, y
912  gRAB_FeedDownSystematicsElossHypothesis->SetPointEXlow( hABbin, ptwidth);
913  gRAB_FeedDownSystematicsElossHypothesis->SetPointEXhigh( hABbin, ptwidth);
914  gRAB_GlobalSystematics->SetPoint( hABbin, pt, RaaCharm ); // i, x, y
915  gRAB_GlobalSystematics->SetPointEXlow(hABbin,0.4); gRAB_GlobalSystematics->SetPointEXhigh(hABbin,0.4);
916  }
917  }
918 
919  //
920  // Filling Eloss scenarii information
921  //
922  // trick in case not fine enough Rb hypothesis to cope with the min/max range
923  // if( RaaCharm>0 && ( (ElossHypo >= MinHypo && ElossHypo <=MaxHypo) || ElossHypo == ElossCentral[ hABbin ] ) && RaaBeauty<=MaxRb ) {
924  // by default better not use it, to monitor when this happens (could affect results)
925  if( RaaCharm>0 && ElossHypo >= MinHypo && ElossHypo <=MaxHypo && RaaBeauty<=MaxRb ) {
926 
927  Double_t Ehigh = ElossMax[ hABbin ] ;
928  Double_t Elow = ElossMin[ hABbin ] ;
929  if ( RaaCharm > Ehigh ) ElossMax[ hABbin ] = RaaCharm ;
930  if ( RaaCharm < Elow ) ElossMin[ hABbin ] = RaaCharm ;
931  if(printout && TMath::Abs(ptprintout-pt)<0.1) {
932  cout<<" Hypothesis " << ElossHypo << " sigma-AB "<< sigmaAB <<", Raa "<< RaaCharm <<", Raa Eloss max "<< ElossMax[hABbin] <<" Raa Eloss min "<< ElossMin[hABbin] << " Rb="<< RaaBeauty <<endl;
933  cout<<" Rb="<< RaaBeauty <<" max "<< RaaBeautyFDhigh <<" min "<< RaaBeautyFDlow <<endl;
934  }
935  Double_t fcEhigh = fcElossMax[ hABbin ] ;
936  Double_t fcElow = fcElossMin[ hABbin ] ;
937  if ( fcAB > fcEhigh ) fcElossMax[ hABbin ] = fcAB ;
938  if ( fcAB < fcElow ) fcElossMin[ hABbin ] = fcAB ;
939  Double_t FDEhigh = FDElossMax[ hABbin ];
940  Double_t FDEmin = FDElossMin[ hABbin ];
941  Double_t RFDhigh = RaaCharmFDhigh>RaaCharmFDlow ? RaaCharmFDhigh : RaaCharmFDlow;
942  Double_t RFDlow = RaaCharmFDlow<RaaCharmFDhigh ? RaaCharmFDlow : RaaCharmFDhigh;
943  if ( RFDhigh > FDEhigh ) FDElossMax[ hABbin ] = RFDhigh ;
944  if ( RFDlow < FDEmin ) FDElossMin[ hABbin ] = RFDlow ;
945  if(printout && TMath::Abs(ptprintout-pt)<0.1)
946  cout<<" Hypothesis " << ElossHypo << " sigma-AB "<< sigmaAB <<", Raa FD-max Eloss max "<< FDElossMax[hABbin] <<" Raa FD-min Eloss min "<< FDElossMin[hABbin] <<endl;
947  }
948 
949 
950  }
951 
952 
953  // Finish filling the y-uncertainties of the Eloss scenarii
954  for (Int_t ibin=0; ibin<=nbins; ibin++){
955  Double_t ipt=0., value =0.;
956  gRAB_ElossHypothesis->GetPoint(ibin,ipt,value);
957  if(ipt<=0) continue;
958  //
959  // Uncertainty on Raa due to the Eloss hypothesis
960  Double_t elossYhigh = TMath::Abs( ElossMax[ibin] - value );
961  Double_t elossYlow = TMath::Abs( value - ElossMin[ibin] );
962  gRAB_ElossHypothesis->SetPointEYhigh(ibin, elossYhigh );
963  gRAB_ElossHypothesis->SetPointEYlow(ibin, elossYlow );
964  gRAB_ElossHypothesis->SetPointEXhigh(ibin, 0.2);
965  gRAB_ElossHypothesis->SetPointEXlow(ibin, 0.2);
966  cout << " pt "<< ipt << " Raa "<< value <<" max "<< ElossMax[ibin] << " min " <<ElossMin[ibin] <<endl;
967  cout<<" Eloss syst +"<< elossYhigh <<" - "<< elossYlow <<endl;
968  // cout << " fc max "<< fcElossMax[ibin] << " fc min " <<fcElossMin[ibin] <<endl;
969  //
970  // Uncertainty on Raa due to the FD unc & Eloss hypothesis
971  Double_t fdElossEYhigh = TMath::Abs( FDElossMax[ibin] - value );
972  Double_t fdElossEYlow = TMath::Abs( value - FDElossMin[ibin] );
973  if(elossFDQuadSum){
974  Double_t fdEYhigh = gRAB_FeedDownSystematics->GetErrorYhigh(ibin);
975  fdElossEYhigh = TMath::Sqrt( elossYhigh*elossYhigh + fdEYhigh*fdEYhigh );
976  Double_t fdEYlow = gRAB_FeedDownSystematics->GetErrorYlow(ibin);
977  fdElossEYlow = TMath::Sqrt( elossYlow*elossYlow + fdEYlow*fdEYlow );
978  }
979  gRAB_FeedDownSystematicsElossHypothesis->SetPointEYhigh(ibin, fdElossEYhigh );
980  gRAB_FeedDownSystematicsElossHypothesis->SetPointEYlow(ibin, fdElossEYlow );
981  gRAB_FeedDownSystematicsElossHypothesis->SetPointEXhigh(ibin, 0.25);
982  gRAB_FeedDownSystematicsElossHypothesis->SetPointEXlow(ibin, 0.25);
983  cout<<" FD & Eloss syst +"<< fdElossEYhigh <<" - "<< fdElossEYlow
984  <<" = + "<< fdElossEYhigh/value <<" - "<< fdElossEYlow/value <<" %" <<endl;
985  //
986  // All uncertainty on Raa (FD unc & Eloss + data)
987  Double_t systdatal = gRAB_DataSystematics->GetErrorYlow(ibin);
988  Double_t systdatah = gRAB_DataSystematics->GetErrorYhigh(ibin);
989  Double_t systgbhUnc = TMath::Sqrt( systdatah*systdatah + fdElossEYhigh*fdElossEYhigh );
990  Double_t systgblUnc = TMath::Sqrt( systdatal*systdatal + fdElossEYlow*fdElossEYlow );
991  gRAB_GlobalSystematics->SetPointEYhigh(ibin,systgbhUnc);
992  gRAB_GlobalSystematics->SetPointEYlow(ibin,systgblUnc);
993  cout<<" Data syst +"<< systdatah <<" - "<< systdatal <<" = + "<< systdatah/value <<" - " << systdatal/value << " % "<<endl;
994  cout<<" Global syst +"<< systgbhUnc <<" - "<< systgblUnc << " = + "<< systgbhUnc/value <<" - "<< systgblUnc/value << " %" <<endl;
995  //
996  }
997 
998  cout<<endl<<" Calculation finished, now drawing"<<endl<<endl;
999 
1000 
1001  gROOT->SetStyle("Plain");
1002  gStyle->SetPalette(1);
1003  gStyle->SetOptStat(0);
1004 
1005 
1006  TCanvas *cRABvsRb = new TCanvas("RABvsRb","RAB vs Rb");
1007  hRABvsRb->Draw("colz");
1008  cRABvsRb->Update();
1009 
1010 // TCanvas *cRABvsRbvsPt = new TCanvas("cRABvsRbvsPt","RAB vs Rb vs pt");
1011 // hRABCharmVsRBeautyVsPt->Draw("lego3z");
1012 // cRABvsRbvsPt->Update();
1013 
1014 
1015  cout<< " Drawing feed-down contribution"<<endl;
1016  TCanvas *cRABvsRbFDl = new TCanvas("RABvsRbFDl","RAB vs Rb (FD low)");
1017  hRABvsRbFDlow->Draw("cont4z");
1018  cRABvsRbFDl->Update();
1019  TCanvas *cRABvsRbFDh = new TCanvas("RABvsRbFDh","RAB vs Rb (FD high)");
1020  hRABvsRbFDhigh->Draw("cont4z");
1021  cRABvsRbFDh->Update();
1022 
1023  TCanvas * cSigmaABptEloss = new TCanvas("cSigmaABptEloss","SigmaAB vs pt, Eloss hypothesis");
1024  TH1D *hSigmaABEloss00= new TH1D("hSigmaABEloss00","hSigmaABEloss00",nbins,limits);
1025  TH1D *hSigmaABEloss05= new TH1D("hSigmaABEloss05","hSigmaABEloss05",nbins,limits);
1026  TH1D *hSigmaABEloss10= new TH1D("hSigmaABEloss10","hSigmaABEloss10",nbins,limits);
1027  TH1D *hSigmaABEloss15= new TH1D("hSigmaABEloss15","hSigmaABEloss15",nbins,limits);
1028  TH1D *hSigmaABEloss20= new TH1D("hSigmaABEloss20","hSigmaABEloss20",nbins,limits);
1029 
1030  delete [] limits;
1031  delete [] binwidths;
1032 
1033  for (Int_t i=0; i<=nSigmaAB->GetEntriesFast(); i++) {
1034  nSigmaAB->GetEntry(i);
1035  if (fdMethod==kfc) {
1036  if( TMath::Abs(Rcb-0.015)<0.009 ) hSigmaABEloss00->Fill(pt,sigmaAB);
1037  if( TMath::Abs(Rcb-0.5)<0.009 ) hSigmaABEloss05->Fill(pt,sigmaAB);
1038  if( TMath::Abs(Rcb-1.0)<0.009 ) hSigmaABEloss10->Fill(pt,sigmaAB);
1039  if( TMath::Abs(Rcb-1.5)<0.009 ) hSigmaABEloss15->Fill(pt,sigmaAB);
1040  if( TMath::Abs(Rcb-2.0)<0.009 ) hSigmaABEloss20->Fill(pt,sigmaAB);
1041  }
1042  else if (fdMethod==kNb) {
1043  if( TMath::Abs(Rb-0.015)<0.009 ) hSigmaABEloss00->Fill(pt,sigmaAB);
1044  if( TMath::Abs(Rb-0.5)<0.009 ) hSigmaABEloss05->Fill(pt,sigmaAB);
1045  if( TMath::Abs(Rb-1.0)<0.009 ) hSigmaABEloss10->Fill(pt,sigmaAB);
1046  if( TMath::Abs(Rb-1.5)<0.009 ) hSigmaABEloss15->Fill(pt,sigmaAB);
1047  if( TMath::Abs(Rb-2.0)<0.009 ) hSigmaABEloss20->Fill(pt,sigmaAB);
1048  }
1049  }
1050  hSigmaABEloss00->SetLineColor(2);
1051  hSigmaABEloss05->SetLineColor(3);
1052  hSigmaABEloss10->SetLineColor(4);
1053  hSigmaABEloss15->SetLineColor(kMagenta+1);
1054  hSigmaABEloss20->SetLineColor(kGreen+2);
1055  hSigmaABEloss00->SetMarkerStyle(22);
1056  hSigmaABEloss05->SetMarkerStyle(26);
1057  hSigmaABEloss10->SetMarkerStyle(20);
1058  hSigmaABEloss15->SetMarkerStyle(25);
1059  hSigmaABEloss20->SetMarkerStyle(21);
1060  if (fdMethod==kNb) {
1061  hSigmaABEloss05->Draw("ph");
1062  hSigmaABEloss10->Draw("phsame");
1063  hSigmaABEloss15->Draw("phsame");
1064  hSigmaABEloss20->Draw("phsame");
1065  }
1066  else {
1067  hSigmaABEloss20->Draw("p");
1068  hSigmaABEloss00->Draw("phsame");
1069  hSigmaABEloss05->Draw("phsame");
1070  hSigmaABEloss10->Draw("phsame");
1071  hSigmaABEloss15->Draw("phsame");
1072  hSigmaABEloss20->Draw("phsame");
1073  }
1074  TLegend *legrcb = new TLegend(0.8,0.8,0.95,0.9);
1075  legrcb->SetFillColor(0);
1076  legrcb->AddEntry(hSigmaABEloss00,"Rc/b=0.0","lp");
1077  legrcb->AddEntry(hSigmaABEloss05,"Rc/b=0.5","lp");
1078  legrcb->AddEntry(hSigmaABEloss10,"Rc/b=1.0","lp");
1079  legrcb->AddEntry(hSigmaABEloss15,"Rc/b=1.5","lp");
1080  legrcb->AddEntry(hSigmaABEloss20,"Rc/b=2.0","lp");
1081  legrcb->Draw();
1082  cSigmaABptEloss->Update();
1083 
1084 
1085  TCanvas * cRABptEloss = new TCanvas("cRABptEloss","RAB vs pt, Eloss hypothesis");
1086  hRABEloss00->SetLineColor(2);
1087  hRABEloss05->SetLineColor(3);
1088  hRABEloss10->SetLineColor(4);
1089  hRABEloss15->SetLineColor(kMagenta+1);
1090  hRABEloss20->SetLineColor(kGreen+2);
1091  hRABEloss00->SetMarkerStyle(22);
1092  hRABEloss05->SetMarkerStyle(26);
1093  hRABEloss10->SetMarkerStyle(20);
1094  hRABEloss15->SetMarkerStyle(25);
1095  hRABEloss20->SetMarkerStyle(21);
1096  if (fdMethod==kNb) {
1097  hRABEloss05->Draw("ph");
1098  hRABEloss10->Draw("phsame");
1099  hRABEloss15->Draw("phsame");
1100  hRABEloss20->Draw("phsame");
1101  }
1102  else {
1103  hRABEloss20->Draw("p");
1104  hRABEloss00->Draw("phsame");
1105  hRABEloss05->Draw("phsame");
1106  hRABEloss10->Draw("phsame");
1107  hRABEloss15->Draw("phsame");
1108  hRABEloss20->Draw("phsame");
1109  }
1110  legrcb = new TLegend(0.8,0.8,0.95,0.9);
1111  legrcb->SetFillColor(0);
1112  if (fdMethod==kfc) {
1113  legrcb->AddEntry(hRABEloss00,"Rc/b=0.0","lp");
1114  legrcb->AddEntry(hRABEloss05,"Rc/b=0.5","lp");
1115  legrcb->AddEntry(hRABEloss10,"Rc/b=1.0","lp");
1116  legrcb->AddEntry(hRABEloss15,"Rc/b=0.5","lp");
1117  legrcb->AddEntry(hRABEloss20,"Rc/b=2.0","lp");
1118  }
1119  else if (fdMethod==kNb) {
1120  legrcb->AddEntry(hRABEloss00,"Rb=0.0","lp");
1121  legrcb->AddEntry(hRABEloss05,"Rb=0.5","lp");
1122  legrcb->AddEntry(hRABEloss10,"Rb=1.0","lp");
1123  legrcb->AddEntry(hRABEloss15,"Rb=0.5","lp");
1124  legrcb->AddEntry(hRABEloss20,"Rb=2.0","lp");
1125  }
1126  legrcb->Draw();
1127  cRABptEloss->Update();
1128 
1129 
1130  cout<< " Drawing summary results"<<endl;
1131  TCanvas * cRABpt = new TCanvas("cRABpt","RAB vs pt, no hypothesis");
1132  hRABEloss10->Draw("");
1133  cRABpt->Update();
1134 
1135  TCanvas * cRABptFDUnc = new TCanvas("cRABptFDUnc","RAB vs pt, FD Uncertainties");
1136  hRABvsRbFDlow_proj->Draw("");
1137  hRABEloss10->Draw("phsame");
1138  hRABvsRbFDhigh_proj->SetLineColor(kMagenta+1);
1139  hRABvsRbFDhigh_proj->Draw("same");
1140  hRABvsRbFDlow_proj->SetLineColor(kGreen+2);
1141  hRABvsRbFDlow_proj->Draw("same");
1142  legrcb = new TLegend(0.8,0.8,0.95,0.9);
1143  legrcb->SetFillColor(0);
1144  legrcb->AddEntry(hRABEloss10,"FD Central","lp");
1145  legrcb->AddEntry(hRABvsRbFDhigh_proj,"FD Upper unc.","l");
1146  legrcb->AddEntry(hRABvsRbFDlow_proj,"FD Lower unc.","l");
1147  legrcb->Draw();
1148  cRABptFDUnc->Update();
1149 
1150  TCanvas *RaaPlot = new TCanvas("RaaPlot","RAB vs pt, plot all");
1151  RaaPlot->SetTopMargin(0.085);
1152  RaaPlot->SetBottomMargin(0.1);
1153  RaaPlot->SetTickx();
1154  RaaPlot->SetTicky();
1155  TH2D *hRaaCanvas = new TH2D("hRaaCanvas"," R_{AB}(c) vs p_{T} (no Eloss hypothesis); p_{t} [GeV/c] ; R_{AA} prompt D",40,0.,40.,100,0.,3.0);
1156  hRaaCanvas->GetXaxis()->SetTitleSize(0.05);
1157  hRaaCanvas->GetXaxis()->SetTitleOffset(0.9);
1158  hRaaCanvas->GetYaxis()->SetTitleSize(0.05);
1159  hRaaCanvas->GetYaxis()->SetTitleOffset(0.9);
1160  hRaaCanvas->Draw();
1161  gRAB_Norm->SetFillStyle(1001);
1162  gRAB_Norm->SetFillColor(kGray+2);
1163  gRAB_Norm->Draw("2");
1164  TLine *line = new TLine(0.0172415,1.0,40.,1.0);
1165  line->SetLineStyle(2);
1166  line->Draw();
1167  hRABvsPt->SetMarkerColor(kBlue);
1168  hRABvsPt->SetMarkerColor(kBlue);
1169  hRABvsPt->SetMarkerStyle(21);
1170  hRABvsPt->SetMarkerSize(1.1);
1171  hRABvsPt->SetLineWidth(2);
1172  hRABvsPt->Draw("psame");
1173  gRAB_DataSystematics->SetLineColor(kBlue);
1174  gRAB_DataSystematics->SetLineWidth(3);
1175  gRAB_DataSystematics->SetLineWidth(2);
1176  gRAB_DataSystematics->SetFillColor(kRed);
1177  gRAB_DataSystematics->SetFillStyle(0);
1178  gRAB_DataSystematics->Draw("2");
1179  gRAB_FeedDownSystematics->SetFillColor(kViolet+1);
1180  gRAB_FeedDownSystematics->SetFillStyle(1001);
1181  gRAB_FeedDownSystematics->Draw("2");
1182  gRAB_ElossHypothesis->SetLineColor(kMagenta-7);
1183  gRAB_ElossHypothesis->SetFillColor(kMagenta-7);
1184  gRAB_ElossHypothesis->SetFillStyle(1001);
1185  gRAB_ElossHypothesis->Draw("2");
1186  hRABvsPt->Draw("psame");
1187  gRAB_DataSystematics->Draw("2");
1188  legrcb = new TLegend(0.5517241,0.6504237,0.8520115,0.8728814,NULL,"brNDC");
1189  legrcb->SetBorderSize(0);
1190  legrcb->SetTextSize(0.03389831);
1191  legrcb->SetLineColor(1);
1192  legrcb->SetLineStyle(1);
1193  legrcb->SetLineWidth(1);
1194  legrcb->SetFillColor(0);
1195  legrcb->SetFillStyle(1001);
1196  if(cc==k020) legrcb->AddEntry(hRABvsPt,"R_{AA} 0-20% CC","pe");
1197  else if(cc==k4080) legrcb->AddEntry(hRABvsPt,"R_{AA} 40-80% CC","pe");
1198  else legrcb->AddEntry(hRABvsPt,"R_{AA} and stat. unc.","pe");
1199  legrcb->AddEntry(gRAB_DataSystematics,"Syst. from data","f");
1200  legrcb->AddEntry(gRAB_ElossHypothesis,"Syst. from R_{AA}(B)","f");
1201  legrcb->AddEntry(gRAB_FeedDownSystematics,"Syst. from B feed-down","f");
1202  legrcb->Draw();
1203  TLatex* tc;
1204  TString system = "Pb-Pb #sqrt{s_{NN}}=2.76 TeV";
1205  if( cc==kpPb0100 || cc==kpPb020 || cc==kpPb1020 || cc==kpPb2040 || cc==kpPb4060 || cc==kpPb60100 ) system = "p-Pb #sqrt{s_{NN}}=5.023 TeV";
1206  if(decay==1) tc =new TLatex(0.18,0.82,Form("D^{0}, %s ",system.Data()));
1207  else if(decay==2) tc =new TLatex(0.18,0.82,Form("D^{+}, %s ",system.Data()));
1208  else if(decay==3) tc =new TLatex(0.18,0.82,Form("D^{*+}, %s ",system.Data()));
1209  else if(decay==4) tc =new TLatex(0.18,0.82,Form("D_{s}^{+}, %s ",system.Data()));
1210  else tc =new TLatex(0.18,0.82,Form("any (?) D meson, %s ",system.Data()));
1211  tc->SetNDC();
1212  tc->SetTextSize(0.038);
1213  tc->SetTextFont(42);
1214  tc->Draw();
1215  RaaPlot->Update();
1216 
1217 
1218  TCanvas *RaaPlotFDEloss = new TCanvas("RaaPlotFDEloss","RAB vs pt, plot FD & ElossUnc");
1219  RaaPlotFDEloss->SetTopMargin(0.085);
1220  RaaPlotFDEloss->SetBottomMargin(0.1);
1221  hRaaCanvas->Draw();
1222  line->Draw();
1223  hRABvsPt->Draw("psame");
1224  gRAB_FeedDownSystematics->SetFillColor(kViolet+1);
1225  gRAB_FeedDownSystematics->SetFillStyle(1001);
1226  gRAB_FeedDownSystematics->Draw("2");
1227  gRAB_ElossHypothesis->SetLineColor(kMagenta-7);
1228  gRAB_ElossHypothesis->SetFillColor(kMagenta-7);
1229  gRAB_ElossHypothesis->SetFillStyle(1001);
1230  gRAB_ElossHypothesis->Draw("2");
1231  gRAB_FeedDownSystematicsElossHypothesis->SetLineColor(kBlack);
1232  gRAB_FeedDownSystematicsElossHypothesis->SetFillStyle(0);
1233  gRAB_FeedDownSystematicsElossHypothesis->SetFillColor(kViolet+1);
1234  gRAB_FeedDownSystematicsElossHypothesis->Draw("2");
1235  hRABvsPt->Draw("psame");
1236  legrcb = new TLegend(0.6,0.6,0.9,0.9);
1237  legrcb->SetBorderSize(0);
1238  legrcb->SetTextSize(0.03389831);
1239  legrcb->SetLineColor(1);
1240  legrcb->SetLineStyle(1);
1241  legrcb->SetLineWidth(1);
1242  legrcb->SetFillColor(0);
1243  legrcb->SetFillStyle(1001);
1244  legrcb->AddEntry(hRABvsPt,"R_{PbPb} and stat. unc.","pe");
1245  legrcb->AddEntry(gRAB_ElossHypothesis,"Energy loss syst.","f");
1246  legrcb->AddEntry(gRAB_FeedDownSystematics,"Feed down syst.","f");
1247  legrcb->AddEntry(gRAB_FeedDownSystematicsElossHypothesis,"Feed down & Eloss syst.","f");
1248  legrcb->Draw();
1249  RaaPlotFDEloss->Update();
1250 
1251 
1252  TCanvas *RaaPlotGlob = new TCanvas("RaaPlotGlob","RAB vs pt, plot Global unc");
1253  RaaPlotGlob->SetTopMargin(0.085);
1254  RaaPlotGlob->SetBottomMargin(0.1);
1255  RaaPlotGlob->SetTickx();
1256  RaaPlotGlob->SetTicky();
1257  hRaaCanvas->Draw();
1258  line->Draw();
1259  hRABvsPt->Draw("psame");
1260  gRAB_DataSystematics->Draw("2");
1261  gRAB_FeedDownSystematicsElossHypothesis->Draw("2");
1262  gRAB_GlobalSystematics->SetLineColor(kRed);
1263  gRAB_GlobalSystematics->SetLineWidth(2);
1264  gRAB_GlobalSystematics->SetFillColor(kRed);
1265  gRAB_GlobalSystematics->SetFillStyle(3002);
1266  gRAB_GlobalSystematics->Draw("2");
1267  hRABvsPt->Draw("psame");
1268  legrcb = new TLegend(0.6,0.6,0.9,0.9);
1269  legrcb->SetBorderSize(0);
1270  legrcb->SetTextSize(0.03389831);
1271  legrcb->SetLineColor(1);
1272  legrcb->SetLineStyle(1);
1273  legrcb->SetLineWidth(1);
1274  legrcb->SetFillColor(0);
1275  legrcb->SetFillStyle(1001);
1276  legrcb->AddEntry(hRABvsPt,"R_{PbPb} and stat. unc.","pe");
1277  legrcb->AddEntry(gRAB_DataSystematics,"Data syst.","f");
1278  legrcb->AddEntry(gRAB_FeedDownSystematicsElossHypothesis,"Feed down & Eloss syst.","f");
1279  legrcb->AddEntry(gRAB_GlobalSystematics,"Global syst.","f");
1280  legrcb->Draw();
1281  RaaPlotGlob->Update();
1282 
1283 
1284 
1285  TCanvas *RaaPlotSimple = new TCanvas("RaaPlotSimple","RAB vs pt, plot Simple unc");
1286  RaaPlotSimple->SetTopMargin(0.085);
1287  RaaPlotSimple->SetBottomMargin(0.1);
1288  RaaPlotSimple->SetTickx();
1289  RaaPlotSimple->SetTicky();
1290  hRaaCanvas->Draw();
1291  line->Draw();
1292  hRABvsPt->Draw("psame");
1293  gRAB_GlobalSystematics->SetLineColor(kBlue);
1294  gRAB_GlobalSystematics->SetLineWidth(2);
1295  gRAB_GlobalSystematics->SetFillStyle(0);
1296  gRAB_GlobalSystematics->Draw("2");
1297  gRAB_Norm->Draw("2");
1298  hRABvsPt->Draw("psame");
1299  legrcb = new TLegend(0.5991379,0.6949153,0.8534483,0.8559322,NULL,"brNDC");
1300  legrcb->SetBorderSize(0);
1301  legrcb->SetTextSize(0.03389831);
1302  legrcb->SetLineColor(1);
1303  legrcb->SetLineStyle(1);
1304  legrcb->SetLineWidth(1);
1305  legrcb->SetFillColor(0);
1306  legrcb->SetFillStyle(1001);
1307  if(cc==k020) legrcb->AddEntry(hRABvsPt,"R_{AA} 0-20% CC","pe");
1308  else if(cc==k4080) legrcb->AddEntry(hRABvsPt,"R_{AA} 40-80% CC","pe");
1309  else legrcb->AddEntry(hRABvsPt,"R_{AA} and stat. unc.","pe");
1310  legrcb->AddEntry(gRAB_GlobalSystematics,"Systematics","f");
1311  legrcb->Draw();
1312  tc->Draw();
1313  RaaPlotSimple->Update();
1314 
1315 
1316  TCanvas *c = new TCanvas("c","");
1317  systematicsAB->DrawErrors();
1318  c->Update();
1319 
1320  TCanvas *cStatUnc = new TCanvas("cStatUnc","stat unc");
1321  cStatUnc->Divide(2,2);
1322  cStatUnc->cd(1);
1323  fhStatUncEffcSigmaAB_Raa->Draw("e");
1324  cStatUnc->cd(2);
1325  fhStatUncEffbSigmaAB_Raa->Draw("e");
1326  cStatUnc->cd(3);
1327  fhStatUncEffcFDAB_Raa->Draw("e");
1328  cStatUnc->cd(4);
1329  fhStatUncEffbFDAB_Raa->Draw("e");
1330  cStatUnc->Update();
1331 
1332  //
1333  // Write the information to a file
1334  //
1335  cout<<endl<< " Save results in the output file"<<endl<<endl;
1336  TFile * out = new TFile(outfile,"recreate");
1337 
1338  ntupleRAB->Write();
1339  hRABvsRcb->Write();
1340  hRABvsRb->Write();
1341  hRABCharmVsRBeautyVsPt->Write();
1342  for(Int_t j=0; j<=nbins; j++) hRCharmVsRBeauty[j]->Write();
1343  // for(Int_t j=0; j<=nbins; j++) hRCharmVsElossHypo[j]->Write();
1344  hRABvsPt->Write();
1345  hRABvsPt_DataSystematics->Write();
1346  gRAB_ElossHypothesis->Write();
1347  gRAB_FeedDownSystematics->Write();
1348  gRAB_fcFeedDownOnly->Write();
1349  gRAB_DataSystematics->Write();
1350  gRAB_DataSystematicsPP->Write();
1351  gSigmaPPSystTheory->Write();
1352  gRAB_DataSystematicsAB->Write();
1353  gRAB_Norm->Write();
1354  gRAB_FeedDownSystematicsElossHypothesis->Write();
1355  gRAB_GlobalSystematics->Write();
1356  if(isScaledAndExtrapRef && hCombinedReferenceFlag) hCombinedReferenceFlag->Write();
1357  systematicsPP->SetName("AliHFSystErrPP");
1358  systematicsPP->Write();
1359  systematicsAB->SetName("AliHFSystErrAA");
1360  systematicsAB->Write();
1361  out->Write();
1362 
1363 }
1364 
1365 //____________________________________________________________
1366 Bool_t PbPbDataSyst(AliHFSystErr *syst, Double_t pt, Int_t cc, Double_t &dataSystUp, Double_t &dataSystDown)
1367 {
1368 
1369  Double_t err=0., errUp=1., errDown=1.;
1370  Bool_t isOk=false;
1371  Double_t pidunc=0.;
1372 
1373  err = syst->GetTotalSystErr(pt)*syst->GetTotalSystErr(pt);
1374  errUp = err ;
1375  errDown = err ;
1376 
1377  // Apply an asymetric PID uncertainty for 2010 PbPb data only
1378  if( syst->GetRunNumber()==10 && syst->GetCollisionType()==1 ) {
1379  if( cc==k07half || cc==k020 || cc==k010 || cc==k1020 || cc==k2040 ) {
1380  if(pt<6) pidunc = 0.15;
1381  else pidunc = 0.05;
1382  errUp = err + pidunc*pidunc - syst->GetPIDEffErr(pt)*syst->GetPIDEffErr(pt);
1383  isOk = true;
1384  }
1385  else if ( cc==k3050 || cc==k4080 || cc==k4060 || cc==k6080 ){
1386  if(pt<3.1) pidunc = 0.10;
1387  else pidunc = 0.05;
1388  errUp = err + pidunc*pidunc - syst->GetPIDEffErr(pt)*syst->GetPIDEffErr(pt);
1389  isOk = true;
1390  }
1391  }
1392  else { isOk = true; }
1393 
1394  dataSystUp = TMath::Sqrt(errUp);
1395  dataSystDown = TMath::Sqrt(errDown);
1396 
1397  return isOk;
1398 }
particularity
Definition: HFPtSpectrum.C:48
const Color_t cc[]
Definition: DrawKs.C:1
double Double_t
Definition: External.C:58
BFDSubtrMethod
Definition: HFPtSpectrum.C:45
energy
Definition: HFPtSpectrum.C:44
centrality
void SetCentrality(TString centrality)
Definition: AliHFSystErr.h:75
TCanvas * c
Definition: TestFitELoss.C:172
void SetIsPbPb2010EnergyScan(Bool_t flag)
Definition: AliHFSystErr.h:104
Bool_t elossFDQuadSum
Int_t GetRunNumber() const
Definition: AliHFSystErr.h:63
void SetIsBDTAnalysis(Bool_t flag)
Definition: AliHFSystErr.h:99
int Int_t
Definition: External.C:63
Definition: External.C:204
float Float_t
Definition: External.C:68
RaavsEP
Definition: HFPtSpectrum.C:46
centestimator
Definition: External.C:252
Definition: External.C:228
Definition: External.C:212
rapidity
Definition: HFPtSpectrum.C:47
void SetIsPass4Analysis(Bool_t flag)
Definition: AliHFSystErr.h:87
void HFPtSpectrumRaa(const char *ppfile="HFPtSpectrum_D0Kpi_method2_rebinnedth_230311_newsigma.root", const char *ABfile="HFPtSpectrum_D0Kpi_PbPbcuts_method2_rebinnedth_230311_newsigma.root", const char *outfile="HFPtSpectrumRaa.root", Int_t decay=1, Double_t sigmaABCINT1B=54.e9, Int_t fdMethod=kNb, Int_t cc=kpp, Int_t Energy=k276, Double_t MinHypo=1./3., Double_t MaxHypo=3.0, Double_t MaxRb=6.0, Bool_t isRbHypo=false, Double_t CentralHypo=1.0, Int_t ccestimator=kV0M, Bool_t isUseTaaForRaa=true, const char *shadRbcFile="", Int_t nSigmaShad=3.0, Int_t isRaavsEP=kPhiIntegrated, Bool_t isScaledAndExtrapRef=kFALSE, Int_t rapiditySlice=kdefault, Int_t analysisSpeciality=kTopological)
Double_t GetTotalSystErr(Double_t pt, Double_t feeddownErr=0) const
void SetIsLowPtAnalysis(Bool_t flag)
Definition: AliHFSystErr.h:83
Double_t ExtractFDSyst(Double_t total, Double_t fd)
Double_t ptprintout
decay
Definition: HFPtSpectrum.C:41
void Init(Int_t decay)
Function to initialize the variables/histograms.
void SetIspPb2011RapidityScan(Bool_t flag)
Definition: AliHFSystErr.h:114
Bool_t printout
void SetCollisionType(Int_t type)
Definition: AliHFSystErr.h:66
Double_t NormABUnc
Bool_t PbPbDataSyst(AliHFSystErr *syst, Double_t pt, Int_t cc, Double_t &dataSystUp, Double_t &dataSystDown)
const Int_t nbins
bool Bool_t
Definition: External.C:53
void DrawErrors(TGraphAsymmErrors *grErrFeeddown=0) const
Int_t FindGraphBin(TGraphAsymmErrors *gr, Double_t pt)
void SetRapidity(TString rapidity)
Settings of rapidity ranges for pPb 0-100% CC.
Definition: AliHFSystErr.h:110
Int_t GetCollisionType() const
Definition: AliHFSystErr.h:72
void SetRunNumber(Int_t number)
Definition: AliHFSystErr.h:59
Double_t NormPPUnc
Double_t GetPIDEffErr(Double_t pt) const