5#include "SimCore/Reweight/GenieReweightProducer.h"
10#include "Framework/EventGen/EventRecord.h"
11#include "Framework/EventGen/HepMC3Converter.h"
12#include "Framework/Interaction/Interaction.h"
13#include "Framework/Messenger/Messenger.h"
14#include "Framework/Utils/AppInit.h"
15#include "Framework/Utils/RunOpt.h"
16#include "HepMC3/GenEvent.h"
17#include "HepMC3/PrintStreams.h"
18#include "RwCalculators/GReWeightAGKY.h"
19#include "RwCalculators/GReWeightDISNuclMod.h"
20#include "RwCalculators/GReWeightDeltaradAngle.h"
21#include "RwCalculators/GReWeightFGM.h"
22#include "RwCalculators/GReWeightFZone.h"
23#include "RwCalculators/GReWeightINuke.h"
24#include "RwCalculators/GReWeightINukeParams.h"
25#include "RwCalculators/GReWeightNonResonanceBkg.h"
26#include "RwCalculators/GReWeightNuXSecCCQE.h"
27#include "RwCalculators/GReWeightNuXSecCCQEaxial.h"
28#include "RwCalculators/GReWeightNuXSecCCQEvec.h"
29#include "RwCalculators/GReWeightNuXSecCCRES.h"
30#include "RwCalculators/GReWeightNuXSecCOH.h"
31#include "RwCalculators/GReWeightNuXSecDIS.h"
32#include "RwCalculators/GReWeightNuXSecNC.h"
33#include "RwCalculators/GReWeightNuXSecNCEL.h"
34#include "RwCalculators/GReWeightNuXSecNCRES.h"
35#include "RwCalculators/GReWeightResonanceDecay.h"
36#include "RwCalculators/GReWeightXSecEmpiricalMEC.h"
37#include "RwCalculators/GReWeightXSecMEC.h"
38#include "RwFramework/GReWeight.h"
39#include "RwFramework/GReWeightI.h"
40#include "RwFramework/GSyst.h"
41#include "RwFramework/GSystSet.h"
42#include "SimCore/Event/EventWeights.h"
43#include "SimCore/Event/HepMC3GenEvent.h"
47GenieReweightProducer::GenieReweightProducer(
const std::string& name,
49 : Producer(name, process) {
50 hep_mc3_converter_ = std::make_unique<genie::HepMC3Converter>();
51 genie_rw_ = std::make_unique<genie::rew::GReWeight>();
55 hepmc3_coll_name_ = ps.
get<std::string>(
"hepmc3CollName");
56 hepmc3_pass_name_ = ps.
get<std::string>(
"hepmc3PassName");
57 event_weights_coll_name_ = ps.
get<std::string>(
"eventWeightsCollName");
58 seed_ = ps.
get<
int>(
"seed");
59 n_weights_ =
static_cast<size_t>(ps.
get<
int>(
"n_weights"));
60 auto var_types_strings = ps.
get<std::vector<std::string> >(
"var_types");
62 message_threshold_file_ = ps.
get<std::string>(
"message_threshold_file");
64 std::default_random_engine generator(seed_);
65 std::normal_distribution<double> normal_distribution;
67 for (
auto const& vt_str : var_types_strings) {
68 auto vtype = ldmx::EventWeights::stringToVariationType(vt_str);
69 for (
size_t i_w = 0; i_w < n_weights_; ++i_w)
70 variation_map_[vtype].push_back(normal_distribution(generator));
74void GenieReweightProducer::reinitializeGenieReweight() {
75 genie_rw_ = std::make_unique<genie::rew::GReWeight>();
78 genie::utils::app_init::MesgThresholds(message_threshold_file_);
80 genie::RunOpt::Instance()->SetTuneName(tune_);
81 if (!genie::RunOpt::Instance()->Tune()) {
82 EXCEPTION_RAISE(
"ConfigurationException",
"No TuneId in RunOption.");
84 genie::RunOpt::Instance()->BuildTune();
86 genie_rw_->AdoptWghtCalc(
"hadro_fzone",
new genie::rew::GReWeightFZone);
87 genie_rw_->AdoptWghtCalc(
"hadro_intranuke",
new genie::rew::GReWeightINuke);
89 auto& syst = genie_rw_->Systematics();
90 for (
auto var : variation_map_)
91 syst.Init(variationTypeToGenieDial(var.first));
95 const std::string genie_tune_par =
"GenieTune";
98 if (par.first.size() >= genie_tune_par.size() &&
99 par.first.compare(par.first.size() - genie_tune_par.size(),
100 genie_tune_par.size(), genie_tune_par) == 0) {
101 new_tune = par.second;
105 if (tune_ != new_tune) {
106 ldmx_log(debug) <<
"Found new tune " << new_tune <<
" (used to be " << tune_
109 reinitializeGenieReweight();
113void GenieReweightProducer::reconfigureGenieReweight(
size_t i_w) {
114 auto& syst = genie_rw_->Systematics();
115 for (
auto var : variation_map_) {
117 ldmx::EventWeights::VariationType::kGENIE_INukeTwkDial_MFP_pi)
118 syst.Set(genie::rew::GSyst::FromString(
"MFP_pi"), var.second[i_w]);
119 else if (var.first ==
120 ldmx::EventWeights::VariationType::kGENIE_INukeTwkDial_MFP_N)
121 syst.Set(genie::rew::GSyst::FromString(
"MFP_N"), var.second[i_w]);
122 else if (var.first ==
123 ldmx::EventWeights::VariationType::kGENIE_INukeTwkDial_FrCEx_pi)
124 syst.Set(genie::rew::GSyst::FromString(
"FrCEx_pi"), var.second[i_w]);
125 else if (var.first ==
126 ldmx::EventWeights::VariationType::kGENIE_INukeTwkDial_FrInel_pi)
127 syst.Set(genie::rew::GSyst::FromString(
"FrInel_pi"), var.second[i_w]);
128 else if (var.first ==
129 ldmx::EventWeights::VariationType::kGENIE_INukeTwkDial_FrAbs_pi)
130 syst.Set(genie::rew::GSyst::FromString(
"FrAbs_pi"), var.second[i_w]);
131 else if (var.first ==
132 ldmx::EventWeights::VariationType::kGENIE_INukeTwkDial_FrPiProd_pi)
133 syst.Set(genie::rew::GSyst::FromString(
"FrPiProd_pi"), var.second[i_w]);
134 else if (var.first ==
135 ldmx::EventWeights::VariationType::kGENIE_INukeTwkDial_FrCEx_N)
136 syst.Set(genie::rew::GSyst::FromString(
"FrCEx_N"), var.second[i_w]);
137 else if (var.first ==
138 ldmx::EventWeights::VariationType::kGENIE_INukeTwkDial_FrInel_N)
139 syst.Set(genie::rew::GSyst::FromString(
"FrInel_N"), var.second[i_w]);
140 else if (var.first ==
141 ldmx::EventWeights::VariationType::kGENIE_INukeTwkDial_FrAbs_N)
142 syst.Set(genie::rew::GSyst::FromString(
"FrAbs_N"), var.second[i_w]);
143 else if (var.first ==
144 ldmx::EventWeights::VariationType::kGENIE_INukeTwkDial_FrPiProd_N)
145 syst.Set(genie::rew::GSyst::FromString(
"FrPiProd_N"), var.second[i_w]);
146 else if (var.first ==
147 ldmx::EventWeights::VariationType::kGENIE_HadrNuclTwkDial_FormZone)
148 syst.Set(genie::rew::GSyst::FromString(
"FormZone"), var.second[i_w]);
150 genie_rw_->Reconfigure();
155 auto hepmc3_col =
event.getObject<std::vector<ldmx::HepMC3GenEvent> >(
156 hepmc3_coll_name_, hepmc3_pass_name_);
161 for (
size_t i_w = 0; i_w < n_weights_; ++i_w) {
162 double running_weight = 1;
164 reconfigureGenieReweight(i_w);
168 for (
size_t i_ev = 0; i_ev < 1; ++i_ev) {
169 auto const& hepmc3_ev = hepmc3_col.at(i_ev);
171 HepMC3::GenEvent hepmc3_genev;
172 hepmc3_genev.read_data(hepmc3_ev);
175 if (i_w == 0) ldmx_log(debug) << hepmc3_genev;
178 auto genie_ev_record_ptr = hep_mc3_converter_->RetrieveGHEP(hepmc3_genev);
181 if (i_w == 0) ldmx_log(debug) << *genie_ev_record_ptr;
184 auto this_weight = genie_rw_->CalcWeight(*genie_ev_record_ptr);
186 running_weight = running_weight * this_weight;
190 ev_weights.addWeight(running_weight);
194 ldmx_log(trace) << ev_weights;
196 event.add(event_weights_coll_name_, ev_weights);
#define DECLARE_PRODUCER(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
Implements an event buffer system for storing event data.
Class which represents the process under execution.
Class encapsulating parameters for configuring a processor.
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Dynamically loadable photonuclear models either from SimCore or external libraries implementing this ...