DAS  3.0
Das Analysis System
ControlPlots

#include <ControlPlots.h>

Public Member Functions

 ControlPlots (TString Name)
 
void operator() (const std::vector< RecJet > &recjets, const double &evW, size_t iJEC=0, size_t iWgt=0)
 
void operator() (const std::vector< GenJet > &genjets, const double &evW, size_t iWgt=0)
 
void Write (TDirectory *D)
 

Public Attributes

const TString name
 
std::optional< TH1F > genpt
 
std::optional< TH1F > recpt
 
std::optional< TH2F > genpt_y
 
std::optional< TH2F > recpt_y
 
std::optional< TH2F > genMjj_y
 
std::optional< TH2F > recMjj_y
 
std::optional< TH3F > genMjj_ysyb
 
std::optional< TH3F > recMjj_ysyb
 
std::optional< TH3F > genptav_ysyb
 
std::optional< TH3F > recptav_ysyb
 
std::optional< TH2F > genpt0_N
 
std::optional< TH2F > recpt0_N
 
std::optional< TH2F > genpt_n
 
std::optional< TH2F > recpt_n
 

Static Public Attributes

static bool isMC = false
 
static bool verbose = false
 

Constructor & Destructor Documentation

◆ ControlPlots()

ControlPlots ( TString  Name)
Todo:
TH3 pt eta rho plots TH2 pt_N TH2 ptmax_Dphi
11  :
12  name(Name)
13 {
14  if (verbose)
15  cout << __func__ << '\t' << Name << endl;
16  recpt = TH1F(Name + "recpt" , ";p_{T};N^{j}_{eff}", nPtBins, pt_edges.data());
17  recpt_y = TH2F(Name + "recpt_y", ";p_{T};y;N^{j}_{eff}", nPtBins, pt_edges.data(), nYbins, y_edges.data());
18  recMjj_y = TH2F(Name + "recMjj_y", ";M_{jj};y;N^{j}_{eff}", nMjjBins, Mjj_edges.data(), nYbins, y_edges.data());
19  recMjj_ysyb = TH3F(Name + "recMjj_ysyb", ";M_{jj};y_{b};y^{*}", nMjjBins, Mjj_edges.data(), nYbins, y_edges.data(), nYbins, y_edges.data());
20  recptav_ysyb = TH3F(Name + "recptav_ysyb", ";p_{T}^{av};y_{b};y^{*}", nPtDijetBins, pt_av_edges.data(), nYbins, y_edges.data(), nYbins, y_edges.data());
21  recpt0_N = TH2F(Name + "recpt0_N", ";p^{leading}_{T};N_{jets};N^{j}_{eff}", nPtBins, pt_edges.data(), maxMult, n_edges.data());
22  recpt_n = TH2F(Name + "recpt_n", ";p_{T};n^{th} jet;N^{j}_{eff}", nPtBins, pt_edges.data(), maxMult, n_edges.data());
23  if (!isMC) return;
24  genpt = TH1F(Name + "genpt" , ";p_{T};N^{j}_{eff}", nPtBins, pt_edges.data());
25  genpt_y = TH2F(Name + "genpt_y", ";p_{T};y;N^{j}_{eff}", nPtBins, pt_edges.data(), nYbins, y_edges.data());
26  genMjj_y = TH2F(Name + "genMjj_y", ";M_{jj};y;N^{j}_{eff}", nMjjBins, Mjj_edges.data(), nYbins, y_edges.data());
27  genMjj_ysyb = TH3F(Name + "genMjj_ysyb", ";M_{jj};y_{b};y^{*}", nMjjBins, Mjj_edges.data(), nYbins, y_edges.data(), nYbins, y_edges.data());
28  genptav_ysyb = TH3F(Name + "genptav_ysyb", ";p_{T}^{av};y_{b};y^{*}", nPtDijetBins, pt_av_edges.data(), nYbins, y_edges.data(), nYbins, y_edges.data());
29  genpt0_N = TH2F(Name + "genpt0_N", ";p^{leading}_{T};N_{jets};N^{j}_{eff}", nPtBins, pt_edges.data(), maxMult, n_edges.data());
30  genpt_n = TH2F(Name + "genpt_n", ";p_{T};n^{th} jet;N^{j}_{eff}", nPtBins, pt_edges.data(), maxMult, n_edges.data());
31 }

Member Function Documentation

◆ operator()() [1/2]

