MARLEY (Model of Argon Reaction Low Energy Yields)
v1.2.0
A Monte Carlo event generator for tens-of-MeV neutrino interactions
|
Abstract base class for objects that deliver output to a file opened by the marley command-line executable. More...
#include <OutputFile.hh>
Public Types | |
enum class | Format { ROOT , HEPEVT , JSON , ASCII } |
Public Member Functions | |
OutputFile (const std::string &name, const std::string &format, const std::string &mode, bool force=false) | |
virtual int_fast64_t | bytes_written ()=0 |
The number of bytes that have been written during the current MARLEY session to this file. | |
virtual void | close (const marley::JSON &json_config, const marley::Generator &gen, const long num_events)=0 |
Close the output file and perform any necessary cleanup. More... | |
bool | mode_is_resume () const |
const std::string & | name () const |
virtual bool | resume (std::unique_ptr< marley::Generator > &gen, long &num_previous_events)=0 |
Load a marley::Generator object whose configuration and state were saved to the metadata in a ROOT or JSON format output file. More... | |
virtual void | write_event (const marley::Event *event)=0 |
Write a new marley::Event to this output file. More... | |
virtual void | write_flux_avg_tot_xsec (double avg_tot_xsec)=0 |
If needed (for the HEPEVT and ASCII formats), write the flux-averaged total cross section to the file. More... | |
Protected Types | |
enum class | Mode { OVERWRITE , APPEND , RESUME } |
Protected Member Functions | |
bool | check_if_file_exists (const std::string &filename) |
Checks that a given file exists (and is readable) | |
virtual void | open ()=0 |
Open the output file for writing. More... | |
std::unique_ptr< marley::Generator > | restore_generator (const marley::JSON &config) |
Helper function for resume() that instantiates a marley::Generator object given the previous JSON configuration object. | |
virtual void | write_generator_state (const marley::JSON &json_config, const marley::Generator &gen, const long num_events)=0 |
For the file formats that support it, write metadata containing the generator state and configuration to the output file. | |
Protected Attributes | |
bool | force_ |
Format | format_ |
Format to use when writing events to the file. | |
Mode | mode_ |
std::string | name_ |
Name of the file to receive output. | |
Abstract base class for objects that deliver output to a file opened by the marley command-line executable.
|
pure virtual |
Close the output file and perform any necessary cleanup.
This function also saves information about the generator configuration and state to those output file formats that support it.
Implemented in marley::TextOutputFile.
|
protectedpure virtual |
Open the output file for writing.
This function will also make whatever preparations are necessary before the first call to write_event()
Implemented in marley::TextOutputFile.
|
pure virtual |
Load a marley::Generator object whose configuration and state were saved to the metadata in a ROOT or JSON format output file.
[out] | num_previous_events | The number of events previously saved to the output file |
Implemented in marley::TextOutputFile.
|
pure virtual |
Write a new marley::Event to this output file.
If a nullptr is passed to this function, then a marley::Error will be thrown
Implemented in marley::TextOutputFile.
|
pure virtual |
If needed (for the HEPEVT and ASCII formats), write the flux-averaged total cross section to the file.
This function is a no-op unless we're at the beginning of the output stream.
Implemented in marley::RootOutputFile, and marley::TextOutputFile.