Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAMoleculeEncounterStepper.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// Author: Mathieu Karamitros
28
29// The code is developed in the framework of the ESA AO7146
30//
31// We would be very happy hearing from you, send us your feedback! :)
32//
33// In order for Geant4-DNA to be maintained and still open-source,
34// article citations are crucial.
35// If you use Geant4-DNA chemistry and you publish papers about your software,
36// in addition to the general paper on Geant4-DNA:
37//
38// Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
39//
40// we would be very happy if you could please also cite the following
41// reference papers on chemistry:
42//
43// J. Comput. Phys. 274 (2014) 841-882
44// Prog. Nucl. Sci. Tec. 2 (2011) 503-508
45
46#pragma once
47
49#include "G4KDTreeResult.hh"
50#include "G4ITTrackHolder.hh"
51#include "G4ITReaction.hh"
52
56
57class G4Molecule;
58
59/**
60 * Given a molecule G4DNAMoleculeEncounterStepper will calculate for its possible reactants
61 * what will be the minimum encounter time and the associated molecules.*
62 *
63 * This model includes dynamical time steps as explained in
64 * "Computer-Aided Stochastic Modeling of the Radiolysis of Liquid Water",
65 * V. Michalik, M. Begusová, E. A. Bigildeev,
66 * Radiation Research, Vol. 149, No. 3 (Mar., 1998), pp. 224-236
67 *
68 */
69
71{
72public:
77
78 virtual void Prepare();
79 virtual G4double CalculateStep(const G4Track&, const G4double&);
81
84
85 void SetVerbose(int);
86 // Final time returned when reaction is available in the reaction table = 1
87 // All details = 2
88
89private:
90 void InitializeForNewTrack();
91
92 class Utils;
93 void CheckAndRecordResults(const Utils&,
94#ifdef G4VERBOSE
95 const G4double reactionRange,
96#endif
98
99 G4bool fHasAlreadyReachedNullTime;
100 const G4DNAMolecularReactionTable*& fMolecularReactionTable;
101 G4VDNAReactionModel* fReactionModel;
102 G4ITTrackHolder* fpTrackContainer;
103 G4ITReactionSet* fReactionSet;
104 G4int fVerbose;
105
106 class Utils
107 {
108 public:
109 Utils(const G4Track& tA, const G4MolecularConfiguration* mB);
110 ~Utils() = default;
111
112 G4double GetConstant() const
113 {
114 return fConstant;
115 }
116
117 const G4Track& fpTrackA;
118 const G4MolecularConfiguration* fpMoleculeB;
119 const G4Molecule* fpMoleculeA;
120 G4double fDA;
121 G4double fDB;
122 G4double fConstant;
123 };
124};
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4DNAMoleculeEncounterStepper(const G4DNAMoleculeEncounterStepper &)=delete
virtual G4double CalculateStep(const G4Track &, const G4double &)
G4DNAMoleculeEncounterStepper & operator=(const G4DNAMoleculeEncounterStepper &)=delete
void SetReactionModel(G4VDNAReactionModel *)
virtual G4double CalculateMinTimeStep(G4double, G4double)