void operator() ( const std::vector< GenJet > &  genjets,
const double &  evW,
size_t  iWgt = 0 
)
34 {
35  assert(isMC);
36 
37  // inclusive jet
38  for (const auto& genjet: genjets) {
39  genpt ->Fill(genjet.p4.Pt(), evW * genjet.weights.at(iWgt));
40  genpt_y->Fill(genjet.p4.Pt(), genjet.AbsRap(), evW * genjet.weights.at(iWgt));
41  }
42 
43  auto inYacceptance = [](const auto& genjet) { return genjet.AbsRap() < ymax; };
44  size_t N = count_if(begin(genjets), end(genjets), inYacceptance);
45 
46  if (N > 0) {
47  auto genjet0 = find_if(begin(genjets), end(genjets), inYacceptance);
48  genpt0_N->Fill(genjet0->p4.Pt(), N, evW * genjet0->weights.front());
49  }
50 
51  for (size_t i = 0; i < genjets.size(); ++i) {
52  const auto& genjet = genjets.at(i);
53  if (genjet.AbsRap() >= ymax) continue;
54  genpt_n->Fill(genjet.p4.Pt(), i+1, evW * genjet.weights.front());
55  }
56 
57  // dijet
58  if (genjets.size() < 2) return;
59  const auto& j0 = genjets.at(0),
60  j1 = genjets.at(1);
61  if ( j0.p4.Pt() < 100. || j1.p4.Pt() < 50. ) return;
62  auto dijet = j0 + j1;
63  auto ymax = max(j0.AbsRap(), j1.AbsRap());
64  auto djW = j0.weights.at(iWgt) * j1.weights.at(iWgt);
65  genMjj_y->Fill(dijet.CorrP4().M(), ymax, evW * djW);
66 
67  auto Y0 = j0.Rapidity(),
68  Y1 = j1.Rapidity();
69  auto ystar = 0.5*abs(Y0-Y1),
70  yboost = 0.5*abs(Y0+Y1);
71  genMjj_ysyb->Fill(dijet.CorrP4().M(), yboost, ystar, evW * djW);
72  auto ptav = 0.5*(j0.p4.Pt() + j1.p4.Pt());
73  genptav_ysyb->Fill(ptav, yboost, ystar, evW * djW);
74 
75 }

◆ operator()() [2/2]

void operator() ( const std::vector< RecJet > &  recjets,
const double &  evW,
size_t  iJEC = 0,
size_t  iWgt = 0 
)
78 {
79  // inclusive jet
80  for (const auto& recjet: recjets) {
81  recpt ->Fill(recjet.CorrPt(iJEC), evW * recjet.weights.at(iWgt));
82  recpt_y->Fill(recjet.CorrPt(iJEC), recjet.AbsRap(), evW * recjet.weights.at(iWgt));
83  }
84 
85  auto inYacceptance = [](const auto& recjet) { return recjet.AbsRap() < ymax; };
86  size_t N = count_if(begin(recjets), end(recjets), inYacceptance);
87 
88  if (N > 0) {
89  auto recjet0 = find_if(begin(recjets), end(recjets), inYacceptance);
90  recpt0_N->Fill(recjet0->p4.Pt(), N, evW * recjet0->weights.front());
91  }
92 
93 
94  for (size_t i = 0; i < recjets.size(); ++i) {
95  const auto& recjet = recjets.at(i);
96  if (recjet.AbsRap() >= ymax) continue;
97  recpt_n->Fill(recjet.p4.Pt(), i+1, evW * recjet.weights.front());
98  }
99  // dijet
100  if (recjets.size() < 2) return;
101  const auto& j0 = recjets.at(0),
102  j1 = recjets.at(1);
103  if ( j0.CorrPt(iJEC) < 100. || j1.CorrPt(iJEC) < 50. ) return;
104  auto ymax = max(j0.AbsRap(), j1.AbsRap());
105  auto djW = j0.weights.at(iWgt) * j1.weights.at(iWgt);
106  auto dijet = j0 + j1;
107  recMjj_y->Fill(dijet.CorrP4().M(), ymax, evW * djW);
108 
109  auto Y0 = j0.Rapidity(),
110  Y1 = j1.Rapidity();
111  auto ystar = 0.5*abs(Y0-Y1),
112  yboost = 0.5*abs(Y0+Y1);
113  recMjj_ysyb->Fill(dijet.CorrP4().M(), yboost, ystar, evW * djW);
114  auto ptav = dijet.HT();
115  recptav_ysyb->Fill(ptav, yboost, ystar, evW * djW);
116 }

◆ Write()

void Write ( TDirectory *  D)
119 {
120  assert(D);
121  auto d = D->mkdir(name);
122  d->cd();
123  vector<TH1*> hs {&*recpt, dynamic_cast<TH1*>(&*recpt_y), dynamic_cast<TH1*>(&*recMjj_y), dynamic_cast<TH1*>(&*recpt0_N), dynamic_cast<TH1*>(&*recpt_n), dynamic_cast<TH1*>(&*recMjj_ysyb), dynamic_cast<TH1*>(&*recptav_ysyb)};
124  if (isMC) {
125  vector<TH1*> hs2 {&*genpt, dynamic_cast<TH1*>(&*genpt_y), dynamic_cast<TH1*>(&*genMjj_y), dynamic_cast<TH1*>(&*genpt0_N), dynamic_cast<TH1*>(&*genpt_n), dynamic_cast<TH1*>(&*genMjj_ysyb), dynamic_cast<TH1*>(&*genptav_ysyb)};
126  hs.insert(hs.end(), hs2.begin(), hs2.end());
127  }
128  for (auto& h: hs) {
129  TString n = h->GetName();
130  n.ReplaceAll(name,"");
131  h->SetDirectory(d);
132  h->Write(n);
133  }
134 }

Member Data Documentation

◆ genMjj_y

std::optional<TH2F> genMjj_y

◆ genMjj_ysyb

std::optional<TH3F> genMjj_ysyb

◆ genpt

std::optional<TH1F> genpt

◆ genpt0_N

std::optional<TH2F> genpt0_N

◆ genpt_n

std::optional<TH2F> genpt_n

◆ genpt_y

std::optional<TH2F> genpt_y

◆ genptav_ysyb

std::optional<TH3F> genptav_ysyb

◆ isMC

bool isMC = false
static

◆ name

const TString name

◆ recMjj_y

std::optional<TH2F> recMjj_y

◆ recMjj_ysyb

std::optional<TH3F> recMjj_ysyb

◆ recpt

std::optional<TH1F> recpt

◆ recpt0_N

std::optional<TH2F> recpt0_N

◆ recpt_n

std::optional<TH2F> recpt_n

◆ recpt_y

std::optional<TH2F> recpt_y

◆ recptav_ysyb

std::optional<TH3F> recptav_ysyb

◆ verbose

bool verbose = false
static

The documentation for this struct was generated from the following files:
DAS::ControlPlots::genMjj_y
std::optional< TH2F > genMjj_y
Definition: ControlPlots.h:26
DAS::y_edges
static const std::vector< double > y_edges
Definition: binnings.h:36
DAS::ControlPlots::recpt
std::optional< TH1F > recpt
Definition: ControlPlots.h:24
DAS::ControlPlots::genptav_ysyb
std::optional< TH3F > genptav_ysyb
Definition: ControlPlots.h:28
DAS::PhysicsObject::p4
FourVector p4
raw four-momentum directly after reconstruction
Definition: PhysicsObject.h:50
DAS::n_edges
static const std::vector< double > n_edges
Definition: binnings.h:37
DAS::ControlPlots::genpt_y
std::optional< TH2F > genpt_y
Definition: ControlPlots.h:25
DAS::pt_edges
static const std::vector< double > pt_edges
Definition: binnings.h:33
DAS::pt_av_edges
static const std::vector< double > pt_av_edges
Definition: binnings.h:35
DAS::nMjjBins
static const int nMjjBins
Definition: binnings.h:40
DAS::maxMult
static const int maxMult
Definition: binnings.h:42
DAS::ControlPlots::genMjj_ysyb
std::optional< TH3F > genMjj_ysyb
Definition: ControlPlots.h:27
DAS::ymax
static const double ymax
Definition: binnings.h:49
DAS::ControlPlots::recpt_y
std::optional< TH2F > recpt_y
Definition: ControlPlots.h:25
DAS::ControlPlots::verbose
static bool verbose
Definition: ControlPlots.h:20
Ntupliser_cfg.genjets
genjets
Definition: Ntupliser_cfg.py:272
DAS::PhysicsObject::AbsRap
float AbsRap(const Uncertainties::Variation &=Uncertainties::nominal) const final
absolute rapidity
Definition: PhysicsObject.h:58
DAS::ControlPlots::genpt
std::optional< TH1F > genpt
Definition: ControlPlots.h:24
recjet
DAS::RecJet recjet
Definition: classes.h:15
DAS::ControlPlots::recpt0_N
std::optional< TH2F > recpt0_N
Definition: ControlPlots.h:29
DAS::ControlPlots::isMC
static bool isMC
Definition: ControlPlots.h:19
DAS::ControlPlots::genpt0_N
std::optional< TH2F > genpt0_N
Definition: ControlPlots.h:29
Ntupliser_cfg.recjets
recjets
Definition: Ntupliser_cfg.py:273
DAS::Mjj_edges
static const std::vector< double > Mjj_edges
Definition: binnings.h:34
DAS::PhysicsObject::CorrPt
float CorrPt(size_t i=0) const
corrected transverse momentum
Definition: PhysicsObject.h:55
DAS::nPtBins
static const int nPtBins
Definition: binnings.h:39
DAS::ControlPlots::recMjj_ysyb
std::optional< TH3F > recMjj_ysyb
Definition: ControlPlots.h:27
DAS::ControlPlots::genpt_n
std::optional< TH2F > genpt_n
Definition: ControlPlots.h:30
DAS::PhysicsObject::weights
Weights weights
object weights
Definition: PhysicsObject.h:52
genjet
DAS::GenJet genjet
Definition: classes.h:14
DAS::nPtDijetBins
static const int nPtDijetBins
Definition: binnings.h:41
DAS::ControlPlots::recptav_ysyb
std::optional< TH3F > recptav_ysyb
Definition: ControlPlots.h:28
DAS::ControlPlots::name
const TString name
Definition: ControlPlots.h:22
DAS::ControlPlots::recpt_n
std::optional< TH2F > recpt_n
Definition: ControlPlots.h:30
DAS::nYbins
static const int nYbins
Definition: binnings.h:38
DAS::ControlPlots::recMjj_y
std::optional< TH2F > recMjj_y
Definition: ControlPlots.h:26