Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id: G4DNASecondOrderReaction.cc 65695 2012-11-27 11:39:12Z gunter $
27//
29#include "G4SystemOfUnits.hh"
30#include "G4Molecule.hh"
33#include "G4DNADamages.hh"
34#include "G4UnitsTable.hh"
35#include "G4ITTrackHolder.hh"
36
37void G4DNASecondOrderReaction::Create()
38{
40 enableAtRestDoIt = false;
41 enableAlongStepDoIt = false;
42 enablePostStepDoIt = true;
43
45
47 // meaning G4DNASecondOrderReaction contains a class inheriting from G4ProcessState
48
49 fIsInitialized = false;
51 fpMaterial = 0;
52 fReactionRate = -1.;
53 fConcentration = -1.;
55 fProposesTimeStep = true;
58
59 verboseLevel = 0;
60}
61
63 G4VITProcess(aName,type),
64 InitProcessState(fpSecondOrderReactionState, fpState)
65{
66 Create();
67}
68
70 G4VITProcess(rhs),
71 InitProcessState(fpSecondOrderReactionState, fpState)
72{
73 Create();
74}
75
77{
78 ;
79}
81{
82 if (this == &rhs) return *this; // handle self assignment
83
84 //assignment operator
85 return *this;
86}
87
89{
91 fIsInGoodMaterial = false;
92}
93
95{
97 fMolarMassOfMaterial = fpMaterial->GetMassOfMolecule()*CLHEP::Avogadro*1e3;
98 fIsInitialized = true;
99}
100
101void
103{
107}
108
109void
111 const G4Material* mat, double reactionRate)
112{
114 {
115 G4ExceptionDescription exceptionDescription ;
116 exceptionDescription << "G4DNASecondOrderReaction was already initialised. ";
117 exceptionDescription << "You cannot set a reaction after initialisation.";
118 G4Exception("G4DNASecondOrderReaction::SetReaction","G4DNASecondOrderReaction001",
119 FatalErrorInArgument,exceptionDescription);
120 }
121 fpMolecularConfiguration = molConf;
122 fpMaterial = mat;
123 fReactionRate = reactionRate;
124}
125
127 G4double /*previousStepSize*/,
128 G4ForceCondition* pForceCond)
129{
130// G4cout << "G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength" << G4endl;
131// G4cout << "For reaction : " << fpMaterial->GetName() << " + " << fpMolecularConfiguration->GetName() << G4endl;
132
133 //_______________________________________________________________________
134 // Check whether the track is in the good material (maybe composite material)
135 const G4Material* material = track.GetMaterial();
136
137 G4Molecule* mol = GetMolecule(track);
138 if(!mol) return DBL_MAX;
140 {
141// G4cout <<"mol->GetMolecularConfiguration() != fpMolecularConfiguration" << G4endl;
142 return DBL_MAX;
143 }
144
145 G4double molDensity = (*fpMoleculeDensity)[material->GetIndex()];
146
147 if(molDensity == 0.0) // ie : not found
148 {
150 {
153 }
154
155// G4cout << " Material " << fpMaterial->GetName() << " not found "
156// <<" | name of current material : " << material->GetName()
157// << G4endl;
158
159 return DBL_MAX; // Becareful return here !!
160 }
161
162// G4cout << " Va calculer le temps d'interaction " << G4endl;
163
165
166// fConcentration = molDensity/fMolarMassOfMaterial;
167 fConcentration = molDensity/CLHEP::Avogadro;
168
169// G4cout << "Concentration : " << fConcentration / (g/mole)<< G4endl;
170
171 //_______________________________________________________________________
172 // Either initialize the lapse of time left
173 // meaning => the track enters for the first time in the material
174 // or substract the previous time step to the previously calculated lapse of time left
175 // meaning => the track has not left this material since the previous call
176
177 G4double previousTimeStep(-1.);
178
179 if(track.GetCurrentStepNumber() > 0)
181
183
184 if ( (previousTimeStep < 0.0) || (fpState->theNumberOfInteractionLengthLeft<=0.0)) {
185 // beggining of tracking (or just after DoIt of this process)
187 } else if ( previousTimeStep > 0.0) {
188 // subtract NumberOfInteractionLengthLeft
189 SubtractNumberOfInteractionLengthLeft(previousTimeStep );
190 } else {
191 // zero step
192 // DO NOTHING
193 }
194
195 // condition is set to "Not Forced"
196 *pForceCond = NotForced;
197
198 // get mean free path
200
201 G4double value;
204 } else {
205 value = DBL_MAX;
206 }
207#ifdef G4VERBOSE
208 if (verboseLevel>2){
209 G4cout << "G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength ";
210 G4cout << "[ " << GetProcessName() << "]" <<G4endl;
211 track.GetDynamicParticle()->DumpInfo();
212 G4cout << " in Material " << track.GetMaterial()->GetName() <<G4endl;
213 G4cout << "InteractionLength= " << value/cm <<"[cm] " <<G4endl;
214 }
215#endif
216
217// G4cout << "currentInteractionLength : " << fpState->currentInteractionLength << G4endl;
218// G4cout << "Returned time : " << G4BestUnit(value,"Time") << G4endl;
219
220 if(value < fReturnedValue)
221 fReturnedValue = value;
222
223 return value*-1;
224 // multiple by -1 to indicate to the tracking system that we are returning a time
225}
226
228{
229 G4Molecule* molecule = GetMolecule(track);
230#ifdef G4VERBOSE
231 if(verboseLevel > 1)
232 {
233 G4cout << "___________" << G4endl;
234 G4cout << ">>> Beginning of G4DNASecondOrderReaction verbose" << G4endl;
235 G4cout << ">>> Returned value : " << G4BestUnit(fReturnedValue,"Time") << G4endl;
236 G4cout << ">>> Time Step : " << G4BestUnit(G4ITTrackHolder::Instance()->GetTimeStep(),"Time") << G4endl;
237 G4cout << ">>> Reaction : " << molecule->GetName() << " + " << fpMaterial->GetName() << G4endl;
238 G4cout << ">>> End of G4DNASecondOrderReaction verbose <<<" << G4endl;
239 }
240#endif
245 return &fParticleChange;
246}
247
@ FatalErrorInArgument
G4ForceCondition
@ NotForced
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:67
G4ProcessType
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
#define InitProcessState(destination, source)
Definition: G4VITProcess.hh:53
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
virtual void AddIndirectDamage(const G4String &baseName, const G4Molecule *molecule, const G4ThreeVector &position, double time)
static G4DNADamages * Instance()
Definition: G4DNADamages.cc:59
static G4DNAMolecularMaterial * Instance()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4DNASecondOrderReaction & operator=(const G4DNASecondOrderReaction &)
const std::vector< double > * fpMoleculeDensity
const G4MolecularConfiguration * fpMolecularConfiguration
G4DNASecondOrderReaction(const G4String &aName="DNASecondOrderReaction", G4ProcessType type=fDecay)
void SetReaction(const G4MolecularConfiguration *, const G4Material *, double)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
SecondOrderReactionState *& fpSecondOrderReactionState
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void DumpInfo(G4int mode=0) const
static G4ITTrackHolder * Instance()
G4double GetMassOfMolecule() const
Definition: G4Material.hh:241
const G4String & GetName() const
Definition: G4Material.hh:177
size_t GetIndex() const
Definition: G4Material.hh:261
G4MolecularConfiguration * GetMolecularConfiguration()
Definition: G4Molecule.hh:273
const G4String & GetName() const
Definition: G4Molecule.cc:259
virtual void Initialize(const G4Track &)
Definition: G4Step.hh:78
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
G4int GetCurrentStepNumber() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void StartTracking(G4Track *)
Definition: G4VITProcess.cc:81
void SetInstantiateProcessState(G4bool flag)
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
Definition: G4VITProcess.cc:96
G4ProcessState * fpState
virtual void ResetNumberOfInteractionLengthLeft()
G4bool fProposesTimeStep
void ProposeTrackStatus(G4TrackStatus status)
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
G4int verboseLevel
Definition: G4VProcess.hh:368
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:351
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:125
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:352
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
#define DBL_MAX
Definition: templates.hh:83