MARLEYGen_module.cc
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////
2 /// \file MARLEYGen_module.cc
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 // standard library includes
10 #include <memory>
11 #include <string>
12 
13 // framework includes
18 #include "art_root_io/TFileService.h"
19 #include "fhiclcpp/types/Table.h"
20 
21 // art extensions
22 #include "nurandom/RandomUtils/NuRandomService.h"
23 
24 // LArSoft includes
28 
31 
32 // ROOT includes
33 #include "TTree.h"
34 
35 namespace evgen {
36  class MarleyGen;
37 }
38 
40 
41  public:
42 
43  using Name = fhicl::Name;
45 
46  /// Ignore the marley_parameters FHiCL table during the validation
47  /// step (MARLEY will take care of that by itself)
48  struct KeysToIgnore {
49  std::set< std::string > operator()()
50  {
51  return { "marley_parameters" };
52  }
53  };
54 
55  /// Collection of configuration parameters for the module
56  struct Config {
57 
59  Name("vertex"),
60  Comment("Configuration for selecting the vertex location(s)")
61  };
62 
63  fhicl::Atom<std::string> module_type_ {
64  Name("module_type"),
65  Comment(""),
66  "MARLEYGen" // default value
67  };
68 
69  }; // struct Config
70 
71  // Type to enable FHiCL parameter validation by art
73 
74  // Configuration-checking constructors
75  explicit MarleyGen(const Parameters& p);
76 
77  virtual void produce(art::Event& e) override;
78  virtual void beginRun(art::Run& run) override;
79 
80  virtual void reconfigure(const Parameters& p);
81 
82  private:
83 
84  // Object that provides an interface to the MARLEY event generator
85  std::unique_ptr<evgen::MARLEYHelper> fMarleyHelper;
86 
87  // Algorithm that allows us to sample vertex locations within the active
88  // volume(s) of the detector
89  std::unique_ptr<evgen::ActiveVolumeVertexSampler> fVertexSampler;
90 
91  // unique_ptr to the current event created by MARLEY
92  std::unique_ptr<marley::Event> fEvent;
93 
94  // the MARLEY event TTree
95  TTree* fEventTree;
96 
97  // Run, subrun, and event numbers from the art::Event being processed
98  uint_fast32_t fRunNumber;
99  uint_fast32_t fSubRunNumber;
100  uint_fast32_t fEventNumber;
101 };
102 
103 //------------------------------------------------------------------------------
105  : EDProducer{ p.get_PSet() },
106  fEvent(new marley::Event), fRunNumber(0), fSubRunNumber(0), fEventNumber(0)
107 {
108  // Configure the module (including MARLEY itself) using the FHiCL parameters
109  this->reconfigure( p );
110 
111  // Create a ROOT TTree using the TFileService that will store the MARLEY
112  // event objects (useful for debugging purposes)
114  fEventTree = tfs->make<TTree>("MARLEY_event_tree",
115  "Neutrino events generated by MARLEY");
116  fEventTree->Branch("event", "marley::Event", fEvent.get());
117 
118  // Add branches that give the art::Event run, subrun, and event numbers for
119  // easy match-ups between the MARLEY and art TTrees. All three are recorded
120  // as 32-bit unsigned integers.
121  fEventTree->Branch("run_number", &fRunNumber, "run_number/i");
122  fEventTree->Branch("subrun_number", &fSubRunNumber, "subrun_number/i");
123  fEventTree->Branch("event_number", &fEventNumber, "event_number/i");
124 
125  produces< std::vector<simb::MCTruth> >();
126  produces< sumdata::RunData, art::InRun >();
127 }
128 
129 //------------------------------------------------------------------------------
131 {
133  run.put(std::make_unique<sumdata::RunData>(geo->DetectorName()));
134 }
135 
136 //------------------------------------------------------------------------------
138 {
139  // Get the run, subrun, and event numbers from the current art::Event
140  fRunNumber = e.run();
141  fSubRunNumber = e.subRun();
142  fEventNumber = e.event();
143 
144  std::unique_ptr< std::vector<simb::MCTruth> >
145  truthcol(new std::vector<simb::MCTruth>);
146 
147  // Get the primary vertex location for this event
149  TLorentzVector vertex_pos = fVertexSampler->sample_vertex_pos(*geo);
150 
151  // Create the MCTruth object, and retrieve the marley::Event object
152  // that was generated as it was created
153  simb::MCTruth truth = fMarleyHelper->create_MCTruth(vertex_pos,
154  fEvent.get());
155 
156  // Write the marley::Event object to the event tree
157  fEventTree->Fill();
158 
159  truthcol->push_back(truth);
160 
161  e.put(std::move(truthcol));
162 }
163 
164 //------------------------------------------------------------------------------
166 {
167  const auto& seed_service = art::ServiceHandle<rndm::NuRandomService>();
168  const auto& geom_service = art::ServiceHandle<geo::Geometry const>();
169 
170  // Create a new evgen::ActiveVolumeVertexSampler object based on the current
171  // configuration
172  fVertexSampler = std::make_unique<evgen::ActiveVolumeVertexSampler>(
173  p().vertex_, *seed_service, *geom_service, "MARLEY_Vertex_Sampler");
174 
175  // Create a new marley::Generator object based on the current configuration
176  fhicl::ParameterSet marley_pset = p.get_PSet().get< fhicl::ParameterSet >(
177  "marley_parameters" );
178  fMarleyHelper = std::make_unique<MARLEYHelper>( marley_pset,
179  *seed_service, "MARLEY" );
180 }
181 
virtual void reconfigure(const Parameters &p)
fhicl::Comment Comment
EventNumber_t event() const
Definition: DataViewImpl.cc:85
LArSoft interface to the MARLEY (Model of Argon Reaction Low Energy Yields) supernova neutrino event ...
virtual void beginRun(art::Run &run) override
EDProducer(fhicl::ParameterSet const &pset)
Definition: EDProducer.h:20
ChannelGroupService::Name Name
art framework interface to geometry description
std::unique_ptr< evgen::MARLEYHelper > fMarleyHelper
Definition: Run.h:17
Algorithm that samples vertex locations uniformly within the active volume of a detector. It is fully experiment-agnostic and multi-TPC aware.
const double e
#define DEFINE_ART_MODULE(klass)
Definition: ModuleMacros.h:67
uint_fast32_t fRunNumber
std::string DetectorName() const
Returns a string with the name of the detector, as configured.
def move(depos, offset)
Definition: depos.py:107
uint_fast32_t fSubRunNumber
std::unique_ptr< evgen::ActiveVolumeVertexSampler > fVertexSampler
Collection of configuration parameters for the module.
p
Definition: test.py:223
SubRunNumber_t subRun() const
Definition: DataViewImpl.cc:78
ProductID put(std::unique_ptr< PROD > &&edp, std::string const &instance={})
Definition: DataViewImpl.h:686
RunNumber_t run() const
Definition: DataViewImpl.cc:71
virtual void produce(art::Event &e) override
auto const & get_PSet() const
Definition: ProducerTable.h:47
#define Comment
uint_fast32_t fEventNumber
Event generator information.
Definition: MCTruth.h:32
MarleyGen(const Parameters &p)
LArSoft geometry interface.
Definition: ChannelGeo.h:16
Event Generation using GENIE, cosmics or single particles.
std::unique_ptr< marley::Event > fEvent
std::set< std::string > operator()()