MARLEY (Model of Argon Reaction Low Energy Yields)  v1.2.0
A Monte Carlo event generator for tens-of-MeV neutrino interactions
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
marley::EventFileReader Class Reference

Object that parses MARLEY output files written in any of the available formats, except for ROOT format. More...

#include <EventFileReader.hh>

Inheritance diagram for marley::EventFileReader:
marley::RootEventFileReader

Public Member Functions

 EventFileReader (const std::string &file_name)
 
double flux_averaged_xsec (bool natural_units=false)
 Returns the flux-averaged total cross section used to produce the events in the file. More...
 
virtual bool next_event (marley::Event &ev)
 Read the next MARLEY event record from the file. More...
 
virtual operator bool () const
 Implicit boolean conversion allows the state of the input stream (or ROOT file) to be tested for readiness to read in another event.
 
EventFileReaderoperator>> (marley::Event &ev)
 Stream operator for reading in the next event.
 

Protected Member Functions

virtual bool deduce_file_format ()
 Helper function that auto-detects which of the available output formats is appropriate for the requested file.
 
void ensure_initialized ()
 This function should be called at the beginning of all public member functions of EventFileReader that interact with data in the file. More...
 
virtual void initialize ()
 Prepares the file for reading the events.
 

Protected Attributes

std::string file_name_
 Name of the file (with any needed path specification) to be read.
 
double flux_avg_tot_xs_ = 0.
 Flux-averaged total cross section (MeV -2) used to produce the events in the file, or zero if that information is not included in a particular format.
 
OutputFile::Format format_
 Format of the output file being read. More...
 
std::ifstream in_
 Input stream used to read from textual output formats.
 
bool initialized_ = false
 Flag that indicates whether initialize() has been called or not. More...
 
marley::JSON json_event_array_
 Used to parse events from JSON-format files.
 
marley::JSON::JSONWrapper< std::deque< marley::JSON > > json_event_array_wrapper_
 Used to iterate over events from JSON-format files.
 
std::deque< marley::JSON >::iterator json_event_iter_
 Iterator to the next JSON event in json_event_array.
 

Detailed Description

Object that parses MARLEY output files written in any of the available formats, except for ROOT format.

For a version of this class that can also handle ROOT-format output files, see the RootEventFileReader class

Member Function Documentation

◆ ensure_initialized()

void marley::EventFileReader::ensure_initialized ( )
protected

This function should be called at the beginning of all public member functions of EventFileReader that interact with data in the file.

It provides necessary initialization as a workaround to calling virtual functions in the constructor

◆ flux_averaged_xsec()

double marley::EventFileReader::flux_averaged_xsec ( bool  natural_units = false)

Returns the flux-averaged total cross section used to produce the events in the file.

For file formats which do not include this information, this function will return zero

Parameters
[in]natural_unitsIf true, then this function will return the flux-averaged total cross section in natural units (MeV -2). If false (default), then 10-42 cm2 will be used.

◆ next_event()

bool marley::EventFileReader::next_event ( marley::Event ev)
virtual

Read the next MARLEY event record from the file.

Parameters
[out]evReference to the object that will be filled with the next event record
Returns
True if reading the next event was successful, or false otherwise. This behavior is designed to be used as a while loop condition for iterating over events in the output file

Reimplemented in marley::RootEventFileReader.

Member Data Documentation

◆ format_

OutputFile::Format marley::EventFileReader::format_
protected

Format of the output file being read.

This format will be determined automatically by deduce_file_format() and does not need to be specified by the user

◆ initialized_

bool marley::EventFileReader::initialized_ = false
protected

Flag that indicates whether initialize() has been called or not.

To avoid problems with using virtual functions in the constructor, we defer nearly all of the initialization to the first call to one of the other public member functions.


The documentation for this class was generated from the following files: