DAS  3.0
Das Analysis System
RivetNtupliser
+ Inheritance diagram for RivetNtupliser:
+ Collaboration diagram for RivetNtupliser:

Public Member Functions

 RIVET_DEFAULT_ANALYSIS_CTOR (RivetNtupliser)
 
void init ()
 
void analyze (const Event &event)
 
void finalize ()
 

Private Attributes

DAS::GenEventevent_
 
vector< DAS::GenJet > * genJets
 
vector< DAS::RecJet > * recJets_vec
 
TFile * file
 
TTree * tree
 
pt::ptree config
 
double AbsEtaMin_V
 
double AbsEtaMax_V
 
double pTJetMin_V
 
double pTJetMax_V
 
bool detector_bool
 

Member Function Documentation

◆ analyze()

void analyze ( const Event &  event)
inline

Get FastJet final states genJets and recoJets, apply selections, fill ntuples with jet observables and weights.

120  {
121  event_->clear();
122  genJets->clear();
123  if (detector_bool) {
124  recJets_vec->clear();
125  }
126  // one hepmc weight per event, in this case the phase space bias weight
127  const auto& weights = event.weights();
128 
129  for (size_t i=0; i < weights.size() ; ++i) {
130  float w = weights[i];
131  event_->weights.push_back(DAS::Weight{w,0});
132  }
133 
134  // genJets Fastjet projections
135  const FastJets& fj = apply<FastJets>(event, "Jets");
136  // const Jets& RivetJets = fj.jets(Cuts::pT > pTJetMin_V*GeV && Cuts::abseta < AbsEtaMax_V);
137  const Jets& RivetJets = fj.jets(Cuts::abseta < AbsEtaMax_V);
138 
139  for (const Jet& j : RivetJets) {
140  DAS::GenJet genJet;
141  genJet.p4.SetPt(j.pT());
142  genJet.p4.SetEta(j.eta());
143  genJet.p4.SetPhi(j.phi());
144  genJet.p4.SetM(j.mass());
145  genJets->push_back(genJet);
146  }
147 
148  // // SmearedJets behaves like a FastJets objects e.g. but they are
149 
150  // recoJets Fastjet projections
151  if (detector_bool) {
152 
153  const SmearedJets& ApplySmearedJets = apply<SmearedJets>(event,"RecoJetDeclaration");
154 
155  const auto& SmearedJetsPostCuts = ApplySmearedJets.jets(Cuts::pT > pTJetMin_V*GeV && Cuts::abseta < AbsEtaMax_V);
156 
157  for (const auto& j : SmearedJetsPostCuts) {
158  DAS::RecJet recJet;
159  recJet.p4.SetPt(j.pT());
160  recJet.p4.SetEta(j.eta());
161  recJet.p4.SetPhi(j.phi());
162  recJet.p4.SetM(j.mass());
163  recJets_vec->push_back(recJet);
164  }
165 
166  }
167 
168  tree->Fill();
169 
170  }

◆ finalize()

void finalize ( )
inline
173  {
174 
175  tree->Write();
176  file->Close();
177 
178  }

◆ init()

void init ( )
inline

Read .json configuration values, set up ttree and branches, declare fastjet finalstates.

Get key-value pairs

49  {
50  const char* fname = getenv("DAS_RIVET_OUTPUT");
51  assert(fname != nullptr);
52  cout << "output filename: " << fname << endl;
53  file = TFile::Open(fname, "RECREATE");
54  tree = new TTree("events","events");
55  // Read .json configurations
56  const char* cfgname = getenv("DAS_RIVET_CONFIG");
57  assert(cfgname != nullptr);
58 
59  fs::path p_cfg = cfgname;
60  pt::ptree config;
61  pt::read_json(p_cfg.string(), config);
62 
63  // Read json Config values
64  float radius = config.get<float>("radius");
65  detector_bool = config.get<bool>("detector");
66 
67  cout << "\nradius: " << radius << endl;
68  cout << "\nreco jets?: " << detector_bool << endl;
69 
70  auto &selection = config.get_child("selection");
71 
73 
74  auto &AbsEtaMin_K = selection.get_child("AbsEtaMin");
75  AbsEtaMin_V = AbsEtaMin_K.get_value<float>();
76 
77  auto &AbsEtaMax_K = selection.get_child("AbsEtaMax");
78  AbsEtaMax_V = AbsEtaMax_K.get_value<double>();
79 
80  auto &pTJetMin_K = selection.get_child("pTJetMin");
81  pTJetMin_V = pTJetMin_K.get_value<float>();
82 
83  auto &pTJetMax_K = selection.get_child("pTJetMax");
84  pTJetMax_V = pTJetMax_K.get_value<float>();
85  // Print out the configurations
86  cout << "\n\tAbsEtaMin: " << AbsEtaMin_V << "\tAbsEtaMax: " << AbsEtaMax_V << "\tpTJetMin: " << pTJetMin_V << "\tpTJetMax: " << pTJetMax_V << endl;
87 
88  // Declare final state FastJets
89  const FinalState finalstate;
90 
91  // Set up ttree branches
92  event_ = new DAS::GenEvent;
93  tree->Branch("genEvent", &event_);
94  genJets = new vector<DAS::GenJet>;
95  genJets->reserve(25);
96  tree->Branch("genJets", &genJets);
97  // Construct recoJets FastJets
98  // if "detector":true in Rivet.json, then save recoJets as smearedjets
99 
100  if (detector_bool) {
101  recJets_vec = new vector<DAS::RecJet>;
102  recJets_vec->reserve(25);
103  tree->Branch("recJets", &recJets_vec);
104  }
105  else {
106  recJets_vec = nullptr;
107  }
108 
109 
110  // Construct genJets FastJets (declare Jets as a type)
111  declare(FastJets(finalstate, FastJets::ANTIKT, radius), "Jets");
112  if (detector_bool) {
113  declare(SmearedJets(FastJets(finalstate, FastJets::ANTIKT, radius), JET_SMEAR_CMS_RUN2), "RecoJetDeclaration");
114  }
115 
116  }

