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

#include <G4ChemReboundTransportation.hh>

Inheritance diagram for G4ChemReboundTransportation:

Classes

struct  G4ITBrownianState

Public Member Functions

 G4ChemReboundTransportation (const G4String &aName="ChemReboundTransportation", const G4DNABoundingBox *=nullptr, G4int verbosityLevel=0)
 ~G4ChemReboundTransportation () override=default
 G4ChemReboundTransportation (const G4ChemReboundTransportation &)=delete
G4ChemReboundTransportationoperator= (const G4ChemReboundTransportation &)=delete
void BuildPhysicsTable (const G4ParticleDefinition &) override
void StartTracking (G4Track *aTrack) override
void ComputeStep (const G4Track &, const G4Step &, G4double, G4double &) override
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *) override
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &) override
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &) override
void SetBoundary (const G4DNABoundingBox *)
G4double GetTimeToBoundary (const G4Track &track)
Public Member Functions inherited from G4ITTransportation
 G4ITTransportation (const G4String &aName="ITTransportation", G4int verbosityLevel=0)
 ~G4ITTransportation () override
 G4ITTransportation (const G4ITTransportation &)
void BuildPhysicsTable (const G4ParticleDefinition &) override
void StartTracking (G4Track *aTrack) override
G4bool IsStepLimitedByGeometry ()
G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *) override
G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &) override
G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *pForceCond) override
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData) override
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &) override
G4PropagatorInFieldGetPropagatorInField ()
void SetPropagatorInField (G4PropagatorInField *pFieldPropagator)
void SetVerboseLevel (G4int verboseLevel)
G4int GetVerboseLevel () const
G4double GetThresholdWarningEnergy () const
G4double GetThresholdImportantEnergy () const
G4int GetThresholdTrials () const
void SetThresholdWarningEnergy (G4double newEnWarn)
void SetThresholdImportantEnergy (G4double newEnImp)
void SetThresholdTrials (G4int newMaxTrials)
G4double GetMaxEnergyKilled () const
G4double GetSumEnergyKilled () const
void ResetKilledStatistics (G4int report=1)
void EnableShortStepOptimisation (G4bool optimise=true)
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 &)

Static Public Member Functions

static G4double calculateNextCoordinate (G4double nextPos, G4double high, G4double low)
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

G4ThreeVector BouncingAction (const G4ThreeVector &nextPosition)
G4double calculateDistanceFromTimeStep (MolConf mol, G4double timeStep)
Protected Member Functions inherited from G4ITTransportation
G4bool DoesGlobalFieldExist ()
void SetInstantiateProcessState (G4bool flag)
G4bool InstantiateProcessState ()
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 ()

Protected Attributes

G4ITReactionSetfReactionSet = G4ITReactionSet::Instance()
G4MaterialfNistWater = nullptr
G4double fMaximumTimeStep = 0
G4double fInternalMinTimeStep
const G4DNABoundingBoxfpBoundingBox = nullptr
Protected Attributes inherited from G4ITTransportation
G4ITNavigator * fLinearNavigator
G4PropagatorInFieldfFieldPropagator
G4ParticleChangeForTransport fParticleChange
G4double fThreshold_Warning_Energy
G4double fThreshold_Important_Energy
G4int fThresholdTrials {10}
G4double fUnimportant_Energy
G4double fSumEnergyKilled {0.0}
G4double fMaxEnergyKilled {0.0}
G4bool fShortStepOptimisation {false}
G4ITSafetyHelperfpSafetyHelper
G4int fVerboseLevel
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

Detailed Description

Definition at line 38 of file G4ChemReboundTransportation.hh.

Constructor & Destructor Documentation

◆ G4ChemReboundTransportation() [1/2]

G4ChemReboundTransportation::G4ChemReboundTransportation ( const G4String & aName = "ChemReboundTransportation",
const G4DNABoundingBox * pB = nullptr,
G4int verbosityLevel = 0 )
explicit

Definition at line 63 of file G4ChemReboundTransportation.cc.

