File indexing completed on 2025-02-22 10:31:24
0001
0002
0003
0004
0005
0006
0007
0008 #pragma once
0009
0010 #include <string>
0011
0012 #include "corecel/Macros.hh"
0013 #include "celeritas/ext/RootUniquePtr.hh"
0014 #include "celeritas/phys/Primary.hh"
0015
0016 #include "EventIOInterface.hh"
0017
0018 namespace celeritas
0019 {
0020 class ParticleParams;
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036 class RootEventReader : public EventReaderInterface
0037 {
0038 public:
0039
0040
0041 using SPConstParticles = std::shared_ptr<ParticleParams const>;
0042 using result_type = std::vector<Primary>;
0043
0044
0045
0046 RootEventReader(std::string const& filename, SPConstParticles params);
0047
0048
0049 CELER_DELETE_COPY_MOVE(RootEventReader);
0050
0051
0052 result_type operator()(EventId event_id);
0053
0054
0055 result_type operator()() final;
0056
0057
0058 size_type num_events() const final { return num_events_; }
0059
0060 private:
0061
0062
0063 SPConstParticles params_;
0064 UPExtern<TFile> tfile_;
0065 UPExtern<TTree> ttree_;
0066 std::size_t num_entries_;
0067 size_type num_events_;
0068 std::size_t entry_count_{0};
0069 EventId expected_event_id_{0};
0070 std::vector<std::size_t> event_to_entry_{0};
0071
0072
0073
0074
0075 char const* tree_name() { return "primaries"; }
0076 };
0077
0078
0079 #if !CELERITAS_USE_ROOT
0080 inline RootEventReader::RootEventReader(std::string const&, SPConstParticles)
0081 {
0082 CELER_DISCARD(params_);
0083 CELER_DISCARD(tfile_);
0084 CELER_DISCARD(ttree_);
0085 CELER_DISCARD(num_entries_);
0086 CELER_DISCARD(num_events_);
0087 CELER_DISCARD(entry_count_);
0088 CELER_DISCARD(expected_event_id_);
0089 CELER_DISCARD(event_to_entry_);
0090 CELER_NOT_CONFIGURED("ROOT");
0091 }
0092
0093 inline RootEventReader::result_type RootEventReader::operator()()
0094 {
0095 CELER_ASSERT_UNREACHABLE();
0096 }
0097 #endif
0098
0099
0100 }