ATLAS Offline Software
Loading...
Searching...
No Matches
DiTauSelectionCuts.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5
6// local include(s)
9
10using namespace TauAnalysisTools;
11
12//______________________________________________________________________________
14 : m_sName(sName)
15 , m_tDTST(tDTST)
16{}
17
18//______________________________________________________________________________
22
23//______________________________________________________________________________
29
30//______________________________________________________________________________
31std::unique_ptr<TH1F> DiTauSelectionCut::CreateControlPlot(const char* sName, const char* sTitle, int iBins, double dXLow, double dXUp)
32{
33 if (m_tDTST->m_bCreateControlPlots)
34 {
35 auto hHist = std::make_unique<TH1F>(sName, sTitle, iBins, dXLow, dXUp);
36 hHist->SetDirectory(0);
37 return hHist;
38 }
39
40 return nullptr;
41}
42
43//______________________________________________________________________________
48
49//______________________________________________________________________________
54
55
56//______________________________________________________________________________
57void DiTauSelectionCut::setProperty(const std::string& name, const std::string& value)
58{
59 std::map<std::string, std::string&>::iterator it = m_mProperties.find(name);
60 if(it == m_mProperties.end() )
61 throw std::runtime_error (("Undeclared property: " + name + "\n").c_str());
62 it->second = value;
63}
64
65//______________________________________________________________________________
66void DiTauSelectionCut::declareProperty(const std::string& name, std::string& loc)
67{
68 std::pair<std::string, std::string&> p(name, loc);
69 m_mProperties.insert(p);
70}
71
72//______________________________________________________________________________
73std::string DiTauSelectionCut::getProperty(const std::string& name)
74{
75 std::map<std::string, std::string&>::iterator it = m_mProperties.find(name);
76 if(it == m_mProperties.end() )
77 throw std::runtime_error (("Undeclared property: " + name + "\n").c_str());
78
79 return it->second;
80}
81
82
83//_______________________________SelectionCutPt_________________________________
84//______________________________________________________________________________
86 : DiTauSelectionCut("CutPt", tDTST)
87{
88 m_hHistCutPre = CreateControlPlot("hPt_pre","Pt_pre;di-#tau-p_{T} [GeV]; events",1000,0,1000);
89 m_hHistCut = CreateControlPlot("hPt_cut","Pt_cut;di-#tau-p_{T} [GeV]; events",1000,0,1000);
90}
91
92//______________________________________________________________________________
93void DiTauSelectionCutPt::fillHistogram(const xAOD::DiTauJet& xDiTau, TH1F& hHist) const
94{
95 hHist.Fill(xDiTau.pt()/1000.);
96}
97
98//______________________________________________________________________________
100{
101 info.addCut( "Pt",
102 "Selection of ditaus according to their transverse momentum" );
103}
104//______________________________________________________________________________
106 asg::AcceptData& acceptData)
107{
108 // save ditau pt in GeV
109 double pt = xDiTau.pt() / 1000.;
110 // in case of only one entry in vector, run for lower limits
111 if (m_tDTST->m_vPtRegion.size() == 1)
112 {
113 if ( pt >= m_tDTST->m_vPtRegion.at(0) )
114 {
115 acceptData.setCutResult( "Pt", true );
116 return true;
117 }
118 }
119 unsigned int iNumPtRegion = m_tDTST->m_vPtRegion.size()/2;
120 for( unsigned int iPtRegion = 0; iPtRegion < iNumPtRegion; iPtRegion++ )
121 {
122 if ( pt >= m_tDTST->m_vPtRegion.at(iPtRegion*2) and pt <= m_tDTST->m_vPtRegion.at(iPtRegion*2+1))
123 {
124 acceptData.setCutResult( "Pt", true );
125 return true;
126 }
127 }
128 m_tDTST->msg() << MSG::VERBOSE << "DiTau failed pt requirement, ditau pt [GeV]: " << pt << endmsg;
129 return false;
130}
131
132//_____________________________SelectionCutAbsEta_______________________________
133//______________________________________________________________________________
135 : DiTauSelectionCut("CutAbsEta", tDTST)
136{
137 m_hHistCutPre = CreateControlPlot("hEta_pre","Eta_pre;di-#tau-#eta; events",100,-3,3);
138 m_hHistCut = CreateControlPlot("hEta_cut","Eta_cut;di-#tau-#eta; events",100,-3,3);
139}
140
141//______________________________________________________________________________
142void DiTauSelectionCutAbsEta::fillHistogram(const xAOD::DiTauJet& xDiTau, TH1F& hHist) const
143{
144 hHist.Fill(xDiTau.eta());
145}
146
147//______________________________________________________________________________
149{
150 info.addCut( "AbsEta",
151 "Selection of ditaus according to their absolute pseudorapidity" );
152}
153//______________________________________________________________________________
155 asg::AcceptData& acceptData)
156{
157 // check regions of eta, if ditau is in one region then return true; false otherwise
158 unsigned int iNumEtaRegion = m_tDTST->m_vAbsEtaRegion.size()/2;
159 for( unsigned int iEtaRegion = 0; iEtaRegion < iNumEtaRegion; iEtaRegion++ )
160 {
161 if ( std::abs( xDiTau.eta() ) >= m_tDTST->m_vAbsEtaRegion.at(iEtaRegion*2) and std::abs( xDiTau.eta() ) <= m_tDTST->m_vAbsEtaRegion.at(iEtaRegion*2+1))
162 {
163 acceptData.setCutResult( "AbsEta", true );
164 return true;
165 }
166 }
167 m_tDTST->msg() << MSG::VERBOSE << "DiTau failed eta requirement, ditau eta: " << xDiTau.eta() << endmsg;
168 return false;
169}
170
171//_______________________________SelectionCutNSubjets_________________________________
172//______________________________________________________________________________
174 : DiTauSelectionCut("CutNSubjets", tDTST)
175{
176 m_hHistCutPre = CreateControlPlot("hNSubjets_pre","NSubjets_pre;di-#tau NSubjets; events",30,0,30);
177 m_hHistCut = CreateControlPlot("hNSubjets_cut","NSubjets_cut;di-#tau NSubjets; events",30,0,30);
178}
179
180//______________________________________________________________________________
181void DiTauSelectionCutNSubjets::fillHistogram(const xAOD::DiTauJet& xDiTau, TH1F& hHist) const
182{
183 hHist.Fill(xDiTau.nSubjets());
184}
185
186//______________________________________________________________________________
188{
189 info.addCut( "NSubjets",
190 "Selection of ditaus according to their number of subjets" );
191}
192//______________________________________________________________________________
194 asg::AcceptData& acceptData)
195{
196 float nsubjets = xDiTau.nSubjets();
197 // in case of only one entry in vector, run for lower limits
198 if (m_tDTST->m_vNSubjetsRegion.size() == 1)
199 {
200 if ( nsubjets >= m_tDTST->m_vNSubjetsRegion.at(0) )
201 {
202 acceptData.setCutResult( "NSubjets", true );
203 return true;
204 }
205 }
206 unsigned int iNumNSubjetsRegion = m_tDTST->m_vNSubjetsRegion.size()/2;
207 for( unsigned int iNSubjetsRegion = 0; iNSubjetsRegion < iNumNSubjetsRegion; iNSubjetsRegion++ )
208 {
209 if ( nsubjets >= m_tDTST->m_vNSubjetsRegion.at(iNSubjetsRegion*2) and nsubjets <= m_tDTST->m_vNSubjetsRegion.at(iNSubjetsRegion*2+1))
210 {
211 acceptData.setCutResult( "NSubjets", true );
212 return true;
213 }
214 }
215 m_tDTST->msg() << MSG::VERBOSE << "DiTau failed NSubjets requirement, ditau number of subjets: " << nsubjets << endmsg;
216 return false;
217}
218
219//____________________________SelectionCutAbsCharge_____________________________
220//______________________________________________________________________________
222 : DiTauSelectionCut("CutAbsCharge", tDTST),
223 m_bDiTauCharge(-1234)
224{
225 m_hHistCutPre = CreateControlPlot("hCharge_pre","Charge_pre;charge; events",7,-3.5,3.5);
226 m_hHistCut = CreateControlPlot("hCharge_cut","Charge_cut;charge; events",7,-3.5,3.5);
227}
228
229//______________________________________________________________________________
230void DiTauSelectionCutAbsCharge::fillHistogram(const xAOD::DiTauJet& /*xTau*/, TH1F& hHist) const
231{
232 hHist.Fill(m_bDiTauCharge);
233}
234
235//______________________________________________________________________________
237{
238 info.addCut( "AbsCharge",
239 "Selection of taus according to their absolute charge" );
240}
241//______________________________________________________________________________
243 asg::AcceptData& acceptData)
244{
245 m_bDiTauCharge = 0;
246 static const SG::ConstAccessor<float> acc_charge ("charge");
247 if ( acc_charge.isAvailable(xTau) ) {
248 m_bDiTauCharge = acc_charge(xTau);
249 } else {
250 for (const auto& xTrack : xTau.trackLinks()) {
251 if (!xTrack.isValid())
252 continue;
253
254 if(xTau.nSubjets() >= 2){
255 for (int i = 0; i < 2; ++i) { // loop over two leading subjets
256 TLorentzVector tlvSubjet = TLorentzVector();
257 tlvSubjet.SetPtEtaPhiE(xTau.subjetPt(i), xTau.subjetEta(i),
258 xTau.subjetPhi(i), xTau.subjetE(i));
259 double dR = tlvSubjet.DeltaR((*xTrack)->p4());
260 if (dR < 0.1) {
261 m_bDiTauCharge += (*xTrack)->charge();
262 break; //prevents double counting of tracks
263 }
264 } // loop over subjets
265 }
266 } // loop over tracks
267 }
268
269 // check charge, if ditau has one of the charges requiered then return true; false otherwise
270 for( unsigned int iCharge = 0; iCharge < m_tDTST->m_vAbsCharges.size(); iCharge++ )
271 {
272 if ( std::abs( m_bDiTauCharge ) == m_tDTST->m_vAbsCharges.at(iCharge) )
273 {
274 acceptData.setCutResult( "AbsCharge", true );
275 return true;
276 }
277 }
278 m_tDTST->msg() << MSG::VERBOSE << "DiTau failed charge requirement, ditau charge: " << m_bDiTauCharge << endmsg;
279 return false;
280
281}
282
283//___________________________SelectionCutOmniScore_____________________________
284//______________________________________________________________________________
286 : DiTauSelectionCut("CutOmniScore", tDTST)
287{
288 m_hHistCutPre = CreateControlPlot("hOmniScore_pre","OmniScore_pre;OmniScore; events",100,0,1);
289 m_hHistCut = CreateControlPlot("hOmniScore_cut","OmniScore_cut;OmniScore; events",100,0,1);
290}
291//______________________________________________________________________________
293{
294 static const SG::ConstAccessor<float> acc_OmniScore("omni_score");
295 hHist.Fill(acc_OmniScore(xTau));
296}
297//______________________________________________________________________________
299{
300 info.addCut( "OmniScore",
301 "Selection of taus according to their OmniScore" );
302}
303//______________________________________________________________________________
305 asg::AcceptData& acceptData)
306{
307 // check OmniScore score, if tau has a OmniScore score in one of the regions requiered then return true; false otherwise
308 static const SG::ConstAccessor<float> acc ("omni_score");
309 float dOmniScore = acc(xTau);
310 unsigned int iNumOmniScoreRegion = m_tDTST->m_vOmniScoreRegion.size()/2;
311 for( unsigned int iOmniScoreRegion = 0; iOmniScoreRegion < iNumOmniScoreRegion; iOmniScoreRegion++ )
312 {
313 if ( dOmniScore >= m_tDTST->m_vOmniScoreRegion.at(iOmniScoreRegion*2) and dOmniScore <= m_tDTST->m_vOmniScoreRegion.at(iOmniScoreRegion*2+1))
314 {
315 acceptData.setCutResult( "OmniScore", true );
316 return true;
317 }
318 }
319 m_tDTST->msg() << MSG::VERBOSE << "Tau failed OmniScore requirement, tau OmniScore: " << dOmniScore << endmsg;
320 return false;
321}
322
323//_____________________________SelectionCutOmniIDWP______________________________
324//_______________________________________________________________________________
326 : DiTauSelectionCut("CutOmniIDWP", tDTST)
327{
328 m_hHistCutPre = CreateControlPlot("hOmniIDWP_pre","OmniIDWP_pre;; events",8,-.5,7.5);
329 m_hHistCut = CreateControlPlot("hOmniIDWP_cut","OmniIDWP_cut;; events",8,-.5,7.5);
330 // only proceed if histograms are defined
331 if (!m_hHistCutPre or !m_hHistCut)
332 return;
333
334 m_hHistCutPre->GetXaxis()->SetBinLabel(1,"!VeryLoose");
335 m_hHistCutPre->GetXaxis()->SetBinLabel(2,"VeryLoose");
336 m_hHistCutPre->GetXaxis()->SetBinLabel(3,"!Loose");
337 m_hHistCutPre->GetXaxis()->SetBinLabel(4,"Loose");
338 m_hHistCutPre->GetXaxis()->SetBinLabel(5,"!Medium");
339 m_hHistCutPre->GetXaxis()->SetBinLabel(6,"Medium");
340 m_hHistCutPre->GetXaxis()->SetBinLabel(7,"!Tight");
341 m_hHistCutPre->GetXaxis()->SetBinLabel(8,"Tight");
342 m_hHistCut->GetXaxis()->SetBinLabel(1,"!VeryLoose");
343 m_hHistCut->GetXaxis()->SetBinLabel(2,"VeryLoose");
344 m_hHistCut->GetXaxis()->SetBinLabel(3,"!Loose");
345 m_hHistCut->GetXaxis()->SetBinLabel(4,"Loose");
346 m_hHistCut->GetXaxis()->SetBinLabel(5,"!Medium");
347 m_hHistCut->GetXaxis()->SetBinLabel(6,"Medium");
348 m_hHistCut->GetXaxis()->SetBinLabel(7,"!Tight");
349 m_hHistCut->GetXaxis()->SetBinLabel(8,"Tight");
350}
351
352//______________________________________________________________________________
354{
355 if(m_tDTST->m_useOmniScore){
356 static const SG::ConstAccessor<char> acc_OmniVeryLoose("omni_score_VL");
357 static const SG::ConstAccessor<char> acc_OmniLoose("omni_score_L");
358 static const SG::ConstAccessor<char> acc_OmniMedium("omni_score_M");
359 static const SG::ConstAccessor<char> acc_OmniTight("omni_score_T");
360 hHist.Fill(acc_OmniVeryLoose(xTau));
361 hHist.Fill(acc_OmniLoose(xTau)+2);
362 hHist.Fill(acc_OmniMedium(xTau)+4);
363 hHist.Fill(acc_OmniTight(xTau)+6);
364 }
365}
366
367//______________________________________________________________________________
369{
370 info.addCut( "OmniIDWP",
371 "Selection of ditaus according to their OmniIDScore" );
372}
373//______________________________________________________________________________
375 asg::AcceptData& acceptData)
376{
377 // check Omni ID working point, if ditau passes OmniID working point then return true; false otherwise
378 bool bPass = false;
379 switch (m_tDTST->m_iOmniIDWP)
380 {
381 case OMNIIDNONE:
382 bPass = true;
383 break;
384 case OMNIIDVERYLOOSE:
385 static const SG::ConstAccessor<char> acc_OmniVeryLoose("omni_score_VL");
386 if (!acc_OmniVeryLoose.isAvailable(xTau)) m_tDTST->msg() << MSG::WARNING << "Omni VeryLoose WP not available" << endmsg;
387 else bPass = acc_OmniVeryLoose(xTau);
388 break;
389 case OMNIIDLOOSE:
390 static const SG::ConstAccessor<char> acc_OmniLoose("omni_score_L");
391 if (!acc_OmniLoose.isAvailable(xTau)) m_tDTST->msg() << MSG::WARNING << "Omni Loose WP not available" << endmsg;
392 else bPass = acc_OmniLoose(xTau);
393 break;
394 case OMNIIDMEDIUM:
395 static const SG::ConstAccessor<char> acc_OmniMedium("omni_score_M");
396 if (!acc_OmniMedium.isAvailable(xTau)) m_tDTST->msg() << MSG::WARNING << "Omni Medium WP not available" << endmsg;
397 else bPass = acc_OmniMedium(xTau);
398 break;
399 case OMNIIDTIGHT:
400 static const SG::ConstAccessor<char> acc_OmniTight("omni_score_T");
401 if (!acc_OmniTight.isAvailable(xTau)) m_tDTST->msg() << MSG::WARNING << "Omni Tight WP not available" << endmsg;
402 else bPass = acc_OmniTight(xTau);
403 break;
404 default:
405 m_tDTST->msg() << MSG::WARNING << "The Omni ID working point with the enum " << m_tDTST->m_iOmniIDWP << " is not available" << endmsg;
406 break;
407 }
408 if (bPass)
409 {
410 acceptData.setCutResult( "OmniIDWP", true );
411 return true;
412 }
413 m_tDTST->msg() << MSG::VERBOSE << "DiTau failed OmniIDWP requirement" << endmsg;
414 return false;
415}
416
417
#define endmsg
Helper class to provide constant type-safe access to aux data.
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
virtual void setAcceptInfo(asg::AcceptInfo &info) const override
virtual void fillHistogram(const xAOD::DiTauJet &xTau, TH1F &hHist) const override
virtual bool accept(const xAOD::DiTauJet &xTau, asg::AcceptData &accept) override
virtual void fillHistogram(const xAOD::DiTauJet &xTau, TH1F &hHist) const override
virtual void setAcceptInfo(asg::AcceptInfo &info) const override
virtual bool accept(const xAOD::DiTauJet &xTau, asg::AcceptData &accept) override
virtual void setAcceptInfo(asg::AcceptInfo &info) const override
virtual void fillHistogram(const xAOD::DiTauJet &xTau, TH1F &hHist) const override
virtual bool accept(const xAOD::DiTauJet &xTau, asg::AcceptData &accept) override
virtual bool accept(const xAOD::DiTauJet &xTau, asg::AcceptData &accept) override
virtual void setAcceptInfo(asg::AcceptInfo &info) const override
virtual void fillHistogram(const xAOD::DiTauJet &xTau, TH1F &hHist) const override
virtual bool accept(const xAOD::DiTauJet &xTau, asg::AcceptData &accept) override
virtual void fillHistogram(const xAOD::DiTauJet &xTau, TH1F &hHist) const override
virtual void setAcceptInfo(asg::AcceptInfo &info) const override
virtual void fillHistogram(const xAOD::DiTauJet &xTau, TH1F &hHist) const override
virtual void setAcceptInfo(asg::AcceptInfo &info) const override
DiTauSelectionCutPt(DiTauSelectionTool *tDTST)
virtual bool accept(const xAOD::DiTauJet &xTau, asg::AcceptData &accept) override
std::unique_ptr< TH1F > CreateControlPlot(const char *sName, const char *sTitle, int iBins, double dXLow, double dXUp)
std::map< std::string, std::string & > m_mProperties
void fillHistogramCut(const xAOD::DiTauJet &xTau)
void fillHistogramCutPre(const xAOD::DiTauJet &xTau)
void setProperty(const std::string &name, const std::string &value)
void declareProperty(const std::string &name, std::string &loc)
DiTauSelectionCut(const std::string &sName, TauAnalysisTools::DiTauSelectionTool *tDTST)
virtual void fillHistogram(const xAOD::DiTauJet &xTau, TH1F &hHist) const =0
std::string getProperty(const std::string &name)
void setCutResult(const std::string &cutName, bool cutResult)
Set the result of a cut, based on the cut name (safer)
Definition AcceptData.h:134
virtual double eta() const
The pseudorapidity ( ) of the particle.
float subjetEta(unsigned int numSubjet) const
virtual double pt() const
The transverse momentum ( ) of the particle.
float subjetE(unsigned int numSubjet) const
float subjetPt(unsigned int numSubjet) const
float subjetPhi(unsigned int numSubjet) const
const TrackParticleLinks_t & trackLinks() const
float nSubjets() const
DiTauJet_v1 DiTauJet
Definition of the current version.
Definition DiTauJet.h:17