Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAScavengerProcess Class Reference

#include <G4DNAScavengerProcess.hh>

+ Inheritance diagram for G4DNAScavengerProcess:

Classes

struct  G4DNAScavengerProcessState
 

Public Types

using MolType = const G4MolecularConfiguration*
 
using Data = const G4DNAMolecularReactionData
 

Public Member Functions

 G4DNAScavengerProcess (const G4String &aName, const G4DNABoundingBox &box, G4ProcessType type=fUserDefined)
 
 ~G4DNAScavengerProcess () override
 
 G4DNAScavengerProcess (const G4DNAScavengerProcess &)=delete
 
G4DNAScavengerProcessoperator= (const G4DNAScavengerProcess &)=delete
 
void StartTracking (G4Track *) override
 
void SetReaction (MolType, Data *pData)
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) override
 
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
- Public Member Functions inherited from G4VITProcess
 G4VITProcess (const G4String &name, G4ProcessType type=fNotDefined)
 
 ~G4VITProcess () override
 
 G4VITProcess (const G4VITProcess &other)
 
G4VITProcessoperator= (const G4VITProcess &other)
 
G4bool operator== (const G4VITProcess &right) const
 
G4bool operator!= (const G4VITProcess &right) const
 
size_t GetProcessID () const
 
G4shared_ptr< G4ProcessState_LockGetProcessState ()
 
void SetProcessState (G4shared_ptr< G4ProcessState_Lock > aProcInfo)
 
void ResetProcessState ()
 
void StartTracking (G4Track *) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4double GetInteractionTimeLeft ()
 
void ResetNumberOfInteractionLengthLeft () override
 
G4bool ProposesTimeStep () const
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual const G4VProcessGetCreatorProcess () const
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Attributes

G4bool fIsInitialized
 
G4double fReturnedValue
 
G4ParticleChange fParticleChange
 
const G4MolecularConfigurationfpMolecularConfiguration
 
std::map< MolType, std::map< MolType, Data * > > fConfMap
 
std::vector< MolTypefpMaterialVector
 
MolType fpMaterialConf
 
const G4DNABoundingBoxfpBoundingBox
 
G4DNAScavengerMaterialfpScavengerMaterial {nullptr}
 
- Protected Attributes inherited from G4VITProcess
G4shared_ptr< G4ProcessStatefpState
 
G4bool fProposesTimeStep
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VITProcess
static const size_t & GetMaxProcessIndex ()
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VITProcess
void RetrieveProcessInfo ()
 
void CreateInfo ()
 
template<typename T >
T * GetState ()
 
virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
virtual void ClearInteractionTimeLeft ()
 
virtual void ClearNumberOfInteractionLengthLeft ()
 
void SetInstantiateProcessState (G4bool flag)
 
G4bool InstantiateProcessState ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 35 of file G4DNAScavengerProcess.hh.

Member Typedef Documentation

◆ Data

◆ MolType

Constructor & Destructor Documentation

◆ G4DNAScavengerProcess() [1/2]

G4DNAScavengerProcess::G4DNAScavengerProcess ( const G4String & aName,
const G4DNABoundingBox & box,
G4ProcessType type = fUserDefined )
explicit

Definition at line 45 of file G4DNAScavengerProcess.cc.

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}
const G4DNABoundingBox * fpBoundingBox
const G4MolecularConfiguration * fpMolecularConfiguration
void SetInstantiateProcessState(G4bool flag)
G4VITProcess(const G4String &name, G4ProcessType type=fNotDefined)
G4bool fProposesTimeStep
G4bool enableAtRestDoIt
G4int verboseLevel
G4bool enableAlongStepDoIt
void SetProcessSubType(G4int)
G4bool enablePostStepDoIt
G4VParticleChange * pParticleChange
#define DBL_MAX
Definition templates.hh:62

◆ ~G4DNAScavengerProcess()

G4DNAScavengerProcess::~G4DNAScavengerProcess ( )
override

Definition at line 65 of file G4DNAScavengerProcess.cc.

66{
67 for(auto& iter : fConfMap)
68 {
69 for(auto& iter2 : iter.second)
70 {
71 delete iter2.second;
72 }
73 }
74}
std::map< MolType, std::map< MolType, Data * > > fConfMap

◆ G4DNAScavengerProcess() [2/2]

