Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNASecondOrderReaction Class Reference

#include <G4DNASecondOrderReaction.hh>

+ Inheritance diagram for G4DNASecondOrderReaction:

Classes

struct  SecondOrderReactionState
 

Public Member Functions

 G4DNASecondOrderReaction (const G4String &aName="DNASecondOrderReaction", G4ProcessType type=fDecay)
 
virtual ~G4DNASecondOrderReaction ()
 
 G4DNASecondOrderReaction (const G4DNASecondOrderReaction &)
 
G4DNASecondOrderReactionoperator= (const G4DNASecondOrderReaction &)
 
void StartTracking (G4Track *)
 
void SetReaction (const G4MolecularConfiguration *, const G4Material *, double)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VITProcess
 G4VITProcess (const G4String &name, G4ProcessType type=fNotDefined)
 
virtual ~G4VITProcess ()
 
 G4VITProcess (const G4VITProcess &other)
 
G4VITProcessoperator= (const G4VITProcess &other)
 
G4int operator== (const G4VITProcess &right) const
 
G4int operator!= (const G4VITProcess &right) const
 
size_t GetProcessID () const
 
G4ProcessState_LockGetProcessState ()
 
void SetProcessState (G4ProcessState_Lock *aProcInfo)
 
virtual void StartTracking (G4Track *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double GetInteractionTimeLeft ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4bool ProposesTimeStep () const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Attributes

G4bool fIsInitialized
 
G4double fReturnedValue
 
const std::vector< double > * fpMoleculeDensity
 
G4double fReactionRate
 
G4double fConcentration
 
G4double fMolarMassOfMaterial
 
SecondOrderReactionState *& fpSecondOrderReactionState
 
G4ParticleChange fParticleChange
 
const G4MolecularConfigurationfpMolecularConfiguration
 
const G4MaterialfpMaterial
 
- Protected Attributes inherited from G4VITProcess
G4ProcessStatefpState
 
G4bool fProposesTimeStep
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VITProcess
static const size_t & GetMaxProcessIndex ()
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VITProcess
void RetrieveProcessInfo ()
 
void CreateInfo ()
 
virtual void ClearInteractionTimeLeft ()
 
virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
virtual void ClearNumberOfInteractionLengthLeft ()
 
void SetInstantiateProcessState (G4bool flag)
 
G4bool InstantiateProcessState ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 46 of file G4DNASecondOrderReaction.hh.

Constructor & Destructor Documentation

◆ G4DNASecondOrderReaction() [1/2]

G4DNASecondOrderReaction::G4DNASecondOrderReaction ( const G4String aName = "DNASecondOrderReaction",
G4ProcessType  type = fDecay 
)

Definition at line 62 of file G4DNASecondOrderReaction.cc.

62 :
63 G4VITProcess(aName,type),
65{
66 Create();
67}
#define InitProcessState(destination, source)
Definition: G4VITProcess.hh:53
SecondOrderReactionState *& fpSecondOrderReactionState
G4ProcessState * fpState

◆ ~G4DNASecondOrderReaction()

G4DNASecondOrderReaction::~G4DNASecondOrderReaction ( )
virtual

Definition at line 76 of file G4DNASecondOrderReaction.cc.

77{
78 ;
79}

◆ G4DNASecondOrderReaction() [2/2]

G4DNASecondOrderReaction::G4DNASecondOrderReaction ( const G4DNASecondOrderReaction rhs)

Definition at line 69 of file G4DNASecondOrderReaction.cc.

69 :
70 G4VITProcess(rhs),
72{
73 Create();
74}

Member Function Documentation

◆ AlongStepDoIt()

virtual G4VParticleChange * G4DNASecondOrderReaction::AlongStepDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 94 of file G4DNASecondOrderReaction.hh.

97 {return 0;}

◆ AlongStepGetPhysicalInteractionLength()

virtual G4double G4DNASecondOrderReaction::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4double  ,
G4double ,
G4GPILSelection  
)
inlinevirtual

Implements G4VProcess.

Definition at line 85 of file G4DNASecondOrderReaction.hh.

91 { return -1.0; }

◆ AtRestDoIt()

virtual G4VParticleChange * G4DNASecondOrderReaction::AtRestDoIt ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VProcess.

Definition at line 79 of file G4DNASecondOrderReaction.hh.

82 {return 0;}

◆ AtRestGetPhysicalInteractionLength()

virtual G4double G4DNASecondOrderReaction::AtRestGetPhysicalInteractionLength ( const G4Track ,
G4ForceCondition  
)
inlinevirtual

Implements G4VProcess.

Definition at line 74 of file G4DNASecondOrderReaction.hh.

77 { return -1.0; }

◆ BuildPhysicsTable()

void G4DNASecondOrderReaction::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VITProcess.

Definition at line 94 of file G4DNASecondOrderReaction.cc.

