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

#include <G4DNABrownianTransportation.hh>

+ Inheritance diagram for G4DNABrownianTransportation:

Classes

struct  G4ITBrownianState
 

Public Member Functions

 G4DNABrownianTransportation (const G4String &aName="DNABrownianTransportation", G4int verbosityLevel=0)
 
virtual ~G4DNABrownianTransportation ()
 
 G4DNABrownianTransportation (const G4DNABrownianTransportation &)=delete
 
G4DNABrownianTransportationoperator= (const G4DNABrownianTransportation &)=delete
 
void SetBrownianAction (G4BrownianAction *)
 
void SetUserBrownianAction (G4VUserBrownianAction *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void StartTracking (G4Track *aTrack)
 
virtual void ComputeStep (const G4Track &, const G4Step &, const G4double, G4double &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &)
 
void UseMaximumTimeBeforeReachingBoundary (bool flag=true)
 
void UseCumulativeDensitFunction (bool flag=true)
 
void UseLimitingTimeSteps (bool flag=true)
 
void SpeedLevel (int level)
 
- Public Member Functions inherited from G4ITTransportation
 G4ITTransportation (const G4String &aName="ITTransportation", G4int verbosityLevel=0)
 
virtual ~G4ITTransportation ()
 
 G4ITTransportation (const G4ITTransportation &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void ComputeStep (const G4Track &, const G4Step &, const double timeStep, double &spaceStep)
 
virtual void StartTracking (G4Track *aTrack)
 
G4bool IsStepLimitedByGeometry ()
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double, G4ForceCondition *pForceCond)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &)
 
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)
 
virtual ~G4VITProcess ()
 
 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 ()
 
