Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAIRTMoleculeEncounterStepper.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/*
28 * G4DNAIRTMoleculeEncounterStepper.hh
29 *
30 * Created on: Jul 23, 2019
31 * Author: W. G. Shin
32 * J. Ramos-Mendez and B. Faddegon
33*/
34
35#pragma once
36
38#include "G4KDTreeResult.hh"
39#include "G4ITReaction.hh"
40#include "G4ITTrackHolder.hh"
41
45
46class G4Molecule;
47
48/**
49 * Given a molecule G4DNAIRTMoleculeEncounterStepper will calculate for its possible reactants
50 * what will be the minimum encounter time and the associated molecules.*
51 *
52 * This model includes dynamical time steps as explained in
53 * "Computer-Aided Stochastic Modeling of the Radiolysis of Liquid Water",
54 * V. Michalik, M. Begusová, E. A. Bigildeev,
55 * Radiation Research, Vol. 149, No. 3 (Mar., 1998), pp. 224-236
56 *
57 */
58
60{
61public:
66
67 void Prepare() override;
68 G4double CalculateStep(const G4Track&, const G4double&) override;
70
73
74 void SetVerbose(int);
75 // Final time returned when reaction is available in the reaction table = 1
76 // All details = 2
77
78private:
79 void InitializeForNewTrack();
80
81 class Utils;
82 void CheckAndRecordResults(const Utils&,
83#ifdef G4VERBOSE
84 const G4double reactionRange,
85#endif
87
88 G4bool fHasAlreadyReachedNullTime{false};
89
90 const G4DNAMolecularReactionTable*& fMolecularReactionTable;
91 G4VDNAReactionModel* fReactionModel{nullptr};
92 G4ITReactionSet* fReactionSet;
93 G4ITTrackHolder* fpTrackContainer;
94 G4int fVerbose{0};
95
96 class Utils
97 {
98 public:
99 Utils(const G4Track& tA, const G4MolecularConfiguration* mB);
100 ~Utils() = default;
101
102 G4double GetConstant() const
103 {
104 return fConstant;
105 }
106
107 const G4Track& fpTrackA;
108 const G4MolecularConfiguration* fpMoleculeB;
109 const G4Molecule* fpMoleculeA;
110 G4double fDA;
111 G4double fDB;
112 G4double fConstant;
113 };
114};
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4DNAIRTMoleculeEncounterStepper & operator=(const G4DNAIRTMoleculeEncounterStepper &)=delete
G4DNAIRTMoleculeEncounterStepper(const G4DNAIRTMoleculeEncounterStepper &)=delete
G4double CalculateMinTimeStep(G4double, G4double) override
G4double CalculateStep(const G4Track &, const G4double &) override