DAS  3.0
Das Analysis System
DAS::MN Namespace Reference

Classes

struct  Hist
 
struct  Obs2Jets
 
struct  ObsMiniJets
 

Functions

void getMNobservables (const vector< fs::path > inputs, const fs::path output, const pt::ptree &config, const int steering, const DT::Slice slice={1, 0})
 
template<typename Jet >
std::optional< std::pair< Jet, Jet > > GetMNJet (std::vector< Jet > jets, std::function< bool(Jet &)> ptcut=[](Jet &jet) {return jet.p4.Pt()< 35 ;})
 
template<typename Jet >
std::vector< Jet > GetMiniJets (std::vector< Jet > jets, const std::pair< Jet, Jet > &MNJets, std::function< bool(Jet &)> ptcut=[](Jet &jet) {return jet.p4.Pt()< 20;})
 

Variables

static const vector< double > binsY = {0., 0.25, .5, 0.75, 1., 1.25, 1.5, 1.75, 2., 2.25, 2.5, 2.75, 3., 3.25, 3.5, 3.75, 4., 4.25, 4.5, 4.75, 5., 5.25, 5.5, 5.75, 6., 6.25, 6.5, 6.75, 7., 7.25, 7.5, 7.75, 8.}
 
static const vector< double > binsPhi = {0.,.1*M_PI, .2*M_PI, .3*M_PI, .4*M_PI, .5*M_PI, .6*M_PI, .7*M_PI, .8*M_PI, .9*M_PI, M_PI}
 
static const double minY = binsY.front()
 
static const double maxY = binsY.back()
 
static const double maxy = 4.7
 
static const double minpt = 28
 
static const double maxpt = 1.e10
 

Detailed Description

The namespace MN_Helper contains the classes that are needed by getMNobservables.cc to build the histograms. Two types of objects are present here: bs2Jets and ObsMiniJets that contains structures to store the observables themselves (they are higher level physics objects) and Hist that contains histograms.

Function Documentation

◆ GetMiniJets()

std::vector< Jet > GetMiniJets ( std::vector< Jet >  jets,
const std::pair< Jet, Jet > &  MNJets,
std::function< bool(Jet &)>  ptcut = [](Jet& jet) {return jet.p4.Pt()<20;} 
)

This function takes as an input a vector of jet, a pair of Mueller-Navelet jets in this event and a selection fuction. All the jet for which the selection function is false are dicarded. Remove all the jets that are at higher-or-equal rapidity than the first member of the pair (considered as the 'leading MN jet') or lower-or-equal rapidity than the second member of the pair (considered a the 'subleading MN jet'). The function is build assuming that the two Mueller-Navelet jets that are provided are also present in the vector of jets of the event. The two Mueller-Navelet jets are removed by the lower-or-equal / higher-or-equal logic of the sample.

Note
there is no consistency check. So if the two MNJets that are given are not present in the vector of jets, it will not throw any error.
if the paire provided is trivial, it will not throw any error. The non-triviallity must be check before hand.
As this function is a template function, it is declared and defined in the common.h. Definitions of template functions are not compiled until they are instantiated with a specific type instance.
Todo:
check the consistency of this condition
Todo:
check the consistency of this condition
129 {
130  jets.erase(std::remove_if(jets.begin(), jets.end(), ptcut), jets.end());
131  jets.erase(std::remove_if(jets.begin(), jets.end(),
132  [MNJets](Jet& jet){return jet.p4.Eta() >= MNJets.first.p4.Eta();}), jets.end());
133 
134  jets.erase(std::remove_if(jets.begin(), jets.end(),
135  [MNJets](Jet& jet){return jet.p4.Eta() <= MNJets.second.p4.Eta();}), jets.end());
136 
137  std::sort( jets.begin(), jets.end() , [](Jet& j1, Jet& j2)
138  {return j1.p4.Eta() > j2.p4.Eta();} );
139 
140  return jets;
141 }

◆ GetMNJet()