66 : G4ITTransportation(aName, verbosity), fpBoundingBox(pB)
67{
68 fVerboseLevel = 0;
69 fpState = std::make_shared<G4ITBrownianState>();
72 fInternalMinTimeStep = 1 * CLHEP::ps;
73}
@ fLowEnergyBrownianTransportation
G4ITTransportation(const G4String &aName="ITTransportation", G4int verbosityLevel=0)
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
G4shared_ptr< G4ProcessState > fpState
void SetProcessSubType(G4int)

Referenced by G4ChemReboundTransportation(), and operator=().

◆ ~G4ChemReboundTransportation()

G4ChemReboundTransportation::~G4ChemReboundTransportation ( )
overridedefault

◆ G4ChemReboundTransportation() [2/2]

G4ChemReboundTransportation::G4ChemReboundTransportation ( const G4ChemReboundTransportation & )
delete

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4ChemReboundTransportation::AlongStepDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Implements G4VProcess.

Definition at line 261 of file G4ChemReboundTransportation.cc.

263{
264 if (GetIT(track)->GetTrackingInfo()->IsLeadingStep()) {
265 // Hoang fixed: - Fixed the state update in AlongStepDoIt (G4ChemReboundTransportation) to
266 // avoid setting incorrect molecule positions (leading tracks)
267 // when the scavenger process is called.
268
269 // G4double spaceStep = DBL_MAX;
270 // const auto molConf = GetMolecule(track)->GetMolecularConfiguration();
271 // spaceStep = calculateDistanceFromTimeStep(molConf, State(theInteractionTimeLeft));
272
273 // State(fGeometryLimitedStep) = false;
274 // State(fTransportEndPosition) =
275 // BouncingAction(track.GetPosition() + spaceStep * G4RandomDirection());
276 // State(fEndPointDistance) = spaceStep;
277 if (fVerboseLevel > 1)
278 //if(GetMolecule(track)->GetName() == "e_aq^-1")
279 {
280 G4cout << "G4DNABrownianTransportation::AlongStepDoIt() : IsLeadingStep " << " trackID : "
281 << track.GetTrackID()
282 << "position : "<<track.GetPosition()
283 << " Molecule name: "
284 << GetMolecule(track)->GetName()<< "State(fTransportEndPosition) : "<<State(fTransportEndPosition) << " State(theInteractionTimeLeft) : "
285 << G4BestUnit(State(theInteractionTimeLeft), "Time")
286 << " Diffusion length : " << G4BestUnit(step.GetStepLength(), "Length")
287 << " within time step : " << G4BestUnit(step.GetDeltaTime(), "Time")
288 << "\t Current global time : " << G4BestUnit(track.GetGlobalTime(), "Time")
289 << " track.GetMomentumDirection() : " << track.GetMomentumDirection()<<" GetIT(track)->GetTrackingInfo()->IsLeadingStep() : "<<GetIT(track)->GetTrackingInfo()->IsLeadingStep() << G4endl;
290 //throw;
291 }
292 }
293
295
296 return &fParticleChange;
297}
#define State(X)
G4IT * GetIT(const G4Track *track)
Definition G4IT.cc:48
G4Molecule * GetMolecule(const G4Track &track)
Definition G4Molecule.cc:75
#define G4BestUnit(a, b)
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForTransport fParticleChange
G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData) override
G4TrackingInformation * GetTrackingInfo()
Definition G4IT.hh:143
const G4String & GetName() const override
G4double GetDeltaTime() const
G4double GetStepLength() const
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
const G4ThreeVector & GetMomentumDirection() const

◆ AlongStepGetPhysicalInteractionLength()

G4double G4ChemReboundTransportation::AlongStepGetPhysicalInteractionLength ( const G4Track & track,
G4double ,
G4double ,
G4double & ,
G4GPILSelection *  )
overridevirtual

Implements G4VProcess.

Definition at line 174 of file G4ChemReboundTransportation.cc.

177{
178 if (!fpBoundingBox->contains(track.GetPosition())) {
180 errMsg << "Point is out of box : " << *fpBoundingBox
181 << " of particle : " << GetIT(track)->GetName() << "(" << track.GetTrackID()
182 << ") : " << track.GetPosition();
184 "ChemReboundTransportation::AlongStepGetPhysicalInteractionLength"
185 "ChemReboundTransportation",
186 "ChemReboundTransportation", FatalErrorInArgument, errMsg);
187 }
188 if (fNistWater != track.GetMaterial()) {
190 errMsg << "This is not water";
192 "ChemReboundTransportation::AlongStepGetPhysicalInteractionLength"
193 "ChemReboundTransportation",
194 "ChemReboundTransportation", FatalErrorInArgument, errMsg);
195 }
196
197 G4double geometryStepLength = DBL_MAX;
198 State(theInteractionTimeLeft) = DBL_MAX;
199
200 auto molConf = GetMolecule(track)->GetMolecularConfiguration();
201 G4ITReactionPerTime& reactionPerTime = fReactionSet->GetReactionsPerTime();
202 auto reaction_i = reactionPerTime.begin();
203 if (reaction_i == reactionPerTime.end()) {
204 State(fGeometryLimitedStep) = false;
205 State(theInteractionTimeLeft) = fMaximumTimeStep;
206 if (fVerboseLevel > 1) {
207 G4cout << "out of reaction " << G4BestUnit(State(theInteractionTimeLeft), "Time") << G4endl;
208 }
209 }
210 else {
211 G4Track* pTrackA = (*reaction_i)->GetReactants().first;
212 G4Track* pTrackB = (*reaction_i)->GetReactant(pTrackA);
213
214 if (&track == pTrackA || &track == pTrackB) {
215 State(theInteractionTimeLeft) = GetTimeToBoundary(track);
216 State(fTimeStepReachedLimit) = false;
217 State(fGeometryLimitedStep) = false;
218
219 if (fVerboseLevel > 1)
220 G4cout << "Molecule A is of type : " << GetMolecule(track)->GetName()
221 << " with trackID : " << track.GetTrackID()
222 << " fMaximumTimeStep : " << G4BestUnit(fMaximumTimeStep, "Time")
223 << " State(theInteractionTimeLeft) : "
224 << G4BestUnit(State(theInteractionTimeLeft), "Time") << G4endl;
225
226 if (State(theInteractionTimeLeft) < fInternalMinTimeStep) {
227 State(fTimeStepReachedLimit) = true;
228 State(theInteractionTimeLeft) = fInternalMinTimeStep;
229 }
230 else if (State(theInteractionTimeLeft) > fMaximumTimeStep) {
231 State(fTimeStepReachedLimit) = true;
232 State(theInteractionTimeLeft) = fMaximumTimeStep;
233 }
234 }
235 else {
236 State(fGeometryLimitedStep) = false;
237 State(theInteractionTimeLeft) = DBL_MAX;
238 }
239 }
240
241 geometryStepLength = calculateDistanceFromTimeStep(molConf, State(theInteractionTimeLeft));
242 State(fTransportEndPosition) =
243 BouncingAction(track.GetPosition() + geometryStepLength * G4RandomDirection());
244
245 State(fTimeStepReachedLimit) = true;
246 State(fCandidateEndGlobalTime) = track.GetGlobalTime() + State(theInteractionTimeLeft);
247 State(fEndGlobalTimeComputed) = true;
248
249#ifdef G4VERBOSE
250 if (fVerboseLevel > 1) {
251 G4cout << "ChemReboundTransportation::AlongStepGetPhysicalInteractionLength = "
252 << G4BestUnit(geometryStepLength, "Length") << " "
253 << G4BestUnit(State(theInteractionTimeLeft), "Time")
254 << " | trackID = " << track.GetTrackID() << G4endl;
255 }
256#endif
257 return geometryStepLength;
258}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
std::multiset< G4ITReactionPtr, compReactionPerTime > G4ITReactionPerTime
G4ThreeVector G4RandomDirection()
double G4double
Definition G4Types.hh:83
G4ThreeVector BouncingAction(const G4ThreeVector &nextPosition)
G4double calculateDistanceFromTimeStep(MolConf mol, G4double timeStep)
G4double GetTimeToBoundary(const G4Track &track)
virtual const G4String & GetName() const =0
const G4MolecularConfiguration * GetMolecularConfiguration() const
G4Material * GetMaterial() const
#define DBL_MAX
Definition templates.hh:62

◆ BouncingAction()

G4ThreeVector G4ChemReboundTransportation::BouncingAction ( const G4ThreeVector & nextPosition)
protected

Definition at line 300 of file G4ChemReboundTransportation.cc.

301{
302 // from Karamitros, Mathieu et al.2020,arXiv:2006.14225 (2020)
303 // https://doi.org/10.48550/arXiv.2006.14225
304
305 G4ThreeVector output;
306
307 G4double RxM = fpBoundingBox->Getxhi();
308 G4double RyM = fpBoundingBox->Getyhi();
309 G4double RzM = fpBoundingBox->Getzhi();
310 G4double Rxm = fpBoundingBox->Getxlo();
311 G4double Rym = fpBoundingBox->Getylo();
312 G4double Rzm = fpBoundingBox->Getzlo();
313
314 G4double x = calculateNextCoordinate(nextPosition.getX(), RxM, Rxm);
315 G4double y = calculateNextCoordinate(nextPosition.getY(), RyM, Rym);
316 G4double z = calculateNextCoordinate(nextPosition.getZ(), RzM, Rzm);
317 output.set(x, y, z);
318 return output;
319}
CLHEP::Hep3Vector G4ThreeVector
double getZ() const
void set(double x, double y, double z)
double getX() const
double getY() const
static G4double calculateNextCoordinate(G4double nextPos, G4double high, G4double low)

Referenced by AlongStepGetPhysicalInteractionLength(), and ComputeStep().

◆ BuildPhysicsTable()

void G4ChemReboundTransportation::BuildPhysicsTable ( const G4ParticleDefinition & particle)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 91 of file G4ChemReboundTransportation.cc.

92{
93 fpSafetyHelper->InitialiseHelper();
95
96 if (fpBoundingBox == nullptr) {
98 errMsg << "fpBoundingBox is nullptr";
100 "ChemReboundTransportation::BuildPhysicsTable"
101 "ChemReboundTransportation",
102 "ChemReboundTransportation", FatalErrorInArgument, errMsg);
103 }
104 G4double halfSize =
105 std::min({fpBoundingBox->halfSideLengthInX(), fpBoundingBox->halfSideLengthInY(),
106 fpBoundingBox->halfSideLengthInZ()});
107 fMaximumTimeStep = (halfSize * halfSize) / (60 * G4H3O::Definition()->GetDiffusionCoefficient());
108}
static G4H3O * Definition()
Definition G4H3O.cc:46
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4ITSafetyHelper * fpSafetyHelper

◆ calculateDistanceFromTimeStep()

G4double G4ChemReboundTransportation::calculateDistanceFromTimeStep ( MolConf mol,
G4double timeStep )
protected

Definition at line 358 of file G4ChemReboundTransportation.cc.

359{
360 G4double diffuCoeff = mol->GetDiffusionCoefficient();
361 if (mol->GetDiffusionCoefficient() < 0) {
362 G4ExceptionDescription exceptionDescription;
363 exceptionDescription << "GetDiffusionCoefficient is negative";
364 G4Exception("ChemReboundTransportation::calculateDistanceFromTimeStep",
365 "ChemReboundTransportation030", FatalErrorInArgument, exceptionDescription);
366 }
367 G4double sqrt_2Dt = sqrt(2 * diffuCoeff * timeStep);
368 G4double x = G4RandGauss::shoot(0, sqrt_2Dt);
369 G4double y = G4RandGauss::shoot(0, sqrt_2Dt);
370 G4double z = G4RandGauss::shoot(0, sqrt_2Dt);
371 return sqrt(x * x + y * y + z * z);
372}

Referenced by AlongStepGetPhysicalInteractionLength(), and ComputeStep().