G4DNAScavengerProcess::G4DNAScavengerProcess ( const G4DNAScavengerProcess & )
delete

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4DNAScavengerProcess::AlongStepDoIt ( const G4Track & ,
const G4Step &  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 78 of file G4DNAScavengerProcess.hh.

79 {
80 return nullptr;
81 }

◆ AlongStepGetPhysicalInteractionLength()

G4double G4DNAScavengerProcess::AlongStepGetPhysicalInteractionLength ( const G4Track & ,
G4double ,
G4double ,
G4double & ,
G4GPILSelection *  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 70 of file G4DNAScavengerProcess.hh.

73 {
74 return -1.0;
75 }

◆ AtRestDoIt()

G4VParticleChange * G4DNAScavengerProcess::AtRestDoIt ( const G4Track & ,
const G4Step &  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 64 of file G4DNAScavengerProcess.hh.

65 {
66 return nullptr;
67 }

◆ AtRestGetPhysicalInteractionLength()

G4double G4DNAScavengerProcess::AtRestGetPhysicalInteractionLength ( const G4Track & ,
G4ForceCondition *  )
inlineoverridevirtual

Implements G4VProcess.

Definition at line 58 of file G4DNAScavengerProcess.hh.

60 {
61 return -1.0;
62 }

◆ BuildPhysicsTable()

void G4DNAScavengerProcess::BuildPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 81 of file G4DNAScavengerProcess.cc.

82{
85 if(fpScavengerMaterial == nullptr)
86 {
87 return;
88 }
89 fIsInitialized = true;
90}
G4DNAScavengerMaterial * fpScavengerMaterial
static G4Scheduler * Instance()
G4VScavengerMaterial * GetScavengerMaterial() const

◆ operator=()

G4DNAScavengerProcess & G4DNAScavengerProcess::operator= ( const G4DNAScavengerProcess & )
delete

◆ PostStepDoIt()

G4VParticleChange * G4DNAScavengerProcess::PostStepDoIt ( const G4Track & track,
const G4Step &  )
overridevirtual

Implements G4VProcess.

Definition at line 247 of file G4DNAScavengerProcess.cc.

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}
#define State(X)
G4IT * GetIT(const G4Track *track)
Definition G4IT.cc:48
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:74
#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
void ReduceNumberMoleculePerVolumeUnitForMaterialConf(MolType, G4double)
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
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
void ProposeTrackStatus(G4TrackStatus status)
static G4VScheduler * Instance()

◆ PostStepGetPhysicalInteractionLength()

G4double G4DNAScavengerProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Implements G4VProcess.

Definition at line 121 of file G4DNAScavengerProcess.cc.

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}
@ NotForced
#define G4UniformRand()
Definition Randomize.hh:52
G4double Volume() const
G4double GetNumberMoleculePerVolumeUnitForMaterialConf(MolType) const
virtual void SubtractNumberOfInteractionLengthLeft(G4double previousStepSize)
G4shared_ptr< G4ProcessState > fpState
void ResetNumberOfInteractionLengthLeft() override

◆ SetReaction()

void G4DNAScavengerProcess::SetReaction ( MolType molConf,
Data * pData )

Definition at line 99 of file G4DNAScavengerProcess.cc.

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}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription

◆ StartTracking()

void G4DNAScavengerProcess::StartTracking ( G4Track * track)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 92 of file G4DNAScavengerProcess.cc.

93{
95 G4VITProcess::fpState = std::make_shared<G4DNAScavengerProcessState>();
97}
void StartTracking(G4Track *) override
virtual void StartTracking(G4Track *)
Definition G4VProcess.cc:87

Member Data Documentation

◆ fConfMap

std::map<MolType , std::map<MolType , Data*> > G4DNAScavengerProcess::fConfMap
protected

◆ fIsInitialized

G4bool G4DNAScavengerProcess::fIsInitialized
protected

Definition at line 93 of file G4DNAScavengerProcess.hh.

Referenced by BuildPhysicsTable(), G4DNAScavengerProcess(), and SetReaction().

◆ fParticleChange

G4ParticleChange G4DNAScavengerProcess::fParticleChange
protected

Definition at line 95 of file G4DNAScavengerProcess.hh.

Referenced by G4DNAScavengerProcess(), and PostStepDoIt().

◆ fpBoundingBox

const G4DNABoundingBox* G4DNAScavengerProcess::fpBoundingBox
protected

Definition at line 101 of file G4DNAScavengerProcess.hh.

Referenced by PostStepGetPhysicalInteractionLength().

◆ fpMaterialConf

MolType G4DNAScavengerProcess::fpMaterialConf
protected

◆ fpMaterialVector

std::vector<MolType> G4DNAScavengerProcess::fpMaterialVector
protected

Definition at line 99 of file G4DNAScavengerProcess.hh.

◆ fpMolecularConfiguration

const G4MolecularConfiguration* G4DNAScavengerProcess::fpMolecularConfiguration
protected

◆ fpScavengerMaterial

G4DNAScavengerMaterial* G4DNAScavengerProcess::fpScavengerMaterial {nullptr}
protected

Definition at line 102 of file G4DNAScavengerProcess.hh.

102{nullptr};

Referenced by BuildPhysicsTable(), PostStepDoIt(), and PostStepGetPhysicalInteractionLength().

◆ fReturnedValue

G4double G4DNAScavengerProcess::fReturnedValue
protected

The documentation for this class was generated from the following files: