ATLAS Offline Software
Loading...
Searching...
No Matches
TauAnalysisToolsExample.cxx File Reference
#include <memory>
#include <cstdlib>
#include <iostream>
#include <TFile.h>
#include <TError.h>
#include <TString.h>
#include "xAODCore/ShallowCopy.h"
#include "xAODEventInfo/EventInfo.h"
#include "xAODTau/TauJetContainer.h"
#include "xAODTau/TauJetAuxContainer.h"
#include "AsgTools/ToolHandle.h"
#include "TauAnalysisTools/TauEfficiencyCorrectionsTool.h"
#include "TauAnalysisTools/TauSelectionTool.h"
#include "TauAnalysisTools/TauSmearingTool.h"
#include "TauAnalysisTools/TauTruthMatchingTool.h"
#include "xAODCore/tools/IOStats.h"
#include "xAODCore/tools/ReadStats.h"

Go to the source code of this file.

Macros

#define CHECK(ARG)
#define RETRIEVE(TYPE, CONTAINER, NAME)

Functions

int main (int argc, char *argv[])

Macro Definition Documentation

◆ CHECK

#define CHECK ( ARG)
Value:
do { \
const bool result = ARG; \
if( ! result ) { \
::Error( "TauAnalysisToolsExample", "Failed to execute: \"%s\"", \
#ARG ); \
return 1; \
} \
} while( false )

Definition at line 42 of file TauAnalysisToolsExample.cxx.

42#define CHECK( ARG ) \
43 do { \
44 const bool result = ARG; \
45 if( ! result ) { \
46 ::Error( "TauAnalysisToolsExample", "Failed to execute: \"%s\"", \
47 #ARG ); \
48 return 1; \
49 } \
50 } while( false )

◆ RETRIEVE

#define RETRIEVE ( TYPE,
CONTAINER,
NAME )
Value:
do { \
if (xEvent.contains<TYPE>(NAME)) \
CHECK( xEvent.retrieve( CONTAINER, NAME ) ); \
else \
Warning("TauAnalysisToolsExample","%s container is not available", NAME); \
} while(false) \
#define TYPE(CODE, TYP, IOTYP)

Definition at line 52 of file TauAnalysisToolsExample.cxx.

52#define RETRIEVE( TYPE, CONTAINER , NAME ) \
53 do { \
54 if (xEvent.contains<TYPE>(NAME)) \
55 CHECK( xEvent.retrieve( CONTAINER, NAME ) ); \
56 else \
57 Warning("TauAnalysisToolsExample","%s container is not available", NAME); \
58 } while(false) \
59

Function Documentation

◆ main()

int main ( int argc,
char * argv[] )

Definition at line 60 of file TauAnalysisToolsExample.cxx.

61{
62 StatusCode::enableFailure();
63
64 // Check if we received a file name:
65 if( argc < 2 )
66 {
67 Error( "TauAnalysisToolsExample", "No file name received!" );
68 Error( "TauAnalysisToolsExample", " Usage: %s [xAOD file name]", "TauAnalysisToolsExample" );
69 return 1;
70 }
71
72 // Initialise the application:
73 CHECK( xAOD::Init( "TauAnalysisToolsExample" ) );
74
75 // Open the input file:
76 const TString sInputFileName = argv[ 1 ];
77 Info( "TauAnalysisToolsExample", "Opening input file: %s", sInputFileName.Data() );
78 std::unique_ptr< TFile > fInputFile( TFile::Open( sInputFileName, "READ" ) );
79 CHECK( fInputFile.get() );
80
81 // Create the output file:
82 TString sOutputFileName = "output.root";
83 if (argc>3)
84 sOutputFileName = TString(argv[3]);
85
86 //Do the trigger efficiency tools, requires correct ilumicalc
87 bool m_doTrigger = bool(argc>4);
88
89 Info( "TauAnalysisToolsExample", "Opening output file: %s", sOutputFileName.Data() );
90 std::unique_ptr< TFile > fOutputFile( TFile::Open( sOutputFileName, "RECREATE" ) );
91 CHECK( fOutputFile.get() );
92
93 // Create a TEvent object:
94 // xAOD::TEvent xEvent( xAOD::TEvent::kClassAccess );
96 CHECK( xEvent.readFrom( fInputFile.get() ) );
97
98 // Connect TEvent with output file :
99 CHECK( xEvent.writeTo( fOutputFile.get() ) );
100
101 Info( "TauAnalysisToolsExample", "Number of events in the file: %i",
102 static_cast< int >( xEvent.getEntries() ) );
103
104 // Decide how many events to run over:
105 Long64_t iEntries = xEvent.getEntries();
106 if( argc > 2 )
107 {
108 const Long64_t iMaxEntries = atoll( argv[ 2 ] );
109 if( iMaxEntries < iEntries )
110 {
111 iEntries = iMaxEntries;
112 }
113 }
114
115 // defining needed Container
116 const xAOD::EventInfo* xEventInfo = 0;
117 const xAOD::TauJetContainer* xTauJetContainer = 0;
118
119 // ===========================================================================
120 // TauSelectionTool
121 // ===========================================================================
122 TauAnalysisTools::TauSelectionTool* TauSelTool = new TauAnalysisTools::TauSelectionTool( "TauSelectionTool" );
123 TauSelTool->msg().setLevel( MSG::DEBUG );
124 // preparation for control hisograms
125 TauSelTool->setOutFile( fOutputFile.get() );
126 CHECK(TauSelTool->setProperty("CreateControlPlots", true ));
127 CHECK(TauSelTool->setProperty("JetIDWP", int(JETIDRNNMEDIUM) ));
128 CHECK(TauSelTool->setProperty("EleIDWP", int(ELEIDRNNLOOSE) ));
129 CHECK(TauSelTool->setProperty("EleIDVersion", 1 ));
130 CHECK(TauSelTool->setProperty("PtMin", 20. ));
131 CHECK(TauSelTool->setProperty("ConfigPath", "" ));
132 CHECK(TauSelTool->setProperty("SelectionCuts", int(CutPt|CutJetIDWP|CutEleIDWP) ));
133 CHECK(TauSelTool->initialize());
134
135 // ===========================================================================
136 // TauSmearingTool
137 // ===========================================================================
138 TauAnalysisTools::TauSmearingTool TauSmeTool( "TauSmearingTool" );
139 TauSmeTool.msg().setLevel( MSG::DEBUG );
140 CHECK(TauSmeTool.setProperty("RecommendationTag","2025-prerec"));
141 CHECK(TauSmeTool.setProperty("Campaign","mc23")); // can be also set to mc20 depending on the mc campaign
142 CHECK(TauSmeTool.initialize());
143
144 // restructure all recommended systematic variations for smearing tool
145 std::vector<CP::SystematicSet> vSmearingSystematicSet;
146 for (auto SystematicsVariation : TauSmeTool.recommendedSystematics())
147 {
148 vSmearingSystematicSet.push_back(CP::SystematicSet());
149 vSmearingSystematicSet.back().insert(SystematicsVariation);
150 }
151
152 // ===========================================================================
153 // TauEfficiencyCorrectionsTool
154 // ===========================================================================
155 TauAnalysisTools::TauEfficiencyCorrectionsTool TauEffCorrTool( "TauEfficiencyCorrectionsTool" );
156 TauEffCorrTool.msg().setLevel( MSG::VERBOSE );
157 CHECK(TauEffCorrTool.setProperty("JetIDLevel", static_cast<int>(TauAnalysisTools::JetID::JETIDRNNMEDIUM)));
158 CHECK(TauEffCorrTool.setProperty("EfficiencyCorrectionTypes", static_cast<int>(TauAnalysisTools::EfficiencyCorrectionType::SFJetIDHadTau)));
159 CHECK(TauEffCorrTool.setProperty("RecommendationTag","2025-prerec"));
160 CHECK(TauEffCorrTool.setProperty("Campaign", "mc23")); // can be also set to mc20 depending on the mc campaign
161 CHECK(TauEffCorrTool.initialize());
162
163 // restructure all recommended systematic variations for efficiency tools
164 std::vector<CP::SystematicSet> vEfficiencyCorrectionsSystematicSet;
165 vEfficiencyCorrectionsSystematicSet.push_back(CP::SystematicSet());
166 for (auto SystematicsVariation : TauEffCorrTool.recommendedSystematics())
167 {
168 vEfficiencyCorrectionsSystematicSet.push_back(CP::SystematicSet());
169 vEfficiencyCorrectionsSystematicSet.back().insert(SystematicsVariation);
170 }
171
172 // ===========================================================================
173 // TauEfficiencyCorrectionsTriggerTool
174 // ===========================================================================
175 TauAnalysisTools::TauEfficiencyCorrectionsTool TauEffTrigTool( "TauEfficiencyCorrectionsTriggerTool" );
176 // restructure all recommended systematic variations for efficiency tools
177 std::vector<CP::SystematicSet> vEfficiencyCorrectionsTriggerSystematicSet;
178 if (m_doTrigger){
179
180 TauEffTrigTool.msg().setLevel( MSG::DEBUG );
181 CHECK(TauEffTrigTool.setProperty("EfficiencyCorrectionTypes", std::vector<int>({SFTriggerHadTau}) ));
182 CHECK(TauEffTrigTool.setProperty("TriggerName", "HLT_tau25_mediumRNN_tracktwoMVA" ));
183 CHECK(TauEffTrigTool.setProperty("JetIDLevel", static_cast<int>(JETIDRNNMEDIUM) ));
184 CHECK(TauEffTrigTool.setProperty("Campaign", "mc23a")); // can be also set to mc23d depending on the mc campaign
185 CHECK(TauEffTrigTool.initialize());
186
187 vEfficiencyCorrectionsTriggerSystematicSet.push_back(CP::SystematicSet());
188 for (auto SystematicsVariation : TauEffTrigTool.recommendedSystematics())
189 {
190 vEfficiencyCorrectionsTriggerSystematicSet.push_back(CP::SystematicSet());
191 vEfficiencyCorrectionsTriggerSystematicSet.back().insert(SystematicsVariation);
192 }
193 }
194 // ===========================================================================
195 // TauTruthMatchingTool
196 // ===========================================================================
197 TauAnalysisTools::TauTruthMatchingTool T2MT( "TauTruthMatchingTool");
198 T2MT.msg().setLevel( MSG::INFO );
199 CHECK(T2MT.setProperty("TruthJetContainerName", "AntiKt4TruthDressedWZJets"));
200 CHECK(T2MT.initialize());
201
202 static const SG::ConstAccessor<char> acc_IsTruthMatched("IsTruthMatched");
203 static const SG::ConstAccessor<char> acc_IsHadronicTau("IsHadronicTau");
204 static const SG::ConstAccessor<size_t> acc_numCharged("numCharged");
205 static const SG::ConstAccessor<ElementLink< xAOD::JetContainer >> acc_truthJetLink("truthJetLink");
206 static const SG::ConstAccessor<double> acc_TauSFRecoHadTau("TauScaleFactorReconstructionHadTau");
207 static const SG::ConstAccessor<double> acc_TauSFJetIDHadTau("TauScaleFactorJetIDHadTau");
208 static const SG::ConstAccessor<double> acc_TauSFEleIDHadTau("TauScaleFactorEleIDHadTau");
209 static const SG::ConstAccessor<double> acc_TauSFEleIDEle("TauScaleFactorEleIDElectron");
210 static const SG::ConstAccessor<double> acc_TauSFTrigHadTau("TauScaleFactorTriggerHadTau");
211
212 // Loop over the events:
213 for( Long64_t iEntry = 0; iEntry < iEntries; ++iEntry )
214 {
215 // Tell the object which entry to look at:
216 xEvent.getEntry( iEntry );
217
218 // Print some event information for fun:
219 RETRIEVE(xAOD::EventInfo, xEventInfo, "EventInfo");
220 if (xEventInfo)
221 Info( "TauAnalysisToolsExample",
222 "===>>> start processing event #%i, "
223 "run #%i %i events processed so far <<<===",
224 static_cast< int >( xEventInfo->eventNumber() ),
225 static_cast< int >( xEventInfo->runNumber() ),
226 static_cast< int >( iEntry ) );
227
228
229 //Check TauJet Container Name
230 const char * m_tauJetContainerName = "TauJets";
231 if (xEvent.contains<xAOD::TauJetContainer>(m_tauJetContainerName)){
232 RETRIEVE(xAOD::TauJetContainer, xTauJetContainer, m_tauJetContainerName);
233 }else{
234 m_tauJetContainerName = "AnalysisTauJets";
235 RETRIEVE(xAOD::TauJetContainer, xTauJetContainer, m_tauJetContainerName);
236 }
237 std::pair< xAOD::TauJetContainer*, xAOD::ShallowAuxContainer* >xTauShallowContainer = xAOD::shallowCopyContainer(*xTauJetContainer);
238 if(iEntry==0){
239 Info( "TauAnalysisToolsExample:: TauJetContainer = ",m_tauJetContainerName);
240 }
241
242 // // copy truth particles to get truthparticle link for truth taus to work
243 if (xEvent.contains<xAOD::TruthParticleContainer>("TruthParticles"))
244 CHECK( xEvent.copy("TruthParticles") );
245
246 // copy taus
247 CHECK( xEvent.copy(m_tauJetContainerName) );
248
249 // copy tracks
250 CHECK( xEvent.copy("InDetTrackParticles") );
251
252 // Print tau properties, using the tools:
253 for ( auto xTau : *xTauShallowContainer.first )
254 {
255 // perform truth matching
256 auto xTruthTau = T2MT.getTruth(*xTau);
257
258 if (static_cast<bool>(acc_IsTruthMatched(*xTau)))
259 {
260 if (xTruthTau->isTau())
261 {
262 if (static_cast<bool>(acc_IsHadronicTau(*xTruthTau)))
263 Info( "TauAnalysisToolsExample",
264 "Tau was matched to a truth hadronic tau, which has %i prongs and a charge of %i",
265 int(acc_numCharged(*xTruthTau)),
266 int(xTruthTau->charge()));
267 else
268 Info( "TauAnalysisToolsExample",
269 "Tau was matched to a truth leptonic tau, which has a charge of %i",
270 int(xTruthTau->charge()));
271 }
272 else if (xTruthTau->isElectron())
273 Info( "TauAnalysisToolsExample",
274 "Tau was matched to a truth electron");
275 else if (xTruthTau->isMuon())
276 Info( "TauAnalysisToolsExample",
277 "Tau was matched to a truth muon");
278 }
279 else
280 Info( "TauAnalysisToolsExample", "Tau was not matched to truth" );
281
282 auto xTruthJetLink = acc_truthJetLink(*xTau);
283 if (xTruthJetLink.isValid())
284 {
285 const xAOD::Jet* xTruthJet = *xTruthJetLink;
286 Info( "TauAnalysisToolsExample",
287 "Tau was matched to a truth jet, which has pt=%g, eta=%g, phi=%g, m=%g",
288 xTruthJet->p4().Pt(),
289 xTruthJet->p4().Eta(),
290 xTruthJet->p4().Phi(),
291 xTruthJet->p4().M());
292 }
293 else
294 Info( "TauAnalysisToolsExample", "Tau was not matched to truth jet" );
295
296 Info( "TauAnalysisToolsExample",
297 "Un-Smeared tau pt: %g ",
298 xTau->ptFinalCalib());
299 CHECK( TauSmeTool.applyCorrection(*xTau) );
300 Info( "TauAnalysisToolsExample",
301 "Smeared tau pt: %g ",
302 xTau->pt());
303
304 for (auto sSystematicSet: vSmearingSystematicSet)
305 {
306 CHECK( TauSmeTool.applySystematicVariation(sSystematicSet)) ;
307 CHECK( TauSmeTool.applyCorrection(*xTau) );
308 //Skip TES uncertainty print out for non-had taus
309 if (static_cast<bool>(acc_IsTruthMatched(*xTau)) && xTruthTau->isTau() && static_cast<bool>(acc_IsHadronicTau(*xTruthTau))){
310 Info( "TauAnalysisToolsExample",
311 "Smeared tau pt: %g for type %s ",
312 xTau->pt(),
313 sSystematicSet.name().c_str());
314 }
315 }
316
317 // Select "good" taus:
318 if( ! TauSelTool->accept( *xTau ) ){
319 Info( "TauAnalysisToolsExample",
320 "Tau does not pass selection tool: pt %g ",
321 xTau->pt());
322 continue;
323 }
324
325 for (auto sSystematicSet: vEfficiencyCorrectionsSystematicSet)
326 {
327 CHECK( TauEffCorrTool.applySystematicVariation(sSystematicSet));
328 CHECK( TauEffCorrTool.applyEfficiencyScaleFactor(*xTau) );
329 Info( "TauAnalysisToolsExample",
330 "SystType %s: RecoSF: %g JetIDSF: %g EleOLRSFHadTau: %g EleRNNSFElectron: %g",
331 sSystematicSet.name().c_str(),
332 acc_TauSFRecoHadTau(*xTau),
333 acc_TauSFJetIDHadTau(*xTau),
334 acc_TauSFEleIDHadTau(*xTau),
335 acc_TauSFEleIDEle(*xTau));
336 }
337
338 for (auto sSystematicSet: vEfficiencyCorrectionsTriggerSystematicSet)
339 {
340 CHECK( TauEffTrigTool.applySystematicVariation(sSystematicSet));
341 CHECK( TauEffTrigTool.applyEfficiencyScaleFactor(*xTau) );
342 Info( "TauAnalysisToolsExample",
343 "SystType %s: Trigger: %g",
344 sSystematicSet.name().c_str(),
345 acc_TauSFTrigHadTau(*xTau));
346 }
347 // print some info about the selected tau:
348 Info( "TauAnalysisToolsExample", "Selected tau: pt = %g MeV, eta = %g, phi = %g, prong = %i, charge = %i",
349 xTau->pt(), xTau->eta(), xTau->phi(), int(xTau->nTracks()), int(xTau->charge()));
350 }
351 if (xTauJetContainer->empty())
352 CHECK (T2MT.retrieveTruthTaus());
353 xEvent.fill();
354 }
355
356 TauSelTool->writeControlHistograms();
357 delete TauSelTool;
358 CHECK( xEvent.finishWritingTo(fOutputFile.get()));
359 fOutputFile.get()->Close();
360 Info( "TauAnalysisToolsExample", "Finished writing to file: %s", sOutputFileName.Data() );
361
362 // smart slimming
364
365 return 0;
366}
#define CHECK(ARG)
#define RETRIEVE(TYPE, CONTAINER, NAME)
MsgStream & msg() const
Class to wrap a set of SystematicVariations.
bool empty() const noexcept
Returns true if the collection is empty.
Helper class to provide constant type-safe access to aux data.
virtual void setOutFile(TFile *fOutFile) override
Set output file for control histograms.
virtual StatusCode initialize() override
Function initialising the tool.
virtual void writeControlHistograms() override
Write control histograms to output file.
virtual asg::AcceptData accept(const xAOD::IParticle *p) const override
Get the decision using a generic IParticle pointer.
uint32_t runNumber() const
The current event's run number.
uint64_t eventNumber() const
The current event's event number.
ReadStats & stats()
Access the object belonging to the current thread.
Definition IOStats.cxx:17
static IOStats & instance()
Singleton object accessor.
Definition IOStats.cxx:11
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition Jet_v1.cxx:71
void printSmartSlimmingBranchList(bool autoIncludeLinks=false) const
Print the accessed variables, formatted for smart slimming.
Tool for accessing xAOD files outside of Athena.
@ kAthenaAccess
Access containers/objects like Athena does.
Error
The different types of error that can be flagged in the L1TopoRDO.
Definition Error.h:16
@ Info
Definition ZDCMsg.h:20
Jet_v1 Jet
Definition of the current "jet version".
StatusCode Init(const char *appname)
Function initialising ROOT/PyROOT for using the ATLAS EDM.
Definition Init.cxx:31
EventInfo_v1 EventInfo
Definition of the latest event info version.
std::pair< std::unique_ptr< T >, std::unique_ptr< ShallowAuxContainer > > shallowCopyContainer(const T &cont, const EventContext &ctx)
Function making a shallow copy of a constant container.
setBGCode setTAP setLVL2ErrorBits bool
TauJetContainer_v3 TauJetContainer
Definition of the current "taujet container version".
TruthParticleContainer_v1 TruthParticleContainer
Declare the latest version of the truth particle container.