◆ calculateNextCoordinate()

G4double G4ChemReboundTransportation::calculateNextCoordinate ( G4double nextPos,
G4double high,
G4double low )
static

Definition at line 322 of file G4ChemReboundTransportation.cc.

324{
325 // from Karamitros, Mathieu et al.2020,arXiv:2006.14225 (2020)
326
327 G4double length = high - low;
328 if (std::abs(length) < 1e-10) {
329 return low;
330 }
331
332 G4double relativePos = std::abs(nextPos - low);
333 if (!std::isfinite(relativePos)) {
334 return low; // no crash
335 }
336
337 G4double n = relativePos / length;
338 if (!std::isfinite(n)) {
339 return low;
340 }
341
342 G4double truncVal = std::floor(n);//n is already positive
343 G4double h = truncVal;
344 if(truncVal > 2.0){
345 h = std::fmod(truncVal, 2.0);
346 }
347
348 G4double mod = relativePos;
349
350 if(relativePos > length) {
351 mod = std::fmod(relativePos, length);
352 }
353
354 return low + h * length + (1 - 2 * h) * std::abs(mod);
355}

Referenced by BouncingAction().

◆ ComputeStep()

void G4ChemReboundTransportation::ComputeStep ( const G4Track & track,
const G4Step & step,
G4double timeStep,
G4double & spaceStep )
overridevirtual

Reimplemented from G4ITTransportation.

Definition at line 111 of file G4ChemReboundTransportation.cc.

113{
114 if (GetIT(track)->GetTrackingInfo()->IsLeadingStep()) {
115 G4ExceptionDescription exceptionDescription;
116 exceptionDescription << "ComputeStep is called while the track has"
117 "the minimum interaction time";
118 exceptionDescription << " so it should not recompute a timeStep ";
119 G4Exception("ChemReboundTransportation::ComputeStep", "ChemReboundTransportation0001",
120 FatalErrorInArgument, exceptionDescription);
121 }
122 State(fGeometryLimitedStep) = false;
123 if (timeStep == 0) {
124 State(fTransportEndPosition) = track.GetPosition();
125 spaceStep = 0.;
126 }
127 else {
128 const auto molConf = GetMolecule(track)->GetMolecularConfiguration();
129 spaceStep = calculateDistanceFromTimeStep(molConf, timeStep);
130 }
131 State(fTransportEndPosition) =
132 BouncingAction(track.GetPosition() + spaceStep * G4RandomDirection());
133 State(fEndPointDistance) = (track.GetPosition() - State(fTransportEndPosition)).mag();
134 if (fVerboseLevel > 1)
135 // if(GetMolecule(track)->GetName() == "e_aq^-1")
136 {
137 G4cout << G4endl;
138 G4cout << "ComputeStep: timeStep : " << G4BestUnit(timeStep, "Time")
139 << " State(theInteractionTimeLeft) : " << State(theInteractionTimeLeft)
140 << " State(fEndPointDistance) : " << G4BestUnit(State(fEndPointDistance), "Length")
141 << " trackID : " << track.GetTrackID()
142 << " Molecule name: " << "track.GetPosition() : " << track.GetPosition()
143 << " State(fTransportEndPosition) : " << State(fTransportEndPosition) << " "
144 << GetMolecule(track)->GetName() << " Diffusion length : " << G4endl;
145 }
146 State(fCandidateEndGlobalTime) = step.GetPreStepPoint()->GetGlobalTime() + timeStep;
147 State(fEndGlobalTimeComputed) = true;
148}
G4double GetGlobalTime() const
G4StepPoint * GetPreStepPoint() const

◆ GetTimeToBoundary()

G4double G4ChemReboundTransportation::GetTimeToBoundary ( const G4Track & track)

Definition at line 375 of file G4ChemReboundTransportation.cc.

376{
377 if (!fpBoundingBox->contains(track.GetPosition())) {
379 errMsg << "Point is out of box : " << *fpBoundingBox
380 << " of particle : " << GetIT(track)->GetName() << "(" << track.GetTrackID()
381 << ") : " << track.GetPosition();
383 "BoundedBrownianAction::GetTimeToBoundary"
384 "BoundedBrownianAction",
385 "BoundedBrownianAction", FatalErrorInArgument, errMsg);
386 }
387 auto diffusionCoefficient = GetMolecule(track)->GetDiffusionCoefficient();
388
389 auto dx = std::min(track.GetPosition().getX() - fpBoundingBox->Getxlo(),
390 fpBoundingBox->Getxhi() - track.GetPosition().getX());
391 auto dy = std::min(track.GetPosition().getY() - fpBoundingBox->Getylo(),
392 fpBoundingBox->Getyhi() - track.GetPosition().getY());
393 auto dz = std::min(track.GetPosition().getZ() - fpBoundingBox->Getzlo(),
394 fpBoundingBox->Getzhi() - track.GetPosition().getZ());
395
396 std::vector<G4double> distanceVector{dx, dy, dz};
397 G4double MinTime = DBL_MAX;
398 for (const auto& it : distanceVector) {
399 G4double distance = it;
400 auto random = G4UniformRand();
401 auto minTime = 1 / (4 * diffusionCoefficient) * pow(distance / InvErfc(random), 2);
402
403 if (MinTime > minTime) {
404 MinTime = minTime;
405 }
406 }
407 return MinTime;
408}
#define G4UniformRand()
Definition Randomize.hh:52
G4double GetDiffusionCoefficient() const

Referenced by AlongStepGetPhysicalInteractionLength().

◆ operator=()

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

◆ PostStepDoIt()

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

Implements G4VProcess.

Definition at line 151 of file G4ChemReboundTransportation.cc.

153{
155
156#ifdef G4VERBOSE
157 // DEBUG
158 if (fVerboseLevel > 1)
159 // if(GetMolecule(track)->GetName() == "e_aq^-1")
160 if (GetIT(track)->GetTrackingInfo()->IsLeadingStep()) {
161 G4cout << "ChemReboundTransportation::PostStepDoIt() :" << " trackID : " << track.GetTrackID()
162 << " Molecule name: " << "prePosition : " << step.GetPreStepPoint()->GetPosition()
163 << " postPostion : " << step.GetPostStepPoint()->GetPosition() << " "
164 << GetMolecule(track)->GetName()
165 << " Diffusion length : " << G4BestUnit(step.GetStepLength(), "Length")
166 << " within time step : " << G4BestUnit(step.GetDeltaTime(), "Time")
167 << "\t Current global time : " << G4BestUnit(track.GetGlobalTime(), "Time")
168 << " track.GetMomentumDirection() : " << track.GetMomentumDirection() << G4endl;
169 }
170#endif
171 return &fParticleChange;
172}
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &) override
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPostStepPoint() const

◆ SetBoundary()

void G4ChemReboundTransportation::SetBoundary ( const G4DNABoundingBox * pBounding)
inline

Definition at line 87 of file G4ChemReboundTransportation.hh.

88{
89 fpBoundingBox = pBounding;
90}

◆ StartTracking()

void G4ChemReboundTransportation::StartTracking ( G4Track * aTrack)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 83 of file G4ChemReboundTransportation.cc.

84{
85 fpState = std::make_shared<G4ITBrownianState>();
88}
void SetInstantiateProcessState(G4bool flag)
void StartTracking(G4Track *aTrack) override

Member Data Documentation

◆ fInternalMinTimeStep

G4double G4ChemReboundTransportation::fInternalMinTimeStep
protected

◆ fMaximumTimeStep

G4double G4ChemReboundTransportation::fMaximumTimeStep = 0
protected

◆ fNistWater

G4Material* G4ChemReboundTransportation::fNistWater = nullptr
protected

◆ fpBoundingBox

const G4DNABoundingBox* G4ChemReboundTransportation::fpBoundingBox = nullptr
protected

◆ fReactionSet

G4ITReactionSet* G4ChemReboundTransportation::fReactionSet = G4ITReactionSet::Instance()
protected

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