ATLAS Offline Software
Loading...
Searching...
No Matches
JetTileCorrectionTester.cxx File Reference

Contains RootCore testing code for the JetTileCorrectionTool. More...

#include <memory>
#include <cstdlib>
#include "TFile.h"
#include "TError.h"
#include "TString.h"
#include "xAODEventInfo/EventInfo.h"
#include "xAODJet/JetContainer.h"
#include "xAODCore/ShallowCopy.h"
#include <AsgTools/StandaloneToolHandle.h>
#include "AsgTools/AsgTool.h"
#include "PATInterfaces/CorrectionCode.h"
#include "JetTileCorrection/JetTileCorrectionTool.h"

Go to the source code of this file.

Macros

#define CHECK(ARG)

Functions

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

Variables

static SG::AuxElement::Accessor< unsigned int > acc_tileok ("TileStatus")
static SG::AuxElement::Accessor< float > acc_ptraw ("Ptraw")

Detailed Description

Contains RootCore testing code for the JetTileCorrectionTool.

Author
Martin Tripiana tripi.nosp@m.ana@.nosp@m.cern..nosp@m.ch

Definition in file JetTileCorrectionTester.cxx.

Macro Definition Documentation

◆ CHECK

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

Definition at line 45 of file JetTileCorrectionTester.cxx.

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

Function Documentation

◆ main()

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

Definition at line 65 of file JetTileCorrectionTester.cxx.

