 |
DAS
3.0
Das Analysis System
|
This functor implements the selection and normalisation of the events in the dataset. The public attributes are the following:
- trigger_lumi -> contains the trigger threshold (integer) and associated lumi
- eff -> contains the efficiency of the trigger at the leadingJet
- ntriggers -> number of trigger used
- inv_total_lumi -> inverse of the total uminosity of the sample
- minpt -> pt turnon of the lowest trigger
- ibit -> index of the highest pt trigger that have fired
- it -> iterator to the
pair<trigger, TriggerLumi>
of the highest pt trigger that have fired
The operator() implements the selection and normalisation itself. A private funtion pointer is used to allow to point to different private methods corresponding to different strategies. This behaviour is switched via 'strategy' argument in the constructor. It also iterate Functor::ibit
down to the highest trigger for which the leading jet pt is above the trigger turnon. This functor returns the leading jet according to the strategy that has been chosen AND multiplies the rec-level event weight to normalise the triggers together and to normalise the full dataset to the correct luminosity, according to the method that has been used. The weights are multiplied by: $\frac{ prescale * \mathcurl{L}_{total}^{-1} }{ \epsilon_{trigger} }$ in method 'presc' and by $ \frac{ \mathcurl{L}_{effective, trigger}^{-1} }{ \epsilon_{trigger} }$ in method 'lumi'. (NOTE: The lumi, trigger thresholds and trigger efficiency curves should be provided)
Two strategies are implemented.
Two methods are implemented:
- 'prescales' (see the private method
Functor::NormalisationCommonPresc()
)
- 'lumi' (see the private method
DataNormalsiation::NormalisationCommonLumi()
) (NOTE: Other methods could be implemented in the future. They should be implemented in this class.)
Technical remark: the switch wrt the different strategies and is done through private methods and private pointer function.
#include <applyDataNormalisation.h>
|
const RecJet & | NormalisationCommonPresc (RecEvent &evnt, const RecJet &leadingJet, const Trigger &trigger) |
|
const RecJet & | NormalisationCommonLumi (RecEvent &evnt, const RecJet &leadingJet, const Trigger &trigger) |
|
const RecJet & | NormalisationStd (RecEvent &evnt, std::vector< RecJet > &recJets, std::vector< FourVector > &hltJets, Trigger &trigger) |
|
const RecJet & | NormalisationFwd (RecEvent &evnt, std::vector< RecJet > &recJets, std::vector< FourVector > &hltJets, Trigger &trigger) |
|
◆ Functor()
Functor |
( |
fs::path |
lumi_file, |
|
|
fs::path |
turnon_file, |
|
|
fs::path |
trigger_curves, |
|
|
const std::string & |
strategy, |
|
|
const std::string & |
method, |
|
|
const int |
year |
|
) |
| |
|
inline |
- Todo:
- document
103 cout <<
"ntriggers = " <<
ntriggers <<
'\n'
105 <<
"minpt = " <<
minpt << endl;
107 if (
method ==
"prescales")
109 else if (
method ==
"lumi")
112 cerr <<
"method ('" <<
method <<
"') is not defined, aborting\n";
116 if (strategy ==
"pt")
118 else if (strategy ==
"eta")
121 cerr <<
"strategy of selection and normalisation ('" << strategy <<
"') is not defined, aborting\n";
◆ NormalisationCommonLumi()
Implements the common part to all the strategies in the 'lumi' method:
- find the luminosity that apply to the event
- set the weight accroding to the effective luminosity and efficiency of the trigger
- remove events with null trigger efficiency (NOTE: The method lumi could be used to force the application of a weight. For instance, inv_effective_prescales * total_inv_lumi)
199 double efficiency =
eff(leadingJet);
200 double inv_eff_lumi = (
it->second).
weight;
206 double norm = inv_eff_lumi / efficiency;
207 evnt.weights *= norm;
◆ NormalisationCommonPresc()
Implements the part that is common to all the strategies in the 'presc' method.
- calculating the prescales
- set the weight according to the prescales, efficiency andi total luminosity
- remove events with null trigger efficiency
- Todo:
- this should be changed there are cases in which preL1min != preL1max
155 double preHLT = trigger.PreHLT[
ibit];
156 double preL1min = trigger.PreL1min[
ibit];
157 double preL1max = trigger.PreL1max[
ibit];
158 if (preL1min != preL1max)
159 cerr <<
"\x1B[31m\e[1m" << preL1min <<
' ' << preL1max <<
"\x1B[30m\e[0m\n";
161 double prescale = preHLT * preL1min;
166 static vector<map<pair<int,int>,
int>> prescales(
ntriggers);
167 pair<int,int> runlum = {evnt.runNo, evnt.lumi};
168 if (prescales.at(
ibit).count(runlum))
169 assert(prescales.at(
ibit).at(runlum) == prescale);
171 prescales.at(
ibit)[runlum] = prescale;
175 double efficiency =
eff(leadingJet);
182 evnt.weights *= norm;
◆ NormalisationFwd()
Implements the normalisation and selection used for Mueller-Navelet jets studies. Return the leadingJet. Calculates the luminosity, efficiency, and prescale and set the wieghts accordingly The selection rules for event are the following: -> Must be the most forward or most backward jet -> Must match an hltJet -> The hltJet and recJet matching together must be above the hlt threshold and turnon of the trigger respectivelly. -> If both the most bwd and most fwd jets satisfies the above conditions, the jets with highest absolute rapidity is defined as the leading jet. If the event do not macth these criteria, the event weights are set to zero and the leadingJet is the most forward jet. The event must be discarted when filling a tree.
298 auto mnorder = [] (
const RecJet& recJet1,
const RecJet& recJet2) {
299 return (recJet1.p4.Eta() > recJet2.p4.Eta());
301 sort(recJets.begin(), recJets.end(), mnorder);
304 auto leadingJet = recJets.begin();
305 auto subleadingJet = prev(recJets.end());
311 return recJets.front();
315 float leading_pt = leadingJet->CorrPt();
316 float subleading_pt = subleadingJet->CorrPt();
318 if (leading_pt <
minpt && subleading_pt <
minpt) {
320 return recJets.front();
326 float leadingHlt_pt = leadingHltTk.Pt();
327 float subleadingHlt_pt = subleadingHltTk.Pt();
334 int recturnon =
it->second.turnon;
335 float hltthreshold =
it->first;
336 if ( (leading_pt >= recturnon) && (leadingHlt_pt >= hltthreshold)
337 && (subleading_pt >= recturnon) && (subleadingHlt_pt >= hltthreshold) ) {
338 if ( std::abs(subleadingJet->p4.Eta()) > std::abs(leadingJet->p4.Eta()) ) {
339 swap(leadingJet, subleadingJet);
340 swap(leading_pt, subleading_pt);
341 swap(leadingHlt_pt, subleadingHlt_pt);
345 if (leading_pt >= recturnon && leadingHlt_pt >= hltthreshold)
break;
346 if (subleading_pt >= recturnon && subleadingHlt_pt >= hltthreshold) {
347 swap(leadingJet, subleadingJet);
348 swap(leading_pt, subleading_pt);
349 swap(leadingHlt_pt, subleadingHlt_pt);
360 return recJets.front();
364 cerr <<
"\x1B[31m\e[1mProblem with `ibit` (leading_pt = " << leading_pt <<
")\x1B[30m\e[0m\n";
366 return recJets.front();
370 if (!trigger.Bit.at(
ibit)) {
372 return recJets.front();
◆ NormalisationStd()
Implements the "standard" normalisation. The jets are ordered in pt. The jets with highest pt within ($|\eta| < 2.5$) is the "leading jet". Find the trigger bit corresponding to: $ max \left( \left{ TriggerThresods | TriggerThresholds < p_{T}^{leading} \right} \right) $ Calculates the efficiency and prescales and reweight the events. (NOTE: this done regardless to the fact that the leadingJet
has matched to an hltJet
)
227 auto leadingInTk =
phaseSel(recJets);
230 if (leadingInTk == recJets.end()) {
232 return recJets.front();
235 auto leading_pt = leadingInTk->CorrPt();
237 if (leading_pt <
minpt) {
239 return recJets.front();
248 if (leading_pt >=
it->second.turnon)
break;
259 cerr <<
"\x1B[31m\e[1mProblem with `trigger_lumi` iterator (leading_pt = " << leading_pt <<
")\x1B[30m\e[0m\n";
261 return recJets.front();
264 cerr <<
"\x1B[31m\e[1mProblem with `ibit` (leading_pt = " << leading_pt <<
")\x1B[30m\e[0m\n";
266 return recJets.front();
270 if (!trigger.Bit.at(
ibit)) {
272 return recJets.front();
◆ operator()()
Implements the selection and normalisation. It sets the RecEvent::weights
so that the different triggers are normalised and that the dataset in normalised to the total luminosity. In this respect it is assumed that the highest pt trigger has prescale 1. And that its lumi is the total luminosity.
Apply some preliminary cuts depending depending on triggers. (For more details see the strategies)
◆ eff
◆ ibit
◆ inv_total_lumi
const double inv_total_lumi |
◆ it
◆ minpt
◆ NormalisationCommonPtr
◆ NormalisationPtr
Pointer functions. NormalisationPtr
is set by the constructor to one of private method implementing a specific normalisation strategy (i.e. normalisation,selection of the leading jet, trigger considered as triggering the event). The NormaloisationCommonPtr
pointer implements the part common to all the strategies but implementing different method to normalise, depending on the argument method
in the constructor. It is used by the private methods pointed by NormalisationPtr
.
◆ ntriggers
◆ triggers_lumi
The documentation for this struct was generated from the following file:
std::vector< DAS::RecJet >::iterator phaseSel(std::vector< DAS::RecJet > &recjets)
Definition: PhaseSelection.h:9
method
Definition: Core-gitclone-lastrun.txt:5
cerr
Definition: Ntupliser_cfg.py:101
const RecJet & NormalisationStd(RecEvent &evnt, std::vector< RecJet > &recJets, std::vector< FourVector > &hltJets, Trigger &trigger)
Definition: applyDataNormalisation.h:219
TriggerEff eff
Definition: applyDataNormalisation.h:68
int year
Definition: Ntupliser_cfg.py:63
const RecJet & NormalisationFwd(RecEvent &evnt, std::vector< RecJet > &recJets, std::vector< FourVector > &hltJets, Trigger &trigger)
Definition: applyDataNormalisation.h:290
const RecJet & NormalisationCommonLumi(RecEvent &evnt, const RecJet &leadingJet, const Trigger &trigger)
Definition: applyDataNormalisation.h:194
void reset(Weights &wgts, const double v=0)
Definition: applyDataNormalisation.h:28
std::map< int, TriggerLumi >::const_reverse_iterator it
Definition: applyDataNormalisation.h:73
const size_t ntriggers
Definition: applyDataNormalisation.h:69
const RecJet &(Functor::* NormalisationPtr)(RecEvent &, std::vector< RecJet > &, std::vector< FourVector > &, Trigger &)
Definition: applyDataNormalisation.h:83
DAS::Weight weight
Definition: classes.h:11
const float minpt
Definition: applyDataNormalisation.h:71
size_t ibit
Definition: applyDataNormalisation.h:72
std::map< int, TriggerLumi > GetLumiFromFiles(const std::filesystem::path &lumi_file, const std::filesystem::path &turnon_file)
path to text file with turn-on points
const double inv_total_lumi
Definition: applyDataNormalisation.h:70
void swap(SharedPtr< T > &a, SharedPtr< T > &b)
Definition: fjcore.hh:412
const std::map< int, TriggerLumi > triggers_lumi
Definition: applyDataNormalisation.h:67
DAS::FourVector match(const DAS::FourVector &jet, const std::vector< DAS::FourVector > *hltJets)
Definition: match.h:7
const RecJet & NormalisationCommonPresc(RecEvent &evnt, const RecJet &leadingJet, const Trigger &trigger)
Definition: applyDataNormalisation.h:148
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< float > > FourVector
Definition: PhysicsObject.h:15
const RecJet &(Functor::* NormalisationCommonPtr)(RecEvent &, const RecJet &, const Trigger &)
Definition: applyDataNormalisation.h:84