Geant4 11.2.2
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)
 
 ~G4DNASecondOrderReaction () override
 
 G4DNASecondOrderReaction (const G4DNASecondOrderReaction &)
 
G4DNASecondOrderReactionoperator= (const G4DNASecondOrderReaction &)
 
void StartTracking (G4Track *) override
 
void SetReaction (const G4MolecularConfiguration *, const G4Material *, double)
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) override
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
- Public Member Functions inherited from G4VITProcess
 G4VITProcess (const G4String &name, G4ProcessType type=fNotDefined)
 
 ~G4VITProcess () override
 
 G4VITProcess (const G4VITProcess &other)
 
G4VITProcessoperator= (const G4VITProcess &other)
 
G4bool operator== (const G4VITProcess &right) const
 
G4bool operator!= (const G4VITProcess &right) const
 
size_t GetProcessID () const
 
G4shared_ptr< G4ProcessState_LockGetProcessState ()
 
void SetProcessState (G4shared_ptr< G4ProcessState_Lock > aProcInfo)
 
void ResetProcessState ()
 
void StartTracking (G4Track *) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4double GetInteractionTimeLeft ()
 
void ResetNumberOfInteractionLengthLeft () override
 
G4bool ProposesTimeStep () const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
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 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 const G4VProcessGetCreatorProcess () const
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Attributes

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

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 ()
 
template<typename T >
T * GetState ()
 
virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
virtual void ClearInteractionTimeLeft ()
 
virtual void ClearNumberOfInteractionLengthLeft ()
 
void SetInstantiateProcessState (G4bool flag)
 
G4bool InstantiateProcessState ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 54 of file G4DNASecondOrderReaction.hh.

Constructor & Destructor Documentation

◆ G4DNASecondOrderReaction() [1/2]

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

Definition at line 71 of file G4DNASecondOrderReaction.cc.

71 :
72 G4VITProcess(aName,type)
73{
74 Create();
75}
G4VITProcess(const G4String &name, G4ProcessType type=fNotDefined)

◆ ~G4DNASecondOrderReaction()

G4DNASecondOrderReaction::~G4DNASecondOrderReaction ( )
overridedefault

◆ G4DNASecondOrderReaction() [2/2]

G4DNASecondOrderReaction::G4DNASecondOrderReaction ( const G4DNASecondOrderReaction & rhs)

Definition at line 77 of file G4DNASecondOrderReaction.cc.

77 :
78 G4VITProcess(rhs)
79{
80 Create();
81}

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4DNASecondOrderReaction::AlongStepDoIt ( const G4Track & ,
const G4Step &  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 102 of file G4DNASecondOrderReaction.hh.

105 {return nullptr;}

◆ AlongStepGetPhysicalInteractionLength()

G4double G4DNASecondOrderReaction::AlongStepGetPhysicalInteractionLength ( const G4Track & ,
G4double ,
G4double ,
G4double & ,
G4GPILSelection *  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 93 of file G4DNASecondOrderReaction.hh.

99 { return -1.0; }

◆ AtRestDoIt()

G4VParticleChange * G4DNASecondOrderReaction::AtRestDoIt ( const G4Track & ,
const G4Step &  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 87 of file G4DNASecondOrderReaction.hh.

90 {return nullptr;}

◆ AtRestGetPhysicalInteractionLength()

G4double G4DNASecondOrderReaction::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition *  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 82 of file G4DNASecondOrderReaction.hh.

85 { return -1.0; }

◆ BuildPhysicsTable()

void G4DNASecondOrderReaction::BuildPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 100 of file G4DNASecondOrderReaction.cc.

101{
103 fMolarMassOfMaterial = fpMaterial->GetMassOfMolecule()*CLHEP::Avogadro*1e3;
104 fIsInitialized = true;
105}
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()
const std::vector< double > * fpMoleculeDensity
G4double GetMassOfMolecule() const

◆ operator=()

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

Definition at line 86 of file G4DNASecondOrderReaction.cc.

87{
88 if (this == &rhs) return *this; // handle self assignment
89
90 //assignment operator
91 return *this;
92}

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 247 of file G4DNASecondOrderReaction.cc.

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}
#define State(X)
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:74
#define G4BestUnit(a, b)
@ fStopAndKill
#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 G4String & GetName() const
const G4String & GetName() const override
void Initialize(const G4Track &) override
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
void ProposeTrackStatus(G4TrackStatus status)
G4int verboseLevel
static G4VScheduler * Instance()
#define DBL_MAX
Definition templates.hh:62

◆ PostStepGetPhysicalInteractionLength()

G4double G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 132 of file G4DNASecondOrderReaction.cc.

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}
@ NotForced
double G4double
Definition G4Types.hh:83
const G4MolecularConfiguration * fpMolecularConfiguration
void DumpInfo(G4int mode=0) const
std::size_t GetIndex() const
const G4MolecularConfiguration * GetMolecularConfiguration() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4shared_ptr< G4ProcessState > fpState
void ResetNumberOfInteractionLengthLeft() override
const G4String & GetProcessName() const

◆ SetReaction()

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

Definition at line 116 of file G4DNASecondOrderReaction.cc.

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}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription

◆ StartTracking()

void G4DNASecondOrderReaction::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 108 of file G4DNASecondOrderReaction.cc.

109{
111 G4VITProcess::fpState = std::make_shared<SecondOrderReactionState>();
113}
void StartTracking(G4Track *) override
virtual void StartTracking(G4Track *)
Definition G4VProcess.cc:87

Member Data Documentation

◆ fConcentration

G4double G4DNASecondOrderReaction::fConcentration
protected

Definition at line 126 of file G4DNASecondOrderReaction.hh.

Referenced by PostStepGetPhysicalInteractionLength().

◆ fIsInitialized

G4bool G4DNASecondOrderReaction::fIsInitialized
protected

Definition at line 120 of file G4DNASecondOrderReaction.hh.

Referenced by BuildPhysicsTable(), and SetReaction().

◆ fMolarMassOfMaterial

G4double G4DNASecondOrderReaction::fMolarMassOfMaterial
protected

Definition at line 127 of file G4DNASecondOrderReaction.hh.

Referenced by BuildPhysicsTable().

◆ fParticleChange

G4ParticleChange G4DNASecondOrderReaction::fParticleChange
protected

Definition at line 128 of file G4DNASecondOrderReaction.hh.

Referenced by PostStepDoIt().

◆ fpMaterial

const G4Material* G4DNASecondOrderReaction::fpMaterial
protected

Definition at line 131 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 124 of file G4DNASecondOrderReaction.hh.

Referenced by BuildPhysicsTable().

◆ fReactionRate

G4double G4DNASecondOrderReaction::fReactionRate
protected

◆ fReturnedValue

G4double G4DNASecondOrderReaction::fReturnedValue
protected

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