Geant4-11
G4DNAGillespieDirectMethod.cc
Go to the documentation of this file.
1// ********************************************************************
2// * License and Disclaimer *
3// * *
4// * The Geant4 software is copyright of the Copyright Holders of *
5// * the Geant4 Collaboration. It is provided under the terms and *
6// * conditions of the Geant4 Software License, included in the file *
7// * LICENSE and available at http://cern.ch/geant4/license . These *
8// * include a list of copyright holders. *
9// * *
10// * Neither the authors of this software system, nor their employing *
11// * institutes,nor the agencies providing financial support for this *
12// * work make any representation or warranty, express or implied, *
13// * regarding this software system or assume any liability for its *
14// * use. Please see the license in the file LICENSE and URL above *
15// * for the full disclaimer and the limitation of liability. *
16// * *
17// * This code implementation is the result of the scientific and *
18// * technical work of the GEANT4 collaboration. *
19// * By using, copying, modifying or distributing the software (or *
20// * any work based on the software) you agree to acknowledge its *
21// * use in resulting scientific publications, and indicate your *
22// * acceptance of all terms of the Geant4 Software license. *
23// ********************************************************************
24//
25
28#include "Randomize.hh"
30#include <memory>
31#include <tuple>
32#include "G4DNAEventSet.hh"
33#include "G4UnitsTable.hh"
35#include "G4Scheduler.hh"
36#include <cassert>
37
39 : fMolecularReactions(G4DNAMolecularReactionTable::Instance())
40 , fpMesh(nullptr)
41 , fTimeStep(0)
42 , fpEventSet(nullptr)
43 , fVerbose(0)
44 , fpScavengerMaterial(nullptr)
45{}
46
48
50{
51 fpEventSet = pEventSet;
52}
53
54//#define DEBUG 1
55
57{
58 auto LengthY = fpMesh->GetBoundingBox(index).Getyhi() -
60 auto LengthX = fpMesh->GetBoundingBox(index).Getxhi() -
62 auto LengthZ = fpMesh->GetBoundingBox(index).Getzhi() -
64 G4double V = LengthY * LengthX * LengthZ;
65 assert(V > 0);
66 return V;
67}
69 MolType moleType)
70{
71 if(moleType->GetDiffusionCoefficient() == 0)
72 {
73 return 0.;
74 }
75 const auto& node = fpMesh->GetVoxelMapList(index);
76 G4double alpha = 0;
77 auto it = node.find(moleType);
78 if(it != node.end())
79 {
80 auto LengthY = fpMesh->GetBoundingBox(index).Getyhi() -
82 G4double d = it->first->GetDiffusionCoefficient() / std::pow(LengthY, 2);
83 alpha = d * it->second;
84
85#ifdef DEBUG
86 G4cout << it->first->GetName() << " " << it->second
87 << " D : " << it->first->GetDiffusionCoefficient()
88 << " LengthY : " << LengthY << " PropensityFunction : " << alpha
89 << G4endl;
90#endif
91 }
92 return alpha;
93}
94
96 ReactionData* data)
97{
98 G4double value;
99 auto ConfA = data->GetReactant1();
100 auto ConfB = data->GetReactant2();
101 G4double scavengerNumber = 0;
102 auto typeANumber = FindScavenging(index, ConfA, scavengerNumber)
103 ? scavengerNumber
104 : ComputeNumberInNode(index, ConfA);
105
106 auto typeBNumber = FindScavenging(index, ConfB, scavengerNumber)
107 ? scavengerNumber
108 : ComputeNumberInNode(index, ConfB);
109
110 if(typeANumber == 0 || typeBNumber == 0)
111 {
112 return 0;
113 }
114
115 auto k =
117 if(ConfA == ConfB)
118 {
119 value = typeANumber * (typeBNumber - 1) * k;
120 }
121 else
122 {
123 value = typeANumber * typeBNumber * k;
124 }
125
126 if(value < 0)
127 {
128 G4cout << "G4DNAGillespieDirectMethod::PropensityFunction for : "
129 << ConfA->GetName() << "(" << typeANumber << ") + "
130 << ConfB->GetName() << "(" << typeBNumber
131 << ") : propensity : " << value
132 << " GetObservedReactionRateConstant : "
134 << " GetEffectiveReactionRadius : "
135 << G4BestUnit(data->GetEffectiveReactionRadius(), "Length")
136 << " k : " << k << " volume : " << VolumeOfNode(index)
137 << " Index : " << index << G4endl;
138 assert(false);
139 }
140
141#ifdef DEBUG
142 if(value > 0)
143 G4cout << "G4DNAGillespieDirectMethod::PropensityFunction for : "
144 << ConfA->GetName() << "(" << typeANumber << ") + "
145 << ConfB->GetName() << "(" << typeBNumber
146 << ") : propensity : " << value
147 << " Time to Reaction : " << G4BestUnit(timeToReaction, "Time")
148 << " k : " << k << " Index : " << index << G4endl;
149#endif
150
151 return value;
152}
153
155{
156 // for Scavenger
159
160 auto begin = fpMesh->begin();
161 auto end = fpMesh->end();
162 for(; begin != end; begin++)
163 {
164 auto key = begin->first;
165#ifdef DEBUG
167#endif
168 CreateEvent(key);
169 }
170}
171
173{
174 fTimeStep = stepTime;
175}
177{
178 G4double r1 = G4UniformRand();
179 G4double r2 = G4UniformRand();
180 auto index = fpMesh->GetIndex(key);
181 G4double dAlpha0 = DiffusiveJumping(index);
182 G4double rAlpha0 = Reaction(index);
183 G4double alphaTotal = dAlpha0 + rAlpha0;
184
185 if(alphaTotal == 0)
186 {
187 return;
188 }
189 auto timeStep = ((1.0 / (alphaTotal)) * std::log(1.0 / r1)) + fTimeStep;
190
191#ifdef DEBUG
192 G4cout << "r2 : " << r2 << " rAlpha0 : " << rAlpha0
193 << " dAlpha0 : " << dAlpha0 << " rAlpha0 / (dAlpha0 + rAlpha0) : "
194 << rAlpha0 / (dAlpha0 + rAlpha0) << G4endl;
195#endif
196 if(r2 < rAlpha0 / alphaTotal)
197 {
198 if(fVerbose > 1)
199 {
200 G4cout << "=>>>>reaction at : " << timeStep << " timeStep : "
201 << G4BestUnit(((1.0 / alphaTotal) * std::log(1.0 / r1)), "Time")
202 << G4endl;
203 }
204 auto rSelectedIter = fReactionDataMap.upper_bound(r2 * alphaTotal);
205 fpEventSet->CreateEvent(timeStep, key, rSelectedIter->second);
206 }
207 else if(dAlpha0 > 0)
208 {
209 if(fVerbose > 1)
210 {
211 G4cout << "=>>>>jumping at : " << timeStep << " timeStep : "
212 << G4BestUnit(((1.0 / alphaTotal) * std::log(1.0 / r1)), "Time")
213 << G4endl;
214 }
215
216 auto dSelectedIter = fJumpingDataMap.upper_bound(r2 * alphaTotal - rAlpha0);
217 auto pDSelected =
218 std::make_unique<std::pair<MolType, Index>>(dSelectedIter->second);
219 fpEventSet->CreateEvent(timeStep, key, std::move(pDSelected));
220 }
221#ifdef DEBUG
222 G4cout << G4endl;
223#endif
224}
225
227{
228 fReactionDataMap.clear();
229 G4double alpha0 = 0;
231 if(dataList.empty())
232 {
233 G4cout << "MolecularReactionTable empty" << G4endl;
234 assert(false);
235 }
236 for(const auto& it : dataList)
237 {
238 auto propensity = PropensityFunction(index, it);
239 if(propensity == 0)
240 {
241 continue;
242 }
243 alpha0 += propensity;
244 fReactionDataMap[alpha0] = it;
245 }
246#ifdef DEBUG
247 G4cout << "Reaction :alpha0 : " << alpha0 << G4endl;
248#endif
249 return alpha0;
250}
251
253{
254 fJumpingDataMap.clear();
255 G4double alpha0 = 0;
256 auto NeighboringVoxels = fpMesh->FindNeighboringVoxels(index);
257 if(NeighboringVoxels.empty())
258 {
259 return 0;
260 }
262 while(iter())
263 {
264 const auto conf = iter.value();
265 auto propensity = PropensityFunction(index, conf);
266 if(propensity == 0)
267 {
268 continue;
269 }
270 for(const auto& it_Neighbor : NeighboringVoxels)
271 {
272 alpha0 += propensity;
273 fJumpingDataMap[alpha0] = std::make_pair(conf, it_Neighbor);
274#ifdef DEBUG
275 G4cout << "mole : " << conf->GetName()
276 << " number : " << ComputeNumberInNode(index, conf)
277 << " propensity : " << propensity << " alpha0 : " << alpha0
278 << G4endl;
279#endif
280 }
281 }
282#ifdef DEBUG
283 G4cout << "DiffusiveJumping :alpha0 : " << alpha0 << G4endl;
284#endif
285 return alpha0;
286}
287
289 const Index& index, MolType type) // depend node ?
290{
291 if(type->GetDiffusionCoefficient() != 0)
292 {
293 const auto& node = fpMesh->GetVoxelMapList(index);
294 const auto& it = node.find(type);
295 return (it != node.end()) ? (it->second) : 0;
296 }
297 else
298 {
299 return 0;
300 }
301}
302
304 MolType moletype,
305 G4double& numberOfScavenger)
306{
307 numberOfScavenger = 0;
308 if(fpScavengerMaterial == nullptr)
309 {
310 return false;
311 }
312 auto volumeOfNode = VolumeOfNode(index);
313 if(G4MoleculeTable::Instance()->GetConfiguration("H2O") == moletype)
314 {
315 auto factor = Avogadro * volumeOfNode;
316 numberOfScavenger = factor;
317 return true;
318 }
319
320 G4double totalNumber =
322 moletype);
323 if(totalNumber == 0)
324 {
325 return false;
326 }
327 else
328 {
329 G4double numberInDouble =
330 volumeOfNode * std::floor(totalNumber) / fpMesh->GetBoundingBox().Volume();
331 auto numberInInterg = (int) (std::floor(numberInDouble));
332 G4double ram = G4UniformRand();
333 G4double change = numberInDouble - numberInInterg;
334 if(ram > change)
335 {
336 numberOfScavenger = numberInInterg;
337 }
338 else
339 {
340 numberOfScavenger = numberInInterg + 1;
341 }
342 return true;
343 }
344}
static const G4double alpha
#define G4BestUnit(a, b)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4double Getxlo() const
G4double Volume() const
G4double Getyhi() const
G4double Getylo() const
G4double Getxhi() const
G4double Getzlo() const
G4double Getzhi() const
void CreateEvent(G4double time, Key index, Event::ReactionData *pReactionData)
std::map< G4double, ReactionData * > fReactionDataMap
G4bool FindScavenging(const Index &index, MolType, G4double &)
G4double VolumeOfNode(const Index &index)
G4double DiffusiveJumping(const Index &index)
G4double PropensityFunction(const Index &index, ReactionData *data)
G4double ComputeNumberInNode(const Index &index, MolType type)
G4DNAMolecularReactionTable * fMolecularReactions
void SetTimeStep(const G4double &stepTime)
G4double Reaction(const Index &index)
std::map< G4double, JumpingData > fJumpingDataMap
G4DNAScavengerMaterial * fpScavengerMaterial
void PrintVoxel(const Index &index)
Definition: G4DNAMesh.cc:98
VoxelMap::iterator end()
Definition: G4DNAMesh.hh:137
std::vector< Index > FindNeighboringVoxels(const Index &index) const
Definition: G4DNAMesh.cc:276
G4Voxel::MapList & GetVoxelMapList(Key key)
Definition: G4DNAMesh.cc:37
const G4DNABoundingBox & GetBoundingBox() const
Definition: G4DNAMesh.cc:211
VoxelMap::iterator begin()
Definition: G4DNAMesh.hh:138
Index GetIndex(Key key) const
Definition: G4DNAMesh.cc:164
G4double GetNumberMoleculePerVolumeUnitForMaterialConf(MolType) const
static G4MoleculeTable * Instance()
G4ConfigurationIterator GetConfigurationIterator()
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:101
G4VScavengerMaterial * GetScavengerMaterial() const
Definition: G4Scheduler.hh:194
float Avogadro
Definition: hepunit.py:252