DAS  3.0
Das Analysis System
DAS::Muons Namespace Reference

Functions

template<class Muon >
bool keepEvent (const vector< Muon > &muons, float pt1, float pt2)
 
void applyDimuonSkim (const vector< fs::path > &inputs, const fs::path &output, const pt::ptree &config, const int steering, const DT::Slice slice={1, 0})
 

Function Documentation

◆ applyDimuonSkim()

void DAS::Muons::applyDimuonSkim ( 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 muons, with configurable pT thresholds.

Parameters
inputsinput ROOT files (n-tuple)
outputoutput ROOT file (n-tuple)
configconfig handled with `Darwin::Tools::Options`
steeringparameters obtained from explicit options
slicenumber and index of slice
44  {1,0}
45  )
46 {
47  cout << __func__ << ' ' << slice << " start" << endl;
48 
49  DT::Flow flow(steering);
50  auto tIn = flow.GetInputTree(inputs, slice);
51  auto [fOut, tOut] = flow.GetOutput(output);
52 
53  DT::MetaInfo metainfo(tOut);
54  metainfo.Check(config);
55  const bool isMC = metainfo.Get<bool>("flags", "isMC");
56 
57  const auto pt1 = config.get<float>("skims.dimuon.pt1");
58  metainfo.Set<float>("skims", "dimuon", "pt1", pt1);
59 
60  const auto pt2 = config.get<float>("skims.dimuon.pt2");
61  metainfo.Set<float>("skims", "dimuon", "pt2", pt2);
62 
63  cout << "Skimming events: pt(mu1) > " << pt1 << "\t pt(mu2) > " << pt2 << endl;
64 
65  // event information
66  auto genEvt = flow.GetBranchReadOnly<GenEvent>("genEvent", DT::facultative);
67 
68  auto recMuons = flow.GetBranchReadOnly<vector<RecMuon>>("recMuons", DT::facultative);
69  RecDimuon * recDimuon = nullptr;
70  if (recMuons != nullptr)
71  recDimuon = flow.GetBranchWriteOnly<RecDimuon>("recDimuon");
72 
73  auto genMuons = flow.GetBranchReadOnly<vector<GenMuon>>("genMuons", DT::facultative);
74  GenDimuon * genDimuon = nullptr;
75  if (genMuons != nullptr)
76  genDimuon = flow.GetBranchWriteOnly<GenDimuon>("genDimuon");
77 
78  if (!recMuons && !genMuons)
79  BOOST_THROW_EXCEPTION( DE::BadInput("No muons in input tree", tIn) );
80 
81  // histogram used in MC for normalisation
82  TH1 * hSumWgt = nullptr;
83  if (isMC) hSumWgt = new TH1F("hSumWgt",";;#sum_{i} w_{i}", 1, -0.5, 0.5);
84 
85  for (DT::Looper looper(tIn); looper(); ++looper) {
86  [[ maybe_unused ]]
87  static auto& cout = steering & DT::verbose ? ::cout : DT::dev_null;
88 
89  bool passesRec = true, passesGen = true;
90 
91  if (recMuons) {
92  passesRec = keepEvent(*recMuons, pt1, pt2);
93  if (passesRec) *recDimuon = recMuons->at(0) + recMuons->at(1);
94  else recDimuon->clear();
95  }
96  if (genMuons) {
97  passesGen = keepEvent(*genMuons, pt1, pt2);
98  if (passesGen) *genDimuon = genMuons->at(0) + genMuons->at(1);
99  else genDimuon->clear();
100  }
101 
102  if ((steering & DT::fill) && (passesRec || passesGen)) tOut->Fill();
103 
104  if (!isMC) continue;
105  auto w = genEvt->weights.front();
106  hSumWgt->Fill(0.0, w);
107  }
108 
109  metainfo.Set<bool>("git", "complete", true);
110  fOut->cd();
111  if (hSumWgt) hSumWgt->Write();
112 
113  cout << __func__ << ' ' << slice << " stop" << endl;
114 }

◆ keepEvent()

bool DAS::Muons::keepEvent ( const vector< Muon > &  muons,
float  pt1,
float  pt2 
)

Checks whether an event should be kept.

32 {
33  return muons.size() >= 2 && muons[0].p4.Pt() > pt1 && muons[1].p4.Pt() > pt2;
34 }
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
Ntupliser_cfg.muons
string muons
Definition: Ntupliser_cfg.py:43
Darwin::Tools::Looper
Facility to loop over a n-tuple, including parallelisation and printing.
Definition: Looper.h:32
DAS::JetEnergy::w
static const float w
Definition: common.h:51
Darwin::Tools::MetaInfo
Generic meta-information for n-tuple (including speficities to Darwin).
Definition: MetaInfo.h:68
DAS::Muons::keepEvent
bool keepEvent(const vector< Muon > &muons, float pt1, float pt2)
Checks whether an event should be kept.
Definition: applyDimuonSkim.cc:31
Ntupliser_cfg.config
config
Definition: Ntupliser_cfg.py:264
Ntupliser_cfg.isMC
string isMC
Definition: Ntupliser_cfg.py:59
jercExample.inputs
def inputs
Definition: jercExample.py:118
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