std::optional< std::pair< Jet, Jet > > GetMNJet ( std::vector< Jet >  jets,
std::function< bool(Jet &)>  ptcut = [](Jet& jet) {return jet.p4.Pt() < 35 ;} 
)

This function takes as an input a vector of jet and a function<bool(Jet&)>. It finds the two jets that are th most bakward and most forward in the event and satisfies the cut implemented in the function<bool(Jet&)> provided as argument. This function is apllied to all the jets. If it return false, the jets are discarded. The out is an optional. So that if there are no jets satifying the cuts in the event, it returns nullopt.

Note
As this function is a template function, it is declared and defined in the common.h. Definitions of template functions are not compiled until they are instantiated with a specific type instance.
Todo:
check the consistency of this condition
Todo:
check the consistency of this condition
99 {
100  jets.erase(remove_if(jets.begin(), jets.end(), ptcut), jets.end());
101  if (jets.size() < 2) return std::nullopt;
102  auto result = minmax_element( jets.begin(), jets.end() , [](Jet& j1, Jet& j2)
103  {return j1.p4.Eta() > j2.p4.Eta();} );
104 
105  if (result.first->p4.Eta() <= 0 || result.second->p4.Eta() >= 0) return std::nullopt;
106  return std::make_optional<std::pair<Jet,Jet>>(*result.first, *result.second);
107 }

◆ getMNobservables()

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

Calculate the observables that can be used to caracterise Mueller-Navelet Jets.

Parameters
inputsinput ROOT files (n-tuples)
outputname of output root file containing the histograms
configconfig file from `DTOptions`
steeringsteering parameters from `DTOptions`
sliceslices for running
184  {1,0}
185  )
186 {
187  cout << __func__ << ' ' << slice << " start" << endl;
188 
189  unique_ptr<TChain> tIn = DT::GetChain(inputs);
190  unique_ptr<TFile> fOut(DT_GetOutput(output));
191  auto tOut = unique_ptr<TTree>(tIn->CloneTree(0));
192 
193  DT::MetaInfo metainfo(tOut);
194  metainfo.Check(config);
195  auto isMC = metainfo.Get<bool>("flags", "isMC");
196 
197  GenEvent * gEv = nullptr;
198  RecEvent * rEv = nullptr;
199  if (isMC)
200  tIn->SetBranchAddress("genEvent", &gEv);
201  tIn->SetBranchAddress("recEvent", &rEv);
202 
203  vector<RecJet> * recJets = nullptr;
204  vector<GenJet> * genJets = nullptr;
205  tIn->SetBranchAddress("recJets", &recJets);
206  if (isMC)
207  tIn->SetBranchAddress("genJets", &genJets);
208 
209  Hist genHist("gen"), // NOTE: one should only declare the object if (isMC)
210  recHist("rec");
211 
212  for (DT::Looper looper(tIn, slice); looper(); ++looper) {
213  [[ maybe_unused ]]
214  static auto& cout = (steering & DT::verbose) == DT::verbose ? ::cout : DT::dev_null;
215 
216  auto MNJets = GetMNJet<RecJet>(*recJets, [](RecJet jet){return jet.CorrPt() < 35;});
217  if (MNJets) {
218  double recy = std::abs( MNJets->first.p4.Eta() - MNJets->second.p4.Eta() );
219  bool LowRecY = recy < minY,
220  HighRecY = recy >= maxY;
221  bool goodRec = (!LowRecY) && (!HighRecY)
222  && MNJets->first.p4.Eta() < maxy && MNJets->second.p4.Eta() < maxy
223  && minpt < MNJets->first.CorrPt() && MNJets->first.CorrPt() < maxpt
224  && minpt < MNJets->second.CorrPt() && MNJets->second.CorrPt() < maxpt;
225 
226  if (goodRec) {
227  double evweight = rEv->weights.front();
228  if (isMC) evweight *= gEv->weights.front();
229  recHist.Fill(*recJets, evweight);
230  Obs2Jets obs2jets(MNJets->first, MNJets->second);
231  auto minijets = GetMiniJets(*recJets, *MNJets, function<bool(RecJet&)>([](RecJet& jet){return jet.CorrPt() < 20;}));
232  ObsMiniJets obsminijets(minijets);
233  recHist.Fill(obs2jets, obsminijets, evweight);
234  }
235  }
236 
237  if (!isMC) continue;
238 
239  auto genMNJets = GetMNJet<GenJet>(*genJets, [](GenJet jet){return jet.p4.Pt() < 35;});
240  if (genMNJets) {
241  double genY = std::abs( genMNJets->first.p4.Eta() - genMNJets->second.p4.Eta() );
242  bool LowGenY = genY < minY,
243  HighGenY = genY >= maxY,
244  goodGen = (!LowGenY) && (!HighGenY)
245  && genMNJets->first.p4.Eta() < maxy && genMNJets->second.p4.Eta() < maxy
246  && minpt < genMNJets->first.p4.Pt() && genMNJets->first.p4.Pt() < maxpt
247  && minpt < genMNJets->second.p4.Pt() && genMNJets->second.p4.Pt() < maxpt;
248 
249  if (goodGen) {
250  double evweight = gEv->weights.front();
251  genHist.Fill(*genJets, evweight);
252  Obs2Jets obs2jets(genMNJets->first, genMNJets->second);
253  auto minijets = GetMiniJets(*genJets, *genMNJets, function<bool(GenJet&)>([](GenJet& jet){return jet.p4.Pt() < 20;}));
254  ObsMiniJets obsminijets(minijets);
255  genHist.Fill(obs2jets, obsminijets, evweight);
256  }
257  }
258  }
259 
260  recHist.Write(fOut.get());
261  if (isMC)
262  genHist.Write(fOut.get());
263 
264  metainfo.Set<bool>("git", "complete", true);
265  tOut->Write();
266 
267  cout << __func__ << ' ' << slice << " end" << endl;
268 }

Variable Documentation

◆ binsPhi

const vector<double> binsPhi = {0.,.1*M_PI, .2*M_PI, .3*M_PI, .4*M_PI, .5*M_PI, .6*M_PI, .7*M_PI, .8*M_PI, .9*M_PI, M_PI}
static

◆ binsY

const vector<double> binsY = {0., 0.25, .5, 0.75, 1., 1.25, 1.5, 1.75, 2., 2.25, 2.5, 2.75, 3., 3.25, 3.5, 3.75, 4., 4.25, 4.5, 4.75, 5., 5.25, 5.5, 5.75, 6., 6.25, 6.5, 6.75, 7., 7.25, 7.5, 7.75, 8.}
static

◆ maxpt

const double maxpt = 1.e10
static

◆ maxy

const double maxy = 4.7
static

◆ maxY

const double maxY = binsY.back()
static

◆ minpt

const double minpt = 28
static

◆ minY

const double minY = binsY.front()
static
Step::verbose
static bool verbose
Definition: Step.h:40
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::MN::maxy
static const double maxy
Definition: getMNobservables.cc:41
DAS::MN::maxY
static const double maxY
Definition: getMNobservables.cc:41
Darwin::Tools::MetaInfo
Generic meta-information for n-tuple (including speficities to Darwin).
Definition: MetaInfo.h:65
DAS::MN::minpt
static const double minpt
Definition: getMNobservables.cc:42
Ntupliser_cfg.jets
string jets
Definition: Ntupliser_cfg.py:41
DAS::MN::minY
static const double minY
Definition: getMNobservables.cc:41
Ntupliser_cfg.config
config
Definition: Ntupliser_cfg.py:263
DAS::MN::maxpt
static const double maxpt
Definition: getMNobservables.cc:42
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
DT_GetOutput
#define DT_GetOutput(output)
Definition: FileUtils.h:96
DAS::MN::GetMiniJets
std::vector< Jet > GetMiniJets(std::vector< Jet > jets, const std::pair< Jet, Jet > &MNJets, std::function< bool(Jet &)> ptcut=[](Jet &jet) {return jet.p4.Pt()< 20;})
Definition: MuellerNavelet.h:128