DAS  3.0
Das Analysis System
DAS::Jets Namespace Reference

Functions

template<class Jet >
bool keepEvent (const vector< Jet > &jets, float pt1, float pt2)
 
void applyDijetSkim (const vector< fs::path > &inputs, const fs::path &output, const pt::ptree &config, const int steering, const DT::Slice slice={1, 0})
 
template<class Jet , class Dimuon >
bool keepEvent (const vector< Jet > &jets, const Dimuon &dimuon, float zpt, float jpt)
 
void applyZJetSkim (const vector< fs::path > &inputs, const fs::path &output, const pt::ptree &config, const int steering, const DT::Slice slice={1, 0})
 

Function Documentation

◆ applyDijetSkim()

void DAS::Jets::applyDijetSkim ( const vector< fs::path > &  inputs,
const fs::path &  output,
const pt::ptree &  config,
const int  steering,
const DT::Slice  slice = {1,0} 
)

Skims the input tree to keep only events with at least two jets, with configurable pT thresholds.

Parameters
inputsinput ROOT files (n-tuples)
outputoutput ROOT file (n-tuple)
configconfig handled with `Darwin::Tools::options`
steeringparameters obtained from explicit options
slicenumber and index of slice
48  {1,0}
49  )
50 {
51  cout << __func__ << ' ' << slice << " start" << endl;
52 
53  DT::Flow flow(steering);
54  auto tIn = flow.GetInputTree(inputs, slice);
55  auto tOut = flow.GetOutputTree(output);
56 
57  DT::MetaInfo metainfo(tOut);
58  metainfo.Check(config);
59 
60  const auto pt1 = config.get<float>("skims.dijet.pt1");
61  metainfo.Set<float>("skims", "dijet", "pt1", pt1);
62 
63  const auto pt2 = config.get<float>("skims.dijet.pt2");
64  metainfo.Set<float>("skims", "dijet", "pt2", pt2);
65 
66  cout << "Skimming events: pt(jet1) > " << pt1 << "\t pt(jet2) > " << pt2 << endl;
67 
68  auto recJets = flow.GetBranchReadOnly<vector<RecJet>>("recJets", DT::facultative);
69  RecDijet * recDijet = nullptr;
70  if (recJets != nullptr)
71  recDijet = flow.GetBranchWriteOnly<RecDijet>("recDijet");
72 
73  auto genJets = flow.GetBranchReadOnly<vector<GenJet>>("genJets", DT::facultative);
74  GenDijet * genDijet = nullptr;
75  if (genJets != nullptr)
76  genDijet = flow.GetBranchWriteOnly<GenDijet>("genDijet");
77 
78  if (!recJets && !genJets)
79  BOOST_THROW_EXCEPTION( DE::BadInput("No jets in input tree", tIn) );
80 
81  for (DT::Looper looper(tIn); looper(); ++looper) {
82  [[ maybe_unused ]]
83  static auto& cout = (steering & DT::verbose) == DT::verbose ? ::cout : DT::dev_null;
84 
85  bool passesRec = true, passesGen = true;
86 
87  if (recJets) {
88  passesRec = keepEvent(*recJets, pt1, pt2);
89  if (passesRec) *recDijet = recJets->at(0) + recJets->at(1);
90  else recDijet->clear();
91  }
92  if (genJets) {
93  passesGen = keepEvent(*genJets, pt1, pt2);
94  if (passesGen) *genDijet = genJets->at(0) + genJets->at(1);
95  else genDijet->clear();
96  }
97 
98  if ((steering & DT::fill) && (passesRec || passesGen)) tOut->Fill();
99  }
100 
101  metainfo.Set<bool>("git", "complete", true);
102 
103  cout << __func__ << ' ' << slice << " stop" << endl;
104 }

◆ applyZJetSkim()

void DAS::Jets::applyZJetSkim ( const vector< fs::path > &  inputs,
const fs::path &  output,
const pt::ptree &  config,
const int  steering,
const DT::Slice  slice = {1,0} 
)

Skims the input tree to keep only events with at least two jets, with configurable pT thresholds.

Parameters
inputsinput ROOT files (n-tuples)
outputoutput ROOT file (n-tuple)
configconfig handled with `Darwin::Tools::options`
steeringparameters obtained from explicit options
slicenumber and index of slice
50  {1,0}
51  )
52 {
53  cout << __func__ << ' ' << slice << " start" << endl;
54 
55  DT::Flow flow(steering);
56  auto tIn = flow.GetInputTree(inputs, slice);
57  auto tOut = flow.GetOutputTree(output);
58 
59  DT::MetaInfo metainfo(tOut);
60  metainfo.Check(config);
61 
62  const auto zpt = config.get<float>("skims.ZJet.zpt");
63  metainfo.Set<float>("skims", "ZJet", "zpt", zpt);
64 
65  const auto jpt = config.get<float>("skims.ZJet.jpt");
66  metainfo.Set<float>("skims", "ZJet", "jpt", jpt);
67 
68  cout << "Skimming events: pt(Z) > " << zpt << "\t pt(jet) > " << jpt << endl;
69 
70  auto recJets = flow.GetBranchReadOnly<vector<RecJet>>("recJets", DT::facultative);
71  auto genJets = flow.GetBranchReadOnly<vector<GenJet>>("genJets", DT::facultative);
72  if (!recJets && !genJets)
73  BOOST_THROW_EXCEPTION( DE::BadInput("No jets in input tree", tIn) );
74 
75  auto recDimuon = flow.GetBranchReadOnly<RecDimuon>("recDimuon", DT::facultative);
76  auto genDimuon = flow.GetBranchReadOnly<GenDimuon>("genDimuon", DT::facultative);
77  if (!recDimuon && !genDimuon)
78  BOOST_THROW_EXCEPTION( DE::BadInput("No dimuon in input tree", tIn) );
79 
80  RecZJet * recZJet = nullptr;
81  if (recJets && recDimuon)
82  recZJet = flow.GetBranchWriteOnly<RecZJet>("recZJet");
83  GenZJet * genZJet = nullptr;
84  if (genJets && genDimuon)
85  genZJet = flow.GetBranchWriteOnly<GenZJet>("genZJet");
86 
87  for (DT::Looper looper(tIn); looper(); ++looper) {
88  [[ maybe_unused ]]
89  static auto& cout = (steering & DT::verbose) == DT::verbose ? ::cout : DT::dev_null;
90 
91  bool passesRec = true, passesGen = true;
92 
93  if (recJets) {
94  passesRec = keepEvent(*recJets, *recDimuon, zpt, jpt);
95  if (passesRec) *recZJet = *recDimuon + recJets->at(0);
96  else recZJet->clear();
97  }
98  if (genJets) {
99  passesGen = keepEvent(*genJets, *genDimuon, zpt, jpt);
100  if (passesGen) *genZJet = *genDimuon + genJets->at(0);
101  else genZJet->clear();
102  }
103 
104  if ((steering & DT::fill) && (passesRec || passesGen)) tOut->Fill();
105  }
106 
107  metainfo.Set<bool>("git", "complete", true);
108 
109  cout << __func__ << ' ' << slice << " stop" << endl;
110 }

◆ keepEvent() [1/2]

bool DAS::Jets::keepEvent ( const vector< Jet > &  jets,
const Dimuon &  dimuon,
float  zpt,
float  jpt 
)

Checks whether an event should be kept.

37 {
38  return jets.size() >= 2 && jets[0].p4.Pt() > jpt &&
39  dimuon && dimuon.CorrP4().Pt() > zpt;
40 }

◆ keepEvent() [2/2]

bool DAS::Jets::keepEvent ( const vector< Jet > &  jets,
float  pt1,
float  pt2 
)

Checks whether an event should be kept.

36 {
37  return jets.size() >= 2 && jets[0].p4.Pt() > pt1 && jets[1].p4.Pt() > pt2;
38 }
Darwin::Tools::fill
@ fill
activate -f to fill the tree
Definition: Options.h:27
Darwin::Tools::Flow
User-friendly handling of input and output n-tuples.
Definition: Flow.h:69
Step::verbose
static bool verbose
Definition: Step.h:40
DAS::GenDijet
Di< GenJet, GenJet > GenDijet
Definition: Di.h:77
Darwin::Tools::Looper
Facility to loop over a n-tuple, including parallelisation and printing.
Definition: Looper.h:32
DAS::GenZJet
Di< GenDimuon, GenJet > GenZJet
Definition: Di.h:79
Darwin::Tools::MetaInfo
Generic meta-information for n-tuple (including speficities to Darwin).
Definition: MetaInfo.h:68
Ntupliser_cfg.jets
string jets
Definition: Ntupliser_cfg.py:41
Ntupliser_cfg.config
config
Definition: Ntupliser_cfg.py:264
DAS::RecZJet
Di< RecDimuon, RecJet > RecZJet
Definition: Di.h:82
DAS::Jets::keepEvent
bool keepEvent(const vector< Jet > &jets, const Dimuon &dimuon, float zpt, float jpt)
Checks whether an event should be kept.
Definition: applyZJetSkim.cc:36
DAS::RecDijet
Di< RecJet, RecJet > RecDijet
Definition: Di.h:80
jercExample.inputs
def inputs
Definition: jercExample.py:118
DAS::Jets::keepEvent
bool keepEvent(const vector< Jet > &jets, float pt1, float pt2)
Checks whether an event should be kept.
Definition: applyDijetSkim.cc:35
Darwin::Tools::dev_null
static std::ostream dev_null(nullptr)
to redirect output stream to nowhere
Darwin::Exceptions::BadInput
Generic exception for ill-defined input (before the event loop).
Definition: exceptions.h:83
DAS::RecDimuon
Di< RecMuon, RecMuon > RecDimuon
Definition: Di.h:81
Darwin::Tools::facultative
@ facultative
mounting branch is facultative
Definition: Flow.h:31
DAS::GenDimuon
Di< GenMuon, GenMuon > GenDimuon
Definition: Di.h:78