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
0028
0029
0030 #ifndef G4DNAIndependentReactionTimeStepper_hh
0031 #define G4DNAIndependentReactionTimeStepper_hh 1
0032
0033 #include "G4VITTimeStepComputer.hh"
0034 #include "G4KDTreeResult.hh"
0035 #include "G4IRTUtils.hh"
0036 #include <memory>
0037 #include <set>
0038 #include "G4ITTrackHolder.hh"
0039 #include "G4ITReaction.hh"
0040 #include "G4ReferenceCast.hh"
0041
0042 class G4VDNAReactionModel;
0043 class G4DNAMolecularReactionTable;
0044 class G4MolecularConfiguration;
0045 class G4Molecule;
0046 class G4ITReactionSet;
0047 class G4ITReactionChange;
0048 class G4VITReactionProcess;
0049 class G4ITTrackHolder;
0050
0051 class G4DNAIndependentReactionTimeStepper : public G4VITTimeStepComputer
0052 {
0053 public:
0054 G4DNAIndependentReactionTimeStepper();
0055 ~G4DNAIndependentReactionTimeStepper() override = default;
0056 G4DNAIndependentReactionTimeStepper(
0057 const G4DNAIndependentReactionTimeStepper&) = delete;
0058 G4DNAIndependentReactionTimeStepper& operator =(
0059 const G4DNAIndependentReactionTimeStepper&) = delete;
0060
0061 void Prepare() override;
0062 G4double CalculateStep(const G4Track&, const G4double&) override;
0063 G4double CalculateMinTimeStep(G4double, G4double) override;
0064
0065 void SetReactionModel(G4VDNAReactionModel*);
0066 G4VDNAReactionModel* GetReactionModel();
0067
0068 std::unique_ptr<G4ITReactionChange> FindReaction(
0069 G4ITReactionSet* pReactionSet, const G4double& currentStepTime = 0,
0070 const G4double& previousStepTime = 0,
0071 const G4bool& reachedUserStepTimeLimit = false);
0072 void SetReactionProcess(G4VITReactionProcess* pReactionProcess);
0073 void SetVerbose(G4int);
0074
0075 private:
0076 void InitializeForNewTrack();
0077 class Utils;
0078 void CheckAndRecordResults(const Utils& utils);
0079
0080 G4double GetTimeToEncounter(const G4Track& trackA, const G4Track& trackB);
0081
0082 G4bool fHasAlreadyReachedNullTime = false;
0083 const G4DNAMolecularReactionTable*& fMolecularReactionTable =
0084 reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable);
0085 G4VDNAReactionModel* fReactionModel = nullptr;
0086 G4ITTrackHolder* fpTrackContainer = G4ITTrackHolder::Instance();
0087 G4ITReactionSet* fReactionSet = G4ITReactionSet::Instance();
0088 G4int fVerbose = 0;
0089 G4double fRCutOff = G4IRTUtils::GetRCutOff();
0090 G4VITReactionProcess* fpReactionProcess = nullptr;
0091 std::map<G4int, G4ThreeVector> fSampledPositions;
0092 std::set<G4int> fCheckedTracks;
0093
0094 class Utils
0095 {
0096 public:
0097 Utils(const G4Track& tA, const G4Track& tB);
0098 ~Utils() = default;
0099 const G4Track& fTrackA;
0100 const G4Track& fTrackB;
0101 const G4Molecule* fpMoleculeA;
0102 const G4Molecule* fpMoleculeB;
0103 };
0104 };
0105 #endif