Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAPolyNucleotideReactionProcess.cc
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/// file: G4DNAPolyNucleotideReactionProcess.cc
27/// brief: This file handls reaction process with DNA geometry.
29
31#include "G4IRTUtils.hh"
32#include "G4ITReactionChange.hh"
34#include "G4Molecule.hh"
35#include "G4VDNAHitModel.hh"
36
37#include <memory>
38
39#ifndef PrepareState
40# define PrepareState() \
41 G4PolyNucleotideReactionState* _state = \
42 this->GetState<G4PolyNucleotideReactionState>()
43#endif
44
45#ifndef State
46# define State(theXInfo) (_state->theXInfo)
47#endif
48
50 const G4String& aName, G4int verbosityLevel)
52 ,
53 fVerbose(verbosityLevel)
54 , fRCutOff(G4IRTUtils::GetDNADistanceCutOff())
55{
57 enableAtRestDoIt = false;
58 enableAlongStepDoIt = false;
59 enablePostStepDoIt = true;
60 fProposesTimeStep = true;
63}
64
66{
67 delete fpDamageModel; // for now, this process handls only one model
68}
69
70G4DNAPolyNucleotideReactionProcess::G4PolyNucleotideReactionState::
71G4PolyNucleotideReactionState()
72{
73 fSampledMinTimeStep = 0;
74 fPreviousTimeAtPreStepPoint = -1;
75}
76
78 const G4Track& trackA, const G4double& /*userTimeStep*/)
79{
81 fHasAlreadyReachedNullTime = false;
82 State(fSampledMinTimeStep) = DBL_MAX;
83 State(theInteractionTimeLeft) = DBL_MAX;
84 State(currentInteractionLength) = -1;
85#ifdef G4VERBOSE
86 if(fVerbose > 1)
87 {
88 auto pMoleculeA = GetMolecule(trackA);
89 G4cout << "________________________________________________________________"
90 "_______"
91 << G4endl;
92 G4cout << "G4DNAPolyNucleotideReactionProcess::CalculateTimleStep"
93 << G4endl;
94 G4cout << "Check done for molecule : " << pMoleculeA->GetName() << " ("
95 << trackA.GetTrackID() << ") " << G4endl;
96 }
97#endif
98 //__________________________________________________________________
99 // Retrieve general informations for making reactions
100 G4double reactionTime =
101 fpDamageModel->CalculateReactionTime(trackA, State(fNodeReactant));
102
103 if(reactionTime < 0)
104 {
105 return DBL_MAX;
106 }
107 State(fSampledMinTimeStep) = reactionTime;
108 State(theInteractionTimeLeft) = State(fSampledMinTimeStep);
109 State(currentInteractionLength) = State(theInteractionTimeLeft);
110
111#ifdef G4VERBOSE
112 if(fVerbose > 1)
113 {
114 G4cout << " theInteractionTimeLeft : " << State(theInteractionTimeLeft)
115 << G4endl;
116 G4cout << " State(fNodeReactant) : " << State(fNodeReactant).index()
117 << G4endl;
118
119 G4cout << "________________________________________________________________"
120 "_______"
121 << G4endl;
122 }
123#endif
124
125 return State(fSampledMinTimeStep);
126}
127
130 const G4Track& track,
131 G4double, // previousStepSize
132 G4ForceCondition* pForceCond)
133{
134 PrepareState();
135 CalculateTimeStep(track);
136 *pForceCond = NotForced;
137
138 G4double previousTimeStep(-1.);
139
140 if(State(fPreviousTimeAtPreStepPoint) != -1)
141 {
142 previousTimeStep =
143 track.GetGlobalTime() - State(fPreviousTimeAtPreStepPoint);
144 }
145
146 State(fPreviousTimeAtPreStepPoint) = track.GetGlobalTime();
147
148 if((fpState->currentInteractionLength <= 0) || (previousTimeStep < 0.0) ||
149 (fpState->theNumberOfInteractionLengthLeft <= 0.0))
150 {
152 }
153 else if(previousTimeStep > 0.0)
154 {
156 previousTimeStep); // TODO need to check this
157 }
158 return State(theInteractionTimeLeft) * -1;
159 // return value * -1;//should regard this !
160}
161
162///////////////////////////////////////////////////////////////////////////////
163
165 const G4Track& track, const G4Step&)
166{
167 PrepareState();
168 auto isReacted = fpDamageModel->DoReaction(
169 track, State(theInteractionTimeLeft), State(fNodeReactant));
170 if(!isReacted)
171 {
172 // no reaction even this process is called
174 return pParticleChange;
175 }
176 fParticleChange.Initialize(track); // To initialise TouchableChange
178 State(fPreviousTimeAtPreStepPoint) = -1;
179 return pParticleChange;
180}
182{
184 G4VITProcess::fpState = std::make_shared<G4PolyNucleotideReactionState>();
186}
187
188#undef State
189#undef PrepareState
G4ForceCondition
@ NotForced
#define State(X)
#define PrepareState()
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:74
@ fUserDefined
@ fStopAndKill
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double, G4ForceCondition *pForceCond) override
G4DNAPolyNucleotideReactionProcess(const G4String &aName="DNAStaticMoleculeReactionProcess", G4int verbosityLevel=0)
G4double CalculateTimeStep(const G4Track &trackA, const G4double &userTimeStep=0)
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &) override
G4int GetTrackID() const
G4double GetGlobalTime() const
virtual G4bool DoReaction(const G4Track &track, const G4double &, const DNANode &)=0
virtual G4double CalculateReactionTime(const G4Track &trackA, DNANode &)=0
void SetInstantiateProcessState(G4bool flag)
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4shared_ptr< G4ProcessState > fpState
void StartTracking(G4Track *) override
G4bool fProposesTimeStep
void ResetNumberOfInteractionLengthLeft() override
void ProposeTrackStatus(G4TrackStatus status)
virtual void Initialize(const G4Track &)
G4bool enableAtRestDoIt
G4bool enableAlongStepDoIt
void SetProcessSubType(G4int)
virtual void StartTracking(G4Track *)
Definition G4VProcess.cc:87
G4bool enablePostStepDoIt
G4VParticleChange * pParticleChange
#define DBL_MAX
Definition templates.hh:62