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

#include <G4DNAPolyNucleotideReactionProcess.hh>

+ Inheritance diagram for G4DNAPolyNucleotideReactionProcess:

Public Member Functions

 G4DNAPolyNucleotideReactionProcess (const G4String &aName="DNAStaticMoleculeReactionProcess", G4int verbosityLevel=0)
 
 ~G4DNAPolyNucleotideReactionProcess () override
 
void SetDNADamageReactionModel (G4VDNAHitModel *pModel)
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
G4double CalculateTimeStep (const G4Track &trackA, const G4double &userTimeStep=0)
 
void StartTracking (G4Track *aTrack) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *pForceCond) override
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &) override
 
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *) override
 
G4DNAPolyNucleotideReactionProcessoperator= (const G4DNAPolyNucleotideReactionProcess &)=delete
 
void SetVerbose (G4int verbose)
 
- Public Member Functions inherited from G4VITDiscreteProcess
 G4VITDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VITDiscreteProcess (G4VITDiscreteProcess &)
 
virtual ~G4VITDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
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)
 
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 ()
 
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 ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool 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 const G4VProcessGetCreatorProcess () const
 
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
 
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

G4VParticleChange fParticleChange
 
- 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)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- 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 37 of file G4DNAPolyNucleotideReactionProcess.hh.

Constructor & Destructor Documentation

◆ G4DNAPolyNucleotideReactionProcess()

G4DNAPolyNucleotideReactionProcess::G4DNAPolyNucleotideReactionProcess ( const G4String aName = "DNAStaticMoleculeReactionProcess",
G4int  verbosityLevel = 0 
)
explicit

Definition at line 46 of file G4DNAPolyNucleotideReactionProcess.cc.

49 , fHasAlreadyReachedNullTime(false)
50 , fVerbose(verbosityLevel)
52 , fpDamageModel(nullptr)
53{
55 enableAtRestDoIt = false;
56 enableAlongStepDoIt = false;
57 enablePostStepDoIt = true;
58 fProposesTimeStep = true;
61}
@ fUserDefined
static G4double GetDNADistanceCutOff()
Definition: G4IRTUtils.cc:62
void SetInstantiateProcessState(G4bool flag)
G4bool fProposesTimeStep
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:363
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:364
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:365
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325

◆ ~G4DNAPolyNucleotideReactionProcess()

G4DNAPolyNucleotideReactionProcess::~G4DNAPolyNucleotideReactionProcess ( )
override

Definition at line 63 of file G4DNAPolyNucleotideReactionProcess.cc.

64{
65 delete fpDamageModel; // for now, this process handls only one model
66}

Member Function Documentation

◆ CalculateTimeStep()

G4double G4DNAPolyNucleotideReactionProcess::CalculateTimeStep ( const G4Track trackA,
const G4double userTimeStep = 0 
)

Definition at line 75 of file G4DNAPolyNucleotideReactionProcess.cc.

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}
#define State(X)
#define PrepareState()
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:74
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4int GetTrackID() const
virtual G4double CalculateReactionTime(const G4Track &trackA, DNANode &)=0
#define DBL_MAX
Definition: templates.hh:62

Referenced by PostStepGetPhysicalInteractionLength().

◆ GetMeanFreePath()

G4double G4DNAPolyNucleotideReactionProcess::GetMeanFreePath ( const G4Track ,
G4double  ,
G4ForceCondition  
)
inlineoverridevirtual

Implements G4VITDiscreteProcess.

Definition at line 60 of file G4DNAPolyNucleotideReactionProcess.hh.

61 {
62 return DBL_MAX;
63 }

◆ IsApplicable()

G4bool G4DNAPolyNucleotideReactionProcess::IsApplicable ( const G4ParticleDefinition )
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 47 of file G4DNAPolyNucleotideReactionProcess.hh.

47{ return true; }

◆ operator=()

G4DNAPolyNucleotideReactionProcess & G4DNAPolyNucleotideReactionProcess::operator= ( const G4DNAPolyNucleotideReactionProcess )
delete

◆ PostStepDoIt()

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

Reimplemented from G4VITDiscreteProcess.

Definition at line 162 of file G4DNAPolyNucleotideReactionProcess.cc.

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}
@ fStopAndKill
virtual G4bool DoReaction(const G4Track &track, const G4double &, const DNANode &)=0
void ProposeTrackStatus(G4TrackStatus status)
virtual void Initialize(const G4Track &)

◆ PostStepGetPhysicalInteractionLength()

G4double G4DNAPolyNucleotideReactionProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  ,
G4ForceCondition pForceCond 
)
overridevirtual

Reimplemented from G4VITDiscreteProcess.

Definition at line 127 of file G4DNAPolyNucleotideReactionProcess.cc.

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}
@ NotForced
G4double CalculateTimeStep(const G4Track &trackA, const G4double &userTimeStep=0)
G4double GetGlobalTime() const
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4shared_ptr< G4ProcessState > fpState
virtual void ResetNumberOfInteractionLengthLeft()

◆ SetDNADamageReactionModel()

void G4DNAPolyNucleotideReactionProcess::SetDNADamageReactionModel ( G4VDNAHitModel pModel)
inline

Definition at line 98 of file G4DNAPolyNucleotideReactionProcess.hh.

100{
101 fpDamageModel = pModel;
102}

◆ SetVerbose()

void G4DNAPolyNucleotideReactionProcess::SetVerbose ( G4int  verbose)
inline

Definition at line 93 of file G4DNAPolyNucleotideReactionProcess.hh.

94{
95 fVerbose = verbose;
96}

◆ StartTracking()

void G4DNAPolyNucleotideReactionProcess::StartTracking ( G4Track aTrack)
overridevirtual

Reimplemented from G4VITProcess.

Definition at line 179 of file G4DNAPolyNucleotideReactionProcess.cc.

180{
182 G4VITProcess::fpState.reset(new G4PolyNucleotideReactionState());
184}
virtual void StartTracking(G4Track *)
Definition: G4VITProcess.cc:85
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:87

Member Data Documentation

◆ fParticleChange

G4VParticleChange G4DNAPolyNucleotideReactionProcess::fParticleChange
protected

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