Geant4  10.01.p01
G4DNAMolecularReaction.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4DNAMolecularReaction.cc 87061 2014-11-24 11:43:34Z gcosmo $
27 //
28 // Author: Mathieu Karamitros (kara@cenbg.in2p3.fr)
29 //
30 // WARNING : This class is released as a prototype.
31 // It might strongly evolve or even disapear in the next releases.
32 //
33 // History:
34 // -----------
35 // 10 Oct 2011 M.Karamitros created
36 //
37 // -------------------------------------------------------------------
38 
41 #include "G4VDNAReactionModel.hh"
42 #include "G4MoleculeFinder.hh"
43 #include "G4UnitsTable.hh"
44 
45 using namespace std;
46 
49  fMolReactionTable(
50  reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
51 {
52  //ctor
53  fVerbose = 0;
54  fReactionModel = 0;
55  fReactionRadius = -1;
56  fDistance = -1;
57 }
58 
60 {
61  //dtor
62  if (fpChanges) delete fpChanges;
63 }
64 
66  G4VITReactionProcess(other),
67  fMolReactionTable(
68  reference_cast<const G4DNAMolecularReactionTable*>(fpReactionTable))
69 {
70  //copy ctor
71  fVerbose = other.fVerbose;
73  fReactionModel = 0;
74  fReactionRadius = -1;
75  fDistance = -1;
76 }
77 
79 {
80  if (this == &rhs) return *this; // handle self assignment
81 
82  fVerbose = rhs.fVerbose;
84  fReactionRadius = -1;
85  fDistance = -1;
86 
87  //assignment operator
88  return *this;
89 }
90 
92  const G4Track& trackB,
93  const double currentStepTime,
94  const double /*previousStepTime*/,
95  bool userStepTimeLimit) /*const*/
96 {
97  G4Molecule* moleculeA = GetMolecule(trackA);
98  G4Molecule* moleculeB = GetMolecule(trackB);
99 
100  if (!fReactionModel)
101  {
102  G4ExceptionDescription exceptionDescription(
103  "You have to give a reaction model to the molecular reaction process");
104  G4Exception("G4DNAMolecularReaction::TestReactibility",
105  "MolecularReaction001", FatalErrorInArgument,
106  exceptionDescription);
107  return false; // makes coverity happy
108  }
109  if (fMolReactionTable == 0)
110  {
111  G4ExceptionDescription exceptionDescription(
112  "You have to give a reaction table to the molecular reaction process");
113  G4Exception("G4DNAMolecularReaction::TestReactibility",
114  "MolecularReaction002", FatalErrorInArgument,
115  exceptionDescription);
116  return false; // makes coverity happy
117  }
118 
119  // Retrieve reaction range
120  fReactionRadius = -1; // reaction Range
121  fReactionRadius = fReactionModel->GetReactionRadius(moleculeA, moleculeB);
122 
123  fDistance = -1; // separation distance
124 
125  if (currentStepTime == 0.)
126  {
127  userStepTimeLimit = false;
128  }
129 
130  G4bool output = fReactionModel->FindReaction(trackA, trackB, fReactionRadius,
131  fDistance, userStepTimeLimit);
132 
133 #ifdef G4VERBOSE
134  // DEBUG
135  if (fVerbose > 1)
136  {
137  G4cout << "\033[1;39;36m" << "G4MolecularReaction " << G4endl;
138  G4cout << "FindReaction returned : " << G4BestUnit(output,"Length") << G4endl;
139  G4cout << " reaction radius : " << G4BestUnit(fReactionRadius,"Length")
140  << " real distance : " << G4BestUnit((trackA.GetPosition() - trackB.GetPosition()).mag(), "Length")
141  << " calculated distance by model (= -1 if real distance > reaction radius and the user limitation step is not reached) : "
142  << G4BestUnit(fDistance,"Length")
143  << G4endl;
144 
145  G4cout << "TrackID A : " << trackA.GetTrackID()
146  << ", TrackID B : " << trackB.GetTrackID()
147  << " | MolA " << moleculeA->GetName()
148  << ", MolB " << moleculeB->GetName()
149  <<"\033[0m\n"
150  << G4endl;
151  G4cout << "--------------------------------------------" << G4endl;
152  }
153 #endif
154  return output;
155 }
156 
158  const G4Track& trackB)
159 {
161  fpChanges->Initialize(trackA, trackB);
162 
163  G4Molecule* moleculeA = GetMolecule(trackA);
164  G4Molecule* moleculeB = GetMolecule(trackB);
165 
166 #ifdef G4VERBOSE
167  // DEBUG
168  if (fVerbose)
169  {
170  G4cout << "G4DNAMolecularReaction::MakeReaction" << G4endl;
171  G4cout<<"TrackA n°" << trackA.GetTrackID()
172  <<"\t | Track B n°" << trackB.GetTrackID() << G4endl;
173 
174  G4cout<<"Track A : Position : " << G4BestUnit(trackA.GetPosition(),"Length")
175  <<"\t Global Time : " << G4BestUnit(trackA.GetGlobalTime(), "Time")<< G4endl;
176 
177  G4cout<<"Track B : Position : " << G4BestUnit(trackB.GetPosition() ,"Length")
178  <<"\t Global Time : " << G4BestUnit(trackB.GetGlobalTime(), "Time")<< G4endl;
179 
180  G4cout<<"Reaction range : " << G4BestUnit(fReactionRadius,"Length")
181  << " \t Separation distance : " << G4BestUnit(fDistance,"Length")<< G4endl;
182  G4cout << "--------------------------------------------" << G4endl;
183  }
184 #endif
185 
186  const G4DNAMolecularReactionData* reactionData = fMolReactionTable
187  ->GetReactionData(moleculeA, moleculeB);
188 
189  G4int nbProducts = reactionData->GetNbProducts();
190 
191  if (nbProducts)
192  {
193  G4double D1 = moleculeA->GetDiffusionCoefficient();
194  G4double D2 = moleculeB->GetDiffusionCoefficient();
195  G4double sqrD1 = sqrt(D1);
196  G4double sqrD2 = sqrt(D2);
197  G4double numerator = sqrD1 + sqrD2;
198  G4ThreeVector reactionSite = sqrD1 / numerator * trackA.GetPosition()
199  + sqrD2 / numerator * trackB.GetPosition();
200 
201  for (G4int j = 0; j < nbProducts; j++)
202  {
203  G4Molecule* product = new G4Molecule(*reactionData->GetProduct(j));
204  G4Track* productTrack = product->BuildTrack(trackA.GetGlobalTime(),
205  reactionSite);
206 
207 // G4cout << ">> G4DNAMolecularReaction::MakeReaction " << G4endl
208 // << "\t track A ("<< trackA.GetTrackID() <<"): " << trackA.GetGlobalTime() << G4endl
209 // << "\t track B ("<< trackB.GetTrackID() <<"): " << trackB.GetGlobalTime() << G4endl;
210 
211  productTrack->SetTrackStatus(fAlive);
212 
213  fpChanges->AddSecondary(productTrack);
214  G4MoleculeFinder::Instance()->Push(productTrack);
215 // G4ITManager<G4Molecule>::Instance()->Push(productTrack);
216  }
217  }
218 
219  fpChanges->KillParents(true);
220 
221  return fpChanges;
222 }
void SetTrackStatus(const G4TrackStatus aTrackStatus)
G4DNAMolecularReaction & operator=(const G4DNAMolecularReaction &other)
Assignment operator.
virtual G4bool TestReactibility(const G4Track &, const G4Track &, const double currentStepTime, const double previousStepTime, bool userStepTimeLimit)
Given two tracks, it tells you whether they can react.
Similar to G4ParticleChange, but deal with two tracks rather than one.
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
CLHEP::Hep3Vector G4ThreeVector
const G4ThreeVector & GetPosition() const
G4DNAMolecularReactionTable sorts out the G4DNAMolecularReactionData for bimolecular reaction...
ReturnType & reference_cast(OriginalType &source)
G4VITReactionProcess defines the reaction between two G4IT.
G4double GetDiffusionCoefficient() const
Returns the diffusion coefficient D.
Definition: G4Molecule.cc:465
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:78
const G4String & GetName() const
Returns the name of the molecule.
Definition: G4Molecule.cc:303
const G4Molecule * GetProduct(G4int i) const
G4GLOB_DLL std::ostream G4cout
virtual ~G4DNAMolecularReaction()
Default destructor.
void Initialize(const G4Track &, const G4Track &, G4VParticleChange *particleChangeA=0, G4VParticleChange *particleChangeB=0)
bool G4bool
Definition: G4Types.hh:79
G4int GetTrackID() const
G4double GetGlobalTime() const
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:67
const G4DNAMolecularReactionTable *& fMolReactionTable
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4DNAMolecularReactionData contains the information relative to a given reaction (eg : °OH + °OH -> H...
G4Track * BuildTrack(G4double globalTime, const G4ThreeVector &Position)
Definition: G4Molecule.cc:328
static G4ITFinder * Instance()
const G4DNAMolecularReactionData * GetReactionData(const G4Molecule *, const G4Molecule *) const
virtual G4bool FindReaction(const G4Track &, const G4Track &, const G4double, G4double &, const G4bool)=0
G4DNAMolecularReaction is the reaction process used in G4DNAMolecularStepByStepModel between two mole...
G4DNAMolecularReaction()
Default constructor.
virtual G4double GetReactionRadius(const G4Molecule *, const G4Molecule *)=0
G4VDNAReactionModel * fReactionModel
G4ITReactionChange * fpChanges
#define G4endl
Definition: G4ios.hh:61
virtual G4ITReactionChange * MakeReaction(const G4Track &, const G4Track &)
Will generate the products of the two given tracks.
Class Description The dynamic molecule holds all the data that change for a molecule It has a pointer...
Definition: G4Molecule.hh:93
double G4double
Definition: G4Types.hh:76
void AddSecondary(G4Track *aSecondary)
virtual void Push(G4Track *track)