Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNASecondOrderReaction.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//
28
29#include "G4DNADamage.hh"
33#include "G4Molecule.hh"
34#include "G4SystemOfUnits.hh"
36#include "G4UnitsTable.hh"
37
38#include <G4VScheduler.hh>
39
40#include <memory>
41
42#ifndef State
43#define State(theXInfo) (GetState<SecondOrderReactionState>()->theXInfo)
44#endif
45
46void G4DNASecondOrderReaction::Create()
47{
49 enableAtRestDoIt = false;
50 enableAlongStepDoIt = false;
51 enablePostStepDoIt = true;
52
54
56 // meaning G4DNASecondOrderReaction contains a class inheriting from G4ProcessState
57
58 fIsInitialized = false;
60 fpMaterial = nullptr;
61 fReactionRate = -1.;
62 fConcentration = -1.;
64 fProposesTimeStep = true;
66 fpMoleculeDensity = nullptr;
67
68 verboseLevel = 0;
69}
70
72 G4VITProcess(aName,type)
73{
74 Create();
75}
76
82
84= default;
85
87{
88 if (this == &rhs) return *this; // handle self assignment
89
90 //assignment operator
91 return *this;
92}
93
99
106
107void
109{
111 G4VITProcess::fpState = std::make_shared<SecondOrderReactionState>();
113}
114
115void
117 const G4Material* mat, double reactionRate)
118{
120 {
121 G4ExceptionDescription exceptionDescription ;
122 exceptionDescription << "G4DNASecondOrderReaction was already initialised. ";
123 exceptionDescription << "You cannot set a reaction after initialisation.";
124 G4Exception("G4DNASecondOrderReaction::SetReaction","G4DNASecondOrderReaction001",
125 FatalErrorInArgument,exceptionDescription);
126 }
127 fpMolecularConfiguration = molConf;
128 fpMaterial = mat;
129 fReactionRate = reactionRate;
130}
131
133 G4double /*previousStepSize*/,
134 G4ForceCondition* pForceCond)
135{
136 // G4cout << "G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength" << G4endl;
137 // G4cout << "For reaction : " << fpMaterial->GetName() << " + " << fpMolecularConfiguration->GetName() << G4endl;
138
139 //_______________________________________________________________________
140 // Check whether the track is in the good material (maybe composite material)
141 const G4Material* material = track.GetMaterial();
142
143 G4Molecule* mol = GetMolecule(track);
144 if(mol == nullptr) return DBL_MAX;
146 {
147 // G4cout <<"mol->GetMolecularConfiguration() != fpMolecularConfiguration" << G4endl;
148 return DBL_MAX;
149 }
150
151 G4double molDensity = (*fpMoleculeDensity)[material->GetIndex()];
152
153 if(molDensity == 0.0) // ie : not found
154 {
155 if(State(fIsInGoodMaterial))
156 {
158 //State(fPreviousTimeAtPreStepPoint) = -1;
159 State(fIsInGoodMaterial) = false;
160 }
161
162 // G4cout << " Material " << fpMaterial->GetName() << " not found "
163 // <<" | name of current material : " << material->GetName()
164 // << G4endl;
165
166 return DBL_MAX; // Becareful return here !!
167 }
168
169 // G4cout << " Va calculer le temps d'interaction " << G4endl;
170
171 State(fIsInGoodMaterial) = true;
172
173 // fConcentration = molDensity/fMolarMassOfMaterial;
174 fConcentration = molDensity/CLHEP::Avogadro;
175 // G4cout << "Concentration : " << fConcentration / (g/mole)<< G4endl;
176
177 //_______________________________________________________________________
178 // Either initialize the lapse of time left
179 // meaning
180 // => the track enters for the first time in the material
181 // or substract the previous time step to the previously calculated lapse
182 // of time left
183 // meaning
184 // => the track has not left this material since the previous call
185 G4double previousTimeStep(-1.);
186
187 if(State(fPreviousTimeAtPreStepPoint) != -1)
188 {
189 previousTimeStep = track.GetGlobalTime() -
190 State(fPreviousTimeAtPreStepPoint) ;
191 }
192
193 State(fPreviousTimeAtPreStepPoint) = track.GetGlobalTime();
194
195 // condition is set to "Not Forced"
196 *pForceCond = NotForced;
197
198 if (
199 (previousTimeStep < 0.0) ||
200 (fpState->theNumberOfInteractionLengthLeft<=0.0)) {
201 // beggining of tracking (or just after DoIt of this process)
203 } else if ( previousTimeStep > 0.0) {
204 // get mean free path
205 // subtract NumberOfInteractionLengthLeft
206 SubtractNumberOfInteractionLengthLeft( previousTimeStep );
207 } else {
208 // zero time step
209 // Force trigerring the process
210 //*pForceCond = Forced;
211 }
212
213 fpState->currentInteractionLength = 1/(fReactionRate*fConcentration);
214
215 // G4cout << "fpState->currentInteractionLength = "
216 // << fpState->currentInteractionLength << G4endl;
217
218 G4double value;
219 if (fpState->currentInteractionLength <DBL_MAX) {
220 value = fpState->theNumberOfInteractionLengthLeft
221 * (fpState->currentInteractionLength);
222 } else {
223 value = DBL_MAX;
224 }
225#ifdef G4VERBOSE
226 if (verboseLevel>2) {
227 G4cout << "G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength ";
228 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
229 track.GetDynamicParticle()->DumpInfo();
230 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
231 G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
232 }
233#endif
234
235// G4cout << "currentInteractionLength : " << fpState->currentInteractionLength << G4endl;
236// G4cout << "Material : " << fpMaterial->GetName()
237// << "ID: " << track.GetTrackID()
238// << " Returned time : " << G4BestUnit(value,"Time") << G4endl;
239
240 if(value < fReturnedValue)
241 fReturnedValue = value;
242
243 return value*-1;
244 // multiple by -1 to indicate to the tracking system that we are returning a time
245}
246
248{
249 G4Molecule* molecule = GetMolecule(track);
250#ifdef G4VERBOSE
251 if(verboseLevel > 1)
252 {
253 G4cout << "___________" << G4endl;
254 G4cout << ">>> Beginning of G4DNASecondOrderReaction verbose" << G4endl;
255 G4cout << ">>> Returned value : " << G4BestUnit(fReturnedValue,"Time") << G4endl;
256 G4cout << ">>> Time Step : " << G4BestUnit(G4VScheduler::Instance()->GetTimeStep(),"Time") << G4endl;
257 G4cout << ">>> Reaction : " << molecule->GetName() << " + " << fpMaterial->GetName() << G4endl;
258 G4cout << ">>> End of G4DNASecondOrderReaction verbose <<<" << G4endl;
259 }
260#endif
265 State(fPreviousTimeAtPreStepPoint) = -1;
266 return &fParticleChange;
267}
268
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ForceCondition
@ NotForced
#define State(X)
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:74
G4ProcessType
#define G4BestUnit(a, b)
@ fStopAndKill
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
virtual void AddIndirectDamage(const G4String &baseName, const G4Molecule *molecule, const G4ThreeVector &position, G4double time)
static G4DNADamage * Instance()
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
G4DNASecondOrderReaction & operator=(const G4DNASecondOrderReaction &)
const std::vector< double > * fpMoleculeDensity
~G4DNASecondOrderReaction() override
const G4MolecularConfiguration * fpMolecularConfiguration
G4DNASecondOrderReaction(const G4String &aName="DNASecondOrderReaction", G4ProcessType type=fDecay)
void SetReaction(const G4MolecularConfiguration *, const G4Material *, double)
void StartTracking(G4Track *) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void DumpInfo(G4int mode=0) const
G4double GetMassOfMolecule() const
std::size_t GetIndex() const
const G4String & GetName() const
const G4MolecularConfiguration * GetMolecularConfiguration() const
const G4String & GetName() const override
void Initialize(const G4Track &) override
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
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)
G4bool enableAtRestDoIt
G4int verboseLevel
G4bool enableAlongStepDoIt
void SetProcessSubType(G4int)
virtual void StartTracking(G4Track *)
Definition G4VProcess.cc:87
G4bool enablePostStepDoIt
G4VParticleChange * pParticleChange
const G4String & GetProcessName() const
static G4VScheduler * Instance()
#define DBL_MAX
Definition templates.hh:62