File indexing completed on 2025-01-18 09:58:41
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
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045 #pragma once
0046
0047 #include "G4VMoleculeCounter.hh"
0048 #include <map>
0049 #include <memory>
0050 #include <set>
0051 #include <vector>
0052
0053
0054
0055 namespace G4 {
0056 namespace MoleculeCounter {
0057 struct TimePrecision
0058 {
0059 bool operator()(const double& a, const double& b) const;
0060 static G4ThreadLocal double fPrecision;
0061 };
0062 }
0063 }
0064 using NbMoleculeAgainstTime = std::map<G4double, G4int, G4::MoleculeCounter::TimePrecision>;
0065 using RecordedTimes = std::unique_ptr<std::set<G4double>>;
0066
0067
0068
0069 class G4MoleculeCounter : public G4VMoleculeCounter
0070 {
0071
0072 public:
0073 using ReactantList = std::vector<Reactant*>;
0074 using CounterMapType = std::map<Reactant*, NbMoleculeAgainstTime>;
0075 using RecordedMolecules = std::unique_ptr<ReactantList>;
0076
0077 static G4MoleculeCounter* Instance();
0078
0079 void Initialize() override;
0080 void ResetCounter() override;
0081
0082
0083 void DontRegister(const G4MoleculeDefinition*) override;
0084 bool IsRegistered(const G4MoleculeDefinition*) override;
0085 void RegisterAll() override;
0086
0087
0088
0089 int GetNMoleculesAtTime(Reactant* molecule, double time);
0090 const NbMoleculeAgainstTime& GetNbMoleculeAgainstTime(Reactant* molecule);
0091
0092 RecordedMolecules GetRecordedMolecules();
0093 RecordedTimes GetRecordedTimes();
0094
0095 void SetVerbose(G4int);
0096 G4int GetVerbose();
0097
0098
0099 static void SetTimeSlice(double);
0100
0101 void Dump();
0102
0103 G4bool IsTimeCheckedForConsistency() const;
0104 void CheckTimeForConsistency(G4bool flag);
0105
0106 #ifdef MOLECULE_COUNTER_TESTING
0107 public:
0108 #else
0109 protected:
0110 #endif
0111 void AddAMoleculeAtTime(Reactant*,
0112 G4double time,
0113 const G4ThreeVector* position = nullptr,
0114 int number = 1) override;
0115 void RemoveAMoleculeAtTime(Reactant*,
0116 G4double time,
0117 const G4ThreeVector* position = nullptr,
0118 int number = 1) override;
0119
0120
0121 protected:
0122 G4bool SearchTimeMap(Reactant* molecule);
0123 int SearchUpperBoundTime(double time, bool sameTypeOfMolecule);
0124
0125 protected:
0126 G4MoleculeCounter();
0127 ~G4MoleculeCounter() override;
0128
0129 CounterMapType fCounterMap;
0130 std::map<const G4MoleculeDefinition*, G4bool> fDontRegister;
0131
0132 G4int fVerbose;
0133 G4bool fCheckTimeIsConsistentWithScheduler;
0134
0135 struct Search
0136 {
0137 Search()
0138 {
0139 fLowerBoundSet = false;
0140 }
0141 CounterMapType::iterator fLastMoleculeSearched;
0142 NbMoleculeAgainstTime::iterator fLowerBoundTime;
0143 bool fLowerBoundSet;
0144 };
0145
0146 std::unique_ptr<Search> fpLastSearch;
0147
0148 friend class G4Molecule;
0149 friend class G4VMoleculeCounter;
0150 };