Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNAMolecularReaction Class Reference

#include <G4DNAMolecularReaction.hh>

Inheritance diagram for G4DNAMolecularReaction:
Collaboration 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
 
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 64 of file G4DNAMolecularReaction.hh.

Constructor & Destructor Documentation

G4DNAMolecularReaction::G4DNAMolecularReaction ( )

Default constructor

Definition at line 48 of file G4DNAMolecularReaction.cc.

48  :
51  reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
52 {
53  //ctor
54  fVerbose = 0;
55  fReactionModel = 0;
56  fReactionRadius = -1;
57  fDistance = -1;
58 }
const G4ITReactionTable * fpReactionTable
const G4DNAMolecularReactionTable *& fMolReactionTable
G4VDNAReactionModel * fReactionModel
G4DNAMolecularReaction::~G4DNAMolecularReaction ( )
virtual

Default destructor

Definition at line 60 of file G4DNAMolecularReaction.cc.

61 {
62  //dtor
63  if (fpChanges) delete fpChanges;
64 }
G4ITReactionChange * fpChanges
G4DNAMolecularReaction::G4DNAMolecularReaction ( const G4DNAMolecularReaction other)

Copy constructor

Parameters
otherObject to copy from

Definition at line 66 of file G4DNAMolecularReaction.cc.

66  :
67  G4VITReactionProcess(other),
69  reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
70 {
71  //copy ctor
72  fVerbose = other.fVerbose;
74  fReactionModel = 0;
75  fReactionRadius = -1;
76  fDistance = -1;
77 }
const G4ITReactionTable * fpReactionTable
const G4DNAMolecularReactionTable *& fMolReactionTable
G4VDNAReactionModel * fReactionModel

Member Function Documentation

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

Will generate the products of the two given tracks

Implements G4VITReactionProcess.

Definition at line 160 of file G4DNAMolecularReaction.cc.

162 {
164  fpChanges->Initialize(trackA, trackB);
165 
166  G4MolecularConfiguration* moleculeA =
168  G4MolecularConfiguration* moleculeB =
170 
171 #ifdef G4VERBOSE
172  // DEBUG
173  if (fVerbose)
174  {
175  G4cout << "G4DNAMolecularReaction::MakeReaction" << G4endl;
176  G4cout<<"TrackA n°" << trackA.GetTrackID()
177  <<"\t | Track B n°" << trackB.GetTrackID() << G4endl;
178 
179  G4cout<<"Track A : Position : " << G4BestUnit(trackA.GetPosition(),"Length")
180  <<"\t Global Time : " << G4BestUnit(trackA.GetGlobalTime(), "Time")<< G4endl;
181 
182  G4cout<<"Track B : Position : " << G4BestUnit(trackB.GetPosition() ,"Length")
183  <<"\t Global Time : " << G4BestUnit(trackB.GetGlobalTime(), "Time")<< G4endl;
184 
185  G4cout<<"Reaction range : " << G4BestUnit(fReactionRadius,"Length")
186  << " \t Separation distance : " << G4BestUnit(fDistance,"Length")<< G4endl;
187  G4cout << "--------------------------------------------" << G4endl;
188  }
189 #endif
190 
191  const G4DNAMolecularReactionData* reactionData = fMolReactionTable
192  ->GetReactionData(moleculeA, moleculeB);
193 
194  G4int nbProducts = reactionData->GetNbProducts();
195 
196  if (nbProducts)
197  {
198  G4double D1 = moleculeA->GetDiffusionCoefficient();
199  G4double D2 = moleculeB->GetDiffusionCoefficient();
200  G4double sqrD1 = sqrt(D1);
201  G4double sqrD2 = sqrt(D2);
202  G4double numerator = sqrD1 + sqrD2;
203  G4ThreeVector reactionSite = sqrD1 / numerator * trackA.GetPosition()
204  + sqrD2 / numerator * trackB.GetPosition();
205 
206  for (G4int j = 0; j < nbProducts; j++)
207  {
208  G4Molecule* product = new G4Molecule(reactionData->GetProduct(j));
209  G4Track* productTrack = product->BuildTrack(trackA.GetGlobalTime(),
210  reactionSite);
211 
212 // G4cout << ">> G4DNAMolecularReaction::MakeReaction " << G4endl
213 // << "\t track A ("<< trackA.GetTrackID() <<"): " << trackA.GetGlobalTime() << G4endl
214 // << "\t track B ("<< trackB.GetTrackID() <<"): " << trackB.GetGlobalTime() << G4endl;
215 
216  productTrack->SetTrackStatus(fAlive);
217 
218  fpChanges->AddSecondary(productTrack);
219  G4MoleculeFinder::Instance()->Push(productTrack);
220  }
221  }
222 
223  fpChanges->KillParents(true);
224 
225  return fpChanges;
226 }
void SetTrackStatus(const G4TrackStatus aTrackStatus)
virtual void Push(G4Track *track)
const G4DNAMolecularReactionData * GetReactionData(G4MolecularConfiguration *, G4MolecularConfiguration *) const
const G4ThreeVector & GetPosition() const
static G4ITFinder * Instance()
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:576
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
void Initialize(const G4Track &, const G4Track &, G4VParticleChange *particleChangeA=0, G4VParticleChange *particleChangeB=0)
G4int GetTrackID() const
G4double GetGlobalTime() const
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
const G4DNAMolecularReactionTable *& fMolReactionTable
G4Track * BuildTrack(G4double globalTime, const G4ThreeVector &Position)
Definition: G4Molecule.cc:391
G4ITReactionChange * fpChanges
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4MolecularConfiguration * GetProduct(G4int i) const
void AddSecondary(G4Track *aSecondary)

Here is the call graph for this function:

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

Assignment operator

Parameters
otherObject to assign from
Returns
A reference to this

Definition at line 79 of file G4DNAMolecularReaction.cc.

80 {
81  if (this == &rhs) return *this; // handle self assignment
82 
83  fVerbose = rhs.fVerbose;
84  fMolReactionTable = rhs.fMolReactionTable;
85  fReactionRadius = -1;
86  fDistance = -1;
87 
88  //assignment operator
89  return *this;
90 }
const G4DNAMolecularReactionTable *& fMolReactionTable
void G4DNAMolecularReaction::SetReactionModel ( G4VDNAReactionModel model)
inline

Definition at line 111 of file G4DNAMolecularReaction.hh.

112 {
114 }
G4VDNAReactionModel * fReactionModel
const XML_Char XML_Content * model
Definition: expat.h:151
void G4DNAMolecularReaction::SetReactionTable ( const G4DNAMolecularReactionTable )
inline
void G4DNAMolecularReaction::SetVerbose ( int  verb)
inline

Definition at line 116 of file G4DNAMolecularReaction.hh.

117 {
118  fVerbose = verb;
119 }
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 92 of file G4DNAMolecularReaction.cc.

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

Here is the call graph for this function:

Member Data Documentation

G4double G4DNAMolecularReaction::fDistance
protected

Definition at line 107 of file G4DNAMolecularReaction.hh.

const G4DNAMolecularReactionTable*& G4DNAMolecularReaction::fMolReactionTable
protected

Definition at line 103 of file G4DNAMolecularReaction.hh.

G4VDNAReactionModel* G4DNAMolecularReaction::fReactionModel
protected

Definition at line 104 of file G4DNAMolecularReaction.hh.

G4double G4DNAMolecularReaction::fReactionRadius
protected

Definition at line 106 of file G4DNAMolecularReaction.hh.

G4int G4DNAMolecularReaction::fVerbose
protected

Definition at line 105 of file G4DNAMolecularReaction.hh.


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