ATLAS Offline Software
Loading...
Searching...
No Matches
Rivet_i.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef RIVET_I_H
6#define RIVET_I_H
7
11#include "GaudiKernel/ToolHandle.h"
13
14#include "Rivet/AnalysisHandler.hh"
15
16#include "AtlasHepMC/GenEvent.h"
17
18#include <vector>
19#include <string>
20#include <memory>
21
22class ISvcLocator;
24//class ITHistSvc;
25
26
31class Rivet_i : public AthAlgorithm {
32public:
33
35 Rivet_i(const std::string &name, ISvcLocator *pSvcLocator);
36
41 virtual StatusCode initialize ATLAS_NOT_THREAD_SAFE () override;
42
44 virtual StatusCode execute() override;
45
48 virtual StatusCode finalize() override;
49
50
51private:
52
54// StatusCode regHist(const AIDA::IDataPointSet& dps, const std::string& path);
55
57// StatusCode regGraph(const AIDA::IDataPointSet& dps, const std::string& path);
58
59 // Check and potentially modify events for correct units, beam particles, ...
60 std::unique_ptr<HepMC::GenEvent> checkEvent(const HepMC::GenEvent& event, const EventContext& ctx);
61
62 // Utility method to convert units of the event
63 void MeV2GeV(HepMC::GenEvent& event);
64
66 //ServiceHandle<ITHistSvc> m_histSvc;
67
69 std::string m_stream;
70
72 std::string m_file;
73
74 //specify a pre-existing yoda file to initialize from
75 std::string m_preload;
76
80 std::string m_anapath;
81
87
92
95
97 ToolHandle<IxAODtoHepMCTool> m_xAODtoHepMCTool{this, "HepMCTool", "xAODtoHepMCTool"};
98
105
107 std::string m_runname;
108
110 std::string m_genEventKey;
111
113 Rivet::AnalysisHandler* m_analysisHandler;
114
116 std::vector<std::string> m_analysisNames;
117
120
123
125 bool m_init;
126
129
131 std::string m_matchWeights;
132
134 std::string m_unmatchWeights;
135
138
141
143 , "EventInfo"
144 , "EventInfo"
145 , "ReadHandleKey for xAOD::EventInfo" };
146};
147
148#endif
Define macros for attributes used to control the static checker.
AthAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
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
virtual StatusCode initialize ATLAS_NOT_THREAD_SAFE() override
Initialise the Rivet interface and Athena services.
std::string m_runname
The name of the run (prepended to plot paths).
Definition Rivet_i.h:107
Property holding a SG store/key/clid from which a ReadHandle is made.
void initialize()