File indexing completed on 2025-01-18 09:58:08
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027 #ifndef G4DNAEventSet_hh
0028 #define G4DNAEventSet_hh 1
0029 #include <list>
0030 #include <map>
0031 #include <G4memory.hh>
0032 #include <set>
0033 #include "G4DNAMesh.hh"
0034 #include <variant>
0035 class G4DNAMolecularReactionTable;
0036 class G4DNAMolecularReactionData;
0037 class Event
0038 {
0039 public:
0040 using Index = G4VDNAMesh::Index;
0041 using MolType = const G4MolecularConfiguration*;
0042 using JumpingData = std::pair<MolType, Index>;
0043 using ReactionData = const G4DNAMolecularReactionData;
0044
0045
0046
0047 using Data = std::pair<std::unique_ptr<JumpingData>, ReactionData*>;
0048
0049 Event(const G4double& time, const Index& index, ReactionData*);
0050 Event(const G4double& time, const Index& index, std::unique_ptr<JumpingData>&&);
0051
0052 virtual ~Event();
0053 inline G4double GetTime() const { return fTimeStep; }
0054 inline Index GetIndex() const { return fIndex; }
0055 void PrintEvent() const;
0056 JumpingData* GetJumpingData() const { return std::get<0>(fData).get(); }
0057 ReactionData* GetReactionData() const { return std::get<1>(fData); }
0058
0059 private:
0060 G4double fTimeStep;
0061 Index fIndex;
0062 Data fData;
0063 };
0064
0065 struct comparatorEventSet
0066 {
0067 G4bool operator()(std::unique_ptr<Event> const& rhs,
0068 std::unique_ptr<Event> const& lhs) const;
0069 };
0070
0071 class IEventSet
0072 {
0073 public:
0074 IEventSet() = default;
0075 ~IEventSet() = default;
0076 };
0077
0078 class G4DNAEventSet : public IEventSet
0079 {
0080 public:
0081 using Index = G4VDNAMesh::Index;
0082 using EventSet = std::set<std::unique_ptr<Event>, comparatorEventSet>;
0083 using EventMap = std::unordered_map<Index, EventSet::iterator,G4VDNAMesh::hashFunc>;
0084 G4DNAEventSet();
0085 virtual ~G4DNAEventSet();
0086
0087 void CreateEvent(const G4double& time, const Index& index,
0088 Event::ReactionData* pReactionData);
0089 void CreateEvent(const G4double& time, const Index& index,
0090 std::unique_ptr<Event::JumpingData> jum);
0091
0092 void AddEvent(std::unique_ptr<Event> pEvent);
0093 void RemoveEventSet()
0094 {
0095 fEventSet.clear();
0096 fEventMap.clear();
0097 }
0098 void RemoveEventOfVoxel(const Index& key);
0099
0100 EventSet::iterator end() { return fEventSet.end(); }
0101 EventSet::iterator begin() { return fEventSet.begin(); }
0102
0103 EventSet::reverse_iterator rend() { return fEventSet.rend(); }
0104 EventSet::reverse_iterator rbegin() { return fEventSet.rbegin(); }
0105 EventSet::const_iterator end() const { return fEventSet.end(); }
0106 EventSet::const_iterator begin() const { return fEventSet.begin(); }
0107 size_t size() { return fEventSet.size(); }
0108 G4bool Empty() { return fEventSet.empty(); }
0109 void RemoveEvent(EventSet::iterator iter);
0110 void PrintEventSet();
0111 private:
0112 EventSet fEventSet;
0113 EventMap fEventMap;
0114 };
0115 #endif