virtual void StartTracking (G4Track *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double GetInteractionTimeLeft ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 BuildPhysicsTable (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 StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
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 Member Functions

G4double ComputeGeomLimit (const G4Track &track, G4double &presafety, G4double limit)
 
void Diffusion (const G4Track &track)
 
- 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

G4bool fUseMaximumTimeBeforeReachingBoundary
 
G4MaterialfNistWater
 
G4bool fUseSchedulerMinTimeSteps
 
G4double fInternalMinTimeStep
 
G4bool fSpeedMeUp
 
const std::vector< G4double > * fpWaterDensity
 
G4BrownianActionfpBrownianAction
 
G4VUserBrownianActionfpUserBrownianAction
 
- Protected Attributes inherited from G4ITTransportation
G4ITNavigator * fLinearNavigator
 
G4PropagatorInFieldfFieldPropagator
 
G4ParticleChangeForTransport fParticleChange
 
G4double fThreshold_Warning_Energy
 
G4double fThreshold_Important_Energy
 
G4int fThresholdTrials
 
G4double fUnimportant_Energy
 
G4double fSumEnergyKilled
 
G4double fMaxEnergyKilled
 
G4bool fShortStepOptimisation
 
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
 

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)
 

Detailed Description

Definition at line 109 of file G4DNABrownianTransportation.hh.

Constructor & Destructor Documentation

◆ G4DNABrownianTransportation() [1/2]

G4DNABrownianTransportation::G4DNABrownianTransportation ( const G4String aName = "DNABrownianTransportation",
G4int  verbosityLevel = 0 
)

Definition at line 120 of file G4DNABrownianTransportation.cc.

121 :
122 G4ITTransportation(aName, verbosity)
123{
124
125 fVerboseLevel = 0;
126
127 fpState.reset(new G4ITBrownianState());
128
129 //ctor
131
133 fpWaterDensity = 0;
134
137 fSpeedMeUp = true;
138
139 fInternalMinTimeStep = 1*picosecond;
141 fpUserBrownianAction = nullptr;
142}
@ fLowEnergyBrownianTransportation
G4VUserBrownianAction * fpUserBrownianAction
const std::vector< G4double > * fpWaterDensity
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
G4shared_ptr< G4ProcessState > fpState
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410

◆ ~G4DNABrownianTransportation()

G4DNABrownianTransportation::~G4DNABrownianTransportation ( )
virtual

Definition at line 144 of file G4DNABrownianTransportation.cc.

145{
147}

◆ G4DNABrownianTransportation() [2/2]

G4DNABrownianTransportation::G4DNABrownianTransportation ( const G4DNABrownianTransportation )
delete

Member Function Documentation

◆ AlongStepDoIt()

G4VParticleChange * G4DNABrownianTransportation::AlongStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Reimplemented from G4ITTransportation.

Definition at line 800 of file G4DNABrownianTransportation.cc.

802{
803#ifdef DEBUG_MEM
804 MemStat mem_first, mem_second, mem_diff;
805#endif
806
807#ifdef DEBUG_MEM
808 mem_first = MemoryUsage();
809#endif
810
811 if (GetIT(track)->GetTrackingInfo()->IsLeadingStep()
812 && State(fComputeLastPosition)
813 && State(fGeometryLimitedStep)//Hoang added 26/8/2019
814 )
815 {
816 //==========================================================================
817 // DEBUG
818 //
819// assert(fabs(State(theInteractionTimeLeft)-
820// G4VScheduler::Instance()->GetTimeStep()) < DBL_EPSILON);
821 //==========================================================================
822
823 G4double spaceStep = DBL_MAX;
824 G4double diffusionCoefficient = GetMolecule(track)->
825 GetDiffusionCoefficient();
826
827 G4double sqrt_2Dt= sqrt(2 * diffusionCoefficient * State(theInteractionTimeLeft));
828 G4double x = G4RandGauss::shoot(0, sqrt_2Dt);
829 G4double y = G4RandGauss::shoot(0, sqrt_2Dt);
830 G4double z = G4RandGauss::shoot(0, sqrt_2Dt);
831
832 if(State(theInteractionTimeLeft) <= fInternalMinTimeStep)
833 {
834 spaceStep = State(fEndPointDistance);
835 State(fGeometryLimitedStep) = true;
836 }
837 else
838 {
839 spaceStep = sqrt(x * x + y * y + z * z);
840
841 if(spaceStep >= State(fEndPointDistance))
842 {
843 State(fGeometryLimitedStep) = true;
844 if (
845 //fSpeedMeUp == false&&
847 && spaceStep >= State(fEndPointDistance))
848 {
849 spaceStep = State(fEndPointDistance);
850 }
851 }
852 else
853 {
854 State(fGeometryLimitedStep) = false;
855 }
856 }
857
858// assert( (spaceStep < State(fEndPointDistance) && State(fGeometryLimitedStep) == false)
859//|| (spaceStep >= State(fEndPointDistance) && State(fGeometryLimitedStep)));
860
861 // Calculate final position
862 //
863 State(fTransportEndPosition) = track.GetPosition()
864 + spaceStep * track.GetMomentumDirection();
865
866 //hoang exp
867 if(fpUserBrownianAction != nullptr)
868 {
869 // Let the user Brownian action class decide what to do
870 G4ThreeVector nextPosition{track.GetPosition().getX() + x,
871 track.GetPosition().getY() + y,
872 track.GetPosition().getZ() + z};
873 fpUserBrownianAction->Transport(nextPosition);
874 State(fTransportEndPosition) = nextPosition;
875 }
876 //hoang exp
877 }
878
879 if(fVerboseLevel)
880 {
882 << "G4DNABrownianTransportation::AlongStepDoIt: "
883 "GeometryLimitedStep = "
884 << State(fGeometryLimitedStep)
885 << RESET_COLOR
886 << G4endl;
887 }
888
889// static G4ThreadLocal G4int noCalls = 0;
890// noCalls++;
891
892// fParticleChange.Initialize(track);
893
895
896#ifdef DEBUG_MEM
897 MemStat mem_intermediaire = MemoryUsage();
898 mem_diff = mem_intermediaire-mem_first;
899 G4cout << "\t\t\t >> || MEM || After calling G4ITTransportation::"
900 "AlongStepDoIt for "<< track.GetTrackID() << ", diff is : "
901 << mem_diff << G4endl;
902#endif
903
904 if(track.GetStepLength() != 0 // && State(fGeometryLimitedStep)
905 //========================================================================
906 // TODO: avoid changing direction after too small time steps
907// && (G4VScheduler::Instance()->GetTimeStep() > fInternalMinTimeStep
908// || fSpeedMeUp == false)
909 //========================================================================
910 )
911 {
912 Diffusion(track);
913 }
914 //else
915 //{
916 // fParticleChange.ProposeMomentumDirection(State(fTransportEndMomentumDir));
917 //}
918/*
919 if (State(fParticleIsLooping))
920 {
921 if ((State(fNoLooperTrials)>= fThresholdTrials))
922 {
923 fParticleChange.ProposeTrackStatus(fStopAndKill);
924 State(fNoLooperTrials) = 0;
925#ifdef G4VERBOSE
926 if ((fVerboseLevel > 1))
927 {
928 G4cout
929 << " G4DNABrownianTransportation is killing track that is looping or stuck "
930 << G4endl;
931 G4cout << " Number of trials = " << State(fNoLooperTrials)
932 << " No of calls to AlongStepDoIt = " << noCalls
933 << G4endl;
934 }
935#endif
936 }
937 else
938 {
939 State(fNoLooperTrials)++;
940 }
941 }
942 else
943 {
944 State(fNoLooperTrials)=0;
945 }
946*/
947#ifdef DEBUG_MEM
948 mem_intermediaire = MemoryUsage();
949 mem_diff = mem_intermediaire-mem_first;
950 G4cout << "\t\t\t >> || MEM || After calling G4DNABrownianTransportation::"
951 "Diffusion for "<< track.GetTrackID() << ", diff is : "
952 << mem_diff << G4endl;
953#endif
954
955 return &fParticleChange;
956}
#define State(X)
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:74
#define GREEN_ON_BLUE
Definition: G4Scheduler.cc:85
#define RESET_COLOR
Definition: G4Scheduler.cc:86
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double getZ() const
double getX() const
double getY() const
void Diffusion(const G4Track &track)
G4ParticleChangeForTransport fParticleChange
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetStepLength() const
virtual void Transport(G4ThreeVector &, G4Track *pTrack=nullptr)=0
MemStat MemoryUsage()
Definition: G4MemStat.cc:55
#define DBL_MAX
Definition: templates.hh:62

◆ AlongStepGetPhysicalInteractionLength()

G4double G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety,
G4GPILSelection selection 
)
virtual

