Geant4-11
G4DNAScavengerProcess.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//
28#include <G4VScheduler.hh>
29#include <memory>
30#include "G4Molecule.hh"
33#include "G4UnitsTable.hh"
37#include "G4DNABoundingBox.hh"
39#include "G4MoleculeFinder.hh"
40#include "G4Scheduler.hh"
41
42#ifndef State
43# define State(theXInfo) (GetState<G4DNAScavengerProcessState>()->theXInfo)
44#endif
46 const G4DNABoundingBox& box,
47 G4ProcessType type)
48 : G4VITProcess(aName, type)
49 , fpBoundingBox(&box)
50 , fpScavengerMaterial(nullptr)
51 , fVerbose(1)
52{
54 enableAtRestDoIt = false;
55 enableAlongStepDoIt = false;
56 enablePostStepDoIt = true;
59 fIsInitialized = false;
61 fpMaterialConf = nullptr;
62 fProposesTimeStep = true;
64 verboseLevel = 0;
65}
66
68{
69 for(auto& iter : fConfMap)
70 {
71 for(auto& iter2 : iter.second)
72 {
73 delete iter2.second;
74 }
75 }
76}
77
80{
82 fIsInGoodMaterial = false;
83}
85{
88 if(fpScavengerMaterial == nullptr)
89 {
90 return;
91 }
92 fIsInitialized = true;
93}
94
96{
98 G4VITProcess::fpState = std::make_shared<G4DNAScavengerProcessState>();
100}
101
103{
105 {
106 G4ExceptionDescription exceptionDescription;
107 exceptionDescription
108 << "G4DNASecondOrderReaction was already initialised. ";
109 exceptionDescription << "You cannot set a reaction after initialisation.";
110 G4Exception("G4DNASecondOrderReaction::SetReaction",
111 "G4DNASecondOrderReaction001", FatalErrorInArgument,
112 exceptionDescription);
113 }
114 auto materialConf = pData->GetReactant1() == molConf ? pData->GetReactant2()
115 : pData->GetReactant1();
116 if(fVerbose > 0)
117 {
118 G4cout << "G4DNAScavengerProcess::SetReaction : " << molConf->GetName()
119 << " materialConf : " << materialConf->GetName() << G4endl;
120 }
121 fConfMap[molConf][materialConf] = pData;
122}
123
125 const G4Track& track, G4double /*previousStepSize*/,
126 G4ForceCondition* pForceCond)
127{
128 G4Molecule* molecule = GetMolecule(track);
129 auto molConf = molecule->GetMolecularConfiguration();
130 // reset
131 fpMolecularConfiguration = nullptr;
132 fpMaterialConf = nullptr;
133
134 // this because process for moleculeDifinition not for configuration
135 // TODO: need change this
136 auto it = fConfMap.find(molConf);
137 if(it == fConfMap.end())
138 {
139 return DBL_MAX;
140 }
141
142 fpMolecularConfiguration = molConf;
143 auto MaterialMap = it->second;
144
145 G4double r1 = G4UniformRand();
146 std::map<G4double /*Propensity*/, std::pair<MolType, G4double>>
147 ReactionDataMap;
148 G4double alpha0 = 0;
149
150 for(const auto& mat_it : MaterialMap)
151 {
152 auto matConf = mat_it.first;
153 G4double numMol =
155 matConf);
156 if(numMol == 0.0) // ie : not found
157 {
158 continue;
159 }
160 if(fVerbose > 1)
161 {
162 G4cout << " Material of " << matConf->GetName() << " : " << numMol
163 << G4endl;
164 }
165 // auto data = fReactionMap[mat_it];
166 auto data = mat_it.second;
167 auto reactionRate = data->GetObservedReactionRateConstant(); //_const
168 G4double propensity =
169 numMol * reactionRate / (fpBoundingBox->Volume() * Avogadro);
170 auto reactionData = std::make_pair(mat_it.first, propensity);
171 if(propensity == 0)
172 {
173 continue;
174 }
175 alpha0 += propensity;
176 ReactionDataMap[alpha0] = reactionData;
177 }
178 if(alpha0 == 0)
179 {
180 if(State(fIsInGoodMaterial))
181 {
183 State(fIsInGoodMaterial) = false;
184 }
185 return DBL_MAX;
186 }
187 auto rSelectedIter = ReactionDataMap.upper_bound(r1 * alpha0);
188
189 fpMaterialConf = rSelectedIter->second.first;
190
191 State(fIsInGoodMaterial) = true;
192 G4double previousTimeStep(-1.);
193
194 if(State(fPreviousTimeAtPreStepPoint) != -1)
195 {
196 previousTimeStep =
197 track.GetGlobalTime() - State(fPreviousTimeAtPreStepPoint);
198 }
199
200 State(fPreviousTimeAtPreStepPoint) = track.GetGlobalTime();
201
202 // condition is set to "Not Forced"
203 *pForceCond = NotForced;
204
205 if((previousTimeStep <= 0.0) ||
206 (fpState->theNumberOfInteractionLengthLeft <= 0.0))
207 {
209 }
210 else if(previousTimeStep > 0.0)
211 {
213 }
214
215 fpState->currentInteractionLength = 1 / (rSelectedIter->second.second);
216 G4double value = DBL_MAX;
217 if(fpState->currentInteractionLength < DBL_MAX)
218 {
219 value = fpState->theNumberOfInteractionLengthLeft *
220 (fpState->currentInteractionLength);
221 }
222#ifdef G4VERBOSE
223 if(verboseLevel > 2)
224 {
225 G4cout << "G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength:: "
226 << molConf->GetName() << G4endl;
227 G4cout << "theNumberOfInteractionLengthLeft : "
228 << fpState->theNumberOfInteractionLengthLeft << G4endl;
229 G4cout << "currentInteractionLength : " << fpState->currentInteractionLength
230 << G4endl;
231 G4cout << "Material : " << fpMaterialConf->GetName()
232 << " ID: " << track.GetTrackID()
233 << " Track Time : " << track.GetGlobalTime()
234 << " name : " << molConf->GetName()
235 << " Track Position : " << track.GetPosition()
236 << " Returned time : " << G4BestUnit(value, "Time") << G4endl;
237 }
238#endif
239
240 if(value < fReturnedValue)
241 {
242 fReturnedValue = value;
243 }
244
245 return value * -1;
246 // multiple by -1 to indicate to the tracking system that we are returning a
247 // time
248}
249
251 const G4Step& /*step*/)
252{
253 G4Molecule* molecule = GetMolecule(track);
254 auto molConf = molecule->GetMolecularConfiguration();
255 if(fpMolecularConfiguration != molConf)
256 {
257 G4cout << "fpMolecularConfiguration : "
259 << " molConf : " << molConf->GetName() << G4endl;
260
261 assert(false);
262 }
263 std::vector<G4Track*> products;
264#ifdef G4VERBOSE
265 if(verboseLevel > 1)
266 {
267 G4cout << "___________" << G4endl;
268 G4cout << ">>> Beginning of G4DNAScavengerProcess verbose" << G4endl;
269 G4cout << ">>> Returned value : " << G4BestUnit(fReturnedValue, "Time")
270 << G4endl;
271 G4cout << ">>> Time Step : "
272 << G4BestUnit(G4VScheduler::Instance()->GetTimeStep(), "Time")
273 << G4endl;
274 G4cout << ">>> Global Time : "
275 << G4BestUnit(G4VScheduler::Instance()->GetGlobalTime(), "Time")
276 << G4endl;
277 G4cout << ">>> Global Time Track : "
278 << G4BestUnit(track.GetGlobalTime(), "Time") << G4endl;
279 G4cout << ">>> Track Position : " << track.GetPosition() << G4endl;
280 G4cout << ">>> Reaction : " << molecule->GetName() << "("
281 << track.GetTrackID() << ") + " << fpMaterialConf->GetName()
282 << G4endl;
283 G4cout << ">>> End of G4DNAScavengerProcess verbose <<<" << G4endl;
284 }
285#endif
286
287 G4double reactionTime = track.GetGlobalTime();
289
290 auto nbSecondaries = data->GetNbProducts();
291
292 for(G4int j = 0; j < nbSecondaries; ++j)
293 {
294 auto pProduct = new G4Molecule(data->GetProduct(j));
295 auto pProductTrack =
296 pProduct->BuildTrack(reactionTime, track.GetPosition());
297 pProductTrack->SetTrackStatus(fAlive);
298 G4ITTrackHolder::Instance()->Push(pProductTrack);
299 G4MoleculeFinder::Instance()->Push(pProductTrack);
300 products.push_back(pProductTrack);
301 }
302
303#ifdef G4VERBOSE
304 if(fVerbose != 0)
305 {
306 G4cout << "At time : " << std::setw(7) << G4BestUnit(reactionTime, "Time")
307 << " Reaction : " << GetIT(track)->GetName() << " ("
308 << track.GetTrackID() << ") + " << fpMaterialConf->GetName() << " ("
309 << "B"
310 << ") -> ";
311 }
312#endif
313 if(nbSecondaries > 0)
314 {
315 for(G4int i = 0; i < nbSecondaries; ++i)
316 {
317#ifdef G4VERBOSE
318 if((fVerbose != 0) && i != 0)
319 {
320 G4cout << " + ";
321 }
322
323 if(fVerbose != 0)
324 {
325 G4cout << GetIT(products.at(i))->GetName() << " ("
326 << products.at(i)->GetTrackID() << ")";
327 }
328#endif
329 }
330 }
331 else
332 {
333#ifdef G4VERBOSE
334 if(fVerbose != 0)
335 {
336 G4cout << "No product";
337 }
338#endif
339 }
340#ifdef G4VERBOSE
341 if(fVerbose != 0)
342 {
343 G4cout << G4endl;
344 }
345#endif
346
351 fpMaterialConf, reactionTime);
352 State(fPreviousTimeAtPreStepPoint) = -1;
353 return &fParticleChange;
354}
#define State(theXInfo)
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ NotForced
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:76
G4ProcessType
#define G4BestUnit(a, b)
@ fAlive
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4double Volume() const
G4double GetNumberMoleculePerVolumeUnitForMaterialConf(MolType) const
void ReduceNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
void StartTracking(G4Track *) override
std::map< MolType, std::map< MolType, Data * > > fConfMap
G4DNAScavengerProcess(const G4String &aName, const G4DNABoundingBox &box, G4ProcessType type=fDecay)
const G4DNABoundingBox * fpBoundingBox
const G4MolecularConfiguration * fpMolecularConfiguration
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4DNAScavengerMaterial * fpScavengerMaterial
void BuildPhysicsTable(const G4ParticleDefinition &) override
void SetReaction(MolType, Data *pData)
G4ParticleChange fParticleChange
virtual void Push(G4Track *track)
static G4ITFinder * Instance()
virtual void Push(G4Track *)
static G4ITTrackHolder * Instance()
virtual const G4String & GetName() const =0
const G4String & GetName() const
const G4MolecularConfiguration * GetMolecularConfiguration() const
Definition: G4Molecule.cc:532
const G4String & GetName() const
Definition: G4Molecule.cc:338
virtual void Initialize(const G4Track &)
static G4Scheduler * Instance()
Definition: G4Scheduler.cc:101
G4VScavengerMaterial * GetScavengerMaterial() const
Definition: G4Scheduler.hh:194
Definition: G4Step.hh:62
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
virtual void StartTracking(G4Track *)
Definition: G4VITProcess.cc:85
void SetInstantiateProcessState(G4bool flag)
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4shared_ptr< G4ProcessState > fpState
virtual void ResetNumberOfInteractionLengthLeft()
G4bool fProposesTimeStep
void ProposeTrackStatus(G4TrackStatus status)
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:359
G4int verboseLevel
Definition: G4VProcess.hh:356
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:360
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
virtual void StartTracking(G4Track *)
Definition: G4VProcess.cc:87
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:361
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
static G4VScheduler * Instance()
Definition: G4VScheduler.cc:48
float Avogadro
Definition: hepunit.py:252
#define DBL_MAX
Definition: templates.hh:62