66{
67
68 // The application's name
69 const char* APP_NAME = argv[ 0 ];
70
71 // Check if we received a file name
72 if(argc < 2) {
73 Error( APP_NAME, "No file name received!" );
74 Error( APP_NAME, " Usage: %s [xAOD file name] [num events]", APP_NAME );
75 return 1;
76 }
77
78 // Initialise the application
80 StatusCode::enableFailure();
81
82 // Open the input file
83 const TString fileName = argv[ 1 ];
84 Info(APP_NAME, "Opening file: %s", fileName.Data());
85 std::unique_ptr<TFile> ifile(TFile::Open(fileName, "READ"));
86 CHECK( ifile.get() );
87
88 // Create a TEvent object
91 CHECK( event.readFrom(ifile.get()) );
92 Info(APP_NAME, "Number of events in the file: %i",
93 static_cast<int>(event.getEntries()));
94
95 // Decide how many events to run over
96 Long64_t entries = event.getEntries();
97 if(argc > 2) {
98 const Long64_t e = atoll(argv[2]);
99 entries = std::min(e, entries);
100 }
101
102 // Initialize the tool
104 tool_jtc.setTypeAndName("CP::JetTileCorrectionTool/JetTileCorrectionTool");
105
106 //Set properties if needed
107 CHECK( tool_jtc.setProperty("CorrectionFileName", "JetTileCorrection/JetTile_pFile_010216.root") ); //default anyway
108 // std::vector<std::string> dead_modules = {"1 04","1 05","1 06","1 07","1 08","1 09",
109 // "2 04","2 05","2 06","2 07","2 08","2 09" }; // LBA/C5-10 : NOT REAL SCENARIO !! just trying to get some magnified effect for testing!
110 // CHECK( tool_jtc.setProperty("UserMaskedRegions", dead_modules));
111
112 CHECK( tool_jtc.retrieve() );
113
114 // Loop over the events
115 std::cout << "Starting loop" << std::endl;
116 for(Long64_t entry = 0; entry < entries; ++entry){
117
118 event.getEntry(entry);
119
120 // Print some event information for fun
121 const xAOD::EventInfo* evtInfo = 0;
122 CHECK( event.retrieve(evtInfo, "EventInfo") );
123 if ( entry%100==0 ){
124 Info(APP_NAME, "===>>> Processing entry %lli, run %u, event %lu <<<===",
125 entry, evtInfo->runNumber(), evtInfo->eventNumber());
126 }
127
128 // Get jets
129 const xAOD::JetContainer* jets = 0;
130 CHECK( event.retrieve(jets, "AntiKt4EMTopoJets") );
131
132 std::pair<xAOD::JetContainer*, xAOD::ShallowAuxContainer*> shallowcopy =
134 std::unique_ptr< xAOD::JetContainer > jets_sc( shallowcopy.first );
135 std::unique_ptr< xAOD::ShallowAuxContainer >
136 jets_scaux( shallowcopy.second );
137
138 for( xAOD::Jet* jet : *jets_sc ){
139
140 if (jet->pt() < 20000. || fabs(jet->eta()) > 2.8) continue;
141
142 //--- apply tile dead module correction
143 const CP::CorrectionCode retCode = tool_jtc->applyCorrection(*jet);
144
146 Warning("JetTileCorrectionTester","No valid pt/eta range. No correction applied."); // It might get too verbosed, just here to illustrate all return values.
147 }
148 else if( retCode != CP::CorrectionCode::Ok ){
149 Error("JetTileCorrectionTester","Failed to apply JetTileCorrection!");
150 //return StatusCode::FAILURE;
151 }
152
153 unsigned int j_status = acc_tileok(*jet);
154
155 std::string str_status="";
156 if(j_status == (unsigned int)JTC::TS::GOOD)
157 str_status = "NotAffected";
158 else if(j_status == (unsigned int)JTC::TS::EDGE)
159 str_status = "EdgeAffected";
160 else if(j_status == (unsigned int)JTC::TS::CORE)
161 str_status = "CoreAffected";
162 else
163 str_status = "Unknown";
164
165 Info(APP_NAME, "Jet status : %s, Pt raw = %.3f GeV, Pt corrected %.3f GeV", str_status.c_str(), acc_ptraw(*jet)*0.001, jet->pt()*0.001);
166 }
167
168
169
170 // //Check status only
171 // tool_jtc->setRJET(0.1); //change jet radius (for tile status checks only!)
172 // for( xAOD::Jet* jet : *jets_sc ){
173
174 // JTC::TS j_status = tool_jtc->getTileStatus(*jet);
175
176 // //or well:
177 // //CHECK( tool_jtc->addTileStatus(*jet) );
178 // //unsigned int j_status = acc_tileok(*jet);
179
180 // std::string str_status="";
181 // if(j_status == JTC::TS::GOOD)
182 // str_status = "NotAffected";
183 // else if(j_status == JTC::TS::EDGE)
184 // str_status = "EdgeAffected";
185 // else if(j_status == JTC::TS::CORE)
186 // str_status = "CoreAffected";
187 // else
188 // str_status = "Unknown";
189
190 // Info(APP_NAME, "Jet status : %s, Pt raw = %.3f GeV, Pt corrected %.3f GeV", str_status.c_str(), acc_ptraw(*jet)*0.001, jet->pt()*0.001);
191
192 // }
193 // //back to default
194 // tool_jtc->setRJET(0.4);
195
196 }
197
198 Info(APP_NAME, "Application finished successfully");
199
200 return 0;
201}
static const SG::Accessor< float > acc_ptraw("Ptraw")
static SG::AuxElement::Accessor< unsigned int > acc_tileok("TileStatus")
#define CHECK(ARG)
Return value from object correction CP tools.
@ OutOfValidityRange
Input object is out of validity range.
@ Ok
The correction was done successfully.
an "initializing" ToolHandle for stand-alone applications
StatusCode setProperty(const std::string &name, T2 &&value)
StatusCode retrieve()
initialize the tool, will succeed if the tool was already initialized
void setTypeAndName(const std::string &typeAndName)
uint32_t runNumber() const
The current event's run number.
uint64_t eventNumber() const
The current event's event number.
Tool for accessing xAOD files outside of Athena.
@ kAthenaAccess
Access containers/objects like Athena does.
A relatively simple transient store for objects created in analysis.
Definition TStore.h:45
double entries
Definition listroot.cxx:49
Error
The different types of error that can be flagged in the L1TopoRDO.
Definition Error.h:16
TestStore store
Definition TestStore.cxx:23
@ 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.
JetContainer_v1 JetContainer
Definition of the current "jet container version".

Variable Documentation

◆ acc_ptraw

SG::AuxElement::Accessor< float > acc_ptraw("Ptraw") ( "Ptraw" )
static

◆ acc_tileok

SG::AuxElement::Accessor< unsigned int > acc_tileok("TileStatus") ( "TileStatus" )
static