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