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  unique_ptr<TChain> tIn = DT::GetChain(inputs);
54  unique_ptr<TFile> fOut(DT_GetOutput(output));
55  auto tOut = unique_ptr<TTree>(tIn->CloneTree(0));
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  vector<RecJet> * recJets = nullptr;
69  RecDijet * recDijet = nullptr;
70  if (branchExists(tIn, "recJets")) {
71  tIn->SetBranchAddress("recJets", &recJets);
72  recDijet = new RecDijet;
73  tOut->Branch("recDijet", &recDijet);
74  }
75 
76  vector<GenJet> * genJets = nullptr;
77  GenDijet * genDijet = nullptr;
78  if (branchExists(tIn, "genJets")) {
79  tIn->SetBranchAddress("genJets", &genJets);
80  genDijet = new GenDijet;
81  tOut->Branch("genDijet", &genDijet);
82  }
83 
84  if (!recJets && !genJets)
85  BOOST_THROW_EXCEPTION( DE::BadInput("No jets in input tree", tIn) );
86 
87  for (DT::Looper looper(tIn, slice); 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, pt1, pt2);
95  if (passesRec) *recDijet = recJets->at(0) + recJets->at(1);
96  else recDijet->clear();
97  }
98  if (genJets) {
99  passesGen = keepEvent(*genJets, pt1, pt2);
100  if (passesGen) *genDijet = genJets->at(0) + genJets->at(1);
101  else genDijet->clear();
102  }
103 
104  if ((steering & DT::fill) == DT::fill && (passesRec || passesGen)) tOut->Fill();
105  }
106 
107  metainfo.Set<bool>("git", "complete", true);
108  tOut->Write();
109 
110  cout << __func__ << ' ' << slice << " stop" << endl;
111 }

◆ 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  unique_ptr<TChain> tIn = DT::GetChain(inputs);
56  unique_ptr<TFile> fOut(DT_GetOutput(output));
57  auto tOut = unique_ptr<TTree>(tIn->CloneTree(0));
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  vector<RecJet> * recJets = nullptr;
71  if (branchExists(tIn, "recJets"))
72  tIn->SetBranchAddress("recJets", &recJets);
73  vector<GenJet> * genJets = nullptr;
74  if (branchExists(tIn, "genJets"))
75  tIn->SetBranchAddress("genJets", &genJets);
76  if (!recJets && !genJets)
77  BOOST_THROW_EXCEPTION( DE::BadInput("No jets in input tree", tIn) );
78 
79  RecDimuon * recDimuon = nullptr;
80  if (branchExists(tIn, "recDimuon"))
81  tIn->SetBranchAddress("recDimuon", &recDimuon);
82  GenDimuon * genDimuon = nullptr;
83  if (branchExists(tIn, "genDimuon"))
84  tIn->SetBranchAddress("genDimuon", &genDimuon);
85  if (!recDimuon && !genDimuon)
86  BOOST_THROW_EXCEPTION( DE::BadInput("No dimuon in input tree", tIn) );
87 
88  RecZJet * recZJet = nullptr;
89  if (recJets && recDimuon) {
90  recZJet = new RecZJet;
91  tOut->Branch("recZJet", &recZJet);
92  }
93  GenZJet * genZJet = nullptr;
94  if (genJets && genDimuon) {
95  genZJet = new GenZJet;
96  tOut->Branch("genZJet", &genZJet);
97  }
98 
99  for (DT::Looper looper(tIn, slice); looper(); ++looper) {
100  [[ maybe_unused ]]
101  static auto& cout = (steering & DT::verbose) == DT::verbose ? ::cout : DT::dev_null;
102 
103  bool passesRec = true, passesGen = true;
104 
105  if (recJets) {
106  passesRec = keepEvent(*recJets, *recDimuon, zpt, jpt);
107  if (passesRec) *recZJet = *recDimuon + recJets->at(0);
108  else recZJet->clear();
109  }
110  if (genJets) {
111  passesGen = keepEvent(*genJets, *genDimuon, zpt, jpt);
112  if (passesGen) *genZJet = *genDimuon + genJets->at(0);
113  else genZJet->clear();
114  }
115 
116  if ((steering & DT::fill) == DT::fill && (passesRec || passesGen)) tOut->Fill();
117  }
118 
119  metainfo.Set<bool>("git", "complete", true);
120  tOut->Write();
121 
122  cout << __func__ << ' ' << slice << " stop" << endl;
123 }

◆ 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 }
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:33
Darwin::Tools::GetChain
std::unique_ptr< TChain > GetChain(std::vector< std::filesystem::path > inputs, const char *name="events")
Load chain from a list of files.
Definition: FileUtils.cc:67
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:65
Ntupliser_cfg.jets
string jets
Definition: Ntupliser_cfg.py:41
Darwin::Tools::fill
@ fill
activate -f to fill the tree
Definition: Options.h:28
Ntupliser_cfg.config
config
Definition: Ntupliser_cfg.py:263
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::branchExists
bool branchExists(const TTreePtr &tree, TString brName)
Definition: toolbox.h:27
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
DT_GetOutput
#define DT_GetOutput(output)
Definition: FileUtils.h:96
Darwin::Exceptions::BadInput
Generic exception for ill-defined input (before the event loop).
Definition: exceptions.h:68
DAS::RecDimuon
Di< RecMuon, RecMuon > RecDimuon
Definition: Di.h:81
DAS::GenDimuon
Di< GenMuon, GenMuon > GenDimuon
Definition: Di.h:78