AliPhysics  068200c (068200c)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliDJetTTreeReader.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 1998-2016, ALICE Experiment at CERN, All rights reserved. *
3  * *
4  * Author: The ALICE Off-line Project. *
5  * Contributors are mentioned in the code where appropriate. *
6  * *
7  * Permission to use, copy, modify and distribute this software and its *
8  * documentation strictly for non-commercial purposes is hereby granted *
9  * without fee, provided that the above copyright notice appears in all *
10  * copies and that both the copyright notice and this permission notice *
11  * appear in the supporting documentation. The authors make no claims *
12  * about the suitability of this software for any purpose. It is *
13  * provided "as is" without express or implied warranty. *
14  **************************************************************************/
15 
16 #include <TFile.h>
17 #include <TH1D.h>
18 #include <TH2D.h>
19 #include <TChain.h>
20 #include <Riostream.h>
21 
23 
24 #include "AliDJetTTreeReader.h"
25 
27 ClassImp(AliDJetTTreeReader);
29 
35  fInputFileNames(),
36  fTreeName(),
37  fDBranchName(),
38  fJetBranchName(),
39  fmassmin(0),
40  fmassmax(0),
41  fmasswidth(0)
42 {
43 }
44 
50  AliDJetVReader(source),
51  fInputFileNames(source.fInputFileNames),
52  fTreeName(source.fTreeName),
53  fDBranchName(source.fDBranchName),
54  fJetBranchName(source.fJetBranchName),
55  fmassmin(source.fmassmin),
56  fmassmax(source.fmassmax),
57  fmasswidth(source.fmasswidth)
58 {
59 }
60 
65 {
66 }
67 
73 {
74  Printf("Generating the chain...");
75  TChain* chain = new TChain(fTreeName);
76  for (std::string fname : fInputFileNames) {
77  chain->Add(fname.c_str());
78  }
79  return chain;
80 }
81 
86 {
87  std::cout << "Extracting input mass plot: " << fpTmin << " to " << fpTmax << std::endl;
88 
89  TTree *tree = GenerateChain();
90  if (!tree) {
91  std::cout << "Error in setting the tree/branch names! Exiting..." << std::endl;
92  return kFALSE;
93  }
94 
97  tree->SetBranchAddress(fDBranchName,&brD);
98  tree->SetBranchAddress(fJetBranchName,&brJet);
99 
100  TString hname = TString::Format("invMass_JetPt_%.0f_%.0f", fpTmin, fpTmax);
102  fMassPlot->Sumw2();
103 
104  for (int k = 0; k < tree->GetEntries(); k++) {
105  tree->GetEntry(k);
106  if (brJet->fEta < -0.5 || brJet->fEta >= 0.5) continue;
107  if (brJet->fPt < fpTmin || brJet->fPt >= fpTmax) continue;
108  for (int j = 0; j < fnDbins; j++) {
109  if (brD->fPt < fDbinpTedges[j] || brD->fPt >= fDbinpTedges[j + 1]) continue;
110  fMassPlot->Fill(brD->fInvMass, 1. / fDEffValues[j]);
111  }//end of D-meson pT bin loop
112  }
113 
114  return kTRUE;
115 }
116 
121 {
122  std::cout << "Extracting input mass plot: " << fpTmin << " to " << fpTmax << std::endl;
123 
124  TTree *tree = GenerateChain();
125  if (!tree) {
126  std::cout << "Error in setting the tree/branch names! Exiting..." << std::endl;
127  return kFALSE;
128  }
129 
132  tree->SetBranchAddress(fDBranchName,&brD);
133  tree->SetBranchAddress(fJetBranchName,&brJet);
134 
135  TString hname = TString::Format("invMass_DPt_%.0f_%.0f", fpTmin, fpTmax);
137  fMassPlot->Sumw2();
138 
139  fMassVsJetPtPlot = new TH2D("hInvMassJetPt", "hInvMassJetPt", (fmassmax-fmassmin) / fmasswidth, fmassmin, fmassmax, fnJetPtbins, fJetPtBinEdges);
140  fMassVsJetPtPlot->Sumw2();
141 
142  fMassVsJetzPlot = new TH2D("hInvMassJetz", "hInvMassJetz", (fmassmax-fmassmin) / fmasswidth, fmassmin, fmassmax, fnJetzbins, fJetzBinEdges);
143  fMassVsJetzPlot->Sumw2();
144 
145  for (int k = 0; k < tree->GetEntries(); k++) {
146  tree->GetEntry(k);
147  if (brJet->fEta < -0.5 || brJet->fEta >= 0.5) continue;
148  if (brJet->fPt < fJetPtBinEdges[0] || brJet->fPt >= fJetPtBinEdges[fnJetPtbins]) continue;
149  if (brD->fPt < fpTmin || brD->fPt >= fpTmax) continue;
150  fMassPlot->Fill(brD->fInvMass);
151  fMassVsJetPtPlot->Fill(brD->fInvMass, brJet->fPt);
152  fMassVsJetzPlot->Fill(brD->fInvMass, brJet->fZ);
153  }
154 
155  return kTRUE;
156 }
Double_t fpTmax
pT upper edge of mass plot to evaluate variations of yields
Double_t * fJetzBinEdges
Jet z bin edges to be used for spectrum.
Lightweight class that encapsulates D meson jets.
std::vector< std::string > fInputFileNames
Name of input file.
Double32_t fInvMass
Invariant mass of the D0 meson candidate in GeV/c2.
Double_t * fDbinpTedges
D-meson pt bin edges values.
Implementation of an abstract class to read the invariant mass histograms used to extract the raw yie...
Int_t fnJetzbins
Number of jet z bins to be used for spectrum.
TH1D * fMassPlot
!Mass spectra to be fitted
Double32_t fPt
Transverse momentum of the jet in GeV/c.
Double_t fmassmin
Mass lower edge of inv.mass plots.
TString fDBranchName
Name of input branch for D meson.
Int_t fnJetPtbins
Number of jet pT bins to be used for spectrum.
Double32_t fPt
Transverse momentum of the D meson in GeV/c.
Double_t * fDEffValues
D-meson efficiency values.
Double_t fmasswidth
Mass plots bin width.
Definition: External.C:228
Definition: External.C:212
Int_t fnDbins
Number of D-meson pT bins (for eff scaling)
TString fTreeName
Name of input TTree.
Declaration of class AliDJetTTreeReader.
Double_t fpTmin
pT lower edge of mass plot to evaluate variations of yields
Bool_t ExtractInputMassPlotEffScale()
UInt_t fMassRebin
Rebin the mass histogram axis.
Implementation of an abstract class to read the invariant mass histograms used to extract the raw yie...
TString fJetBranchName
Name of input branch for jet.
Double_t fmassmax
Mass upper edge of inv.mass plots.
TH2D * fMassVsJetPtPlot
!Mass vs jet pt (SB method)
bool Bool_t
Definition: External.C:53
Double_t * fJetPtBinEdges
Jet pT bin edges to be used for spectrum.
Lightweight class that encapsulates D0.
TH2D * fMassVsJetzPlot
!Mass vs jet z (SB method)