◆ RIVET_DEFAULT_ANALYSIS_CTOR()

RIVET_DEFAULT_ANALYSIS_CTOR ( RivetNtupliser  )

Constructor.

Member Data Documentation

◆ AbsEtaMax_V

double AbsEtaMax_V
private

◆ AbsEtaMin_V

double AbsEtaMin_V
private

◆ config

pt::ptree config
private

◆ detector_bool

bool detector_bool
private

◆ event_

DAS::GenEvent* event_
private

◆ file

TFile* file
private

output file, tree, json config, and parameters

output file for ntuple

◆ genJets

vector<DAS::GenJet>* genJets
private

◆ pTJetMax_V

double pTJetMax_V
private

◆ pTJetMin_V

double pTJetMin_V
private

◆ recJets_vec

vector<DAS::RecJet>* recJets_vec
private

◆ tree

TTree* tree
private

The documentation for this class was generated from the following file:
selection
Di< const Jet, const Jet > selection(const TTreeReaderArray< Jet > &jets, const Uncertainties::Variation &v=Uncertainties::nominal)
Mueller-Navelet jet selection.
Definition: MNjets.cc:24
jercExample.fname
fname
Definition: jercExample.py:88
DAS::PhysicsObject::p4
FourVector p4
raw four-momentum directly after reconstruction
Definition: PhysicsObject.h:47
Rivet::RivetNtupliser::detector_bool
bool detector_bool
Definition: RivetNtupliser.cc:39
Rivet::RivetNtupliser::pTJetMax_V
double pTJetMax_V
Definition: RivetNtupliser.cc:38
Rivet::RivetNtupliser::file
TFile * file
output file, tree, json config, and parameters
Definition: RivetNtupliser.cc:34
DAS::RecJet
Definition: Jet.h:37
DAS::AbstractEvent::weights
Weights weights
e.g. cross section normalisation
Definition: Event.h:23
Rivet::RivetNtupliser::event_
DAS::GenEvent * event_
Definition: RivetNtupliser.cc:29
DAS::JetEnergy::w
static const float w
Definition: common.h:51
Rivet::RivetNtupliser::pTJetMin_V
double pTJetMin_V
Definition: RivetNtupliser.cc:38
DAS::GenJet
class GenJet
Definition: Jet.h:9
Rivet::RivetNtupliser::tree
TTree * tree
Definition: RivetNtupliser.cc:35
Rivet::RivetNtupliser::AbsEtaMin_V
double AbsEtaMin_V
Definition: RivetNtupliser.cc:38
Rivet::RivetNtupliser::genJets
vector< DAS::GenJet > * genJets
Definition: RivetNtupliser.cc:30
DAS::Weight
Definition: Weight.h:16
Rivet::RivetNtupliser::config
pt::ptree config
Definition: RivetNtupliser.cc:37
DAS::GenEvent::clear
void clear()
to clear for each new event in n-tupliser
Definition: Event.cc:9
Rivet::RivetNtupliser::AbsEtaMax_V
double AbsEtaMax_V
Definition: RivetNtupliser.cc:38
weights
DAS::Weights weights
Definition: classes.h:12
Ntupliser_cfg.radius
dictionary radius
Definition: Ntupliser_cfg.py:35
Rivet::RivetNtupliser::recJets_vec
vector< DAS::RecJet > * recJets_vec
Definition: RivetNtupliser.cc:31
DAS::GenEvent
Definition: Event.h:38