ATLAS Offline Software
Loading...
Searching...
No Matches
Rivet_i.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// Implementation file for Athena-Rivet interface
6
7#include "Rivet_i.h"
8#include "LogLevels.h"
9
10#include "AtlasHepMC/GenEvent.h"
12
16
17#include "GaudiKernel/IAppMgrUI.h"
18#include "GaudiKernel/Bootstrap.h"
19#include "GaudiKernel/ITHistSvc.h"
20
21#include "Rivet/Rivet.hh"
22#include "Rivet/Analysis.hh"
23#include "Rivet/Config/RivetConfig.hh"
24#include "Rivet/Tools/RivetYODA.hh"
25
26
27#include <cstdlib>
28#include <cstdio>
29#include <memory>
30#include <regex>
31
32
34
35Rivet_i::Rivet_i(const std::string& name, ISvcLocator* pSvcLocator) :
36 AthAlgorithm(name, pSvcLocator),
37 m_needsConversion(false),
39 m_init(false)
40{
41 // Options
42 declareProperty("McEventKey", m_genEventKey="GEN_EVENT");
44 declareProperty("CrossSection", m_crossSection=0.0);
45 declareProperty("CrossSectionUncertainty", m_crossSection_uncert=0.0);
46 declareProperty("Stream", m_stream="/Rivet");
47 declareProperty("RunName", m_runname="");
48 declareProperty("HistoFile", m_file="Rivet.yoda.gz");
49 declareProperty("HistoPreload", m_preload="");
50 declareProperty("AnalysisPath", m_anapath="");
51 declareProperty("IgnoreBeamCheck", m_ignorebeams=false);
52 declareProperty("DoRootHistos", m_doRootHistos=false);
53 declareProperty("SkipWeights", m_skipweights=false);
54 declareProperty("MatchWeights", m_matchWeights="");
55 declareProperty("UnmatchWeights", m_unmatchWeights="");
56 declareProperty("NominalWeightName", m_nominalWeightName="");
57 declareProperty("WeightCap", m_weightcap=-1.0);
58 declareProperty("AddMissingBeamParticles",m_patchBeams=false);
59}
60
61std::string getenv_str(const std::string& key) {
62 char const* val = getenv(key.c_str());
63 return val == NULL ? std::string() : std::string(val);
64}
65
66
67StatusCode Rivet_i::initialize ATLAS_NOT_THREAD_SAFE () {
68 // ^ use of setenv
69 ATH_MSG_INFO("Using Rivet version " << Rivet::version());
70
71 ATH_CHECK(m_evtInfoKey.initialize());
72
73 // Get tool for xAOD::Truth to HepMC::GenEvent conversion
74 ATH_CHECK(m_xAODtoHepMCTool.retrieve());
75
76
77 // Set RIVET_ANALYSIS_PATH based on alg setup
78
79 // First set (overwrite, if necessary) the RIVET_ANALYSIS_PATH variable
80 std::string env_rap(getenv_str("RIVET_ANALYSIS_PATH"));
81 if (!m_anapath.empty()) {
82 ATH_MSG_INFO("Setting Rivet plugin analyses loader env path: " << m_anapath);
83 if (!env_rap.empty()) ATH_MSG_INFO("Overwriting environment's RIVET_ANALYSIS_PATH = " << env_rap << "!");
84 setenv("RIVET_ANALYSIS_PATH", m_anapath.c_str(), 1);
85 }
86
87 // Get the final form of RIVET_ANALYSIS_PATH and talk about it
88 env_rap = getenv_str("RIVET_ANALYSIS_PATH");
89 if (!env_rap.empty()) ATH_MSG_DEBUG("Loading Rivet plugin analyses from env path: " << env_rap);
90
91 // Set up analysis handler
92 #if RIVET_VERSION_CODE >= 40100
93 m_analysisHandler = new Rivet::AnalysisHandler();
94 #else
95 m_analysisHandler = new Rivet::AnalysisHandler(m_runname);
96 #endif
97 assert(m_analysisHandler);
98
99 #if RIVET_VERSION_CODE >= 40000
100 m_analysisHandler->setCheckBeams(!m_ignorebeams); //< Whether to do beam ID/energy consistency checks
101 m_analysisHandler->matchWeightNames(m_matchWeights); //< Only run on a subset of the multi-weights
102 m_analysisHandler->unmatchWeightNames(m_unmatchWeights); //< Veto a subset of the multi-weights
103 #else
104 m_analysisHandler->setIgnoreBeams(m_ignorebeams); //< Whether to do beam ID/energy consistency checks
105 m_analysisHandler->selectMultiWeights(m_matchWeights); //< Only run on a subset of the multi-weights
106 m_analysisHandler->deselectMultiWeights(m_unmatchWeights); //< Veto a subset of the multi-weights
107 #endif
108 m_analysisHandler->skipMultiWeights(m_skipweights); //< Only run on the nominal weight
109 m_analysisHandler->setNominalWeightName(m_nominalWeightName);
110 if (m_weightcap>0) m_analysisHandler->setWeightCap(m_weightcap);
111
112 // Set Rivet native log level to match Athena
113 Rivet::Log::setLevel("Rivet", rivetLevel(msg().level()));
114
115 // Get all available analysis names
116 if (msgLvl(MSG::VERBOSE)) {
117 std::vector<std::string> analysisNames = Rivet::AnalysisLoader::analysisNames();
118 ATH_MSG_VERBOSE("List of available Rivet analyses:");
119 for (const std::string& a : analysisNames) ATH_MSG_VERBOSE(" " + a);
120 }
121
122 // Add analyses
123 for (const std::string& a : m_analysisNames) {
124 ATH_MSG_INFO("Loading Rivet analysis " << a);
125 m_analysisHandler->addAnalysis(a);
126 Rivet::Log::setLevel("Rivet.Analysis."+a, rivetLevel(msg().level()));
127 }
128
129 // Initialise Rivet
130 // m_analysisHandler->init();
131
132 //load a pre-existing yoda file to initialize histograms
133 if (m_preload!= "") {
134 m_analysisHandler->readData(m_preload);
135 }
136
137 return StatusCode::SUCCESS;
138}
139
140
141StatusCode Rivet_i::execute() {
142 ATH_MSG_DEBUG("Rivet_i execute");
143
144 const EventContext& ctx = getContext();
145
147 ATH_MSG_DEBUG("Rivet_i needs xAOD::Truth to HepMC::GenEvent conversion? " << m_needsConversion);
148
149 std::unique_ptr<HepMC::GenEvent> checkedEvent;
150 if (m_needsConversion) {
151 const xAOD::TruthEventContainer* truthCollection = nullptr;
152 if (evtStore()->retrieve(truthCollection, "TruthEvents").isFailure()) {
153 ATH_MSG_ERROR("Could not retrieve TruthEvents collection, aborting.");
154 return StatusCode::FAILURE;
155 }
156
158 ATH_MSG_DEBUG("Attempt xAOD::Truth to HepMC::GenEvent conversion.");
159 std::vector<HepMC::GenEvent>&& hepmc_evts = m_xAODtoHepMCTool->getHepMCEvents(truthCollection, eventInfo.cptr());
160
161 if (hepmc_evts.empty()) {
162 ATH_MSG_ERROR("Conversion didn't yield HepMC::GenEvent, aborting.");
163 return StatusCode::FAILURE;
164 }
165 checkedEvent = checkEvent(std::move(hepmc_evts[0]), ctx);
166 }
167 else {
168 // Get the event collection
170 const McEventCollection* eventCollection = nullptr;
171 StatusCode sc = evtStore()->retrieve(eventCollection, m_genEventKey);
172 if (sc.isFailure() || eventCollection == 0) {
173 ATH_MSG_ERROR("Unable to retrieve event collection from StoreGate with key " << m_genEventKey);
174 return StatusCode::FAILURE;
175 } else {
176 ATH_MSG_DEBUG("Retrieved event collection from StoreGate with key " << m_genEventKey);
177 }
178
179 // Get the first event in the event collection
182 const HepMC::GenEvent* event = eventCollection->front();
183 if (event == nullptr) {
184 ATH_MSG_ERROR("Rivet_i received a null HepMC event");
185 return StatusCode::FAILURE;
186 }
187
188 // Modify the event units etc. if necessary
189 checkedEvent = checkEvent(*event, ctx);
190 }
191
192 if(!checkedEvent) {
193 ATH_MSG_ERROR("Check on HepMC event failed!");
194 return StatusCode::FAILURE;
195 }
196
197 // Initialize Rivet (on the first event only)
198 if (!m_init) {
199 m_analysisHandler->init(*checkedEvent);
200 m_init = true;
201 }
202
203 // Analyse the event
204 m_analysisHandler->analyze(*checkedEvent);
205
206 return StatusCode::SUCCESS;
207}
208
209
210StatusCode Rivet_i::finalize() {
211 ATH_MSG_INFO("Rivet_i finalizing");
212
213 // Setting cross-section in Rivet
214 // If no user-specified cross-section available,
215 // take whatever is in the GenEvent
217
218 // Call Rivet finalize
219 m_analysisHandler->finalize();
220
221
222
223 // Write out YODA file (add .yoda suffix if missing)
224 if (m_file.find(".yoda") == std::string::npos) m_file += ".yoda.gz";
225 auto pos = m_file.find(".gz.");
226 if (pos == std::string::npos) {
227 m_analysisHandler->writeData(m_file);
228 }
229 else {//filename is e.g. something.yoda.gz.1
230 //write to output file name without the ".1"
231 m_analysisHandler->writeData(m_file.substr(0, pos)+".gz");
232 // then rename it as the requested output filename
233 if (std::rename((m_file.substr(0, pos)+".gz").c_str(), m_file.c_str()) !=0) {
234 std::string error_msg = "Impossible to rename ";
235 error_msg += m_file.substr(0, pos)+".gz";
236 error_msg += " as ";
237 error_msg += m_file;
238 ATH_MSG_ERROR(error_msg.c_str());
239 return StatusCode::FAILURE;
240 }
241 }
242
243 return StatusCode::SUCCESS;
244}
245
246
250 return a->momentum().e() > b->momentum().e();
251}
252
253inline std::vector<std::string> split(const std::string& input, const std::string& regex) {
254 // passing -1 as the submatch index parameter performs splitting
255 std::regex re(regex);
256 std::sregex_token_iterator
257 first{input.begin(), input.end(), re, -1},
258 last;
259 return {first, last};
260}
261
262std::unique_ptr<HepMC::GenEvent> Rivet_i::checkEvent(const HepMC::GenEvent& event, const EventContext& ctx) {
263 auto modEvent = std::make_unique<HepMC::GenEvent>(event);
264
265 if (!m_needsConversion) {
266 // overwrite the HEPMC dummy event number with the proper ATLAS event number
268 modEvent->set_event_number(static_cast<int>(evtInfo->eventNumber()));
269 }
270
271 // weight-name cleaning
272#ifdef HEPMC3
273 std::shared_ptr<HepMC3::GenRunInfo> modRunInfo;
274 if (event.run_info()) {
275 modRunInfo = std::make_shared<HepMC3::GenRunInfo>(*(event.run_info().get()));
276 if (event.run_info()->weight_names().empty() && event.weights().size() == 1) {
277 modRunInfo->set_weight_names({"Default"});
278 }
279 }
280 else {
281 ATH_MSG_DEBUG("No run info, event weights size is " << event.weights().size() );
282 modRunInfo = std::make_shared<HepMC3::GenRunInfo>();
283 std::vector<std::string> w_names;
284 for (size_t i = 0; i < event.weights().size(); i++) { w_names.push_back(std::string("badweight") + std::to_string(i)); }
285 modRunInfo->set_weight_names(w_names);
286 }
287 modEvent->set_run_info(std::move(modRunInfo));
288 std::vector<std::string> w_names = modEvent->weight_names();
289 if (w_names.size()) {
290 std::vector<std::pair<std::string,std::string> > w_subs = {
291 {" nominal ",""},
292 {" set = ","_"},
293 {" = ","_"},
294 {"=",""},
295 {",",""},
296 {".",""},
297 {":",""},
298 {" ","_"},
299 {"#","num"},
300 {"\n","_"},
301 {"/","over"}
302 };
303 for (std::string& wname : w_names) {
304 for (const auto& sub : w_subs) {
305 size_t start_pos = wname.find(sub.first);
306 while (start_pos != std::string::npos) {
307 wname.replace(start_pos, sub.first.length(), sub.second);
308 start_pos = wname.find(sub.first);
309 }
310 }
311 }
312 modEvent->run_info()->set_weight_names(w_names);
313 }
314#else
315 const HepMC::WeightContainer& old_wc = event.weights();
316 std::vector<std::string> old_wnames = old_wc.weight_names();
317 if (old_wnames.size()) {
318 HepMC::WeightContainer& new_wc = modEvent->weights();
319 new_wc.clear();
320 std::vector<std::pair<std::string,std::string> > w_subs = {
321 {" nominal ",""},
322 {" set = ","_"},
323 {" = ","_"},
324 {"=",""},
325 {",",""},
326 {".",""},
327 {":",""},
328 {" ","_"},
329 {"#","num"},
330 {"\n","_"},
331 {"/","over"}
332 };
333 std::map<std::string, double> new_name_to_value;
334 std::map<std::string, std::string> old_name_to_new_name;
335 for (const std::string& old_name : old_wnames) {
336 std::string wname = std::string(old_name);
337 double value = old_wc[old_name];
338 for (const auto& sub : w_subs) {
339 size_t start_pos = wname.find(sub.first);
340 while (start_pos != std::string::npos) {
341 wname.replace(start_pos, sub.first.length(), sub.second);
342 start_pos = wname.find(sub.first);
343 }
344 }
345 new_name_to_value[wname] = value;
346 old_name_to_new_name[old_name] = wname;
347 }
348 auto itEnd = old_name_to_new_name.end();
349 for (const std::string& old_name : old_wnames) {
350 if (old_name_to_new_name.find(old_name) == itEnd) continue;
351 const std::string& new_name = old_name_to_new_name[old_name];
352 new_wc[ new_name ] = new_name_to_value[new_name];
353 }
354 // end of weight-name cleaning
355 }
356#endif
357
358#ifdef HEPMC3
359 modEvent->set_units(HepMC3::Units::GEV, HepMC3::Units::MM);
360 if (modEvent->particles().size() == 1) modEvent->set_beam_particles(modEvent->particles().front(), modEvent->particles().front());
361 if (m_patchBeams) {
362 // workaround in case beam particles were missing in xAOD::Truth
363 // first work out centre-of-mass energy
364 double etotal = 0;
365 for (const auto& p : modEvent->particles()) {
366 if (!MC::isStable(p)) continue;
367 etotal += p->momentum().e();
368 }
369 const double ebeam = 0.5*etotal;
370 // create dummy beam particles
371 HepMC::GenParticlePtr b1 = std::make_shared<HepMC::GenParticle>(HepMC::FourVector(0,0,ebeam,ebeam),2212,4);
372 HepMC::GenParticlePtr b2 = std::make_shared<HepMC::GenParticle>(HepMC::FourVector(0,0,-1*ebeam,ebeam),2212,4);
373 // remove particles incorrectly labelled as beam particles
374 std::vector<HepMC::GenParticlePtr> notBeams;
375 for (const HepMC::GenParticlePtr& p : modEvent->beams()) {
376 if (!MC::isBeam(p)) notBeams.push_back(p);
377 }
378 //AV: the loop is over shared pointers! Should be const auto&
379 for (const auto& bp : notBeams) bp->production_vertex()->remove_particle_out(bp);
380 // add dummy beam particles
381 modEvent->set_beam_particles(std::move(b1), std::move(b2));
382 }
383 if (modEvent->beams().front()->momentum().e() > 50000.0) {
384 MeV2GeV(*modEvent);
385 }
386#else
387 modEvent->use_units(HepMC::Units::GEV, HepMC::Units::MM);
388 if (modEvent->particles_size() == 1) modEvent->set_beam_particles(*modEvent->particles_begin(), *modEvent->particles_begin());
389 if (modEvent->beam_particles().first->momentum().e() > 50000.0) {
390 MeV2GeV(*modEvent);
391 }
392#endif
393
394 return modEvent;
395}
396
397void Rivet_i::MeV2GeV(HepMC::GenEvent& evt) {
398#ifdef HEPMC3
399 for (auto& p: evt.particles()) {
400 p->set_momentum(p->momentum()*0.001);
401#else
402 for (HepMC::GenParticlePtr p: evt) {
403 const HepMC::FourVector& mom = p->momentum();
404 p->set_momentum(HepMC::FourVector (mom.px()*0.001,
405 mom.py()*0.001,
406 mom.pz()*0.001,
407 mom.e()*0.001));
408#endif
409 p->set_generated_mass(p->generated_mass()*0.001);
410 }
411}
const boost::regex re(r_e)
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_DEBUG(x)
Helpers for checking error return status codes and reporting errors.
ATLAS-specific HepMC functions.
static Double_t a
static Double_t sc
Rivet::Log::Level rivetLevel(MSG::Level gaudiLevel)
Definition LogLevels.h:21
bool cmpGenParticleByEDesc(HepMC::ConstGenParticlePtr a, HepMC::ConstGenParticlePtr b)
Helper function to sort GenParticles by descending energy.
Definition Rivet_i.cxx:249
std::string getenv_str(const std::string &key)
Definition Rivet_i.cxx:61
StatusCode Rivet_i::initialize ATLAS_NOT_THREAD_SAFE()
Install fatal handler with default options.
Definition Rivet_i.cxx:67
std::vector< std::string > split(const std::string &input, const std::string &regex)
Definition Rivet_i.cxx:253
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
const T * front() const
Access the first element in the collection as an rvalue.
This defines the McEventCollection, which is really just an ObjectVector of McEvent objectsFile: Gene...
bool m_skipweights
Skip variation weights and only run nominal.
Definition Rivet_i.h:128
double m_weightcap
Weight cap to set allowed maximum for weights.
Definition Rivet_i.h:140
std::vector< std::string > m_analysisNames
A list of names of the analyses to run (set from the job properties)
Definition Rivet_i.h:116
bool m_ignorebeams
Whether to avoid the beam consistency checks.
Definition Rivet_i.h:86
double m_crossSection
The cross section for this run of events, set from the job properties.
Definition Rivet_i.h:119
bool m_patchBeams
Flag to insert beam protons when they are unavailable in the event.
Definition Rivet_i.h:94
virtual StatusCode finalize() override
Finalise each analysis and commit the plots to an AIDA tree and the THistSvc ROOT tree.
Definition Rivet_i.cxx:210
std::string m_anapath
The analysis plugin search path.
Definition Rivet_i.h:80
bool m_init
Flag to determine whether Rivet init has already happened (in execute())
Definition Rivet_i.h:125
ToolHandle< IxAODtoHepMCTool > m_xAODtoHepMCTool
A tool to convert xAOD::Truth to HepMC::GenEvent.
Definition Rivet_i.h:97
std::string m_file
The base file name to write results to.
Definition Rivet_i.h:72
std::string m_unmatchWeights
String of weight names (or regex) to veto multiweights.
Definition Rivet_i.h:134
std::string m_matchWeights
String of weight names (or regex) to select multiweights.
Definition Rivet_i.h:131
virtual StatusCode execute() override
Run the Rivet analyses on one event, which is retrieved from StoreGate.
Definition Rivet_i.cxx:141
SG::ReadHandleKey< xAOD::EventInfo > m_evtInfoKey
Definition Rivet_i.h:142
double m_crossSection_uncert
The uncertainity of the cross section for this run of events, set from the job properties.
Definition Rivet_i.h:122
bool m_needsConversion
Do we need to convert xAOD::Truth back to HepMC::GenEvemt?
Definition Rivet_i.h:91
Rivet::AnalysisHandler * m_analysisHandler
A Rivet analysis handler.
Definition Rivet_i.h:113
bool m_doRootHistos
Will we convert Rivet's internal histo format into a ROOT histo for streaming with THistSvc?
Definition Rivet_i.h:104
std::string m_stream
A pointer to the THistSvc.
Definition Rivet_i.h:69
std::string m_preload
Definition Rivet_i.h:75
Rivet_i(const std::string &name, ISvcLocator *pSvcLocator)
Standard algorithm constructor.
Definition Rivet_i.cxx:35
std::unique_ptr< HepMC::GenEvent > checkEvent(const HepMC::GenEvent &event, const EventContext &ctx)
Book an AIDA::IDataPointSet into the THistSvc as a TH1D at path.
Definition Rivet_i.cxx:262
std::string m_genEventKey
The GenEvent StoreGate key (default "GEN_EVENT")
Definition Rivet_i.h:110
std::string m_nominalWeightName
String to specify non-standard nominal weight.
Definition Rivet_i.h:137
void MeV2GeV(HepMC::GenEvent &event)
Definition Rivet_i.cxx:397
std::string m_runname
The name of the run (prepended to plot paths).
Definition Rivet_i.h:107
const_pointer_type cptr()
Dereference the pointer.
GenParticle * GenParticlePtr
Definition GenParticle.h:37
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
bool isStable(const T &p)
Identify if the particle is stable, i.e. has not decayed.
bool isBeam(const T &p)
Identify if the particle is beam particle.
TruthEventContainer_v1 TruthEventContainer
Declare the latest version of the truth event container.
MsgStream & msg
Definition testRead.cxx:32