Reimplemented from G4ITTransportation.

Definition at line 569 of file G4DNABrownianTransportation.cc.

574{
575#ifdef G4VERBOSE
576 if(fVerboseLevel)
577 {
578 G4cout << " G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength - track ID: "
579 << track.GetTrackID() << G4endl;
580 G4cout << "In volume : " << track.GetVolume()->GetName()
581 << " position : " << G4BestUnit(track.GetPosition() , "Length") << G4endl;
582 }
583#endif
584
585 G4double geometryStepLength =
587 track, previousStepSize, currentMinimumStep, currentSafety,
588 selection);
589
590 if(geometryStepLength==0)
591 {
592// G4cout << "geometryStepLength==0" << G4endl;
593 if(State(fGeometryLimitedStep))
594 {
595// G4cout << "if(State(fGeometryLimitedStep))" << G4endl;
596 G4TouchableHandle newTouchable = new G4TouchableHistory;
597
598 newTouchable->UpdateYourself(State(fCurrentTouchableHandle)->GetVolume(),
599 State(fCurrentTouchableHandle)->GetHistory());
600
601 fLinearNavigator->SetGeometricallyLimitedStep();
602 fLinearNavigator->LocateGlobalPointAndUpdateTouchableHandle(
603 track.GetPosition(), track.GetMomentumDirection(),
604 newTouchable, true);
605
606 if(newTouchable->GetVolume() == 0)
607 {
608 return 0;
609 }
610
611 State(fCurrentTouchableHandle) = newTouchable;
612
613 //=======================================================================
614 // TODO: speed up navigation update
615// geometryStepLength = fLinearNavigator->ComputeStep(track.GetPosition(),
616// track.GetMomentumDirection(),
617// currentMinimumStep,
618// currentSafety);
619 //=======================================================================
620
621
622 //=======================================
623 // Longer but safer ...
624 geometryStepLength =
626 track, previousStepSize, currentMinimumStep, currentSafety,
627 selection);
628
629 }
630 }
631
632 //============================================================================
633 // DEBUG
634 // G4cout << "geometryStepLength: " << G4BestUnit(geometryStepLength, "Length")
635 // << " | trackID: " << track.GetTrackID()
636 // << G4endl;
637 //============================================================================
638//Hoang exp
639 if(fpUserBrownianAction != nullptr)
640 {
641 // Let the user Brownian action class decide what to do
642 G4double distanceToBoundary = fpUserBrownianAction->GetDistanceToBoundary(track);
643 geometryStepLength = distanceToBoundary;
644 }
645
646//Hoang exp
647
648
649 G4double diffusionCoefficient = 0;
650
651 /*
652 G4Material* material = track.GetMaterial();
653 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
654
655 if (waterDensity == 0.0)
656 {
657 if(fpBrownianAction)
658 {
659 diffusionCoefficient = fpBrownianAction->GetDiffusionCoefficient(material,
660 GetMolecule(track));
661 }
662
663 if(diffusionCoefficient <= 0)
664 {
665 State(fGeometryLimitedStep) = false;
666 State(theInteractionTimeLeft) = 0;
667 State(fTransportEndPosition) = track.GetPosition();
668 return 0;
669 }
670
671 }
672 else
673 {
674 diffusionCoefficient = GetMolecule(track)->GetDiffusionCoefficient();
675 }
676 */
677
678 diffusionCoefficient = GetMolecule(track)->GetDiffusionCoefficient();
679
680 // To avoid divide by zero of diffusionCoefficient
681 if(diffusionCoefficient <= 0)
682 {
683 State(fGeometryLimitedStep) = false;
684 State(theInteractionTimeLeft) = DBL_MAX;
685 State(fTransportEndPosition) = track.GetPosition();
686 return 0;
687 }
688
689
690 State(fComputeLastPosition) = false;
691 State(fTimeStepReachedLimit) = false;
692
693 if (State(fGeometryLimitedStep))
694 {
695 // 95 % of the space step distribution is lower than
696 // d_95 = 2 * sqrt(2*D*t)
697 // where t is the corresponding time step
698 // so by inversion :
700 {
701 if(fSpeedMeUp)
702 {
703 State(theInteractionTimeLeft) = (geometryStepLength * geometryStepLength)
704 / (diffusionCoefficient); // d_50 - use straight line
705 }
706 else
707 {
708 State(theInteractionTimeLeft) = (currentSafety * currentSafety)
709 / (diffusionCoefficient); // d_50 - use safety
710
711 //=====================================================================
712 // State(theInteractionTimeLeft) = (currentSafety * currentSafety)
713 // / (8 * diffusionCoefficient); // d_95
714 //=====================================================================
715 }
716 State(fComputeLastPosition) = true;
717 }
718 else
719 // Will use a random time - this is precise but long to compute in certain
720 // circumstances (many particles - small volumes)
721 {
722 State(fRandomNumber) = G4UniformRand();
723 State(theInteractionTimeLeft) = 1 / (4 * diffusionCoefficient)
724 * pow(geometryStepLength / InvErfc(State(fRandomNumber)),2);
725
726 State(fTransportEndPosition) = geometryStepLength*
727 track.GetMomentumDirection() + track.GetPosition();
728 }
729
731 {
732 G4double minTimeStepAllowed = G4VScheduler::Instance()->GetLimitingTimeStep();
733 //========================================================================
734 // TODO
735// double currentMinTimeStep = G4VScheduler::Instance()->GetTimeStep();
736 //========================================================================
737
738 if (State(theInteractionTimeLeft) < minTimeStepAllowed)
739 {
740 State(theInteractionTimeLeft) = minTimeStepAllowed;
741 State(fTimeStepReachedLimit) = true;
742 State(fComputeLastPosition) = true;
743 }
744 }
745 else if(State(theInteractionTimeLeft) < fInternalMinTimeStep)
746 // TODO: find a better way when fForceLimitOnMinTimeSteps is not used
747 {
748 State(fTimeStepReachedLimit) = true;
749 State(theInteractionTimeLeft) = fInternalMinTimeStep;
751 {
752 State(fComputeLastPosition) = true;
753 }
754 }
755
756 State(fCandidateEndGlobalTime) =
757 track.GetGlobalTime() + State(theInteractionTimeLeft);
758
759 State(fEndGlobalTimeComputed) = true; // MK: ADDED ON 20/11/2014
760
761 State(fPathLengthWasCorrected) = false;
762 }
763 else
764 {
765 // Transform geometrical step
766 geometryStepLength = 2
767 * sqrt(diffusionCoefficient * State(theInteractionTimeLeft))
768 * InvErf(G4UniformRand());
769 State(fPathLengthWasCorrected) = true;
770 //State(fEndPointDistance) = geometryStepLength;
771 State(fTransportEndPosition) = geometryStepLength*
772 track.GetMomentumDirection() + track.GetPosition();
773 }
774
775#ifdef G4VERBOSE
776 // DEBUG
777 if (fVerboseLevel > 1)
778 {
780 << "G4DNABrownianTransportation::AlongStepGetPhysicalInteractionLength = "
781 << G4BestUnit(geometryStepLength, "Length")
782 << " | trackID = "
783 << track.GetTrackID()
784 << RESET_COLOR
785 << G4endl;
786 }
787#endif
788
789// assert(geometryStepLength < State(fEndPointDistance)
790// || (geometryStepLength >= State(fEndPointDistance) && State(fGeometryLimitedStep)));
791
792 return geometryStepLength;
793}
#define G4UniformRand()
Definition: Randomize.hh:52
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4ITNavigator * fLinearNavigator
G4double GetDiffusionCoefficient() const
Definition: G4Molecule.cc:514
G4VPhysicalVolume * GetVolume() const
G4double GetGlobalTime() const
const G4String & GetName() const
static G4VScheduler * Instance()
Definition: G4VScheduler.cc:48
virtual double GetLimitingTimeStep() const
Definition: G4VScheduler.hh:84
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34
virtual void UpdateYourself(G4VPhysicalVolume *pPhysVol, const G4NavigationHistory *history=nullptr)
Definition: G4VTouchable.cc:66
virtual G4double GetDistanceToBoundary(const G4Track &)=0

◆ BuildPhysicsTable()

void G4DNABrownianTransportation::BuildPhysicsTable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4ITTransportation.

Definition at line 167 of file G4DNABrownianTransportation.cc.

168{
169 if(verboseLevel > 0)
170 {
171 G4cout << G4endl<< GetProcessName() << ": for "
172 << setw(24) << particle.GetParticleName()
173 << "\tSubType= " << GetProcessSubType() << G4endl;
174 }
175 // Initialize water density pointer
177 GetDensityTableFor(G4Material::GetMaterial("G4_WATER"));
178
181}
static G4DNAMolecularMaterial * Instance()
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4ITSafetyHelper * fpSafetyHelper
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:691
const G4String & GetParticleName() const
G4int verboseLevel
Definition: G4VProcess.hh:360
G4int GetProcessSubType() const
Definition: G4VProcess.hh:404
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

◆ ComputeGeomLimit()

G4double G4DNABrownianTransportation::ComputeGeomLimit ( const G4Track track,
G4double presafety,
G4double  limit 
)
protected

Definition at line 550 of file G4DNABrownianTransportation.cc.

553{
554 G4double res = DBL_MAX;
555 if(track.GetVolume() != fpSafetyHelper->GetWorldVolume())
556 {
557 G4TrackStateManager& trackStateMan = GetIT(track)->GetTrackingInfo()
559 fpSafetyHelper->LoadTrackState(trackStateMan);
561 track.GetStep()->GetPreStepPoint()->GetPosition(),
562 track.GetMomentumDirection(),
563 limit, presafety);
565 }
566 return res;
567}
G4VPhysicalVolume * GetWorldVolume()
G4double CheckNextStep(const G4ThreeVector &position, const G4ThreeVector &direction, const G4double currentMaxStep, G4double &newSafety)
G4TrackingInformation * GetTrackingInfo()
Definition: G4IT.hh:143
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPreStepPoint() const
virtual void LoadTrackState(G4TrackStateManager &manager)
virtual void ResetTrackState()
const G4Step * GetStep() const
G4TrackStateManager & GetTrackStateManager()