95{
97 fMolarMassOfMaterial = fpMaterial->GetMassOfMolecule()*CLHEP::Avogadro*1e3;
98 fIsInitialized = true;
99}
static G4DNAMolecularMaterial * Instance()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
const std::vector< double > * fpMoleculeDensity
G4double GetMassOfMolecule() const
Definition: G4Material.hh:241

◆ operator=()

G4DNASecondOrderReaction & G4DNASecondOrderReaction::operator= ( const G4DNASecondOrderReaction rhs)

Definition at line 80 of file G4DNASecondOrderReaction.cc.

81{
82 if (this == &rhs) return *this; // handle self assignment
83
84 //assignment operator
85 return *this;
86}

◆ PostStepDoIt()

G4VParticleChange * G4DNASecondOrderReaction::PostStepDoIt ( const G4Track track,
const G4Step  
)
virtual

Implements G4VProcess.

Definition at line 227 of file G4DNASecondOrderReaction.cc.

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}
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:67
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
@ fStopAndKill
#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 G4ITTrackHolder * Instance()
const G4String & GetName() const
Definition: G4Material.hh:177
const G4String & GetName() const
Definition: G4Molecule.cc:259
virtual void Initialize(const G4Track &)
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
void ProposeTrackStatus(G4TrackStatus status)
G4int verboseLevel
Definition: G4VProcess.hh:368
#define DBL_MAX
Definition: templates.hh:83

◆ PostStepGetPhysicalInteractionLength()

G4double G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VProcess.

Definition at line 126 of file G4DNASecondOrderReaction.cc.

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}
@ NotForced
double G4double
Definition: G4Types.hh:64
const G4MolecularConfiguration * fpMolecularConfiguration
void DumpInfo(G4int mode=0) const
size_t GetIndex() const
Definition: G4Material.hh:261
G4MolecularConfiguration * GetMolecularConfiguration()
Definition: G4Molecule.hh:273
G4int GetCurrentStepNumber() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
Definition: G4VITProcess.cc:96
virtual void ResetNumberOfInteractionLengthLeft()
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ SetReaction()

void G4DNASecondOrderReaction::SetReaction ( const G4MolecularConfiguration molConf,
const G4Material mat,
double  reactionRate 
)

Definition at line 110 of file G4DNASecondOrderReaction.cc.

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}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

◆ StartTracking()

void G4DNASecondOrderReaction::StartTracking ( G4Track track)
virtual

Reimplemented from G4VITProcess.

Definition at line 102 of file G4DNASecondOrderReaction.cc.

103{
105 G4VITProcess::fpState = new SecondOrderReactionState();
107}
virtual void StartTracking(G4Track *)
Definition: G4VITProcess.cc:81
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:125

Member Data Documentation

◆ fConcentration

G4double G4DNASecondOrderReaction::fConcentration
protected

Definition at line 118 of file G4DNASecondOrderReaction.hh.

Referenced by PostStepGetPhysicalInteractionLength().

◆ fIsInitialized

G4bool G4DNASecondOrderReaction::fIsInitialized
protected

Definition at line 112 of file G4DNASecondOrderReaction.hh.

Referenced by BuildPhysicsTable(), and SetReaction().

◆ fMolarMassOfMaterial

G4double G4DNASecondOrderReaction::fMolarMassOfMaterial
protected

Definition at line 119 of file G4DNASecondOrderReaction.hh.

Referenced by BuildPhysicsTable().

◆ fParticleChange

G4ParticleChange G4DNASecondOrderReaction::fParticleChange
protected

Definition at line 121 of file G4DNASecondOrderReaction.hh.

Referenced by PostStepDoIt().

◆ fpMaterial

const G4Material* G4DNASecondOrderReaction::fpMaterial
protected

Definition at line 124 of file G4DNASecondOrderReaction.hh.

Referenced by BuildPhysicsTable(), PostStepDoIt(), and SetReaction().

◆ fpMolecularConfiguration

const G4MolecularConfiguration* G4DNASecondOrderReaction::fpMolecularConfiguration
protected

◆ fpMoleculeDensity

const std::vector<double>* G4DNASecondOrderReaction::fpMoleculeDensity
protected

Definition at line 116 of file G4DNASecondOrderReaction.hh.

Referenced by BuildPhysicsTable().

◆ fpSecondOrderReactionState

SecondOrderReactionState*& G4DNASecondOrderReaction::fpSecondOrderReactionState
protected

Definition at line 120 of file G4DNASecondOrderReaction.hh.

Referenced by PostStepGetPhysicalInteractionLength().

◆ fReactionRate

G4double G4DNASecondOrderReaction::fReactionRate
protected

◆ fReturnedValue

G4double G4DNASecondOrderReaction::fReturnedValue
protected

The documentation for this class was generated from the following files: