|
DAS
3.0
Das Analysis System
|
Functor to apply weights from pre-firing maps
Source: following link + private chats with Laurent Thomas: https://twiki.cern.ch/twiki/bin/viewauth/CMS/L1ECALPrefiringWeightRecipe
Note: we ignore photons
In contrast to the R-scan analysis (AN-18-044), we take the value from the map itself (i.e. we do not fit or smooth anything).
#include <applyPrefiringWeights.h>
|
using | Operation = std::tuple< float, float, float >(Functor::*)(int, const RecJet &) const |
|
|
| Functor (const fs::path &file, const std::string name, const std::string type, bool applySyst) |
|
void | operator() (RecEvent &evnt, const std::vector< RecJet > &recjets) |
|
void | Write (TDirectory *dir) |
|
◆ Operation
using Operation = std::tuple<float,float,float> (Functor::*)(int,const RecJet&) const |
◆ Functor()
Functor |
( |
const fs::path & |
file, |
|
|
const std::string |
name, |
|
|
const std::string |
type, |
|
|
bool |
applySyst |
|
) |
| |
|
inline |
Functor to apply prefiring correction.
- Parameters
-
file | file name containing the map |
name | name of the object to fetch in the file |
type | prefiring option, either "smooth" or "binned" |
141 if (!fs::exists(file))
142 BOOST_THROW_EXCEPTION( fs::filesystem_error(
"File does not exists.",
143 file, make_error_code(errc::no_such_file_or_directory)) );
145 auto f = make_unique<TFile>(file.c_str(),
"READ");
146 string hname = file.stem().string();
147 auto h = unique_ptr<TH2>(
f->Get<TH2>(hname.c_str()));
148 if (h.get() ==
nullptr) {
149 string what =
"Map " + hname +
" can't be found";
155 if (type ==
"binned")
157 else if (type ==
"smooth")
160 BOOST_THROW_EXCEPTION( invalid_argument(
"Types may only be 'smooth' and 'binned'.") );
165 cout <<
"Prefiring weights are ready" << endl;
◆ GetBinnedWeights()
std::tuple<float,float,float> GetBinnedWeights |
( |
int |
run, |
|
|
const RecJet & |
recjet |
|
) |
| const |
|
inlineprivate |
Subroutine to get weights for a give event directly from the maps.
184 if (std::abs(
eta) < 2.0 || std::abs(
eta) >= 3.0)
return { 1, 1, 1};
188 int etabin = MapNo->GetXaxis()->FindBin(
eta),
189 ptbin = min(MapNo->GetYaxis()->FindBin(
pt ), MapNo->GetNbinsY());
191 auto corr = MapNo->GetBinContent(etabin, ptbin),
192 varUp = MapUp->GetBinContent(etabin, ptbin),
193 varDn = MapDn->GetBinContent(etabin, ptbin);
197 return {corr, varUp, varDn};
◆ GetMultiplicityPlots()
static std::vector<TH2 *> GetMultiplicityPlots |
( |
TH2 * |
templ, |
|
|
int |
n = 7 |
|
) |
| |
|
inlinestatic |
Multiplicity plots: not used anywhere, just to compare how the different samples are populating the forward region...
119 for (
int i = 1; i <= n; ++i) {
120 auto h =
dynamic_cast<TH2*
>(templ->Clone(Form(
"multiplicity_nbin%d", i)));
◆ GetSmoothWeights()
std::tuple<float, float, float> GetSmoothWeights |
( |
int |
run, |
|
|
const RecJet & |
recjet |
|
) |
| const |
|
inlineprivate |
Subroutine to get weights for a give event from the smooth functions.
210 int etabin = 1 + abs(
eta)*2;
213 if (etabin < 5 || etabin > 6)
return { 1, 1, 1};
218 auto eval =
f->Eval(
pt);
220 float corr = 1- eval,
226 return {corr, varUp, varDn};
◆ operator()()
void operator() |
( |
RecEvent & |
evnt, |
|
|
const std::vector< RecJet > & |
recjets |
|
) |
| |
|
inline |
Operator overloading for event-by-event call.
- Todo:
- decorrelations
238 tuple<float,float,float> wAllJets { 1, 1, 1 };
243 get<0>(wAllJets) *= get<0>(wJet);
244 get<1>(wAllJets) *= get<1>(wJet);
245 get<2>(wAllJets) *= get<2>(wJet);
253 wEvt0 = evnt.weights.front(),
258 if (i >= n)
continue;
267 float rw0 = evnt.weights.front();
268 assert(get<0>(wAllJets) > 0);
269 evnt.weights /= get<0>(wAllJets);
272 evnt.weights.push_back(Weight{rw0 / get<1>(wAllJets)});
273 evnt.weights.push_back(Weight{rw0 / get<2>(wAllJets)});
◆ PrepareBinnedWeights()
static std::tuple<TH2*, TH2*, TH2*> PrepareBinnedWeights |
( |
const std::unique_ptr< TH2 > & |
hIn | ) |
|
|
inlinestatic |
Prepare weights including variations.
The original maps contain the probability for one jet to cause the prefiring
- we convert it to a weight = 1 - proba
- we derive the uncertainty in addition, following the official recommendation: "the uncertainty is taken as the maximum between 20 percents of the object prefiring probability
and the statistical uncertainty associated to the considered bin in the prefiring map."
- Todo:
- change prescription to better account for decorrelations
63 auto MapNo =
dynamic_cast<TH2*
>(hIn->Clone(
"nominal_map")),
64 MapUp =
dynamic_cast<TH2*
>(hIn->Clone(
"upper_map" )),
65 MapDn =
dynamic_cast<TH2*
>(hIn->Clone(
"lower_map" ));
66 MapNo->SetDirectory(0);
67 MapUp->SetDirectory(0);
68 MapDn->SetDirectory(0);
69 for (
int etabin = 0; etabin <= hIn->GetNbinsX()+1; ++etabin)
70 for (
int ptbin = 0; ptbin <= hIn->GetNbinsY()+1; ++ ptbin) {
71 auto proba = hIn->GetBinContent(etabin, ptbin),
72 error = hIn->GetBinError (etabin, ptbin);
75 bool Forward = (proba >
feps);
76 auto unc = Forward ? max(proba * 0.2
f, error) : 0.
f;
82 MapNo->SetBinContent(etabin, ptbin,
weight);
83 MapUp->SetBinContent(etabin, ptbin,
weight+
unc);
84 MapDn->SetBinContent(etabin, ptbin,
weight-
unc);
90 return { MapNo, MapUp, MapDn };
◆ PrepareSmoothWeights()
static std::tuple<TF1*, TF1*> PrepareSmoothWeights |
( |
TString |
hname, |
|
|
std::unique_ptr< TFile > |
f |
|
) |
| |
|
inlinestatic |
Prepare weights including variations.
The smooth function contain the probability for one jet to cause the prefiring
- as it is a function, the conversion to 1-weight will be done while running
- for the uncertainty, we will just vary the weight while applying it with 20%
- Todo:
- Assertion to exception
- Todo:
- Assertion to exception
103 cout << hname << endl;
104 auto f5 =
f->Get<TF1>(hname +
"_ybin5"),
105 f6 =
f->Get<TF1>(hname +
"_ybin6");
106 assert(f5 !=
nullptr);
107 assert(f6 !=
nullptr);
◆ Write()
void Write |
( |
TDirectory * |
dir | ) |
|
|
inline |
Write maps and control plots to output root file.
287 for (TH2 * h: {MapNo, MapUp, MapDn}) {
288 if (h ==
nullptr)
continue;
289 h->SetDirectory(dir);
290 TString
name = h->GetName();
296 if (f5 !=
nullptr) f5->Write(
"ybin5");
297 if (f6 !=
nullptr) f6->Write(
"ybin6");
300 m->SetDirectory(dir);
◆ GetWeights
◆ multiplicities
std::vector<TH2 *> multiplicities |
jets in the forward region in the same binning as the maps (before correction)
◆ prefBinnedWeights
std::tuple<TH2*, TH2*, TH2*> prefBinnedWeights |
weights extracted from the maps, indiced with run number (three variations)
◆ prefSmoothWeights
std::tuple<TF1 *, TF1*> prefSmoothWeights |
weights extracted from smooth functions, indiced with run number (two rapidity bins)
◆ syst
The documentation for this struct was generated from the following file:
name
Definition: DYToLL_M-50_13TeV_pythia8_cff_GEN_SIM_RECOBEFMIX_DIGI_L1_DIGI2RAW_L1Reco_RECO.py:48
pt
Definition: jmarExample.py:19
FourVector p4
raw four-momentum directly after reconstruction
Definition: PhysicsObject.h:47
std::tuple< float, float, float > GetBinnedWeights(int run, const RecJet &recjet) const
Subroutine to get weights for a give event directly from the maps.
Definition: applyPrefiringWeights.h:171
f
Definition: Ntupliser_cfg.py:252
std::vector< float > scales
energy scale corrections and variations
Definition: PhysicsObject.h:48
const Operation GetWeights
Definition: applyPrefiringWeights.h:45
static const auto feps
Definition: applyPrefiringWeights.h:25
std::tuple< TH2 *, TH2 *, TH2 * > prefBinnedWeights
weights extracted from the maps, indiced with run number (three variations)
Definition: applyPrefiringWeights.h:41
DAS::RecJet recjet
Definition: classes.h:15
DAS::Weight weight
Definition: classes.h:11
string unc
Definition: jercExample.py:70
static std::tuple< TH2 *, TH2 *, TH2 * > PrepareBinnedWeights(const std::unique_ptr< TH2 > &hIn)
Definition: applyPrefiringWeights.h:59
static std::tuple< TF1 *, TF1 * > PrepareSmoothWeights(TString hname, std::unique_ptr< TFile > f)
Definition: applyPrefiringWeights.h:99
recjets
Definition: Ntupliser_cfg.py:269
float CorrPt(size_t i=0) const
corrected transverse momentum
Definition: PhysicsObject.h:52
static std::vector< TH2 * > GetMultiplicityPlots(TH2 *templ, int n=7)
Definition: applyPrefiringWeights.h:114
std::vector< TH2 * > multiplicities
jets in the forward region in the same binning as the maps (before correction)
Definition: applyPrefiringWeights.h:49
const bool syst
Definition: applyPrefiringWeights.h:47
Weights weights
object weights
Definition: PhysicsObject.h:49
std::tuple< TF1 *, TF1 * > prefSmoothWeights
weights extracted from smooth functions, indiced with run number (two rapidity bins)
Definition: applyPrefiringWeights.h:42
std::tuple< float, float, float > GetSmoothWeights(int run, const RecJet &recjet) const
Subroutine to get weights for a give event from the smooth functions.
Definition: applyPrefiringWeights.h:202
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< float > > FourVector
Definition: PhysicsObject.h:15
eta
DeepAK8/ParticleNet tagging.
Definition: jmarExample.py:19