10 #include "cetlib_except/exception.h" 16 #include "nurandom/RandomUtils/NuRandomService.h" 19 #include "TInterpreter.h" 23 #include "marley/Event.hh" 24 #include "marley/RootJSONConfig.hh" 29 constexpr
double MeV_to_GeV = 1
e-3;
34 rndm::NuRandomService& rand_service,
const std::string& helper_name )
35 : fHelperName( helper_name )
47 rndm::NuRandomService::seed_t marley_seed = rand_service.registerEngine(
48 [
this](rndm::NuRandomService::EngineId
const& ,
49 rndm::NuRandomService::seed_t lar_seed) ->
void 52 auto seed =
static_cast<uint_fast64_t
>( lar_seed );
67 uint_fast64_t marley_cast_seed =
static_cast<uint_fast64_t
>( marley_seed );
82 const std::vector<marley::Particle*>& particles,
83 const TLorentzVector& vtx_pos,
bool track )
87 for (
const marley::Particle*
p : particles ) {
92 int pdg =
p->pdg_code();
93 double mass =
p->mass() * MeV_to_GeV;
94 double px =
p->px() * MeV_to_GeV;
95 double py =
p->py() * MeV_to_GeV;
96 double pz =
p->pz() * MeV_to_GeV;
97 double E =
p->total_energy() * MeV_to_GeV;
98 TLorentzVector mom( px, py, pz, E );
101 if ( track ) status = 1;
104 "MARLEY", -1 , mass, status );
113 const TLorentzVector& vtx_pos, marley::Event* marley_event )
126 const marley::Particle& nu =
event.projectile();
127 const marley::Particle& lep =
event.ejectile();
128 double qx = nu.px() - lep.px();
129 double qy = nu.py() - lep.py();
130 double qz = nu.pz() - lep.pz();
131 double Enu = nu.total_energy();
132 double Elep = lep.total_energy();
133 double Q2 = qx*qx + qy*qy + qz*qz -
std::pow(Enu - Elep, 2);
137 double bjorken_x = Q2 / (2 *
event.target().mass() * (Enu - Elep));
138 double inelasticity_y = 1. - Elep / Enu;
144 const marley::Particle& res =
event.residue();
145 double hadronic_mass_W = res.mass() +
event.Ex();
152 marley_utils::get_nucleus_pid(18, 40),
153 marley_utils::NEUTRON,
155 hadronic_mass_W * MeV_to_GeV,
158 Q2 * std::pow(MeV_to_GeV, 2)
161 if ( marley_event ) *marley_event =
event;
183 searchPath.
find_file( fileName, fullName );
185 if ( fullName.empty() )
187 <<
"Cannot find MARLEY " << fileType <<
" data file '" 197 if ( json.has_key(key) ) {
199 marley::JSON&
value = json.at(key);
201 if ( value.is_array() ) {
204 for (
auto& element : value.array_range() ) {
212 <<
"Missing \"" << key <<
"\" key in the MARLEY parameters.";
222 marley::JSON& json = mpsw.
get_json();
231 if ( json.has_key(
"source") ) {
232 marley::JSON& source_object = json.at(
"source" );
234 if ( source_object.has_key(
"tfile") ) {
241 " the JSON configuration\n" << json.dump_string() <<
'\n';
242 marley::RootJSONConfig
config( json );
246 config.create_generator() );
252 static bool already_loaded_marley_dict =
false;
254 if ( already_loaded_marley_dict )
return;
265 if ( gROOT->GetVersionInt() >= 60000 ) {
267 <<
" detected. Loading class information\nfrom headers" 268 <<
" \"marley/Particle.hh\" and \"marley/Event.hh\"";
269 TInterpreter::EErrorCode* ec =
new TInterpreter::EErrorCode();
270 gInterpreter->ProcessLine(
"#include \"marley/Particle.hh\"", ec );
272 <<
"Error loading MARLEY header Particle.hh. For MARLEY headers stored" 273 <<
" in /path/to/include/marley/, please add /path/to/include" 274 <<
" to your ROOT_INCLUDE_PATH environment variable and" 276 gInterpreter->ProcessLine(
"#include \"marley/Event.hh\"" );
277 if ( *ec != 0 )
throw cet::exception(
"MARLEYHelper" ) <<
"Error loading" 278 <<
" MARLEY header Event.hh. For MARLEY headers stored in" 279 <<
" /path/to/include/marley/, please add /path/to/include" 280 <<
" to your ROOT_INCLUDE_PATH environment variable and" 287 already_loaded_marley_dict =
true;
LArSoft interface to the MARLEY (Model of Argon Reaction Low Energy Yields) supernova neutrino event ...
void load_marley_dictionaries()
std::unique_ptr< marley::Generator > fMarleyGenerator
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
void add_marley_particles(simb::MCTruth &truth, const std::vector< marley::Particle * > &particles, const TLorentzVector &vtx_pos, bool track)
double Q2(const Interaction *const i)
void SetOrigin(simb::Origin_t origin)
MARLEYHelper(const fhicl::ParameterSet &pset, rndm::NuRandomService &rand_service, const std::string &generator_name)
simb::MCTruth create_MCTruth(const TLorentzVector &vtx_pos, marley::Event *marley_event=nullptr)
Concrete fhicl::ParameterSetWalker that converts a fhicl::ParameterSet into a marley::JSON object...
void SetNeutrino(int CCNC, int mode, int interactionType, int target, int nucleon, int quark, double w, double x, double y, double qsqr)
const marley::JSON & get_json() const
std::stringstream fMarleyLogStream
std::string find_file(const std::string &fileName, const std::string &fileType)
void reconfigure(const fhicl::ParameterSet &pset)
void walk(ParameterSetWalker &psw) const
#define MF_LOG_INFO(category)
void Add(simb::MCParticle const &part)
n.b.: this group is similar but not quite, entirely unlike GENIE ScatteringType convention ...
void load_full_paths_into_json(marley::JSON &json, const std::string &array_name, bool missing_ok=false)
void line(double t, double *p, double &x, double &y, double &z)
std::string find_file(std::string const &filename) const
Event generator information.
cet::coded_exception< error, detail::translate > exception
Event finding and building.