MARLEY (Model of Argon Reaction Low Energy Yields)  v1.2.0
A Monte Carlo event generator for tens-of-MeV neutrino interactions
RootOutputFile.hh
1 //
5 // This file is part of MARLEY (Model of Argon Reaction Low Energy Yields)
6 //
7 // MARLEY is free software: you can redistribute it and/or modify it under the
8 // terms of version 3 of the GNU General Public License as published by the
9 // Free Software Foundation.
10 //
11 // For the full text of the license please see COPYING or
12 // visit http://opensource.org/licenses/GPL-3.0
13 //
14 // Please respect the MCnet academic usage guidelines. See GUIDELINES
15 // or visit https://www.montecarlonet.org/GUIDELINES for details.
16 
17 #pragma once
18 
19 // standard library includes
20 #include "marley/OutputFile.hh"
21 
22 // ROOT includes
23 #include "TFile.h"
24 #include "TTree.h"
25 
26 namespace marley {
27 
28  class RootOutputFile : public OutputFile {
29  public:
30 
31  RootOutputFile(const std::string& name, const std::string& format,
32  const std::string& mode, bool force = false);
33 
34  virtual ~RootOutputFile() = default;
35 
36  inline int_fast64_t bytes_written() override {
37  return file_->GetBytesWritten();
38  }
39 
40  // This function is a no-op for the ROOT format (we will take
41  // care of saving this information when we write the generator
42  // state variables)
43  inline virtual void write_flux_avg_tot_xsec(double /*avg_tot_xsec*/)
44  override {}
45 
46  private:
47 
48  virtual void open() override;
49 
50  void load_marley_headers();
51 
52  // Write the internal state string of the random number generator to
53  // disk, as well as the generator's JSON configuration. This will allow
54  // MARLEY to resume event generation from where it left off with no
55  // loss of consistency. This trick is based on
56  // http://tinyurl.com/hb7rqsj
57  void write_generator_state(const marley::JSON& json_config,
58  const marley::Generator& gen, const long /*num_events*/) override;
59 
60  // Clean up once we're done with the ROOT file. Save some information
61  // about the current generator configuration as we clean things up.
62  virtual void close(const marley::JSON& json_config,
63  const marley::Generator& gen, const long dummy) override;
64 
65  virtual bool resume(std::unique_ptr<marley::Generator>& gen,
66  long& num_previous_events) override;
67 
68  virtual void write_event(const marley::Event* event) override;
69 
70  // TFile object that will be used to access the ROOT file
71  std::unique_ptr<TFile> file_ = nullptr;
72 
73  // This is a bare pointer, but ROOT will associate it with file_, so we
74  // don't want to delete it ourselves or let a smart pointer do it.
75  TTree* tree_ = nullptr;
76  };
77 
78 }
Container for ingoing and outgoing momentum 4-vectors from a reaction.
Definition: Event.hh:66
The MARLEY Event generator.
Definition: Generator.hh:42
Definition: JSON.hh:62
Abstract base class for objects that deliver output to a file opened by the marley command-line execu...
Definition: OutputFile.hh:29
Definition: RootOutputFile.hh:28
int_fast64_t bytes_written() override
The number of bytes that have been written during the current MARLEY session to this file.
Definition: RootOutputFile.hh:36
virtual void write_flux_avg_tot_xsec(double) override
If needed (for the HEPEVT and ASCII formats), write the flux-averaged total cross section to the file...
Definition: RootOutputFile.hh:43