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

#include <G4DNAMolecularReaction.hh>

+ Inheritance diagram for G4DNAMolecularReaction:

Public Member Functions

 G4DNAMolecularReaction ()
 
virtual ~G4DNAMolecularReaction ()
 
 G4DNAMolecularReaction (const G4DNAMolecularReaction &other)
 
G4DNAMolecularReactionoperator= (const G4DNAMolecularReaction &other)
 
virtual G4bool TestReactibility (const G4Track &, const G4Track &, const double currentStepTime, const double previousStepTime, bool userStepTimeLimit)
 
virtual G4ITReactionChangeMakeReaction (const G4Track &, const G4Track &)
 
void SetReactionModel (G4VDNAReactionModel *)
 
void SetReactionTable (const G4DNAMolecularReactionTable *)
 
void SetVerbose (int)
 
- Public Member Functions inherited from G4VITReactionProcess
 G4VITReactionProcess ()
 
virtual ~G4VITReactionProcess ()
 
 G4VITReactionProcess (const G4VITReactionProcess &other)
 
G4VITReactionProcessoperator= (const G4VITReactionProcess &other)
 
virtual void Initialize ()
 
virtual G4bool IsApplicable (G4ITType, G4ITType) const
 
virtual G4bool TestReactibility (const G4Track &, const G4Track &, const double, const double, bool)=0
 
virtual G4ITReactionChangeMakeReaction (const G4Track &, const G4Track &)=0
 
void SetReactionTable (const G4ITReactionTable *)
 
void ResetChanges ()
 

Protected Attributes

const G4DNAMolecularReactionTable *& fMolReactionTable
 
G4VDNAReactionModelfReactionModel
 
G4int fVerbose
 
G4double fReactionRadius
 
G4double fDistance
 
- Protected Attributes inherited from G4VITReactionProcess
const G4ITReactionTablefpReactionTable
 
G4ITReactionChangefpChanges
 
G4String fName
 

Detailed Description

G4DNAMolecularReaction is the reaction process used in G4DNAMolecularStepByStepModel between two molecules. After the global track steps, it test whether the molecules can react. If so, the reaction is made.

Definition at line 55 of file G4DNAMolecularReaction.hh.

Constructor & Destructor Documentation

◆ G4DNAMolecularReaction() [1/2]

G4DNAMolecularReaction::G4DNAMolecularReaction ( )

Default constructor

Definition at line 47 of file G4DNAMolecularReaction.cc.

48 fMolReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
49{
50 //ctor
51 fVerbose = 0;
53 fReactionRadius = -1;
54 fDistance = -1;
55}
const G4DNAMolecularReactionTable *& fMolReactionTable
G4VDNAReactionModel * fReactionModel
const G4ITReactionTable * fpReactionTable

◆ ~G4DNAMolecularReaction()

G4DNAMolecularReaction::~G4DNAMolecularReaction ( )
virtual

Default destructor

Definition at line 57 of file G4DNAMolecularReaction.cc.

58{
59 //dtor
60 if(fpChanges) delete fpChanges;
61}
G4ITReactionChange * fpChanges

◆ G4DNAMolecularReaction() [2/2]

G4DNAMolecularReaction::G4DNAMolecularReaction ( const G4DNAMolecularReaction other)

Copy constructor

Parameters
otherObject to copy from

Definition at line 63 of file G4DNAMolecularReaction.cc.

