All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MARLEYTimeGen_module.cc
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////
2 /// @file MARLEYTimeGen_module.cc
3 /// @brief Module that allows for sampling neutrino energies and times from
4 /// time-dependent supernova spectra. This module uses MARLEY to help generate
5 /// events.
6 ///
7 /// @author Steven Gardiner <sjgardiner@ucdavis.edu>
8 //////////////////////////////////////////////////////////////////////////////
9 
10 // standard library includes
11 #include <fstream>
12 #include <limits>
13 #include <memory>
14 #include <regex>
15 #include <string>
16 
17 // framework includes
18 #include "art/Framework/Core/EDProducer.h"
19 #include "art/Framework/Core/ModuleMacros.h"
20 #include "art/Framework/Principal/Event.h"
21 #include "art/Framework/Principal/Run.h"
22 #include "art/Framework/Services/Registry/ServiceHandle.h"
23 #include "art_root_io/TFileService.h"
24 #include "art_root_io/TFileDirectory.h"
25 #include "canvas/Persistency/Common/Assns.h"
26 #include "cetlib_except/exception.h"
27 #include "fhiclcpp/ParameterSet.h"
28 #include "fhiclcpp/types/OptionalAtom.h"
29 #include "fhiclcpp/types/Table.h"
30 #include "messagefacility/MessageLogger/MessageLogger.h"
31 #include "art/Persistency/Common/PtrMaker.h"
32 
33 // art extensions
34 #include "nurandom/RandomUtils/NuRandomService.h"
35 
36 // LArSoft includes
42 #include "nusimdata/SimulationBase/MCTruth.h"
43 
44 // ROOT includes
45 #include "TFile.h"
46 #include "TH1D.h"
47 #include "TH2D.h"
48 #include "TTree.h"
49 
50 // MARLEY includes
51 #include "marley/marley_root.hh"
52 #include "marley/marley_utils.hh"
53 #include "marley/Integrator.hh"
54 
55 // Anonymous namespace for definitions local to this source file
56 namespace {
57 
58  // The number of times to attempt to sample an energy uniformly
59  // before giving up
60  constexpr int MAX_UNIFORM_ENERGY_ITERATIONS = 1000;
61 
62  // Neutrino vertices generated using unbiased sampling are assigned
63  // unit weight
64  constexpr double ONE = 1.;
65 
66  // Number of sampling points to use for numerical integration
67  // (via the Clenshaw-Curtis method)
68  constexpr int NUM_INTEGRATION_SAMPLING_POINTS = 100;
69 
70  // Multiply an energy in MeV by this factor to convert to erg
71  // (based on 2017 PDG "Physical Constants" article)
72  constexpr double MeV_to_erg = 1.6021766208e-6; // erg / MeV
73 
74  // The PDG codes that represent one of the varieties of "nu_x"
75  constexpr std::array<int, 4> NU_X_PDG_CODES
76  {
77  marley_utils::MUON_NEUTRINO, marley_utils::MUON_ANTINEUTRINO,
78  marley_utils::TAU_NEUTRINO, marley_utils::TAU_ANTINEUTRINO,
79  };
80 
81  // Helper function that tests whether a given PDG code represents
82  // a "nu_x"
83  bool is_nux(int pdg_code) {
84  return std::find(std::begin(NU_X_PDG_CODES),
85  std::end(NU_X_PDG_CODES), pdg_code) != std::end(NU_X_PDG_CODES);
86  }
87 
88  // Matches comment lines and empty lines in a "fit"-format spectrum file
89  const std::regex rx_comment_or_empty = std::regex("\\s*(#.*)?");
90 
91  // Compute the flux-averaged total cross section (MeV^[-2]) for all
92  // defined reactions and for a particular neutrino source
93  double flux_averaged_total_xs(marley::NeutrinoSource& nu_source,
94  marley::Generator& gen);
95 
96  // Overloaded version that allows access to the flux integral
97  // (which will be loaded into source_integ) and the xs * flux integral
98  // (which will be loaded into tot_xs_integ)
99  double flux_averaged_total_xs(marley::NeutrinoSource& nu_source,
100  marley::Generator& gen, double& source_integ, double& tot_xs_integ);
101 
102  // Returns a numerical integral of the function f on the interval
103  // [x_min, x_max]
104  double integrate(const std::function<double(double)>& f, double x_min,
105  double x_max);
106 }
107 
108 namespace evgen {
109  class MarleyTimeGen;
110 }
111 
112 class evgen::MarleyTimeGen : public art::EDProducer {
113 
114  public:
115 
116  using Name = fhicl::Name;
117  using Comment = fhicl::Comment;
118 
119  /// Ignore the marley_parameters FHiCL table during the validation
120  /// step (MARLEY will take care of that by itself)
121  struct KeysToIgnore {
122  std::set< std::string > operator()()
123  {
124  return { "marley_parameters" };
125  }
126  };
127 
128  /// Collection of configuration parameters for the module
129  struct Config {
130 
131  fhicl::Table<evgen::ActiveVolumeVertexSampler::Config> vertex_ {
132  Name("vertex"),
133  Comment("Configuration for selecting the vertex location(s)")
134  };
135 
136  fhicl::Atom<std::string> module_type_ {
137  Name("module_type"),
138  Comment(""),
139  "MARLEYTimeGen" // default value
140  };
141 
142  fhicl::Atom<std::string> sampling_mode_ {
143  Name("sampling_mode"),
144  Comment("Technique to use when sampling times and energies. Valid"
145  " options are \"histogram\", \"uniform time\","
146  " and \"uniform energy\""),
147  std::string("histogram") // default value
148  };
149 
150  fhicl::Atom<unsigned int> nu_per_event_ {
151  Name("nu_per_event"),
152  Comment("The number of neutrino vertices to generate in"
153  " each art::Event"),
154  1u // default value
155  };
156 
157  fhicl::Atom<std::string> spectrum_file_format_ {
158  Name("spectrum_file_format"),
159  Comment("Format to assume for the neutrino spectrum file."
160  " Valid options are \"th2d\" (a ROOT file containing a "
161  " TH2D object) and \"fit\" (an ASCII text file containing"
162  " fit parameters for each time bin). This parameter is not"
163  " case-sensitive."),
164  "th2d" // default value
165  };
166 
167  fhicl::Atom<std::string> spectrum_file_ {
168  Name("spectrum_file"),
169  Comment("Name of a file that contains a representation"
170  " of the incident (not cross section weighted) neutrino spectrum"
171  " as a function of time and energy.")
172  };
173 
174  fhicl::OptionalAtom<std::string> pinching_parameter_type_ {
175  Name("pinching_parameter_type"),
176  Comment("Type of pinching parameter to assume when parsing"
177  " the time-dependent fit parameters for the incident neutrino"
178  " spectrum. Valid options are \"alpha\" and \"beta\". This"
179  " parameter is not case-sensitive."),
180  [this]() -> bool {
181  auto spectrum_file_format
182  = marley_utils::to_lowercase(spectrum_file_format_());
183  return (spectrum_file_format == "fit");
184  }
185  };
186 
187  fhicl::OptionalAtom<std::string> namecycle_ {
188  Name("namecycle"),
189  Comment("Name of the TH2D object to use to represent the"
190  " incident neutrino spectrum. This value should match the"
191  " name of the TH2D as given in the ROOT file specified"
192  " in the \"spectrum_file\" parameter. The TH2D should use "
193  " time bins on the X axis (seconds) and energy bins on the "
194  " Y axis (MeV)."),
195  [this]() -> bool {
196  auto spectrum_file_format
197  = marley_utils::to_lowercase(spectrum_file_format_());
198  return (spectrum_file_format == "th2d");
199  }
200  };
201 
202  fhicl::OptionalAtom<double> fit_Emin_ {
203  Name("fit_Emin"),
204  Comment("Minimum allowed neutrino energy (MeV) for a \"fit\" format"
205  " spectrum file"),
206  [this]() -> bool {
207  auto spectrum_file_format
208  = marley_utils::to_lowercase(spectrum_file_format_());
209  return (spectrum_file_format == "fit");
210  }
211  };
212 
213  fhicl::OptionalAtom<double> fit_Emax_ {
214  Name("fit_Emax"),
215  Comment("Maximum allowed neutrino energy (MeV) for a \"fit\" format"
216  " spectrum file"),
217  [this]() -> bool {
218  auto spectrum_file_format
219  = marley_utils::to_lowercase(spectrum_file_format_());
220  return (spectrum_file_format == "fit");
221  }
222  };
223 
224  }; // struct Config
225 
226  // Type to enable FHiCL parameter validation by art
227  using Parameters = art::EDProducer::Table<Config, KeysToIgnore>;
228 
229  // @brief Configuration-checking constructor
230  explicit MarleyTimeGen(const Parameters& p);
231 
232  virtual void produce(art::Event& e) override;
233  virtual void beginRun(art::Run& run) override;
234 
235  virtual void reconfigure(const Parameters& p);
236 
237  protected:
238 
239  /// @brief Stores parsed fit parameters from a single time bin and neutrino
240  /// type in a "fit"-format spectrum file
242  public:
243 
244  FitParameters(double Emean, double alpha, double luminosity)
245  : fEmean(Emean), fAlpha(alpha), fLuminosity(luminosity) {}
246 
247  /// @brief Mean neutrino energy (MeV)
248  double Emean() const { return fEmean; }
249  /// @brief Pinching parameter (dimensionless)
250  double Alpha() const { return fAlpha; }
251  /// @brief Luminosity (erg / s)
252  double Luminosity() const { return fLuminosity; }
253 
254  /// @brief Set the mean neutrino energy (MeV)
255  void set_Emean(double Emean) { fEmean = Emean; }
256  /// @brief Set the pinching parameter
257  void set_Alpha(double alpha) { fAlpha = alpha; }
258  /// @brief Set the luminosity (erg / s)
259  void set_Luminosity(double lum) { fLuminosity = lum; }
260 
261  /// @brief Converts an iterator that points to a FitParameters
262  /// object into an iterator that points to that object's
263  /// fLuminosity member.
264  /// @details This function helps us to be able to sample time bins
265  /// with a std::discrete_distribution using the bin luminosities
266  /// without redundnant storage.
267  template<typename It> static marley::IteratorToMember<It, double>
269  {
270  return marley::IteratorToMember<It, double>(
272  }
273 
274  protected:
275  double fEmean; ///< Mean neutrino energy (MeV)
276  double fAlpha; ///< Pinching parameter
277  double fLuminosity; ///< Luminosity (erg / s)
278  };
279 
280  /// @brief Stores fitting parameters for all neutrino types from a single
281  /// time bin in a "fit"-format spectrum file
282  class TimeFit {
283  public:
284  TimeFit(double time, double Emean_nue, double alpha_nue,
285  double lum_nue, double Emean_nuebar, double alpha_nuebar,
286  double lum_nuebar, double Emean_nux, double alpha_nux,
287  double lum_nux) : fTime(time),
288  fNueFitParams(Emean_nue, alpha_nue, lum_nue),
289  fNuebarFitParams(Emean_nuebar, alpha_nuebar, lum_nuebar),
290  fNuxFitParams(Emean_nux, alpha_nux, lum_nux) {}
291 
292  /// @brief Returns the low edge (in s) of the time bin that uses
293  /// these fit parameters
294  /// @details Time zero is defined to be the core bounce
295  double Time() const { return fTime; }
296 
297  protected:
298  /// @brief Helper function that returns a pointer-to-member
299  /// for the FitParameters object appropriate for a given
300  /// neutrino type
301  /// @param pdg_code PDG code for the neutrino type of interest
303  int pdg_code)
304  {
305  if (pdg_code == marley_utils::ELECTRON_NEUTRINO) {
306  return &TimeFit::fNueFitParams;
307  }
308  else if (pdg_code == marley_utils::ELECTRON_ANTINEUTRINO) {
310  }
311  else if ( is_nux(pdg_code) )
312  {
313  // The PDG code represents one of the varieties of nu_x
314  return &TimeFit::fNuxFitParams;
315  }
316  else throw cet::exception("MARLEYTimeGen") << "Invalid neutrino"
317  << " PDG code " << pdg_code << " encountered in MARLEYTimeGen"
318  << "::TimeFit::GetFitParametersMemberPointer()";
319  return nullptr;
320  }
321 
322  public:
323 
324  /// @brief Retrieves fit parameters for a specific neutrino type
325  /// for this time bin
326  /// @param pdg_code The PDG code for the desired neutrino type
327  const FitParameters& GetFitParameters(int pdg_code) const {
328  return this->*GetFitParametersMemberPointer(pdg_code);
329  }
330 
331  /// @brief Converts an iterator that points to a TimeFit
332  /// object into an iterator that points to the appropriate
333  /// (based on the neutrino PDG code) FitParameters member
334  /// owned by the TimeFit object.
335  /// @details This function helps us to be able to sample time bins
336  /// with a std::discrete_distribution using the bin luminosities
337  /// without redundnant storage.
338  /// @param pdg_code PDG code for the neutrino type of interest
339  /// @param iterator An iterator to a TimeFit object that will be
340  /// converted
341  template<typename It> static marley::IteratorToMember<It,
343  It iterator)
344  {
345  return marley::IteratorToMember<It, FitParameters>(
346  iterator, GetFitParametersMemberPointer(pdg_code) );
347  }
348 
349  protected:
350  /// @brief Time bin low edge (s)
351  double fTime;
352 
353  /// @brief Fitting parameters for electron neutrinos in this time bin
355 
356  /// @brief Fitting parameters for electron antineutrinos in this time
357  /// bin
359 
360  /// @brief Fitting parameters for non-electron-flavor neutrinos in this
361  /// time bin
363  };
364 
365  /// @brief Creates a simb::MCTruth object using a uniformly-sampled
366  /// neutrino energy
367  /// @details This function samples a reacting neutrino energy uniformly
368  /// from the full range of energies allowed by the incident spectrum and
369  /// the currently defined reactions.
370  simb::MCTruth make_uniform_energy_mctruth(double E_min, double E_max,
371  double& E_nu, const TLorentzVector& vertex_pos);
372 
373  /// @brief Create a MARLEY neutrino source object using a set of fit
374  /// parameters for a particular time bin
375  std::unique_ptr<marley::NeutrinoSource> source_from_time_fit(
376  const TimeFit& fit);
377 
378  /// @brief Create simb::MCTruth and sim::SupernovaTruth objects using
379  /// spectrum information from a ROOT TH2D
380  void create_truths_th2d(simb::MCTruth& mc_truth,
381  sim::SupernovaTruth& sn_truth, const TLorentzVector& vertex_pos);
382 
383  /// @brief Create simb::MCTruth and sim::SupernovaTruth objects using a
384  /// neutrino spectrum described by a previously-parsed "fit"-format file
385  void create_truths_time_fit(simb::MCTruth& mc_truth,
386  sim::SupernovaTruth& sn_truth, const TLorentzVector& vertex_pos);
387 
388  /// @brief Helper function that makes a final dummy TimeFit object so that
389  /// the final real time bin can have a right edge
390  void make_final_timefit(double time);
391 
392  /// @brief Makes ROOT histograms showing the emitted neutrinos in each time
393  /// bin when using a "fit"-format spectrum file
394  void make_nu_emission_histograms() const;
395 
396  /// @brief Object that provides an interface to the MARLEY event generator
397  std::unique_ptr<evgen::MARLEYHelper> fMarleyHelper;
398 
399  /// @brief Algorithm that allows us to sample vertex locations within the
400  /// active volume(s) of the detector
401  std::unique_ptr<evgen::ActiveVolumeVertexSampler> fVertexSampler;
402 
403  /// @brief unique_ptr to the current event created by MARLEY
404  std::unique_ptr<marley::Event> fEvent;
405 
406  /// @brief ROOT TH2D that contains the time-dependent spectrum to use when
407  /// sampling neutrino times and energies.
408  /// @details This member is only used when reading the spectrum from a ROOT
409  /// file.
410  std::unique_ptr<TH2D> fSpectrumHist;
411 
412  /// @brief Vector that contains the fit parameter information for each time
413  /// bin when using a "fit"-format spectrum file.
414  /// @details This member is unused when the spectrum is read from a ROOT
415  /// file.
416  std::vector<TimeFit> fTimeFits;
417 
418  /// @enum TimeGenSamplingMode
419  /// @brief Enumerated type that defines the allowed ways that a neutrino's
420  /// energy and arrival time may be sampled
421  ///
422  /// @var TimeGenSamplingMode::HISTOGRAM
423  /// @brief Sample from the input spectrum histogram directly (after cross
424  /// section weighting), without any biasing
425  ///
426  /// @var TimeGenSamplingMode::UNIFORM_TIME
427  /// @brief Sample energies without biasing, sample times uniformly
428  /// over the entire allowed range
429  ///
430  /// @var TimeGenSamplingMode::UNIFORM_ENERGY
431  /// @brief Sample energies uniformly over the entire allowed range,
432  /// sample times without biasing
434 
435  /// @brief Represents the sampling mode to use when selecting neutrino
436  /// times and energies
438 
439  /// @enum PinchParamType
440  /// @brief Enumerated type that defines the pinching parameter conventions
441  /// that are understood by this module
442  ///
443  /// @var PinchParamType::ALPHA
444  /// @brief Use the fitting convention in which the pinching parameter
445  /// is called @f$\alpha@f$
446  /// @details In the "alpha" convention, the probability density
447  /// function @f$p(E_\nu)@f$ describing the distribution of neutrino
448  /// energies @f$E_\nu@f$ is given by
449  /// @f$p(E_\nu) = @f$\left(\frac{\alpha + 1}{\left<E_\nu\right>}
450  /// \right)^{\alpha + 1}\frac{E_\nu^\alpha}{\Gamma(\alpha + 1)}
451  /// \exp\left(\frac{-[\alpha + 1]E_\nu}{\left<E_\nu\right>}\right)@f$
452  ///
453  /// @var PinchParamType::BETA
454  /// @brief Use the fitting convention in which the pinching parameter
455  /// is called @f$\beta@f$
456  /// @details In the "beta" convention, the probability density
457  /// function @f$p(E_\nu)@f$ describing the distribution of neutrino
458  /// energies @f$E_\nu@f$ is given by
459  /// @f$p(E_\nu) = @f$\left(\frac{\beta}{\left<E_\nu\right>}
460  /// \right)^{\beta}\frac{E_\nu^(\beta - 1)}{\Gamma(\beta)}
461  /// \exp\left(\frac{-\beta E_\nu}{\left<E_\nu\right>}\right)@f$
462  /// Note that this differs from the "alpha" convention via
463  /// @f$\beta = \alpha + 1@f$.
464  enum class PinchParamType { ALPHA, BETA };
465 
466  /// @brief The pinching parameter convention to use when interpreting the
467  /// time-dependent fits
469 
470  /// @enum SpectrumFileFormat
471  /// @brief Enumerated type that defines the allowed neutrino spectrum input
472  /// file formats
473  ///
474  /// @var SpectrumFileFormat::RootTH2D
475  /// @brief The incident neutrino spectrum is described by a TH2D object
476  /// stored in a ROOT file
477  /// @details The X axis of the TH2D should be the time axis (s) and the
478  /// Y axis should be the energy axis (MeV).
479  // TODO: add information about the bin value units
480  ///
481  /// @var SpectrumFileFormat::FIT
482  /// @brief The incident neutrino spectrum is described by a set of fitting
483  /// parameters for each time bin
484  // TODO: add a "fit"-format description here
485  enum class SpectrumFileFormat { RootTH2D, FIT };
486 
487  /// @brief Format to assume for the neutrino spectrum input file
489 
490  /// @brief The event TTree created by MARLEY
491  /// @details This tree will be saved to the "hist" output file for
492  /// validation purposes. The tree contains the same information as the
493  /// generated simb::MCTruth objects, but in MARLEY's internal format
494  TTree* fEventTree;
495 
496  /// @brief Run number from the art::Event being processed
497  uint_fast32_t fRunNumber;
498  /// @brief Subrun number from the art::Event being processed
499  uint_fast32_t fSubRunNumber;
500  /// @brief Event number from the art::Event being processed
501  uint_fast32_t fEventNumber;
502 
503  /// @brief Time since the supernova core bounce for the current MARLEY
504  /// neutrino vertex
505  double fTNu;
506 
507  /// @brief Statistical weight for the current MARLEY neutrino vertex
508  double fWeight;
509 
510  /// @brief Flux-averaged total cross section (fm<sup>2</sup>, average is
511  /// taken over all energies and times for all defined reactions) used by
512  /// MARLEY to generate neutrino vertices
514 
515  /// @brief The number of MARLEY neutrino vertices to generate in each
516  /// art::Event
517  unsigned int fNeutrinosPerEvent;
518 
519  /// @brief Minimum neutrino energy to consider when using a "fit"-format
520  /// spectrum file
521  double fFitEmin;
522 
523  /// @brief Maximum neutrino energy to consider when using a "fit"-format
524  /// spectrum file
525  double fFitEmax;
526 };
527 
528 //------------------------------------------------------------------------------
530  : EDProducer{ p.get_PSet() }, fEvent(new marley::Event), fRunNumber(0),
531  fSubRunNumber(0), fEventNumber(0), fTNu(0.), fFluxAveragedCrossSection(0.)
532 {
533  // Configure the module (including MARLEY itself) using the FHiCL parameters
534  this->reconfigure(p);
535 
536  // Create a ROOT TTree using the TFileService that will store the MARLEY
537  // event objects (useful for debugging purposes)
538  art::ServiceHandle<art::TFileService const> tfs;
539  fEventTree = tfs->make<TTree>("MARLEY_event_tree",
540  "Neutrino events generated by MARLEY");
541  fEventTree->Branch("event", "marley::Event", fEvent.get());
542 
543  // Add branches that give the art::Event run, subrun, and event numbers for
544  // easy match-ups between the MARLEY and art TTrees. All three are recorded
545  // as 32-bit unsigned integers.
546  fEventTree->Branch("run_number", &fRunNumber, "run_number/i");
547  fEventTree->Branch("subrun_number", &fSubRunNumber, "subrun_number/i");
548  fEventTree->Branch("event_number", &fEventNumber, "event_number/i");
549  fEventTree->Branch("tSN", &fTNu, "tSN/D");
550  fEventTree->Branch("weight", &fWeight, "weight/D");
551 
552  produces< std::vector<simb::MCTruth> >();
553  produces< std::vector<sim::SupernovaTruth> >();
554  produces< art::Assns<simb::MCTruth, sim::SupernovaTruth> >();
556 }
557 
558 //------------------------------------------------------------------------------
559 void evgen::MarleyTimeGen::beginRun(art::Run& run) {
560  art::ServiceHandle<geo::Geometry const> geo;
561  run.put(std::make_unique<sumdata::RunData>(geo->DetectorName()));
562 }
563 
564 //------------------------------------------------------------------------------
565 void evgen::MarleyTimeGen::create_truths_th2d(simb::MCTruth& mc_truth,
566  sim::SupernovaTruth& sn_truth, const TLorentzVector& vertex_pos)
567 {
568  // Get a reference to the generator object created by MARLEY (we'll need
569  // to do a few fancy things with it other than just creating events)
570  marley::Generator& gen = fMarleyHelper->get_generator();
571 
572  if (fSamplingMode == TimeGenSamplingMode::HISTOGRAM)
573  {
574  // Generate a MARLEY event using the time-integrated spectrum
575  // (the generator was already configured to use it by reconfigure())
576  mc_truth = fMarleyHelper->create_MCTruth(vertex_pos,
577  fEvent.get());
578 
579  // Find the time distribution corresponding to the selected energy bin
580  double E_nu = fEvent->projectile().total_energy();
581  int E_bin_index = fSpectrumHist->GetYaxis()->FindBin(E_nu);
582  TH1D* t_hist = fSpectrumHist->ProjectionX("dummy_time_hist", E_bin_index,
583  E_bin_index);
584  double* time_bin_weights = t_hist->GetArray();
585 
586  // Sample a time bin from the distribution
587  std::discrete_distribution<int> time_dist;
588  std::discrete_distribution<int>::param_type time_params(time_bin_weights,
589  &(time_bin_weights[t_hist->GetNbinsX() + 1]));
590  int time_bin_index = gen.sample_from_distribution(time_dist, time_params);
591 
592  // Sample a time uniformly from within the selected time bin
593  double t_min = t_hist->GetBinLowEdge(time_bin_index);
594  double t_max = t_min + t_hist->GetBinWidth(time_bin_index);
595  // sample a time on [ t_min, t_max )
596  fTNu = gen.uniform_random_double(t_min, t_max, false);
597  // Unbiased sampling was used, so assign this neutrino vertex a
598  // unit statistical weight
599  fWeight = ONE;
600 
601  sn_truth = sim::SupernovaTruth(fTNu, fWeight, fFluxAveragedCrossSection,
603  }
604 
605  else if (fSamplingMode == TimeGenSamplingMode::UNIFORM_TIME)
606  {
607  // Generate a MARLEY event using the time-integrated spectrum
608  // (the generator was already configured to use it by reconfigure())
609  mc_truth = fMarleyHelper->create_MCTruth(vertex_pos,
610  fEvent.get());
611 
612  // Sample a time uniformly
613  TAxis* time_axis = fSpectrumHist->GetXaxis();
614  // underflow bin has index zero
615  double t_min = time_axis->GetBinLowEdge(1);
616  double t_max = time_axis->GetBinLowEdge(time_axis->GetNbins() + 1);
617  // sample a time on [ t_min, t_max )
618  fTNu = gen.uniform_random_double(t_min, t_max, false);
619 
620  // Get the value of the true dependent probability density (probability
621  // of the sampled time given the sampled energy) to use as a biasing
622  // correction in the neutrino vertex weight.
623  double E_nu = fEvent->projectile().total_energy();
624  int E_bin_index = fSpectrumHist->GetYaxis()->FindBin(E_nu);
625  TH1D* t_hist = fSpectrumHist->ProjectionX("dummy_time_hist", E_bin_index,
626  E_bin_index);
627  int t_bin_index = t_hist->FindBin(fTNu);
628  double weight_bias = t_hist->GetBinContent(t_bin_index) * (t_max - t_min)
629  / ( t_hist->Integral() * t_hist->GetBinWidth(t_bin_index) );
630 
631  fWeight = weight_bias;
632 
633  sn_truth = sim::SupernovaTruth(fTNu, fWeight, fFluxAveragedCrossSection,
635  }
636 
637  else if (fSamplingMode == TimeGenSamplingMode::UNIFORM_ENERGY)
638  {
639  // Select a time bin using the energy-integrated spectrum
640  TH1D* t_hist = fSpectrumHist->ProjectionX("dummy_time_hist");
641  double* time_bin_weights = t_hist->GetArray();
642 
643  // Sample a time bin from the distribution
644  std::discrete_distribution<int> time_dist;
645  std::discrete_distribution<int>::param_type time_params(time_bin_weights,
646  &(time_bin_weights[t_hist->GetNbinsX() + 1]));
647  int time_bin_index = gen.sample_from_distribution(time_dist, time_params);
648 
649  // Sample a time uniformly from within the selected time bin
650  double t_min = t_hist->GetBinLowEdge(time_bin_index);
651  double t_max = t_min + t_hist->GetBinWidth(time_bin_index);
652  // sample a time on [ t_min, t_max )
653  fTNu = gen.uniform_random_double(t_min, t_max, false);
654 
655  // Sample an energy uniformly over the entire allowed range
656  // underflow bin has index zero
657  TAxis* energy_axis = fSpectrumHist->GetYaxis();
658  double E_min = energy_axis->GetBinLowEdge(1);
659  double E_max = energy_axis->GetBinLowEdge(energy_axis->GetNbins() + 1);
660 
661  double E_nu = std::numeric_limits<double>::lowest();
662 
663  // Generate a MARLEY event using a uniformly sampled energy
664  mc_truth = make_uniform_energy_mctruth(E_min, E_max, E_nu, vertex_pos);
665 
666  // Get the value of the true dependent probability density (probability
667  // of the sampled energy given the sampled time) to use as a biasing
668  // correction in the neutrino vertex weight.
669  //
670  // Get a 1D projection of the energy spectrum for the sampled time bin
671  TH1D* energy_spect = fSpectrumHist->ProjectionY("energy_spect",
672  time_bin_index, time_bin_index);
673 
674  // Create a new MARLEY neutrino source object using this projection (this
675  // will create a normalized probability density that we can use) and load
676  // it into the generator.
677  auto nu_source = marley_root::make_root_neutrino_source(
678  marley_utils::ELECTRON_NEUTRINO, energy_spect);
679  double new_source_E_min = nu_source->get_Emin();
680  double new_source_E_max = nu_source->get_Emax();
681  gen.set_source(std::move(nu_source));
682  // NOTE: The marley::Generator object normalizes the E_pdf to unity
683  // automatically, but just in case, we redo it here.
684  double E_pdf_integ = integrate([&gen](double E_nu)
685  -> double { return gen.E_pdf(E_nu); }, new_source_E_min,
686  new_source_E_max);
687 
688  // Compute the likelihood ratio that we need to bias the neutrino vertex
689  // weight
690  double weight_bias = (gen.E_pdf(E_nu) / E_pdf_integ) * (E_max - E_min);
691 
692  fWeight = weight_bias;
693 
694  sn_truth = sim::SupernovaTruth(fTNu, fWeight, fFluxAveragedCrossSection,
696  }
697 
698  else {
699  throw cet::exception("MARLEYTimeGen") << "Unrecognized sampling mode"
700  << " encountered in evgen::MarleyTimeGen::produce()";
701  }
702 
703 }
704 
705 //------------------------------------------------------------------------------
707 {
708  art::ServiceHandle<geo::Geometry const> geo;
709 
710  // Get the run, subrun, and event numbers from the current art::Event
711  fRunNumber = e.run();
712  fSubRunNumber = e.subRun();
713  fEventNumber = e.event();
714 
715  // Prepare associations and vectors of truth objects that will be produced
716  // and loaded into the current art::Event
717  std::unique_ptr< std::vector<simb::MCTruth> >
718  truthcol(new std::vector<simb::MCTruth>);
719 
720  std::unique_ptr< std::vector<sim::SupernovaTruth> >
721  sn_truthcol(new std::vector<sim::SupernovaTruth>);
722 
723  std::unique_ptr< art::Assns<simb::MCTruth, sim::SupernovaTruth> >
724  truth_assns(new art::Assns<simb::MCTruth, sim::SupernovaTruth>);
725 
726  // Create temporary truth objects that we will use to load the event
727  simb::MCTruth truth;
728  sim::SupernovaTruth sn_truth;
729 
730  for (unsigned int n = 0; n < fNeutrinosPerEvent; ++n) {
731 
732  // Sample a primary vertex location for this event
733  TLorentzVector vertex_pos = fVertexSampler->sample_vertex_pos(*geo);
734 
735  // Reset the neutrino's time-since-supernova to a bogus value (for now)
736  fTNu = std::numeric_limits<double>::lowest();
737 
738  if (fSpectrumFileFormat == SpectrumFileFormat::RootTH2D) {
739  create_truths_th2d(truth, sn_truth, vertex_pos);
740  }
741  else if (fSpectrumFileFormat == SpectrumFileFormat::FIT) {
742  create_truths_time_fit(truth, sn_truth, vertex_pos);
743  }
744  else {
745  throw cet::exception("MARLEYTimeGen") << "Invalid spectrum file"
746  << " format encountered in evgen::MarleyTimeGen::produce()";
747  }
748 
749  // Write the marley::Event object to the event tree
750  fEventTree->Fill();
751 
752  // Add the truth objects to the appropriate vectors
753  truthcol->push_back(truth);
754 
755  sn_truthcol->push_back(sn_truth);
756 
757  // Associate the last entries in each of the truth object vectors (the
758  // truth objects that were just created for the current neutrino vertex)
759  // with each other
760  truth_assns->addSingle(art::PtrMaker<simb::MCTruth>{e}(truthcol->size() - 1),
761  art::PtrMaker<sim::SupernovaTruth>{e}(sn_truthcol->size() - 1));
762  }
763 
764  // Load the completed truth object vectors and associations into the event
765  e.put(std::move(truthcol));
766 
767  e.put(std::move(sn_truthcol));
768 
769  e.put(std::move(truth_assns));
770 }
771 
772 //------------------------------------------------------------------------------
774 {
775  const auto& seed_service = art::ServiceHandle<rndm::NuRandomService>();
776  const auto& geom_service = art::ServiceHandle<geo::Geometry const>();
777 
778  // Create a new evgen::ActiveVolumeVertexSampler object based on the current
779  // configuration
780  fVertexSampler = std::make_unique<evgen::ActiveVolumeVertexSampler>(
781  p().vertex_, *seed_service, *geom_service, "MARLEY_Vertex_Sampler");
782 
783  // Create a new marley::Generator object based on the current configuration
784  const fhicl::ParameterSet marley_pset = p.get_PSet()
785  .get< fhicl::ParameterSet >( "marley_parameters" );
786  fMarleyHelper = std::make_unique<MARLEYHelper>( marley_pset,
787  *seed_service, "MARLEY" );
788 
789  // Get the number of neutrino vertices per event from the FHiCL parameters
790  fNeutrinosPerEvent = p().nu_per_event_();
791 
792  // Determine the current sampling mode from the FHiCL parameters
793  const std::string& samp_mode_str = p().sampling_mode_();
794  if (samp_mode_str == "histogram")
795  fSamplingMode = TimeGenSamplingMode::HISTOGRAM;
796  else if (samp_mode_str == "uniform time")
797  fSamplingMode = TimeGenSamplingMode::UNIFORM_TIME;
798  else if (samp_mode_str == "uniform energy")
799  fSamplingMode = TimeGenSamplingMode::UNIFORM_ENERGY;
800  else throw cet::exception("MARLEYTimeGen")
801  << "Invalid sampling mode \"" << samp_mode_str << "\""
802  << " specified for the MARLEYTimeGen module.";
803 
804  MF_LOG_INFO("MARLEYTimeGen") << fNeutrinosPerEvent << " neutrino vertices"
805  << " will be generated for each art::Event using the \"" << samp_mode_str
806  << "\" sampling mode.";
807 
808  // Retrieve the time-dependent neutrino spectrum from the spectrum file.
809  // Use different methods depending on the file's format.
810  std::string spectrum_file_format = marley_utils::to_lowercase(
811  p().spectrum_file_format_());
812 
813  if (spectrum_file_format == "th2d")
814  fSpectrumFileFormat = SpectrumFileFormat::RootTH2D;
815  else if (spectrum_file_format == "fit") {
816  fSpectrumFileFormat = SpectrumFileFormat::FIT;
817 
818  std::string pinch_type;
819  if ( !p().pinching_parameter_type_(pinch_type) ) {
820  throw cet::exception("MARLEYTimeGen") << "Missing pinching parameter"
821  << " type for a \"fit\" format spectrum file";
822  }
823 
824  marley_utils::to_lowercase_inplace(pinch_type);
825  if (pinch_type == "alpha") fPinchType = PinchParamType::ALPHA;
826  else if (pinch_type == "beta") fPinchType = PinchParamType::BETA;
827  else throw cet::exception("MARLEYTimeGen")
828  << "Invalid pinching parameter type \"" << pinch_type
829  << "\" specified for the MARLEYTimeGen module.";
830 
831  if ( !p().fit_Emin_(fFitEmin) ) throw cet::exception("MARLEYTimeGen")
832  << "Missing minimum energy for a \"fit\" format spectrum"
833  << " used by the MARLEYTimeGen module.";
834 
835  if ( !p().fit_Emax_(fFitEmax) ) throw cet::exception("MARLEYTimeGen")
836  << "Missing maximum energy for a \"fit\" format spectrum"
837  << " used by the MARLEYTimeGen module.";
838 
839  if (fFitEmax < fFitEmin) throw cet::exception("MARLEYTimeGen")
840  << "Maximum energy is less than the minimum energy for"
841  << " a \"fit\" format spectrum used by the MARLEYTimeGen module.";
842  }
843  else throw cet::exception("MARLEYTimeGen")
844  << "Invalid spectrum file format \"" << p().spectrum_file_format_()
845  << "\" specified for the MARLEYTimeGen module.";
846 
847  // Determine the full file name (including path) of the spectrum file
848  std::string full_spectrum_file_name
849  = fMarleyHelper->find_file(p().spectrum_file_(), "spectrum");
850 
851  marley::Generator& gen = fMarleyHelper->get_generator();
852 
853  if (fSpectrumFileFormat == SpectrumFileFormat::RootTH2D) {
854 
855  // Retrieve the time-dependent neutrino flux from a ROOT file
856  std::unique_ptr<TFile> spectrum_file
857  = std::make_unique<TFile>(full_spectrum_file_name.c_str(), "read");
858  TH2D* temp_h2 = nullptr;
859  std::string namecycle;
860  if ( !p().namecycle_(namecycle) ) {
861  throw cet::exception("MARLEYTimeGen") << "Missing namecycle for"
862  << " a TH2D spectrum file";
863  }
864 
865  spectrum_file->GetObject(namecycle.c_str(), temp_h2);
866  fSpectrumHist.reset(temp_h2);
867 
868  // Disassociate the TH2D from its parent TFile. If we fail to do this,
869  // then ROOT will auto-delete the TH2D when the TFile goes out of scope.
870  fSpectrumHist->SetDirectory(nullptr);
871 
872  // Compute the flux-averaged total cross section using MARLEY. This will be
873  // used to compute neutrino vertex weights for the sim::SupernovaTruth
874  // objects.
875 
876  // Get a 1D projection of the energy spectrum (integrated over time)
877  TH1D* energy_spect = fSpectrumHist->ProjectionY("energy_spect");
878 
879  // Create a new MARLEY neutrino source object using this projection
880  // TODO: replace the hard-coded electron neutrino PDG code here (and in
881  // several other places in this source file) when you're ready to use
882  // MARLEY with multiple neutrino flavors
883  std::unique_ptr<marley::NeutrinoSource> nu_source
884  = marley_root::make_root_neutrino_source(marley_utils::ELECTRON_NEUTRINO,
885  energy_spect);
886 
887  // Factor of hbar_c^2 converts from MeV^(-2) to fm^2
888  fFluxAveragedCrossSection = marley_utils::hbar_c2
889  * flux_averaged_total_xs(*nu_source, gen);
890 
891  // For speed, sample energies first whenever possible (and then sample from
892  // an energy-dependent timing distribution). This avoids unnecessary calls
893  // to MARLEY to change the energy spectrum.
894  if (fSamplingMode == TimeGenSamplingMode::HISTOGRAM ||
895  fSamplingMode == TimeGenSamplingMode::UNIFORM_TIME)
896  {
897  gen.set_source(std::move(nu_source));
898  }
899 
900  } // spectrum_file_format == "th2d"
901 
902  else if (fSpectrumFileFormat == SpectrumFileFormat::FIT) {
903 
904  // Clear out the old parameterized spectrum, if one exists
905  fTimeFits.clear();
906 
907  std::ifstream fit_file(full_spectrum_file_name);
908  std::string line;
909 
910  bool found_end = false;
911 
912  // current line number
913  int line_num = 0;
914  // number of lines checked in last call to marley_utils::get_next_line()
915  int lines_checked = 0;
916 
917  double old_time = std::numeric_limits<double>::lowest();
918 
919  while ( line = marley_utils::get_next_line(fit_file, rx_comment_or_empty,
920  false, lines_checked), line_num += lines_checked, !line.empty() )
921  {
922  if (found_end) {
923  MF_LOG_WARNING("MARLEYTimeGen") << "Trailing content after last time"
924  << " bin found on line " << line_num << " of the spectrum file "
925  << full_spectrum_file_name;
926  }
927 
928  double time;
929  double Emean_nue, alpha_nue, luminosity_nue;
930  double Emean_nuebar, alpha_nuebar, luminosity_nuebar;
931  double Emean_nux, alpha_nux, luminosity_nux;
932  std::istringstream iss(line);
933  bool ok_first = static_cast<bool>( iss >> time );
934 
935  if (time <= old_time) throw cet::exception("MARLEYTimeGen")
936  << "Time bin left edges given in the spectrum file must be"
937  << " strictly increasing. Invalid time bin value found on line "
938  << line_num << " of the spectrum file " << full_spectrum_file_name;
939  else old_time = time;
940 
941  bool ok_rest = static_cast<bool>( iss >> Emean_nue >> alpha_nue
942  >> luminosity_nue >> Emean_nuebar >> alpha_nuebar >> luminosity_nuebar
943  >> Emean_nux >> alpha_nux >> luminosity_nux
944  );
945 
946  if (ok_first) {
947  // We haven't reached the final bin, so add another time bin
948  // in the typical way.
949  if (ok_rest) {
950  fTimeFits.emplace_back(time, Emean_nue, alpha_nue, luminosity_nue,
951  Emean_nuebar, alpha_nuebar, luminosity_nuebar, Emean_nux,
952  alpha_nux, luminosity_nux);
953  }
954  else {
955  make_final_timefit(time);
956  found_end = true;
957  }
958  }
959  else throw cet::exception("MARLEYTimeGen") << "Parse error on line "
960  << line_num << " of the spectrum file " << full_spectrum_file_name;
961  }
962 
963  if (!found_end) {
964 
965  size_t num_time_bins = fTimeFits.size();
966  if (num_time_bins < 2) throw cet::exception("MARLEYTimeGen")
967  << "Missing right edge for the final time bin in the spectrum file "
968  << full_spectrum_file_name << ". Unable to guess a bin width for the "
969  << " final bin.";
970 
971  // Guess that the width of the penultimate time bin and the width of
972  // the final time bin are the same
973  double delta_t_bin = fTimeFits.back().Time()
974  - fTimeFits.at(num_time_bins - 2).Time();
975 
976  double last_bin_right_edge = fTimeFits.back().Time() + delta_t_bin;
977 
978  make_final_timefit(last_bin_right_edge);
979 
980  MF_LOG_WARNING("MARLEYTimeGen") << "Missing right"
981  << " edge for the final time bin in the spectrum file "
982  << full_spectrum_file_name << ". Assuming a width of "
983  << delta_t_bin << " s for the final bin.";
984  }
985 
986 
987  // Compute the flux-averaged total cross section for the fitted spectrum.
988  // We will need this to compute neutrino vertex weights.
989  std::vector< std::unique_ptr< marley::NeutrinoSource > > fit_sources;
990  for (const auto& fit : fTimeFits) {
991  fit_sources.emplace_back( source_from_time_fit(fit) );
992  }
993 
994  // TODO: add handling for non-nu_e neutrino types when suitable data become
995  // available in MARLEY
996  auto temp_source = std::make_unique<marley::FunctionNeutrinoSource>(
997  marley_utils::ELECTRON_NEUTRINO, fFitEmin, fFitEmax,
998  [&fit_sources, this](double E_nu) -> double {
999  double flux = 0.;
1000  for (size_t s = 0; s < fit_sources.size(); ++s) {
1001 
1002  const TimeFit& current_fit = this->fTimeFits.at(s);
1003  const auto& current_source = fit_sources.at(s);
1004  const FitParameters& params = current_fit.GetFitParameters(
1005  current_source->get_pid() );
1006 
1007  double lum = params.Luminosity();
1008 
1009  // Skip entries with zero luminosity, since they won't contribute
1010  // anything to the overall integral. Skip negative luminosity ones as
1011  // well, just in case.
1012  if (lum <= 0.) continue;
1013 
1014  flux += lum * current_source->pdf(E_nu);
1015  }
1016  return flux;
1017  }
1018  );
1019 
1020  double flux_integ = 0.;
1021  double tot_xs_integ = 0.;
1022  flux_averaged_total_xs(*temp_source, gen, flux_integ, tot_xs_integ);
1023 
1024  // Factor of hbar_c^2 converts from MeV^(-2) to fm^2
1025  fFluxAveragedCrossSection = marley_utils::hbar_c2
1026  * tot_xs_integ / flux_integ;
1027 
1028  make_nu_emission_histograms();
1029 
1030  } // spectrum_file_format == "fit"
1031 
1032  else {
1033  throw cet::exception("MARLEYTimeGen") << "Unrecognized neutrino spectrum"
1034  << " file format \"" << p().spectrum_file_format_() << "\" encountered"
1035  << " in evgen::MarleyTimeGen::reconfigure()";
1036  }
1037 
1038  MF_LOG_INFO("MARLEYTimeGen") << "The flux-averaged total cross section"
1039  << " predicted by MARLEY for the current supernova spectrum is "
1040  << fFluxAveragedCrossSection << " fm^2";
1041 
1042 }
1043 
1044 //------------------------------------------------------------------------------
1046  sim::SupernovaTruth& sn_truth, const TLorentzVector& vertex_pos)
1047 {
1048  // Get a reference to the generator object created by MARLEY (we'll need
1049  // to do a few fancy things with it other than just creating events)
1050  marley::Generator& gen = fMarleyHelper->get_generator();
1051 
1052  // Initialize the time bin index to something absurdly large. This will help
1053  // us detect strange bugs that arise when it is sampled incorrectly.
1054  size_t time_bin_index = std::numeric_limits<size_t>::max();
1055 
1056  // Create an object to represent the discrete time distribution given in
1057  // the spectrum file. Use the luminosity for each bin as its sampling weight.
1058  // This distribution will actually be used to sample a neutrino arrival time
1059  // unless we're using the uniform time sampling mode, in which case it will
1060  // be used to help calculate the neutrino vertex weight.
1061  int source_pdg_code = gen.get_source().get_pid();
1062 
1063  const auto fit_params_begin = TimeFit::make_FitParameters_iterator(
1064  source_pdg_code, fTimeFits.begin() );
1065  const auto fit_params_end = TimeFit::make_FitParameters_iterator(
1066  source_pdg_code, fTimeFits.end() );
1067 
1068  const auto lum_begin = FitParameters::make_luminosity_iterator(
1069  fit_params_begin);
1070  const auto lum_end = FitParameters::make_luminosity_iterator(
1071  fit_params_end);
1072 
1073  std::discrete_distribution<size_t> time_dist(lum_begin, lum_end);
1074 
1075  if (fSamplingMode == TimeGenSamplingMode::HISTOGRAM
1076  || fSamplingMode == TimeGenSamplingMode::UNIFORM_ENERGY)
1077  {
1078  time_bin_index = gen.sample_from_distribution(time_dist);
1079  }
1080 
1081  else if (fSamplingMode == TimeGenSamplingMode::UNIFORM_TIME) {
1082  int last_time_index = 0;
1083  if (fTimeFits.size() > 0) last_time_index = fTimeFits.size() - 1;
1084  std::uniform_int_distribution<size_t> uid(0, last_time_index);
1085  // for c2: time_bin_index has already been declared
1086  time_bin_index = gen.sample_from_distribution(uid);
1087  }
1088 
1089  else {
1090  throw cet::exception("MARLEYTimeGen") << "Unrecognized sampling mode"
1091  << " encountered in evgen::MarleyTimeGen::produce()";
1092  }
1093 
1094  // Sample a time uniformly from within the selected time bin. Note that
1095  // the entries in fTimeFits use the time_ member to store the bin left
1096  // edges. The module creates fTimeFits in such a way that its last element
1097  // will always have luminosity_ == 0. (zero sampling weight), so we may
1098  // always add one to the sampled bin index without worrying about going
1099  // off the edge.
1100  double t_min = fTimeFits.at(time_bin_index).Time();
1101  double t_max = fTimeFits.at(time_bin_index + 1).Time();
1102  // sample a time on [ t_min, t_max )
1103  fTNu = gen.uniform_random_double(t_min, t_max, false);
1104 
1105  // Create a "beta-fit" neutrino source using the correct parameters for the
1106  // sampled time bin. This will be used to sample a neutrino energy unless
1107  // we're using the uniform time sampling mode. For uniform time sampling,
1108  // it will be used to determine the neutrino event weight.
1109  const auto& fit = fTimeFits.at(time_bin_index);
1110  std::unique_ptr<marley::NeutrinoSource> nu_source = source_from_time_fit(fit);
1111 
1112  if (fSamplingMode == TimeGenSamplingMode::HISTOGRAM
1113  || fSamplingMode == TimeGenSamplingMode::UNIFORM_TIME)
1114  {
1115  // Replace the generator's old source with the new one for the current
1116  // time bin
1117  gen.set_source(std::move(nu_source));
1118 
1119  // Generate a MARLEY event using the updated source
1120  mc_truth = fMarleyHelper->create_MCTruth(vertex_pos, fEvent.get());
1121 
1122  if (fSamplingMode == TimeGenSamplingMode::HISTOGRAM) {
1123  // Unbiased sampling creates neutrino vertices with unit weight
1124  fWeight = ONE;
1125  sn_truth = sim::SupernovaTruth(fTNu, fWeight, fFluxAveragedCrossSection,
1126  sim::kUnbiased);
1127  }
1128  else {
1129  // fSamplingMode == TimeGenSamplingMode::UNIFORM_TIME
1130 
1131  // Multiply by the likelihood ratio in order to correct for uniform
1132  // time sampling if we're using that biasing method
1133  double weight_bias = time_dist.probabilities().at(time_bin_index)
1134  / (t_max - t_min) * ( fTimeFits.back().Time()
1135  - fTimeFits.front().Time() );
1136 
1137  fWeight = weight_bias;
1138  sn_truth = sim::SupernovaTruth(fTNu, fWeight, fFluxAveragedCrossSection,
1140  }
1141  }
1142 
1143  else if (fSamplingMode == TimeGenSamplingMode::UNIFORM_ENERGY)
1144  {
1145  double E_nu = std::numeric_limits<double>::lowest();
1146  mc_truth = make_uniform_energy_mctruth(fFitEmin, fFitEmax, E_nu,
1147  vertex_pos);
1148 
1149  // Get the value of the true dependent probability density (probability
1150  // of the sampled energy given the sampled time) to use as a biasing
1151  // correction in the neutrino vertex weight.
1152 
1153  // Load the generator with the neutrino source that represents the
1154  // true (i.e., unbiased) energy probability distribution. This will
1155  // create a normalized probability density that we can use to determine
1156  // the neutrino vertex weight.
1157  double nu_source_E_min = nu_source->get_Emin();
1158  double nu_source_E_max = nu_source->get_Emax();
1159  gen.set_source(std::move(nu_source));
1160 
1161  // NOTE: The marley::Generator object normalizes the E_pdf to unity
1162  // automatically, but just in case, we redo it here.
1163  double E_pdf_integ = integrate([&gen](double E_nu)
1164  -> double { return gen.E_pdf(E_nu); }, nu_source_E_min,
1165  nu_source_E_max);
1166 
1167  // Compute the likelihood ratio that we need to bias the neutrino vertex
1168  // weight
1169  double weight_bias = (gen.E_pdf(E_nu) / E_pdf_integ)
1170  * (fFitEmax - fFitEmin);
1171 
1172  fWeight = weight_bias;
1173  sn_truth = sim::SupernovaTruth(fTNu, fWeight, fFluxAveragedCrossSection,
1175  }
1176 
1177  else {
1178  throw cet::exception("MARLEYTimeGen") << "Unrecognized sampling mode"
1179  << " encountered in evgen::MarleyTimeGen::produce()";
1180  }
1181 
1182 }
1183 
1184 //------------------------------------------------------------------------------
1185 std::unique_ptr<marley::NeutrinoSource>
1187 {
1188  // Create a "beta-fit" neutrino source using the given fit parameters.
1189 
1190  // The two common fitting schemes (alpha and beta) differ in their
1191  // definitions by \beta = \alpha + 1.
1192  // TODO: add handling for non-nu_e neutrino types
1193  const FitParameters& nue_parameters
1194  = fit.GetFitParameters(marley_utils::ELECTRON_NEUTRINO);
1195 
1196  double beta = nue_parameters.Alpha();
1197  if (fPinchType == PinchParamType::ALPHA) beta += 1.;
1198  else if (fPinchType != PinchParamType::BETA) {
1199  throw cet::exception("MARLEYTimeGen") << "Unreognized pinching parameter"
1200  << " type encountered in evgen::MarleyTimeGen::source_from_time_fit()";
1201  }
1202 
1203  // Create the new source
1204  std::unique_ptr<marley::NeutrinoSource> nu_source
1205  = std::make_unique<marley::BetaFitNeutrinoSource>(
1206  marley_utils::ELECTRON_NEUTRINO, fFitEmin, fFitEmax,
1207  nue_parameters.Emean(), beta);
1208 
1209  return nu_source;
1210 }
1211 
1212 //------------------------------------------------------------------------------
1214  double E_max, double& E_nu, const TLorentzVector& vertex_pos)
1215 {
1216  marley::Generator& gen = fMarleyHelper->get_generator();
1217 
1218  // Sample an energy uniformly over the entire allowed range
1219  double total_xs;
1220  int j = 0;
1221  do {
1222  // We have to check that the cross section is nonzero for the sampled
1223  // energy (otherwise we'll generate an unphysical event). However, if the
1224  // user has given us a histogram that is below threshold, the
1225  // program could get stuck here endlessly, sampling rejected energy
1226  // after rejected energy. Just in case, we cap the total number of tries
1227  // and quit if things don't work out.
1228  if (j >= MAX_UNIFORM_ENERGY_ITERATIONS) {
1229  throw cet::exception("MARLEYTimeGen") << "Exceeded the maximum of "
1230  << MAX_UNIFORM_ENERGY_ITERATIONS << " while attempting to sample"
1231  << " a neutrino energy uniformly.";
1232  }
1233  // Sample an energy uniformly on [ E_min, E_max )
1234  E_nu = gen.uniform_random_double(E_min, E_max, false);
1235  total_xs = 0.;
1236  // Check that at least one defined reaction has a non-zero total
1237  // cross section at the sampled energy. If this is not the case, try
1238  // again.
1239  for (const auto& react : gen.get_reactions()) {
1240  total_xs += react->total_xs(marley_utils::ELECTRON_NEUTRINO, E_nu);
1241  }
1242 
1243  ++j;
1244  } while (total_xs <= 0.);
1245 
1246  // Replace the existing neutrino source with a monoenergetic one at the
1247  // neutrino energy that was just sampled above
1248  std::unique_ptr<marley::NeutrinoSource> nu_source
1249  = std::make_unique<marley::MonoNeutrinoSource>(
1250  marley_utils::ELECTRON_NEUTRINO, E_nu);
1251  gen.set_source(std::move(nu_source));
1252 
1253  // Generate a MARLEY event using the new monoenergetic source
1254  auto mc_truth = fMarleyHelper->create_MCTruth(vertex_pos, fEvent.get());
1255 
1256  return mc_truth;
1257 }
1258 
1259 //------------------------------------------------------------------------------
1261 {
1262  // The final time bin has zero luminosity, and therefore zero sampling
1263  // weight. We need it to be present so that the last nonzero weight bin
1264  // has a right edge.
1265  fTimeFits.emplace_back(time, 0., 0., 0., 0., 0., 0., 0., 0., 0.);
1266 }
1267 
1268 //------------------------------------------------------------------------------
1270 {
1271  // To make a histogram with variable size bins, ROOT needs the bin
1272  // low edges to be contiguous in memory. This is not true of the
1273  // stored times in fTimeFits, so we'll need to make a temporary copy
1274  // of them.
1275  std::vector<double> temp_times;
1276  for (const auto& fit : fTimeFits) temp_times.push_back(fit.Time());
1277 
1278  // If, for some reason, there are fewer than two time bins, just return
1279  // without making the histograms.
1280  // TODO: consider throwing an exception here
1281  if (temp_times.size() < 2) return;
1282 
1283  // Get the number of time bins
1284  int num_bins = temp_times.size() - 1;
1285 
1286  // Create some ROOT TH1D objects using the TFileService. These will store
1287  // the number of emitted neutrinos of each type in each time bin
1288  art::ServiceHandle<art::TFileService const> tfs;
1289  TH1D* nue_hist = tfs->make<TH1D>("NueEmission",
1290  "Number of emitted #nu_{e}; time (s); number of neutrinos in time bin",
1291  num_bins, temp_times.data());
1292  TH1D* nuebar_hist = tfs->make<TH1D>("NuebarEmission",
1293  "Number of emitted #bar{#nu}_{e}; time (s); number of neutrinos in"
1294  " time bin", num_bins, temp_times.data());
1295  TH1D* nux_hist = tfs->make<TH1D>("NuxEmission",
1296  "Number of emitted #nu_{x}; time (s); number of neutrinos in time bin",
1297  num_bins, temp_times.data());
1298 
1299  // Load the histograms with the emitted neutrino counts
1300  for (size_t b = 1; b < temp_times.size(); ++b) {
1301  const TimeFit& current_fit = fTimeFits.at(b - 1);
1302  const TimeFit& next_fit = fTimeFits.at(b);
1303  double bin_deltaT = next_fit.Time() - current_fit.Time();
1304 
1305  const auto& nue_params = current_fit.GetFitParameters(
1306  marley_utils::ELECTRON_NEUTRINO);
1307  const auto& nuebar_params = current_fit.GetFitParameters(
1308  marley_utils::ELECTRON_ANTINEUTRINO);
1309  const auto& nux_params = current_fit.GetFitParameters(
1310  marley_utils::MUON_NEUTRINO);
1311 
1312  // Convert from bin luminosity to number of neutrinos by
1313  // multiplying by the bin time interval and dividing by the
1314  // mean neutrino energy
1315  double num_nue = 0.;
1316  double num_nuebar = 0.;
1317  double num_nux = 0.;
1318  if (nue_params.Emean() != 0.) num_nue = nue_params.Luminosity()
1319  * bin_deltaT / (nue_params.Emean() * MeV_to_erg);
1320  if (nuebar_params.Emean() != 0.) num_nuebar = nuebar_params.Luminosity()
1321  * bin_deltaT / (nuebar_params.Emean() * MeV_to_erg);
1322  if (nux_params.Emean() != 0.) num_nux = nux_params.Luminosity()
1323  * bin_deltaT / (nux_params.Emean() * MeV_to_erg);
1324 
1325  nue_hist->SetBinContent(b, num_nue);
1326  nuebar_hist->SetBinContent(b, num_nuebar);
1327  nux_hist->SetBinContent(b, num_nux);
1328  }
1329 
1330 }
1331 
1332 
1333 //------------------------------------------------------------------------------
1334 // Anonymous namespace function definitions
1335 namespace {
1336 
1337  double flux_averaged_total_xs(marley::NeutrinoSource& nu_source,
1338  marley::Generator& gen, double& source_integ, double& tot_xs_integ)
1339  {
1340  // Get an integral of the source PDF (in case it isn't normalized to 1)
1341  source_integ = integrate(
1342  [&nu_source](double E_nu) -> double { return nu_source.pdf(E_nu); },
1343  nu_source.get_Emin(), nu_source.get_Emax()
1344  );
1345 
1346  tot_xs_integ = integrate(
1347  [&nu_source, &gen](double E_nu) -> double
1348  {
1349  double xs = 0.;
1350  for (const auto& react : gen.get_reactions()) {
1351  xs += react->total_xs(marley_utils::ELECTRON_NEUTRINO, E_nu);
1352  }
1353  return xs * nu_source.pdf(E_nu);
1354  }, nu_source.get_Emin(), nu_source.get_Emax()
1355  );
1356 
1357  return tot_xs_integ / source_integ;
1358  }
1359 
1360  double flux_averaged_total_xs(marley::NeutrinoSource& nu_source,
1361  marley::Generator& gen)
1362  {
1363  double dummy1, dummy2;
1364  return flux_averaged_total_xs(nu_source, gen, dummy1, dummy2);
1365  }
1366 
1367  double integrate(const std::function<double(double)>& f, double x_min,
1368  double x_max)
1369  {
1370  static marley::Integrator integrator(NUM_INTEGRATION_SAMPLING_POINTS);
1371  return integrator.num_integrate(f, x_min, x_max);
1372  }
1373 
1374 }
1375 
1376 DEFINE_ART_MODULE(evgen::MarleyTimeGen)
double fLuminosity
Luminosity (erg / s)
SpectrumFileFormat
Enumerated type that defines the allowed neutrino spectrum input file formats.
virtual void produce(art::Event &e) override
LArSoft interface to the MARLEY (Model of Argon Reaction Low Energy Yields) supernova neutrino event ...
FitParameters fNuebarFitParams
Fitting parameters for electron antineutrinos in this time bin.
SpectrumFileFormat fSpectrumFileFormat
Format to assume for the neutrino spectrum input file.
uint_fast32_t fEventNumber
Event number from the art::Event being processed.
double Alpha() const
Pinching parameter (dimensionless)
Stores fitting parameters for all neutrino types from a single time bin in a &quot;fit&quot;-format spectrum fi...
TimeFit(double time, double Emean_nue, double alpha_nue, double lum_nue, double Emean_nuebar, double alpha_nuebar, double lum_nuebar, double Emean_nux, double alpha_nux, double lum_nux)
double fTime
Time bin low edge (s)
void create_truths_time_fit(simb::MCTruth &mc_truth, sim::SupernovaTruth &sn_truth, const TLorentzVector &vertex_pos)
Create simb::MCTruth and sim::SupernovaTruth objects using a neutrino spectrum described by a previou...
void set_Alpha(double alpha)
Set the pinching parameter.
fhicl::Atom< std::string > sampling_mode_
virtual void reconfigure(const Parameters &p)
pdgs p
Definition: selectors.fcl:22
PinchParamType fPinchType
The pinching parameter convention to use when interpreting the time-dependent fits.
double fEmean
Mean neutrino energy (MeV)
produces< sumdata::RunData, art::InRun >()
FitParameters fNueFitParams
Fitting parameters for electron neutrinos in this time bin.
for pfile in ack l reconfigure(.*) override"` do echo "checking $
Stores parsed fit parameters from a single time bin and neutrino type in a &quot;fit&quot;-format spectrum file...
std::unique_ptr< evgen::ActiveVolumeVertexSampler > fVertexSampler
Algorithm that allows us to sample vertex locations within the active volume(s) of the detector...
fhicl::Atom< std::string > module_type_
std::unique_ptr< evgen::MARLEYHelper > fMarleyHelper
Object that provides an interface to the MARLEY event generator.
double fFluxAveragedCrossSection
Flux-averaged total cross section (fm2, average is taken over all energies and times for all defined ...
double fWeight
Statistical weight for the current MARLEY neutrino vertex.
TTree * fEventTree
The event TTree created by MARLEY.
double Time() const
Returns the low edge (in s) of the time bin that uses these fit parameters.
void set_Luminosity(double lum)
Set the luminosity (erg / s)
uint_fast32_t fSubRunNumber
Subrun number from the art::Event being processed.
Arrival times were sampled uniformly.
Algorithm that samples vertex locations uniformly within the active volume of a detector. It is fully experiment-agnostic and multi-TPC aware.
uint_fast32_t fRunNumber
Run number from the art::Event being processed.
std::unique_ptr< marley::NeutrinoSource > source_from_time_fit(const TimeFit &fit)
Create a MARLEY neutrino source object using a set of fit parameters for a particular time bin...
FitParameters fNuxFitParams
Fitting parameters for non-electron-flavor neutrinos in this time bin.
TimeGenSamplingMode
Enumerated type that defines the allowed ways that a neutrino&#39;s energy and arrival time may be sample...
static marley::IteratorToMember< It, double > make_luminosity_iterator(It it)
Converts an iterator that points to a FitParameters object into an iterator that points to that objec...
art::EDProducer::Table< Config, KeysToIgnore > Parameters
MF_LOG_WARNING("SimWire")<< "SimWire is an example module that works for the "<< "MicroBooNE detector. Each experiment should implement "<< "its own version of this module to simulate electronics "<< "response."
void create_truths_th2d(simb::MCTruth &mc_truth, sim::SupernovaTruth &sn_truth, const TLorentzVector &vertex_pos)
Create simb::MCTruth and sim::SupernovaTruth objects using spectrum information from a ROOT TH2D...
double fTNu
Time since the supernova core bounce for the current MARLEY neutrino vertex.
auto end(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:585
PinchParamType
Enumerated type that defines the pinching parameter conventions that are understood by this module...
void set_Emean(double Emean)
Set the mean neutrino energy (MeV)
virtual void beginRun(art::Run &run) override
static FitParameters TimeFit::* GetFitParametersMemberPointer(int pdg_code)
Helper function that returns a pointer-to-member for the FitParameters object appropriate for a given...
simb::MCTruth make_uniform_energy_mctruth(double E_min, double E_max, double &E_nu, const TLorentzVector &vertex_pos)
Creates a simb::MCTruth object using a uniformly-sampled neutrino energy.
BEGIN_PROLOG vertical distance to the surface Name
double fFitEmax
Maximum neutrino energy to consider when using a &quot;fit&quot;-format spectrum file.
fhicl::OptionalAtom< double > fit_Emax_
static marley::IteratorToMember< It, FitParameters > make_FitParameters_iterator(int pdg_code, It iterator)
Converts an iterator that points to a TimeFit object into an iterator that points to the appropriate ...
MarleyTimeGen(const Parameters &p)
const FitParameters & GetFitParameters(int pdg_code) const
Retrieves fit parameters for a specific neutrino type for this time bin.
FitParameters(double Emean, double alpha, double luminosity)
fhicl::OptionalAtom< std::string > namecycle_
double Luminosity() const
Luminosity (erg / s)
std::unique_ptr< TH2D > fSpectrumHist
ROOT TH2D that contains the time-dependent spectrum to use when sampling neutrino times and energies...
fhicl::Atom< std::string > spectrum_file_format_
auto begin(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:573
std::unique_ptr< marley::Event > fEvent
unique_ptr to the current event created by MARLEY
unsigned int fNeutrinosPerEvent
The number of MARLEY neutrino vertices to generate in each art::Event.
then echo File list $list not found else cat $list while read file do echo $file sed s
Definition: file_to_url.sh:60
Stores extra MC truth information that is recorded when generating events using a time-dependent supe...
Collection of configuration parameters for the module.
Sample directly from cross-section weighted spectrum.
TimeGenSamplingMode fSamplingMode
Represents the sampling mode to use when selecting neutrino times and energies.
double fFitEmin
Minimum neutrino energy to consider when using a &quot;fit&quot;-format spectrum file.
do i e
void make_final_timefit(double time)
Helper function that makes a final dummy TimeFit object so that the final real time bin can have a ri...
fhicl::Atom< unsigned int > nu_per_event_
fhicl::OptionalAtom< std::string > pinching_parameter_type_
std::vector< TimeFit > fTimeFits
Vector that contains the fit parameter information for each time bin when using a &quot;fit&quot;-format spectr...
art::ServiceHandle< art::TFileService > tfs
fhicl::Table< evgen::ActiveVolumeVertexSampler::Config > vertex_
void make_nu_emission_histograms() const
Makes ROOT histograms showing the emitted neutrinos in each time bin when using a &quot;fit&quot;-format spectr...
Neutrino energies were sampled uniformly.
art framework interface to geometry description
fhicl::Atom< std::string > spectrum_file_
fhicl::OptionalAtom< double > fit_Emin_
double Emean() const
Mean neutrino energy (MeV)