Geant4 11.1.1
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.
30#include "G4Molecule.hh"
31#include "G4ITReactionChange.hh"
32#include "G4IRTUtils.hh"
33#include "G4VDNAHitModel.hh"
35
36#ifndef PrepareState
37# define PrepareState() \
38 G4PolyNucleotideReactionState* _state = \
39 this->GetState<G4PolyNucleotideReactionState>()
40#endif
41
42#ifndef State
43# define State(theXInfo) (_state->theXInfo)
44#endif
45
47 const G4String& aName, G4int verbosityLevel)
49 , fHasAlreadyReachedNullTime(false)
50 , fVerbose(verbosityLevel)
51 , fRCutOff(G4IRTUtils::GetDNADistanceCutOff())
52 , fpDamageModel(nullptr)
53{
55 enableAtRestDoIt = false;
56 enableAlongStepDoIt = false;
57 enablePostStepDoIt = true;
58 fProposesTimeStep = true;
61}
62
64{
65 delete fpDamageModel; // for now, this process handls only one model
66}
67
68G4DNAPolyNucleotideReactionProcess::G4PolyNucleotideReactionState::
69G4PolyNucleotideReactionState()
70{
71 fSampledMinTimeStep = 0;
72 fPreviousTimeAtPreStepPoint = -1;
73}
74
76 const G4Track& trackA, const G4double& /*userTimeStep*/)
77{
79 fHasAlreadyReachedNullTime = false;
80 State(fSampledMinTimeStep) = DBL_MAX;
81 State(theInteractionTimeLeft) = DBL_MAX;
82 State(currentInteractionLength) = -1;
83#ifdef G4VERBOSE
84 if(fVerbose > 1)
85 {
86 auto pMoleculeA = GetMolecule(trackA);
87 G4cout << "________________________________________________________________"
88 "_______"
89 << G4endl;
90 G4cout << "G4DNAPolyNucleotideReactionProcess::CalculateTimleStep"
91 << G4endl;
92 G4cout << "Check done for molecule : " << pMoleculeA->GetName() << " ("
93 << trackA.GetTrackID() << ") " << G4endl;
94 }
95#endif
96 //__________________________________________________________________
97 // Retrieve general informations for making reactions
98 G4double reactionTime =
99 fpDamageModel->CalculateReactionTime(trackA, State(fNodeReactant));
100
101 if(reactionTime < 0)
102 {
103 return DBL_MAX;
104 }
105 State(fSampledMinTimeStep) = reactionTime;
106 State(theInteractionTimeLeft) = State(fSampledMinTimeStep);
107 State(currentInteractionLength) = State(theInteractionTimeLeft);
108
109#ifdef G4VERBOSE
110 if(fVerbose > 1)
111 {
112 G4cout << " theInteractionTimeLeft : " << State(theInteractionTimeLeft)
113 << G4endl;
114 G4cout << " State(fNodeReactant) : " << State(fNodeReactant).index()
115 << G4endl;
116
117 G4cout << "________________________________________________________________"
118 "_______"
119 << G4endl;
120 }
121#endif
122
123 return State(fSampledMinTimeStep);
124}
125
128 const G4Track& track,
129 G4double, // previousStepSize
130 G4ForceCondition* pForceCond)
131{
132 PrepareState();
133 CalculateTimeStep(track);
134 *pForceCond = NotForced;
135
136 G4double previousTimeStep(-1.);
137
138 if(State(fPreviousTimeAtPreStepPoint) != -1)
139 {
140 previousTimeStep =
141 track.GetGlobalTime() - State(fPreviousTimeAtPreStepPoint);
142 }
143
144 State(fPreviousTimeAtPreStepPoint) = track.GetGlobalTime();
145
146 if((fpState->currentInteractionLength <= 0) || (previousTimeStep < 0.0) ||
147 (fpState->theNumberOfInteractionLengthLeft <= 0.0))
148 {
150 }
151 else if(previousTimeStep > 0.0)
152 {
154 previousTimeStep); // TODO need to check this
155 }
156 return State(theInteractionTimeLeft) * -1;
157 // return value * -1;//should regard this !
158}
159
160///////////////////////////////////////////////////////////////////////////////
161
163 const G4Track& track, const G4Step&)
164{
165 PrepareState();
166 auto isReacted = fpDamageModel->DoReaction(
167 track, State(theInteractionTimeLeft), State(fNodeReactant));
168 if(!isReacted)
169 {
170 // no reaction even this process is called
172 return pParticleChange;
173 }
174 fParticleChange.Initialize(track); // To initialise TouchableChange
176 State(fPreviousTimeAtPreStepPoint) = -1;
177 return pParticleChange;
178}
180{
182 G4VITProcess::fpState.reset(new G4PolyNucleotideReactionState());
184}
185
186#undef State
187#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:57
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
Definition: G4Step.hh:62
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
virtual void StartTracking(G4Track *)
Definition: G4VITProcess.cc:85
void SetInstantiateProcessState(G4bool flag)
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4shared_ptr< G4ProcessState > fpState
virtual void ResetNumberOfInteractionLengthLeft()
G4bool fProposesTimeStep
void ProposeTrackStatus(G4TrackStatus status)
virtual void Initialize(const G4Track &)
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:363
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:364
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:87
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:365
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325
#define DBL_MAX
Definition: templates.hh:62