Geant4 11.2.2
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
 
 G4DNAPolyNucleotideReactionProcess (const G4DNAPolyNucleotideReactionProcess &)=delete
 
G4DNAPolyNucleotideReactionProcessoperator= (const G4DNAPolyNucleotideReactionProcess &)=delete
 
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
 
void SetVerbose (G4int verbose)
 
- Public Member Functions inherited from G4VITDiscreteProcess
 G4VITDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VITDiscreteProcess (G4VITDiscreteProcess &)
 
 ~G4VITDiscreteProcess () override
 
G4VITDiscreteProcessoperator= (const G4VITDiscreteProcess &right)=delete
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) override
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &) 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 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

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)
 
- Protected Member Functions inherited from G4VITDiscreteProcess
- 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() [1/2]

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

Definition at line 49 of file G4DNAPolyNucleotideReactionProcess.cc.

52 ,
53 fVerbose(verbosityLevel)
55{
57 enableAtRestDoIt = false;
58 enableAlongStepDoIt = false;
59 enablePostStepDoIt = true;
60 fProposesTimeStep = true;
63}
@ fUserDefined
static G4double GetDNADistanceCutOff()
Definition G4IRTUtils.cc:62
void SetInstantiateProcessState(G4bool flag)
G4bool fProposesTimeStep
G4bool enableAtRestDoIt
G4bool enableAlongStepDoIt
void SetProcessSubType(G4int)
G4bool enablePostStepDoIt
G4VParticleChange * pParticleChange

◆ ~G4DNAPolyNucleotideReactionProcess()

G4DNAPolyNucleotideReactionProcess::~G4DNAPolyNucleotideReactionProcess ( )
override

Definition at line 65 of file G4DNAPolyNucleotideReactionProcess.cc.

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

◆ G4DNAPolyNucleotideReactionProcess() [2/2]

G4DNAPolyNucleotideReactionProcess::G4DNAPolyNucleotideReactionProcess ( const G4DNAPolyNucleotideReactionProcess & )
delete

Member Function Documentation

◆ CalculateTimeStep()

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

Definition at line 77 of file G4DNAPolyNucleotideReactionProcess.cc.

79{
81 fHasAlreadyReachedNullTime = false;
82 State(fSampledMinTimeStep) = DBL_MAX;
83 State(theInteractionTimeLeft) = DBL_MAX;
84 State(currentInteractionLength) = -1;
85#ifdef G4VERBOSE
86 if(fVerbose > 1)
87 {
88 auto pMoleculeA = GetMolecule(trackA);
89 G4cout << "________________________________________________________________"
90 "_______"
91 << G4endl;
92 G4cout << "G4DNAPolyNucleotideReactionProcess::CalculateTimleStep"
93 << G4endl;
94 G4cout << "Check done for molecule : " << pMoleculeA->GetName() << " ("
95 << trackA.GetTrackID() << ") " << G4endl;
96 }
97#endif
98 //__________________________________________________________________
99 // Retrieve general informations for making reactions
100 G4double reactionTime =
101 fpDamageModel->CalculateReactionTime(trackA, State(fNodeReactant));
102
103 if(reactionTime < 0)
104 {
105 return DBL_MAX;
106 }
107 State(fSampledMinTimeStep) = reactionTime;
108 State(theInteractionTimeLeft) = State(fSampledMinTimeStep);
109 State(currentInteractionLength) = State(theInteractionTimeLeft);
110
111#ifdef G4VERBOSE
112 if(fVerbose > 1)
113 {
114 G4cout << " theInteractionTimeLeft : " << State(theInteractionTimeLeft)
115 << G4endl;
116 G4cout << " State(fNodeReactant) : " << State(fNodeReactant).index()
117 << G4endl;
118
119 G4cout << "________________________________________________________________"
120 "_______"
121 << G4endl;
122 }
123#endif
124
125 return State(fSampledMinTimeStep);
126}
#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:67
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 64 of file G4DNAPolyNucleotideReactionProcess.hh.

65 {
66 return DBL_MAX;
67 }

◆ IsApplicable()

G4bool G4DNAPolyNucleotideReactionProcess::IsApplicable ( const G4ParticleDefinition & )
inlineoverridevirtual

Reimplemented from G4VProcess.

Definition at line 51 of file G4DNAPolyNucleotideReactionProcess.hh.

51{ return true; }

◆ operator=()

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

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 164 of file G4DNAPolyNucleotideReactionProcess.cc.

166{
167 PrepareState();
168 auto isReacted = fpDamageModel->DoReaction(
169 track, State(theInteractionTimeLeft), State(fNodeReactant));
170 if(!isReacted)
171 {
172 // no reaction even this process is called
174 return pParticleChange;
175 }
176 fParticleChange.Initialize(track); // To initialise TouchableChange
178 State(fPreviousTimeAtPreStepPoint) = -1;
179 return pParticleChange;
180}
@ 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

Implements G4VProcess.

Definition at line 129 of file G4DNAPolyNucleotideReactionProcess.cc.

133{
134 PrepareState();
135 CalculateTimeStep(track);
136 *pForceCond = NotForced;
137
138 G4double previousTimeStep(-1.);
139
140 if(State(fPreviousTimeAtPreStepPoint) != -1)
141 {
142 previousTimeStep =
143 track.GetGlobalTime() - State(fPreviousTimeAtPreStepPoint);
144 }
145
146 State(fPreviousTimeAtPreStepPoint) = track.GetGlobalTime();
147
148 if((fpState->currentInteractionLength <= 0) || (previousTimeStep < 0.0) ||
149 (fpState->theNumberOfInteractionLengthLeft <= 0.0))
150 {
152 }
153 else if(previousTimeStep > 0.0)
154 {
156 previousTimeStep); // TODO need to check this
157 }
158 return State(theInteractionTimeLeft) * -1;
159 // return value * -1;//should regard this !
160}
@ NotForced
G4double CalculateTimeStep(const G4Track &trackA, const G4double &userTimeStep=0)
G4double GetGlobalTime() const
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4shared_ptr< G4ProcessState > fpState
void ResetNumberOfInteractionLengthLeft() override

◆ SetDNADamageReactionModel()

void G4DNAPolyNucleotideReactionProcess::SetDNADamageReactionModel ( G4VDNAHitModel * pModel)
inline

Definition at line 99 of file G4DNAPolyNucleotideReactionProcess.hh.

101{
102 fpDamageModel = pModel;
103}

◆ SetVerbose()

void G4DNAPolyNucleotideReactionProcess::SetVerbose ( G4int verbose)
inline

Definition at line 94 of file G4DNAPolyNucleotideReactionProcess.hh.

95{
96 fVerbose = verbose;
97}

◆ StartTracking()

void G4DNAPolyNucleotideReactionProcess::StartTracking ( G4Track * aTrack)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 181 of file G4DNAPolyNucleotideReactionProcess.cc.

182{
184 G4VITProcess::fpState = std::make_shared<G4PolyNucleotideReactionState>();
186}
void StartTracking(G4Track *) override
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: