remage
Simulation framework for HPGe-based experiments
 
Loading...
Searching...
No Matches
RMGGrabmayrGCReader.hh
1// Copyright (C) 2024 Moritz Neuberger <https://orcid.org/0009-0001-8471-9076>
2// Copyright (C) 2024 Eric Esch <https://orcid.org/0009-0000-4920-9313>
3//
4// This program is free software: you can redistribute it and/or modify it under
5// the terms of the GNU Lesser General Public License as published by the Free
6// Software Foundation, either version 3 of the License, or (at your option) any
7// later version.
8//
9// This program is distributed in the hope that it will be useful, but WITHOUT
10// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more
12// details.
13//
14// You should have received a copy of the GNU Lesser General Public License
15// along with this program. If not, see <https://www.gnu.org/licenses/>.
16
17#ifndef _RMG_GRABMAYR_GC_READER_HH_
18#define _RMG_GRABMAYR_GC_READER_HH_
19
20#include <cmath>
21#include <fstream>
22#include <iostream>
23#include <memory>
24#include <string>
25#include <vector>
26
27#include "G4GenericMessenger.hh"
28#include "globals.hh"
29
30// Modified from WLGDPetersGammaCascadeReader originally contributed by Moritz Neuberger
39 G4int en;
40 G4int ex;
41 G4int m;
42 G4int em;
43 std::vector<G4int> eg;
44};
45
46
54class RMGGrabmayrGCReader {
55 public:
56
58 static RMGGrabmayrGCReader* GetInstance();
59 ~RMGGrabmayrGCReader();
60 // RMGGrabmayrGCReader& operator=(const RMGGrabmayrGCReader&) = delete;
61
63 G4bool IsApplicable(G4int z, G4int a);
65 void CloseFiles();
66
71 GammaCascadeLine GetNextEntry(G4int z, G4int a);
72
73 private:
74
75 static G4ThreadLocal RMGGrabmayrGCReader* instance;
76 RMGGrabmayrGCReader();
77 // std::vector<std::unique_ptr<std::ifstream>> files;
78 // map holding the corresponding file for each isotope
79 std::map<std::pair<G4int, G4int>, std::unique_ptr<std::ifstream>> fCascadeFiles;
80 std::unique_ptr<G4GenericMessenger> fGenericMessenger;
81 G4int fGammaCascadeRandomStartLocation = 0;
82
83 void SetGammaCascadeFile(G4int z, G4int a, G4String file_name);
84 void SetGammaCascadeRandomStartLocation(int answer);
85 void SetStartLocation(std::ifstream& file) const;
86
87 void RandomizeFiles();
88 void DefineCommands();
89
90 // Nested messenger class
91 class GCMessenger : public G4UImessenger {
92 public:
93
94 GCMessenger(RMGGrabmayrGCReader* reader);
95 ~GCMessenger();
96
97 void SetNewValue(G4UIcommand* command, G4String newValues) override;
98
99 private:
100
101 RMGGrabmayrGCReader* fReader;
102 G4UIcommand* fGammaFileCmd;
103
104 void GammaFileCmd(const std::string& parameters);
105 };
106
107 std::unique_ptr<GCMessenger> fUIMessenger;
108};
109
110
111#endif
static RMGGrabmayrGCReader * GetInstance()
Thread-local singleton accessor.
Definition RMGGrabmayrGCReader.cc:28
void CloseFiles()
Close all open cascade files held by this thread.
Definition RMGGrabmayrGCReader.cc:37
GammaCascadeLine GetNextEntry(G4int z, G4int a)
Read and return the next cascade entry for the (Z, A) isotope.
Definition RMGGrabmayrGCReader.cc:52
G4bool IsApplicable(G4int z, G4int a)
Whether a cascade file has been registered for the (Z, A) isotope.
Definition RMGGrabmayrGCReader.cc:45
One pre-computed neutron-capture gamma cascade.
Definition RMGGrabmayrGCReader.hh:38
G4int m
Cascade multiplicity (number of eg entries).
Definition RMGGrabmayrGCReader.hh:41
std::vector< G4int > eg
Photon energies of the cascade [keV].
Definition RMGGrabmayrGCReader.hh:43
G4int en
Neutron kinetic energy bin [keV].
Definition RMGGrabmayrGCReader.hh:39
G4int em
Missing energy not carried by the listed photons [keV].
Definition RMGGrabmayrGCReader.hh:42
G4int ex
Capture-state excitation energy [keV].
Definition RMGGrabmayrGCReader.hh:40