Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAMakeReaction.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
27
28#include "G4DNAMakeReaction.hh"
31#include "G4Molecule.hh"
32#include "G4MoleculeFinder.hh"
33#include "G4ITReactionChange.hh"
34#include "Randomize.hh"
35#include "G4SystemOfUnits.hh"
36#include "G4ITReaction.hh"
38#include "G4Scheduler.hh"
39#include "G4UnitsTable.hh"
42
53
59
64
66 const G4Track& /*trackB*/,
67 G4double currentStepTime,
68 G4bool /*userStepTimeLimit*/) /*const*/
69{
70 fTimeStep = currentStepTime;
71 return true;
72}
73
74std::unique_ptr<G4ITReactionChange>
76 const G4Track &trackB)
77{
78 auto & tA = const_cast<G4Track&>(trackA);
79 auto & tB = const_cast<G4Track&>(trackB);
80 UpdatePositionForReaction( tA , tB );//TODO: should change it
81
82 std::unique_ptr<G4ITReactionChange> pChanges(new G4ITReactionChange());
83 pChanges->Initialize(trackA, trackB);
84
85 const auto pMoleculeA = GetMolecule(trackA)->GetMolecularConfiguration();
86 const auto pMoleculeB = GetMolecule(trackB)->GetMolecularConfiguration();
87
88 const auto pReactionData = fMolReactionTable->GetReactionData(pMoleculeA, pMoleculeB);
89 const G4int nbProducts = pReactionData->GetNbProducts();
90 //add Equilibrium process for particle-based models.
91 if(fpScavengerMaterial != nullptr) { fpScavengerMaterial->SetEquilibrium(pReactionData, trackA.GetGlobalTime()); }
92
93 // Notify molecule (reaction) counter
94 if (G4MoleculeCounterManager::Instance()->GetIsActive()) {
96 }
97 if (nbProducts != 0)
98 {
99 const G4double D1 = pMoleculeA->GetDiffusionCoefficient();
100 const G4double D2 = pMoleculeB->GetDiffusionCoefficient();
101 const G4double sqrD1 = D1 == 0. ? 0. : std::sqrt(D1);
102 const G4double sqrD2 = D2 == 0. ? 0. : std::sqrt(D2);
103 const G4double inv_numerator = 1./(sqrD1 + sqrD2);
104 const G4ThreeVector reactionSite = sqrD2 * inv_numerator * tA.GetPosition()
105 + sqrD1 * inv_numerator * tB.GetPosition();
106
108 auto randP = (1-u) * tA.GetPosition() + u * tB.GetPosition();
109
110 for (G4int j = 0; j < nbProducts; ++j)
111 {
112 auto product = pReactionData->GetProduct(j);
113
114 if(fpScavengerMaterial != nullptr) {
115 auto isScavenger = fpScavengerMaterial->find(product);
116 if (isScavenger) {
117 fpScavengerMaterial->AddNumberMoleculePerVolumeUnitForMaterialConf(
118 product, trackA.GetGlobalTime());
119 continue;
120 }
121 }
122
123 auto pProduct = new G4Molecule(product);
124 auto pProductTrack = pProduct->BuildTrack(trackA.GetGlobalTime(), (reactionSite + randP)/2);
125 pProductTrack->SetTrackStatus(fAlive);
126 G4ITTrackHolder::Instance()->Push(pProductTrack);
127 pChanges->AddSecondary(pProductTrack);
128 }
129 }
130 pChanges->KillParents(true);
131 return pChanges;
132}
133
135{
136 fpReactionModel = pReactionModel;
137}
138
140 G4Track& trackB)
141{
142 const auto pMoleculeA = GetMolecule(trackA)->GetMolecularConfiguration();
143 const auto pMoleculeB = GetMolecule(trackB)->GetMolecularConfiguration();
144 G4double D1 = pMoleculeA->GetDiffusionCoefficient();
145 G4double D2 = pMoleculeB->GetDiffusionCoefficient();
146
147 G4double reactionRadius = fpReactionModel->GetReactionRadius( pMoleculeA, pMoleculeB );
148 G4ThreeVector p1 = trackA.GetPosition();
149 G4ThreeVector p2 = trackB.GetPosition();
150
151 G4ThreeVector S1 = p1 - p2;
152 G4double distance = S1.mag();
153
154 if(D1 == 0)
155 {
156 trackB.SetPosition(p1);
157 return;
158 }
159 if(D2 == 0)
160 {
161 trackA.SetPosition(p2);
162 return;
163 }
164
165 if(distance == 0)
166 {
167 G4ExceptionDescription exceptionDescription;
168 exceptionDescription << "Two particles are overlap: "
169 <<GetMolecule(trackA)->GetName()
170 <<" and "<<GetMolecule(trackB)->GetName()
171 <<" at "<<trackA.GetPosition();
172 G4Exception("G4DNAMakeReaction::PrepareForReaction()",
173 "G4DNAMakeReaction003",
174 FatalErrorInArgument,exceptionDescription);
175 }
176 S1.setMag(reactionRadius);
177
178 const G4double dt = fTimeStep;//irt - actualize molecule time
179
180 if(dt > 0)// irt > 0
181 {
182 G4double s12 = 2.0 * D1 * dt;
183 G4double s22 = 2.0 * D2 * dt;
184 G4double sigma = s12 + ( s12 * s12 ) / s22;
185 G4double alpha = reactionRadius * distance / (2 * (D1 + D2) * dt );
186
187 G4ThreeVector S2 = (p1 + ( s12 / s22 ) * p2) +
188 G4ThreeVector(G4RandGauss::shoot(0.0, sigma),
189 G4RandGauss::shoot(0.0, sigma),
190 G4RandGauss::shoot(0.0, sigma));
191
192 S1.setPhi(rad * G4UniformRand() * 2.0 * CLHEP::pi);
193
194 S1.setTheta(rad * std::acos( 1.0 + (1. / alpha) *
195 std::log(1.0 - G4UniformRand() *
196 (1.-std::exp(-2.0 * alpha)))));
197
198 const G4ThreeVector R1 = (D1 * S1 + D2 * S2) / (D1 + D2);
199 const G4ThreeVector R2 = D2 * (S2 - S1) / (D1 + D2);
200
201 trackA.SetPosition(R1);
202 trackB.SetPosition(R2);
203 }
204}
205
206
207std::vector<std::unique_ptr<G4ITReactionChange>>
209 const G4double currentStepTime,
210 const G4double globalTime,
211 const G4bool /*reachedUserStepTimeLimit*/)
212{
213 std::vector<std::unique_ptr<G4ITReactionChange>> ReactionInfo;
214 auto stepper = dynamic_cast<G4DNAIndependentReactionTimeStepper*>(fpTimeStepper);
215 if (stepper == nullptr) {
216 return ReactionInfo;
217 }else {
218 G4double StepTime = 0;
219 do {
220 auto pReactionChange = stepper->FindReaction(pReactionSet, StepTime, globalTime);
221 if (pReactionChange != nullptr) {
222// G4cout<<" time : "<<globalTime<<" "<<pReactionChange->GetTrackA()->GetTrackID()
223// <<" + "<<pReactionChange->GetTrackB()->GetTrackID()<<G4endl;
224 ReactionInfo.push_back(std::move(pReactionChange));
225 }
226 else{
227 break;
228 }
229 }while(StepTime == currentStepTime);
230}
231 return ReactionInfo;
232}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:75
ReturnType & reference_cast(OriginalType &source)
CLHEP::Hep3Vector G4ThreeVector
@ fAlive
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
void setTheta(double)
double mag() const
void setMag(double)
void setPhi(double)
std::unique_ptr< G4ITReactionChange > MakeReaction(const G4Track &, const G4Track &) override
G4VITTimeStepComputer * fpTimeStepper
std::vector< std::unique_ptr< G4ITReactionChange > > FindReaction(G4ITReactionSet *, const G4double, const G4double, const G4bool) override
G4DNAScavengerMaterial * fpScavengerMaterial
G4VDNAReactionModel * fpReactionModel
void SetReactionModel(G4VDNAReactionModel *)
void SetTimeStepComputer(G4VITTimeStepComputer *)
void UpdatePositionForReaction(G4Track &, G4Track &)
G4bool TestReactibility(const G4Track &, const G4Track &, G4double currentStepTime, G4bool userStepTimeLimit) override
const G4DNAMolecularReactionTable *& fMolReactionTable
void Push(G4Track *) override
static G4ITTrackHolder * Instance()
void RecordReaction(const G4DNAMolecularReactionData *, G4double, G4int=1)
static G4MoleculeCounterManager * Instance()
const G4MolecularConfiguration * GetMolecularConfiguration() const
const G4String & GetName() const override
static G4Scheduler * Instance()
G4VScavengerMaterial * GetScavengerMaterial() const
void SetPosition(const G4ThreeVector &aValue)
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
const G4ITReactionTable * fpReactionTable