64 fMolReactionTable(reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
65{
66 //copy ctor
67 fVerbose = other.fVerbose ;
70 fReactionRadius = -1;
71 fDistance = -1;
72}

Member Function Documentation

◆ MakeReaction()

G4ITReactionChange * G4DNAMolecularReaction::MakeReaction ( const G4Track trackA,
const G4Track trackB 
)
virtual

Will generate the products of the two given tracks

Implements G4VITReactionProcess.

Definition at line 149 of file G4DNAMolecularReaction.cc.

150{
152 fpChanges->Initialize(trackA, trackB);
153
154 G4Molecule* moleculeA = GetMolecule(trackA);
155 G4Molecule* moleculeB = GetMolecule(trackB);
156
157#ifdef G4VERBOSE
158 // DEBUG
159 if(fVerbose)
160 {
161 G4cout << "G4DNAMolecularReaction::MakeReaction" << G4endl;
162 G4cout<<"TrackA n°" << trackA.GetTrackID()
163 <<"\t | Track B n°" << trackB.GetTrackID() << G4endl;
164
165 G4cout<<"Track A : Position : " << G4BestUnit(trackA.GetPosition(),"Length")
166 <<"\t Global Time : " << G4BestUnit(trackA.GetGlobalTime(), "Time")<< G4endl;
167
168 G4cout<<"Track B : Position : " << G4BestUnit(trackB.GetPosition() ,"Length")
169 <<"\t Global Time : " << G4BestUnit(trackB.GetGlobalTime(), "Time")<< G4endl;
170
171 G4cout<<"Reaction range : " << G4BestUnit(fReactionRadius,"Length")
172 << " \t Separation distance : " << G4BestUnit(fDistance,"Length")<< G4endl;
173 G4cout << "--------------------------------------------" << G4endl;
174 }
175#endif
176
177 const G4DNAMolecularReactionData* reactionData = fMolReactionTable->GetReactionData(moleculeA, moleculeB);
178
179 G4int nbProducts = reactionData->GetNbProducts();
180
181 if (nbProducts)
182 {
183 G4double D1 = moleculeA->GetDiffusionCoefficient();
184 G4double D2 = moleculeB->GetDiffusionCoefficient();
185 G4double sqrD1 = sqrt(D1);
186 G4double sqrD2 = sqrt(D2);
187 G4double numerator = sqrD1 + sqrD2;
188 G4ThreeVector reactionSite = sqrD1/numerator * trackA.GetPosition()
189 + sqrD2/numerator * trackB.GetPosition() ;
190
191 for (G4int j=0 ; j < nbProducts; j++)
192 {
193 G4Molecule* product = new G4Molecule(*reactionData->GetProduct(j));
194 G4Track* productTrack = product->BuildTrack(trackA.GetGlobalTime(),
195 reactionSite);
196
197 productTrack->SetTrackStatus(fAlive);
198
199 fpChanges->AddSecondary(productTrack);
200 G4ITManager<G4Molecule>::Instance()->Push(productTrack);
201 }
202 }
203
204 fpChanges->KillParents(true);
205
206 return fpChanges;
207}
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:67
@ fAlive
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
const G4Molecule * GetProduct(G4int i) const
const G4DNAMolecularReactionData * GetReactionData(const G4Molecule *, const G4Molecule *) const
static G4ITManager< T > * Instance()
void AddSecondary(G4Track *aSecondary)
void Initialize(const G4Track &, const G4Track &, G4VParticleChange *particleChangeA=0, G4VParticleChange *particleChangeB=0)
G4Track * BuildTrack(G4double globalTime, const G4ThreeVector &Position)
Definition: G4Molecule.cc:279
G4double GetDiffusionCoefficient() const
Definition: G4Molecule.cc:405
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const

◆ operator=()

G4DNAMolecularReaction & G4DNAMolecularReaction::operator= ( const G4DNAMolecularReaction other)

Assignment operator

Parameters
otherObject to assign from
Returns
A reference to this

Definition at line 74 of file G4DNAMolecularReaction.cc.

75{
76 if (this == &rhs) return *this; // handle self assignment
77
78 fVerbose = rhs.fVerbose ;
79 fMolReactionTable = rhs.fMolReactionTable;
80 fReactionRadius = -1;
81 fDistance = -1;
82
83 //assignment operator
84 return *this;
85}

◆ SetReactionModel()

void G4DNAMolecularReaction::SetReactionModel ( G4VDNAReactionModel model)
inline

Definition at line 102 of file G4DNAMolecularReaction.hh.

103{
104 fReactionModel = model;
105}

◆ SetReactionTable()

void G4DNAMolecularReaction::SetReactionTable ( const G4DNAMolecularReactionTable )
inline

◆ SetVerbose()

void G4DNAMolecularReaction::SetVerbose ( int  verb)
inline

Definition at line 107 of file G4DNAMolecularReaction.hh.

108{
109 fVerbose = verb;
110}

◆ TestReactibility()

G4bool G4DNAMolecularReaction::TestReactibility ( const G4Track trackA,
const G4Track trackB,
const double  currentStepTime,
const double  previousStepTime,
bool  userStepTimeLimit 
)
virtual

Given two tracks, it tells you whether they can react

Implements G4VITReactionProcess.

Definition at line 87 of file G4DNAMolecularReaction.cc.

92{
93 G4Molecule* moleculeA = GetMolecule(trackA);
94 G4Molecule* moleculeB = GetMolecule(trackB);
95
97 {
98 G4ExceptionDescription exceptionDescription ("You have to give a reaction model to the molecular reaction process");
99 G4Exception("G4DNAMolecularReaction::TestReactibility","MolecularReaction001",
100 FatalErrorInArgument,exceptionDescription);
101 return false; // makes coverity happy
102 }
103 if(fMolReactionTable==0)
104 {
105 G4ExceptionDescription exceptionDescription ("You have to give a reaction table to the molecular reaction process");
106 G4Exception("G4DNAMolecularReaction::TestReactibility","MolecularReaction002",
107 FatalErrorInArgument,exceptionDescription);
108 return false; // makes coverity happy
109 }
110
111 // Retrieve reaction range
112 fReactionRadius = -1 ; // reaction Range
113 fReactionRadius = fReactionModel -> GetReactionRadius(moleculeA, moleculeB);
114
115 fDistance = -1 ; // separation distance
116
117 if(currentStepTime == 0.)
118 {
119 userStepTimeLimit = false;
120 }
121
122 G4bool output = fReactionModel->FindReaction(trackA, trackB, fReactionRadius,fDistance, userStepTimeLimit);
123
124#ifdef G4VERBOSE
125 // DEBUG
126 if(fVerbose > 1)
127 {
128 G4cout << "\033[1;39;36m" << "G4MolecularReaction "<< G4endl;
129 G4cout << "FindReaction returned : " << G4BestUnit(output,"Length") << G4endl;
130 G4cout << " reaction radius : " << G4BestUnit(fReactionRadius,"Length")
131 << " real distance : " << G4BestUnit((trackA.GetPosition() - trackB.GetPosition()).mag(), "Length")
132 << " calculated distance by model (= -1 if real distance > reaction radius and the user limitation step is not reached) : "
133 << G4BestUnit(fDistance,"Length")
134 << G4endl;
135
136 G4cout << "TrackID A : " << trackA.GetTrackID()
137 << ", TrackID B : " << trackB.GetTrackID()
138 << " | MolA " << moleculeA->GetName()
139 << ", MolB " << moleculeB->GetName()
140 <<"\033[0m\n"
141 << G4endl;
142 G4cout << "--------------------------------------------" << G4endl;
143 }
144#endif
145 return output ;
146}
@ FatalErrorInArgument
bool G4bool
Definition: G4Types.hh:67
const G4String & GetName() const
Definition: G4Molecule.cc:259
virtual G4bool FindReaction(const G4Track &, const G4Track &, const G4double, G4double &, const G4bool)=0
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

Member Data Documentation

◆ fDistance

G4double G4DNAMolecularReaction::fDistance
protected

◆ fMolReactionTable

const G4DNAMolecularReactionTable*& G4DNAMolecularReaction::fMolReactionTable
protected

◆ fReactionModel

G4VDNAReactionModel* G4DNAMolecularReaction::fReactionModel
protected

◆ fReactionRadius

G4double G4DNAMolecularReaction::fReactionRadius
protected

◆ fVerbose

G4int G4DNAMolecularReaction::fVerbose
protected

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