Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAIndependentReactionTimeStepper.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// 20/2/2019
28// Author: HoangTRAN
29
30#ifndef G4DNAIndependentReactionTimeStepper_hh
31#define G4DNAIndependentReactionTimeStepper_hh 1
32
34#include "G4KDTreeResult.hh"
35#include "G4IRTUtils.hh"
36#include <memory>
37#include <set>
38#include "G4ITTrackHolder.hh"
39#include "G4ITReaction.hh"
40#include "G4ReferenceCast.hh"
41#include <unordered_map>
42
46class G4Molecule;
47class G4ITReactionSet;
50class G4ITTrackHolder;
51
53{
54 public:
61
62 void Prepare() override;
63 G4double CalculateStep(const G4Track&, const G4double&) override;
65
68
69 std::unique_ptr<G4ITReactionChange> FindReaction(
70 G4ITReactionSet* pReactionSet,
71 G4double& currentStepTime,
72 const G4double globalTime);
73 void SetReactionProcess(G4VITReactionProcess* pReactionProcess);
74 void SetVerbose(G4int);
75
76 private:
77 void InitializeForNewTrack();
78 class Utils;
79 void CheckAndRecordResults(G4double reactionTime, const Utils& utils);
80
81 G4double GetTimeToEncounter(const G4Track& trackA, const G4Track& trackB);
82
83 const G4DNAMolecularReactionTable*& fMolecularReactionTable =
85 G4VDNAReactionModel* fReactionModel = nullptr;
86 G4ITTrackHolder* fpTrackContainer = G4ITTrackHolder::Instance();
88 G4int fVerbose = 0;
90 G4VITReactionProcess* fpReactionProcess = nullptr;
91 std::vector<const G4Track *> fSecondaries;
92 std::unordered_map<G4int, G4ThreeVector> fSampledPositions;
93 std::set<G4int> fCheckedTracks;
94 void InitializeReactions(G4double currentGlobalTime);
95 G4bool fIsInitialized = false;
96 G4double GetNextReactionTime();
97 const G4ITReaction* GetNextReaction();
98
99 class Utils
100 {
101 public:
102 Utils(const G4Track& tA, const G4Track& tB);
103 ~Utils() = default;
104
105 G4Track* fpTrackA{nullptr};
106 G4Track* fpTrackB{nullptr};
107 const G4Molecule* fpMoleculeA{nullptr};
108 const G4Molecule* fpMoleculeB{nullptr};
109 };
110};
111#endif
ReturnType & reference_cast(OriginalType &source)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4double CalculateStep(const G4Track &, const G4double &) override
std::unique_ptr< G4ITReactionChange > FindReaction(G4ITReactionSet *pReactionSet, G4double &currentStepTime, const G4double globalTime)
G4DNAIndependentReactionTimeStepper(const G4DNAIndependentReactionTimeStepper &)=delete
G4double CalculateMinTimeStep(G4double, G4double) override
G4DNAIndependentReactionTimeStepper & operator=(const G4DNAIndependentReactionTimeStepper &)=delete
~G4DNAIndependentReactionTimeStepper() override=default
void SetReactionProcess(G4VITReactionProcess *pReactionProcess)
static G4double GetRCutOff()
Definition G4IRTUtils.cc:39
static G4ITReactionSet * Instance()
static G4ITTrackHolder * Instance()
const G4ITReactionTable * fpReactionTable