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