◆ ComputeStep()

void G4DNABrownianTransportation::ComputeStep ( const G4Track track,
const G4Step step,
const G4double  timeStep,
G4double spaceStep 
)
virtual

Reimplemented from G4ITTransportation.

Definition at line 183 of file G4DNABrownianTransportation.cc.

187{
188 // G4cout << "G4ITBrownianTransportation::ComputeStep" << G4endl;
189
190 /* If this method is called, this step
191 * cannot be geometry limited.
192 * In case the step is limited by the geometry,
193 * this method should not be called.
194 * The fTransportEndPosition calculated in
195 * the method AlongStepIL should be taken
196 * into account.
197 * In order to do so, the flag IsLeadingStep
198 * is on. Meaning : this track has the minimum
199 * interaction length over all others.
200 */
201 if(GetIT(track)->GetTrackingInfo()->IsLeadingStep())
202 {
203 const G4VITProcess* ITProc = ((const G4VITProcess*) step.GetPostStepPoint()
205 G4bool makeException = true;
206
207 if(ITProc && ITProc->ProposesTimeStep()) makeException = false;
208
209 if(makeException)
210 {
211
212 G4ExceptionDescription exceptionDescription;
213 exceptionDescription << "ComputeStep is called while the track has"
214 "the minimum interaction time";
215 exceptionDescription << " so it should not recompute a timeStep ";
216 G4Exception("G4DNABrownianTransportation::ComputeStep",
217 "G4DNABrownianTransportation001", FatalErrorInArgument,
218 exceptionDescription);
219 }
220 }
221
222 State(fGeometryLimitedStep) = false;
223
224 G4Molecule* molecule = GetMolecule(track);
225
226 if(timeStep > 0)
227 {
228 spaceStep = DBL_MAX;
229
230 // TODO : generalize this process to all kind of Brownian objects
231 G4double diffCoeff = molecule->GetDiffusionCoefficient(track.GetMaterial(),
232 track.GetMaterial()->GetTemperature());
233
234 static G4double sqrt_2 = sqrt(2.);
235 G4double sqrt_Dt = sqrt(diffCoeff*timeStep);
236 G4double sqrt_2Dt = sqrt_2*sqrt_Dt;
237 G4double x = G4RandGauss::shoot(0,sqrt_2Dt);
238 G4double y = G4RandGauss::shoot(0,sqrt_2Dt);
239 G4double z = G4RandGauss::shoot(0,sqrt_2Dt);
240
241 if(State(fTimeStepReachedLimit)== true)
242 {
243 //========================================================================
244 State(fGeometryLimitedStep) = true;// important
245 //========================================================================
246 spaceStep = State(fEndPointDistance);
247 // G4cout << "State(fTimeStepReachedLimit)== true" << G4endl;
248 }
249 else
250 {
251 spaceStep = sqrt(x*x + y*y + z*z);
252
253 if(spaceStep >= State(fEndPointDistance))
254 {
255 //G4cout << "spaceStep >= State(fEndPointDistance)" << G4endl;
256 //======================================================================
257 State(fGeometryLimitedStep) = true;// important
258 //======================================================================
259/*
260 if(fSpeedMeUp)
261 {
262 G4cout << "SpeedMeUp" << G4endl;
263 }
264 else
265*/
266 if(fUseSchedulerMinTimeSteps == false)// jump over barrier NOT used
267 {
268#ifdef G4VERBOSE
269 if (fVerboseLevel > 1)
270 {
272 << "G4ITBrownianTransportation::ComputeStep() : "
273 << "Step was limited to boundary"
274 << RESET_COLOR
275 << G4endl;
276 }
277#endif
278 //TODO
279 if(State(fRandomNumber)>=0) // CDF is used
280 {
281 /*
282 //=================================================================
283 State(fGeometryLimitedStep) = true;// important
284 //=================================================================
285 spaceStep = State(fEndPointDistance);
286 */
287
288 //==================================================================
289 // BE AWARE THAT THE TECHNIQUE USED BELOW IS A 1D APPROXIMATION
290 // Cumulative density function for the 3D case is not yet
291 // implemented
292 //==================================================================
293// G4cout << GREEN_ON_BLUE
294// << "G4ITBrownianTransportation::ComputeStep() : "
295// << "A random number was selected"
296// << RESET_COLOR
297// << G4endl;
298 G4double value = State(fRandomNumber)+(1-State(fRandomNumber))*G4UniformRand();
299 G4double invErfc = InvErfc(value);
300 spaceStep = invErfc*2*sqrt_Dt;
301
302 if(State(fTimeStepReachedLimit)== false)
303 {
304 //================================================================
305 State(fGeometryLimitedStep) = false;// important
306 //================================================================
307 }
308 //==================================================================
309 // DEBUG
310// if(spaceStep > State(fEndPointDistance))
311// {
312// G4cout << "value = " << value << G4endl;
313// G4cout << "invErfc = " << invErfc << G4endl;
314// G4cout << "spaceStep = " << G4BestUnit(spaceStep, "Length")
315// << G4endl;
316// G4cout << "end point distance= " << G4BestUnit(State(fEndPointDistance), "Length")
317// << G4endl;
318// }
319//
320// assert(spaceStep <= State(fEndPointDistance));
321 //==================================================================
322
323 }
324 else if(fUseMaximumTimeBeforeReachingBoundary == false) // CDF is used
325 {
326 G4double min_randomNumber = Erfc(State(fEndPointDistance)/2*sqrt_Dt);
327 G4double value = min_randomNumber+(1-min_randomNumber)*G4UniformRand();
328 G4double invErfc = InvErfc(value);
329 spaceStep = invErfc*2*sqrt_Dt;
330 if(spaceStep >= State(fEndPointDistance))
331 {
332 //================================================================
333 State(fGeometryLimitedStep) = true;// important
334 //================================================================
335 }
336 else if(State(fTimeStepReachedLimit)== false)
337 {
338 //================================================================
339 State(fGeometryLimitedStep) = false;// important
340 //================================================================
341 }
342 }
343 else // CDF is NOT used
344 {
345 //==================================================================
346 State(fGeometryLimitedStep) = true;// important
347 //==================================================================
348 spaceStep = State(fEndPointDistance);
349 //TODO
350
351 /*
352 //==================================================================
353 // 1D approximation to place the brownian between its starting point
354 // and the geometry boundary
355 //==================================================================
356 double min_randomNumber = Erfc(State(fEndPointDistance)/2*sqrt_Dt);
357 double value = State(fRandomNumber)+(1-State(fRandomNumber))*G4UniformRand();
358 double invErfc = InvErfc(value*G4UniformRand());
359 spaceStep = invErfc*2*sqrt_Dt;
360 State(fGeometryLimitedStep) = false;
361 */
362 }
363 }
364
365 State(fTransportEndPosition)= spaceStep*
366// step.GetPostStepPoint()->GetMomentumDirection()
368 + track.GetPosition();
369 }
370 else
371 {
372 //======================================================================
373 State(fGeometryLimitedStep) = false;// important
374 //======================================================================
375 State(fTransportEndPosition)= spaceStep*step.GetPostStepPoint()->
376 GetMomentumDirection() + track.GetPosition();
377 }
378 }
379 if(fpUserBrownianAction != nullptr)
380 {
381 // Let the user Brownian action class decide what to do:
382 G4ThreeVector nextPosition{track.GetPosition().getX() + x,
383 track.GetPosition().getY() + y,
384 track.GetPosition().getZ() + z};
385 fpUserBrownianAction->Transport(nextPosition);
386 State(fTransportEndPosition) = nextPosition;
387 }
388 }
389 else
390 {
391 spaceStep = 0.;
392 State(fTransportEndPosition) = track.GetPosition();
393 State(fGeometryLimitedStep) = false;
394 }
395
396 State(fCandidateEndGlobalTime) = step.GetPreStepPoint()->GetGlobalTime()
397 + timeStep;
398 State(fEndGlobalTimeComputed) = true;
399
400#ifdef G4VERBOSE
401 // DEBUG
402 if (fVerboseLevel > 1)
403 {
405 << "G4ITBrownianTransportation::ComputeStep() : "
406 << " trackID : " << track.GetTrackID() << " : Molecule name: "
407 << molecule->GetName() << G4endl
408 << "Initial position:" << G4BestUnit(track.GetPosition(), "Length")
409 << G4endl
410 << "Initial direction:" << track.GetMomentumDirection() << G4endl
411 << "Final position:" << G4BestUnit(State(fTransportEndPosition), "Length")
412 << G4endl
413 << "Initial magnitude:" << G4BestUnit(track.GetPosition().mag(), "Length")
414 << G4endl
415 << "Final magnitude:" << G4BestUnit(State(fTransportEndPosition).mag(), "Length")
416 << G4endl
417 << "Diffusion length : "
418 << G4BestUnit(spaceStep, "Length")
419 << " within time step : " << G4BestUnit(timeStep,"Time")
420 << G4endl
421 << "State(fTimeStepReachedLimit)= " << State(fTimeStepReachedLimit) << G4endl
422 << "State(fGeometryLimitedStep)=" << State(fGeometryLimitedStep) << G4endl
423 << "End point distance was: " << G4BestUnit(State(fEndPointDistance), "Length")
424 << G4endl
425 << RESET_COLOR
426 << G4endl<< G4endl;
427 }
428#endif
429
430//==============================================================================
431// DEBUG
432//assert(spaceStep < State(fEndPointDistance)
433// || (spaceStep >= State(fEndPointDistance) && State(fGeometryLimitedStep)));
434//assert(track.GetMomentumDirection() == State(fTransportEndMomentumDir));
435//==============================================================================
436}
@ 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
bool G4bool
Definition: G4Types.hh:86
double mag() const
G4double GetTemperature() const
Definition: G4Material.hh:177
const G4String & GetName() const
Definition: G4Molecule.cc:336
G4double GetGlobalTime() const
const G4VProcess * GetProcessDefinedStep() const
G4StepPoint * GetPostStepPoint() const
G4Material * GetMaterial() const
G4bool ProposesTimeStep() const

◆ Diffusion()

void G4DNABrownianTransportation::Diffusion ( const G4Track track)
protected

Definition at line 462 of file G4DNABrownianTransportation.cc.

463{
464
465#ifdef DEBUG_MEM
466 MemStat mem_first = MemoryUsage();
467#endif
468
469#ifdef G4VERBOSE
470 // DEBUG
471 if (fVerboseLevel > 1)
472 {
473 G4cout << GREEN_ON_BLUE << setw(18)
474 << "G4DNABrownianTransportation::Diffusion :" << setw(8)
475 << GetIT(track)->GetName() << "\t trackID:" << track.GetTrackID()
476 << "\t" << " Global Time = "
477 << G4BestUnit(track.GetGlobalTime(), "Time")
478 << RESET_COLOR
479 << G4endl
480 << G4endl;
481 }
482#endif
483
484/*
485 fParticleChange.ProposePosition(State(fTransportEndPosition));
486 //fParticleChange.ProposeEnergy(State(fTransportEndKineticEnergy));
487 fParticleChange.SetMomentumChanged(State(fMomentumChanged));
488
489 fParticleChange.ProposeGlobalTime(State(fCandidateEndGlobalTime));
490 fParticleChange.ProposeLocalTime(State(fCandidateEndGlobalTime));
491 fParticleChange.ProposeTrueStepLength(track.GetStepLength());
492*/
493 const G4Material* material = track.GetMaterial();
494
495 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
496
497 if (waterDensity == 0.0)
498 {
500 {
501 // Let the user Brownian action class decide what to do
504 return;
505 }
506 else
507 {
508#ifdef G4VERBOSE
509 if(fVerboseLevel)
510 {
511 G4cout << "A track is outside water material : trackID = "
512 << track.GetTrackID() << " (" << GetMolecule(track)->GetName() <<")"
513 << G4endl;
514 G4cout << "Local Time : " << G4BestUnit(track.GetGlobalTime(), "Time")
515 << G4endl;
516 G4cout << "Step Number :" << track.GetCurrentStepNumber() << G4endl;
517 }
518#endif
521 return;// &fParticleChange is the final returned object
522 }
523 }
524
525
526 #ifdef DEBUG_MEM
527 MemStat mem_intermediaire = MemoryUsage();
528 mem_diff = mem_intermediaire-mem_first;
529 G4cout << "\t\t\t >> || MEM || In G4DNABrownianTransportation::Diffusion "
530 "after dealing with waterDensity for "<< track.GetTrackID()
531 << ", diff is : " << mem_diff << G4endl;
532 #endif
533
535 State(fMomentumChanged) = true;
537 //
538 #ifdef DEBUG_MEM
539 mem_intermediaire = MemoryUsage();
540 mem_diff = mem_intermediaire-mem_first;
541 G4cout << "\t\t\t >> || MEM || In G4DNABrownianTransportation::"
542 "After proposing new direction to fParticleChange for "
543 << track.GetTrackID() << ", diff is : " << mem_diff << G4endl;
544 #endif
545
546 return;// &fParticleChange is the final returned object
547}
G4ThreeVector G4RandomDirection()
@ fStopAndKill
virtual void Transport(const G4Track &, G4ParticleChangeForTransport &)=0
virtual const G4String & GetName() const =0
size_t GetIndex() const
Definition: G4Material.hh:255
void SetMomentumChanged(G4bool b)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4int GetCurrentStepNumber() const
void ProposeTrackStatus(G4TrackStatus status)

Referenced by AlongStepDoIt().

◆ operator=()

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

◆ PostStepDoIt()

G4VParticleChange * G4DNABrownianTransportation::PostStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Reimplemented from G4ITTransportation.

Definition at line 438 of file G4DNABrownianTransportation.cc.

440{
442
443#ifdef G4VERBOSE
444 // DEBUG
445 if (fVerboseLevel > 1)
446 {
447 G4cout << GREEN_ON_BLUE << "G4ITBrownianTransportation::PostStepDoIt() :"
448 << " trackID : " << track.GetTrackID() << " Molecule name: "
449 << GetMolecule(track)->GetName() << G4endl;
450 G4cout << "Diffusion length : "
451 << G4BestUnit(step.GetStepLength(), "Length")
452 <<" within time step : " << G4BestUnit(step.GetDeltaTime(),"Time")
453 << "\t Current global time : "
454 << G4BestUnit(track.GetGlobalTime(),"Time")
455 << RESET_COLOR
456 << G4endl<< G4endl;
457 }
458#endif
459 return &fParticleChange;
460}
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &)
G4double GetDeltaTime() const
G4double GetStepLength() const

