MARLEYHelper.cxx
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////
2 /// \file MARLEYHelper.cxx
3 /// \brief LArSoft interface to the MARLEY (Model of Argon Reaction Low Energy
4 /// Yields) supernova neutrino event generator
5 ///
6 /// \author Steven Gardiner <sjgardiner@ucdavis.edu>
7 //////////////////////////////////////////////////////////////////////////////
8 
9 // framework includes
10 #include "cetlib_except/exception.h"
12 
13 // LArSoft includes
16 #include "nurandom/RandomUtils/NuRandomService.h"
17 
18 // ROOT includes
19 #include "TInterpreter.h"
20 #include "TROOT.h"
21 
22 // MARLEY includes
23 #include "marley/Event.hh"
24 #include "marley/RootJSONConfig.hh"
25 
26 namespace {
27  // We need to convert from MARLEY's energy units (MeV) to LArSoft's
28  // (GeV) using this conversion factor
29  constexpr double MeV_to_GeV = 1e-3;
30 }
31 
32 //------------------------------------------------------------------------------
34  rndm::NuRandomService& rand_service, const std::string& helper_name )
35  : fHelperName( helper_name )
36 {
37  // Configure MARLEY using the FHiCL parameters
38  this->reconfigure( pset );
39 
40  // Register this MARLEY generator with the NuRandomService. For simplicity,
41  // we use a lambda as the seeder function (see NuRandomService.h for
42  // details). This allows the SeedService to automatically re-seed MARLEY
43  // whenever necessary. The user can set an explicit seed for MARLEY in the
44  // FHiCL configuration using the "seed" parameter. If you need to get the
45  // seed for MARLEY from the SeedService, note that we're using use the value
46  // of the input variable helper_name as its generator instance name.
47  rndm::NuRandomService::seed_t marley_seed = rand_service.registerEngine(
48  [this](rndm::NuRandomService::EngineId const& /* unused */,
49  rndm::NuRandomService::seed_t lar_seed) -> void
50  {
51  if ( fMarleyGenerator && fMarleyGenerator.get() ) {
52  auto seed = static_cast<uint_fast64_t>( lar_seed );
53  fMarleyGenerator->reseed( seed );
54  }
55  },
56  fHelperName, pset, { "seed" }
57  );
58 
59  // Unless I'm mistaken, the call to registerEngine should seed the generator
60  // with the seed from the FHiCL configuration file if one is included, but it
61  // doesn't appear to do so (as of 16 Aug 2016, larsoft v06_03_00). As a
62  // workaround, I manually reseed the generator (if needed) here using the
63  // result of the call to registerEngine, which will be the seed from the
64  // FHiCL file if one was given.
65  // TODO: figure out what's going on here, and remove this workaround as
66  // needed
67  uint_fast64_t marley_cast_seed = static_cast<uint_fast64_t>( marley_seed );
68  if ( marley_cast_seed != fMarleyGenerator->get_seed() ) {
69  fMarleyGenerator->reseed( marley_cast_seed );
70  }
71 
72  // Log initialization information from the MARLEY generator
74  fMarleyLogStream = std::stringstream();
75 
76  // Do any needed setup of the MARLEY class dictionaries
78 }
79 
80 //------------------------------------------------------------------------------
82  const std::vector<marley::Particle*>& particles,
83  const TLorentzVector& vtx_pos, bool track )
84 {
85  // Loop over the vector of MARLEY particles and add simb::MCParticle
86  // versions of each of them to the MCTruth object.
87  for ( const marley::Particle* p : particles ) {
88  // Treat all of these particles as primaries, which have negative
89  // track IDs by convention
90  int trackID = -1 * ( truth.NParticles() + 1 );
91 
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 );
99 
100  int status = 0; // don't track the particles in LArG4 by default
101  if ( track ) status = 1;
102 
103  simb::MCParticle part( trackID /* trackID to use in Geant4 */, pdg,
104  "MARLEY", -1 /* primary particle */, mass, status );
105 
106  part.AddTrajectoryPoint( vtx_pos, mom );
107  truth.Add( part );
108  }
109 }
110 
111 //------------------------------------------------------------------------------
113  const TLorentzVector& vtx_pos, marley::Event* marley_event )
114 {
115  simb::MCTruth truth;
116 
118 
119  marley::Event event = fMarleyGenerator->create_event();
120 
121  // Add the initial and final state particles to the MCTruth object.
122  add_marley_particles( truth, event.get_initial_particles(), vtx_pos, false );
123  add_marley_particles( truth, event.get_final_particles(), vtx_pos, true );
124 
125  // calculate a few parameters for the call to SetNeutrino
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);
134 
135  // For definitions of Bjorken x, etc., a good reference is Mark Thomson's
136  // set of slides on deep inelastic scattering (http://tinyurl.com/hcn5n6l)
137  double bjorken_x = Q2 / (2 * event.target().mass() * (Enu - Elep));
138  double inelasticity_y = 1. - Elep / Enu;
139 
140  // Include the initial excitation energy of the final-state nucleus when
141  // calculating W (the final-state invariant mass of the hadronic system)
142  // since the other parameters (X, Y) also take into account the 2-2
143  // scattering reaction only.
144  const marley::Particle& res = event.residue();
145  double hadronic_mass_W = res.mass() + event.Ex();
146 
147  // TODO: do a more careful job of setting the parameters here
148  truth.SetNeutrino(
149  simb::kCC, // change when MARLEY can handle NC
150  simb::kUnknownInteraction, // not sure what the mode should be
151  simb::kUnknownInteraction, // not sure what the interaction type should be
152  marley_utils::get_nucleus_pid(18, 40), // Ar-40 PDG code
153  marley_utils::NEUTRON, // nucleon PDG
154  0, // MARLEY handles low enough energies that we shouldn't need HitQuark
155  hadronic_mass_W * MeV_to_GeV,
156  bjorken_x, // dimensionless
157  inelasticity_y, // dimensionless
158  Q2 * std::pow(MeV_to_GeV, 2)
159  );
160 
161  if ( marley_event ) *marley_event = event;
162 
163  // Process the MARLEY logging messages (if any) captured by our
164  // stringstream and forward them to the messagefacility logger
166  while( std::getline(fMarleyLogStream, line) ) {
168  }
169 
170  // Reset the MARLEY log stream
171  fMarleyLogStream = std::stringstream();
172 
173  return truth;
174 }
175 
176 //------------------------------------------------------------------------------
178  const std::string& fileType )
179 {
180  cet::search_path searchPath( "FW_SEARCH_PATH" );
181 
182  std::string fullName;
183  searchPath.find_file( fileName, fullName );
184 
185  if ( fullName.empty() )
186  throw cet::exception( "MARLEYHelper" )
187  << "Cannot find MARLEY " << fileType << " data file '"
188  << fileName << '\'';
189 
190  return fullName;
191 }
192 
193 //------------------------------------------------------------------------------
195  marley::JSON& json, const std::string& key, bool missing_ok )
196 {
197  if ( json.has_key(key) ) {
198 
199  marley::JSON& value = json.at(key);
200 
201  if ( value.is_array() ) {
202  // Replace each file name (which may appear in the FHiCL configuration
203  // without a full path) with the full path found using cetlib
204  for ( auto& element : value.array_range() ) {
205  element = find_file( element.to_string(), key );
206  }
207  }
208 
209  else value = find_file(value.to_string(), key);
210  }
211  else if ( !missing_ok ) throw cet::exception("MARLEYHelper")
212  << "Missing \"" << key << "\" key in the MARLEY parameters.";
213 }
214 
215 //------------------------------------------------------------------------------
217 {
218  // Convert the FHiCL parameters into a JSON object that MARLEY can understand
220  pset.walk( mpsw );
221 
222  marley::JSON& json = mpsw.get_json();
223 
224  // Update the reaction and structure data file names to the full paths
225  // using cetlib to search for them
226  load_full_paths_into_json( json, "reactions", false );
227  load_full_paths_into_json( json, "structure", true );
228 
229  // Also update the path for a neutrino source spectrum given in a ROOT
230  // TFile
231  if ( json.has_key("source") ) {
232  marley::JSON& source_object = json.at( "source" );
233 
234  if ( source_object.has_key("tfile") ) {
235  load_full_paths_into_json( source_object, "tfile" );
236  }
237  }
238 
239  // Create a new MARLEY configuration based on the JSON parameters
240  MF_LOG_INFO( "MARLEYHelper " + fHelperName ) << "MARLEY will now use"
241  " the JSON configuration\n" << json.dump_string() << '\n';
242  marley::RootJSONConfig config( json );
243 
244  // Create a new marley::Generator object based on the current configuration
245  fMarleyGenerator = std::make_unique<marley::Generator>(
246  config.create_generator() );
247 }
248 
249 //------------------------------------------------------------------------------
251 {
252  static bool already_loaded_marley_dict = false;
253 
254  if ( already_loaded_marley_dict ) return;
255 
256  // Current (24 July 2016) versions of ROOT 6 require runtime
257  // loading of headers for custom classes in order to use
258  // dictionaries correctly. If we're running ROOT 6+, do the
259  // loading here, and give the user guidance if there are any
260  // problems.
261  //
262  // This is the same technique used in the MARLEY source code
263  // for the executable (src/marley.cc). If you change how this
264  // code works, please sync changes with the executable as well.
265  if ( gROOT->GetVersionInt() >= 60000 ) {
266  MF_LOG_INFO( "MARLEYHelper " + fHelperName ) << "ROOT 6 or greater"
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 );
271  if ( *ec != 0 ) throw cet::exception( "MARLEYHelper " + fHelperName )
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"
275  << " try again.";
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"
281  << " try again.";
282  }
283 
284  // No further action is required for ROOT 5 because the compiled
285  // dictionaries (which are linked to this algorithm) contain all of
286  // the needed information
287  already_loaded_marley_dict = true;
288 }
LArSoft interface to the MARLEY (Model of Argon Reaction Low Energy Yields) supernova neutrino event ...
std::unique_ptr< marley::Generator > fMarleyGenerator
Definition: MARLEYHelper.h:75
void AddTrajectoryPoint(TLorentzVector const &position, TLorentzVector const &momentum)
Definition: MCParticle.h:257
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)
Definition: KineUtils.cxx:1064
void SetOrigin(simb::Origin_t origin)
Definition: MCTruth.h:82
MARLEYHelper(const fhicl::ParameterSet &pset, rndm::NuRandomService &rand_service, const std::string &generator_name)
std::string string
Definition: nybbler.cc:12
simb::MCTruth create_MCTruth(const TLorentzVector &vtx_pos, marley::Event *marley_event=nullptr)
constexpr T pow(T x)
Definition: pow.h:72
Concrete fhicl::ParameterSetWalker that converts a fhicl::ParameterSet into a marley::JSON object...
int NParticles() const
Definition: MCTruth.h:75
void SetNeutrino(int CCNC, int mode, int interactionType, int target, int nucleon, int quark, double w, double x, double y, double qsqr)
Definition: MCTruth.cxx:30
const marley::JSON & get_json() const
std::stringstream fMarleyLogStream
Definition: MARLEYHelper.h:82
std::string fHelperName
Definition: MARLEYHelper.h:78
std::string find_file(const std::string &fileName, const std::string &fileType)
const double e
fileName
Definition: dumpTree.py:9
def key(type, name=None)
Definition: graph.py:13
static Config * config
Definition: config.cpp:1054
void reconfigure(const fhicl::ParameterSet &pset)
p
Definition: test.py:223
void walk(ParameterSetWalker &psw) const
#define MF_LOG_INFO(category)
void Add(simb::MCParticle const &part)
Definition: MCTruth.h:80
n.b.: this group is similar but not quite, entirely unlike GENIE ScatteringType convention ...
Definition: MCNeutrino.h:80
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)
Supernova neutrinos.
Definition: MCTruth.h:25
std::string find_file(std::string const &filename) const
Definition: search_path.cc:96
Event generator information.
Definition: MCTruth.h:32
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
Event finding and building.