◆ SetBrownianAction()

void G4DNABrownianTransportation::SetBrownianAction ( G4BrownianAction brownianAction)
inline

Definition at line 245 of file G4DNABrownianTransportation.hh.

246{
247 fpBrownianAction = brownianAction;
248}

◆ SetUserBrownianAction()

void G4DNABrownianTransportation::SetUserBrownianAction ( G4VUserBrownianAction brownianAction)
inline

Definition at line 250 of file G4DNABrownianTransportation.hh.

251{
252 fpUserBrownianAction = brownianAction;
253}

◆ SpeedLevel()

void G4DNABrownianTransportation::SpeedLevel ( int  level)
inline

Definition at line 172 of file G4DNABrownianTransportation.hh.

173 {
174 if(level < 0) level =0;
175 else if(level > 2) level = 2;
176
177 switch(level)
178 {
179 case 0:
180 fSpeedMeUp = false;
182 return;
183
184 case 1:
185 fSpeedMeUp = true;
187 return;
188
189 case 2:
190 //======================================================================
191 // NB: BE AWARE THAT IF NO MIN TIME STEPS NO TIME STEPS HAVE BEEN
192 // PROVIDED TO G4Scheduler THIS LEVEL MIGHT BE SLOWER THAN LEVEL 1
193 //======================================================================
194 fSpeedMeUp = true;
196 return;
197 }
198 }

◆ StartTracking()

void G4DNABrownianTransportation::StartTracking ( G4Track aTrack)
virtual

Reimplemented from G4ITTransportation.

Definition at line 158 of file G4DNABrownianTransportation.cc.

159{
160 fpState.reset(new G4ITBrownianState());
161// G4cout << "G4DNABrownianTransportation::StartTracking : "
162// "Initialised track State" << G4endl;
165}
void SetInstantiateProcessState(G4bool flag)
virtual void StartTracking(G4Track *aTrack)

◆ UseCumulativeDensitFunction()

void G4DNABrownianTransportation::UseCumulativeDensitFunction ( bool  flag = true)
inline

Definition at line 156 of file G4DNABrownianTransportation.hh.

157 {
158 if(flag == true)
159 {
161 return;
162 }
164 }

◆ UseLimitingTimeSteps()

void G4DNABrownianTransportation::UseLimitingTimeSteps ( bool  flag = true)
inline

Definition at line 167 of file G4DNABrownianTransportation.hh.

168 {
170 }

◆ UseMaximumTimeBeforeReachingBoundary()

void G4DNABrownianTransportation::UseMaximumTimeBeforeReachingBoundary ( bool  flag = true)
inline

Definition at line 147 of file G4DNABrownianTransportation.hh.

148 {
150 }

Member Data Documentation

◆ fInternalMinTimeStep

G4double G4DNABrownianTransportation::fInternalMinTimeStep
protected

◆ fNistWater

G4Material* G4DNABrownianTransportation::fNistWater
protected

Definition at line 230 of file G4DNABrownianTransportation.hh.

Referenced by G4DNABrownianTransportation().

◆ fpBrownianAction

G4BrownianAction* G4DNABrownianTransportation::fpBrownianAction
protected

◆ fpUserBrownianAction

◆ fpWaterDensity

const std::vector<G4double>* G4DNABrownianTransportation::fpWaterDensity
protected

◆ fSpeedMeUp

G4bool G4DNABrownianTransportation::fSpeedMeUp
protected

◆ fUseMaximumTimeBeforeReachingBoundary

G4bool G4DNABrownianTransportation::fUseMaximumTimeBeforeReachingBoundary
protected

◆ fUseSchedulerMinTimeSteps

G4bool G4DNABrownianTransportation::fUseSchedulerMinTimeSteps
protected

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