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

{ The transportation method implemented is the one from Ermak-McCammon : J. Chem. Phys. 69, 1352 (1978)} More...

#include <G4DNABrownianTransportation.hh>

+ Inheritance diagram for G4DNABrownianTransportation:

Classes

struct  G4ITBrownianState
 

Public Member Functions

 G4DNABrownianTransportation (const G4String &aName="DNABrownianTransportation", G4int verbosityLevel=1)
 
virtual ~G4DNABrownianTransportation ()
 
 G4DNABrownianTransportation (const G4DNABrownianTransportation &other)
 
G4DNABrownianTransportationoperator= (const G4DNABrownianTransportation &other)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void StartTracking (G4Track *aTrack)
 
virtual void ComputeStep (const G4Track &, const G4Step &, const double, double &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &)
 
- Public Member Functions inherited from G4ITTransportation
 G4ITTransportation (const G4String &aName="ITTransportation", G4int verbosityLevel=1)
 
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)
 
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)
 
G4int operator== (const G4VITProcess &right) const
 
G4int operator!= (const G4VITProcess &right) const
 
size_t GetProcessID () const
 
G4ProcessState_LockGetProcessState ()
 
void SetProcessState (G4ProcessState_Lock *aProcInfo)
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

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 ()
 
virtual void ClearInteractionTimeLeft ()
 
virtual void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
virtual void ClearNumberOfInteractionLengthLeft ()
 
void SetInstantiateProcessState (G4bool flag)
 
G4bool InstantiateProcessState ()
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4ITBrownianState *constfpBrownianState
 
G4bool fUseMaximumTimeBeforeReachingBoundary
 
G4MaterialfNistWater
 
const std::vector< G4double > * fpWaterDensity
 
- Protected Attributes inherited from G4ITTransportation
G4ITTransportationState *constfTransportationState
 
G4ITNavigatorfLinearNavigator
 
G4PropagatorInFieldfFieldPropagator
 
G4ParticleChangeForTransport fParticleChange
 
G4double fThreshold_Warning_Energy
 
G4double fThreshold_Important_Energy
 
G4int fThresholdTrials
 
G4double fUnimportant_Energy
 
G4double fSumEnergyKilled
 
G4double fMaxEnergyKilled
 
G4bool fShortStepOptimisation
 
G4SafetyHelperfpSafetyHelper
 
G4int fVerboseLevel
 
- Protected Attributes inherited from G4VITProcess
G4ProcessStatefpState
 
G4bool fProposesTimeStep
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

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

{ The transportation method implemented is the one from Ermak-McCammon : J. Chem. Phys. 69, 1352 (1978)}

Definition at line 49 of file G4DNABrownianTransportation.hh.

Constructor & Destructor Documentation

◆ G4DNABrownianTransportation() [1/2]

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

Definition at line 83 of file G4DNABrownianTransportation.cc.

83 :
84 G4ITTransportation(aName, verbosity),
86{
87 //ctor
89 verboseLevel = 1;
93}
#define InitProcessState(destination, source)
Definition: G4VITProcess.hh:53
const std::vector< G4double > * fpWaterDensity
G4ITBrownianState *const & fpBrownianState
G4ITTransportationState *const & fTransportationState
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403

◆ ~G4DNABrownianTransportation()

G4DNABrownianTransportation::~G4DNABrownianTransportation ( )
virtual

Definition at line 95 of file G4DNABrownianTransportation.cc.

96{;}

◆ G4DNABrownianTransportation() [2/2]

G4DNABrownianTransportation::G4DNABrownianTransportation ( const G4DNABrownianTransportation other)

Definition at line 98 of file G4DNABrownianTransportation.cc.

98 :
99 G4ITTransportation(right),
101{
102 //copy ctor
104 verboseLevel = right.verboseLevel;
105 fUseMaximumTimeBeforeReachingBoundary = right.fUseMaximumTimeBeforeReachingBoundary;
106 fNistWater = right.fNistWater;
107 fpWaterDensity = right.fpWaterDensity ;
108}

Member Function Documentation

◆ AlongStepDoIt()

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

Reimplemented from G4ITTransportation.

Definition at line 362 of file G4DNABrownianTransportation.cc.

364{
366 Diffusion(track);
367 return &fParticleChange;
368}
void Diffusion(const G4Track &track)
G4ParticleChangeForTransport fParticleChange
virtual G4VParticleChange * AlongStepDoIt(const G4Track &track, const G4Step &stepData)

◆ AlongStepGetPhysicalInteractionLength()

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

Reimplemented from G4ITTransportation.

Definition at line 314 of file G4DNABrownianTransportation.cc.

320{
322 previousStepSize,
323 currentMinimumStep,
324 currentSafety,
325 selection);
326
327 G4double diffusionCoefficient = GetMolecule(track)->GetDiffusionCoefficient();
328 // G4double diffusionCoefficient = GetITBrownianObject(track)->GetDiffusionCoefficient(track.GetMaterial());
329
330 if(State(fGeometryLimitedStep))
331 {
332 // 99 % of the space step distribution is lower than
333 // d_99 = 8 * sqrt(D*t)
334 // where t is the corresponding time step
335 // so by inversion :
337 {
338 State(theInteractionTimeLeft) = (geometryStepLength*geometryStepLength)/(64 * diffusionCoefficient);
339 }
340 else // Will use a random time
341 {
342 State(theInteractionTimeLeft) = 1/(4*diffusionCoefficient) * pow(geometryStepLength/InvErfc(G4UniformRand()),2);
343 }
344
345 State(fCandidateEndGlobalTime) = track.GetGlobalTime() + State(theInteractionTimeLeft);
346 State(fPathLengthWasCorrected) = false;
347 }
348 else
349 {
350 geometryStepLength = 2*sqrt(diffusionCoefficient*State(theInteractionTimeLeft) ) *InvErf(G4UniformRand());
351 State(fPathLengthWasCorrected) = true;
352 State(endpointDistance) = geometryStepLength;
353 }
354
355 return geometryStepLength ;
356}
#define State(theXInfo)
G4Molecule * GetMolecule(const G4Track &track)
Definition: G4Molecule.cc:67
double G4double
Definition: G4Types.hh:64
#define G4UniformRand()
Definition: Randomize.hh:53
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4double GetDiffusionCoefficient() const
Definition: G4Molecule.cc:405
G4double GetGlobalTime() const

◆ BuildPhysicsTable()

void G4DNABrownianTransportation::BuildPhysicsTable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4ITTransportation.

Definition at line 129 of file G4DNABrownianTransportation.cc.

130{
131 if(verboseLevel > 0)
132 {
133 G4cout << G4endl << GetProcessName() << ": for "
134 << setw(24) << particle.GetParticleName()
135 << "\tSubType= " << GetProcessSubType() << G4endl;
136 }
137
138 // Initialize water density pointer
140}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
const std::vector< double > * GetDensityTableFor(const G4Material *) const
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
const G4String & GetParticleName() const
G4int GetProcessSubType() const
Definition: G4VProcess.hh:397
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

◆ ComputeStep()

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

Reimplemented from G4ITTransportation.

Definition at line 142 of file G4DNABrownianTransportation.cc.

146{
147 // G4cout << "G4ITBrownianTransportation::ComputeStep" << G4endl;
148
149 // If this method is called, this step
150 // cannot be geometry limited.
151 // In case the step is limited by the geometry,
152 // this method should not be called.
153 // The fTransportEndPosition calculated in
154 // the method AlongStepIL should be taken
155 // into account.
156 // In order to do so, the flag IsLeadingStep
157 // is on. Meaning : this track has the minimum
158 // interaction length over all others.
159 if(GetIT(track)->GetTrackingInfo()->IsLeadingStep())
160 {
161 const G4VITProcess* ITProc = ((const G4VITProcess*) step.GetPostStepPoint()->GetProcessDefinedStep());
162 bool makeException = true;
163
164 if(ITProc && ITProc->ProposesTimeStep()) makeException = false;
165
166 if(makeException)
167 {
168
169 G4ExceptionDescription exceptionDescription ;
170 exceptionDescription << "ComputeStep is called while the track has the minimum interaction time";
171 exceptionDescription << " so it should not recompute a timeStep ";
172 G4Exception("G4DNABrownianTransportation::ComputeStep","G4DNABrownianTransportation001",
173 FatalErrorInArgument,exceptionDescription);
174 }
175 }
176
177 State(fGeometryLimitedStep) = false;
178 // TODO : generalize this process to all kind of brownian objects
179 // G4ITBrownianObject* ITBrown = GetITBrownianObject(track) ;
180 // G4double diffCoeff = ITBrown->GetDiffusionCoefficient(track.GetMaterial());
181 G4Molecule* molecule = GetMolecule(track) ;
182 G4double diffCoeff = molecule->GetDiffusionCoefficient();
183
184 if(timeStep > 0)
185 {
186 spaceStep = DBL_MAX;
187
188 while(spaceStep > State(endpointDistance))
189 // Probably inefficient when the track is close close to boundaries
190 // it goes with fUserMaximumTimeBeforeReachingBoundary == false
191 // fUserMaximumTimeBeforeReachingBoundary == true, it should never loop
192 {
193 G4double x = G4RandGauss::shoot(0,sqrt(2*diffCoeff*timeStep));
194 G4double y = G4RandGauss::shoot(0,sqrt(2*diffCoeff*timeStep));
195 G4double z = G4RandGauss::shoot(0,sqrt(2*diffCoeff*timeStep));
196
197 spaceStep = sqrt(x*x + y*y + z*z);
198 }
199 // State(fTransportEndPosition).set(x + track.GetPosition().x(),
200 // y + track.GetPosition().y(),
201 // z + track.GetPosition().z());
202
203 State(fTransportEndPosition)= spaceStep*step.GetPostStepPoint()->GetMomentumDirection() + track.GetPosition();
204 }
205 else
206 {
207 spaceStep = 0. ;
208 State(fTransportEndPosition) = track.GetPosition() ;
209 }
210
211 State(fCandidateEndGlobalTime) = step.GetPreStepPoint()->GetGlobalTime() + timeStep ;
212 State(fEndGlobalTimeComputed) = true ;
213
214#ifdef G4VERBOSE
215 // DEBUG
216 if(fVerboseLevel>1)
217 {
219 << "G4ITBrownianTransportation::ComputeStep() : "
220 << " trackID : " << track.GetTrackID()
221 << " : Molecule name: " << molecule-> GetName()
222 << G4endl
223 << "Diffusion length : " << G4BestUnit(spaceStep, "Length")
224 << " within time step : " << G4BestUnit(timeStep,"Time")
225 << RESET
226 << G4endl<< G4endl;
227 }
228#endif
229}
#define GREEN_ON_BLUE
@ FatalErrorInArgument
G4IT * GetIT(const G4Track *track)
Definition: G4IT.cc:48
G4double GetGlobalTime() const
const G4VProcess * GetProcessDefinedStep() const
const G4ThreeVector & GetMomentumDirection() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4bool ProposesTimeStep() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
#define DBL_MAX
Definition: templates.hh:83

◆ Diffusion()

void G4DNABrownianTransportation::Diffusion ( const G4Track track)
protected

Definition at line 253 of file G4DNABrownianTransportation.cc.

255{
256
257#ifdef G4VERBOSE
258 // DEBUG
259 if (fVerboseLevel>1)
260 {
262 << setw(18)<< "G4DNABrownianTransportation::Diffusion :"
263 << setw(8) << GetIT(track)->GetName()
264 << "\t trackID:" << track.GetTrackID() <<"\t"
265 << " Global Time = " << G4BestUnit(track.GetGlobalTime(),"Time")
266 << RESET
267 << G4endl<< G4endl;
268 }
269#endif
270
271 G4Material* material = track.GetMaterial();
272// if (material != fNistWater && material->GetBaseMaterial() != fNistWater)
273
274 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
275
276 if(waterDensity == 0.0)
277// if (material == nistwater || material->GetBaseMaterial() == nistwater)
278 {
279 G4cout << "A track is outside water material : trackID"<< track.GetTrackID() << " (" << GetMolecule(track)->GetName() <<")" << G4endl;
280 G4cout << "Local Time : "<< (track.GetLocalTime()) /s<<G4endl;
281 G4cout<< "Step Number :" << track.GetCurrentStepNumber() <<G4endl;
282
285
286 // Got pb with :
287 // fParticleChange.ProposeTrackStatus(fStopAndKill);
288 // It makes the tracks going straight without killing them
289
290 return ; // &fParticleChange is the final returned object
291 }
292
293 G4double costheta = (2*G4UniformRand()-1);
294 G4double theta = acos (costheta);
295 G4double phi = 2*pi*G4UniformRand();
296
297 G4double xMomentum = cos(phi)* sin(theta);
298 G4double yMomentum = sin(theta)*sin(phi);
299 G4double zMomentum = costheta;
300
301 fParticleChange.ProposeMomentumDirection(xMomentum, yMomentum, zMomentum);
302 State(fMomentumChanged) = true;
304
305 // G4cout << "BROWN : Propose new direction :" << G4ThreeVector(xMomentum, yMomentum, zMomentum) << G4endl;
306
307 // Alternative
308 //fParticleChange.ProposeMomentumDirection(G4RandomDirection());
309
310 return; // &fParticleChange is the final returned object
311}
@ fStopButAlive
virtual const G4String & GetName() const =0
size_t GetIndex() const
Definition: G4Material.hh:261
const G4String & GetName() const
Definition: G4Molecule.cc:259
void SetMomentumChanged(G4bool b)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4int GetCurrentStepNumber() const
G4double GetLocalTime() const
G4Material * GetMaterial() const
void ProposeTrackStatus(G4TrackStatus status)
const G4double pi

Referenced by AlongStepDoIt().

◆ operator=()

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

Definition at line 110 of file G4DNABrownianTransportation.cc.

111{
112 if (this == &rhs) return *this; // handle self assignment
113 //assignment operator
114 return *this;
115}

◆ PostStepDoIt()

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

Reimplemented from G4ITTransportation.

Definition at line 231 of file G4DNABrownianTransportation.cc.

232{
234
235#ifdef G4VERBOSE
236 // DEBUG
237 if(fVerboseLevel>1)
238 {
240 << "G4ITBrownianTransportation::PostStepDoIt() :"
241 << " trackID : " << track.GetTrackID()
242 << " Molecule name: " << GetMolecule(track)-> GetName() << G4endl;
243 G4cout<< "Diffusion length : "<< G4BestUnit(step.GetStepLength(),"Length") <<" within time step : "
244 << G4BestUnit(step.GetDeltaTime(),"Time") << "\t"
245 << " Current global time : " << G4BestUnit(track.GetGlobalTime(),"Time")
246 << RESET
247 << G4endl<< G4endl;
248 }
249#endif
250 return &fParticleChange ;
251}
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &)
G4double GetDeltaTime() const
G4double GetStepLength() const

◆ StartTracking()

void G4DNABrownianTransportation::StartTracking ( G4Track aTrack)
virtual

Reimplemented from G4ITTransportation.

Definition at line 122 of file G4DNABrownianTransportation.cc.

123{
124 fpState = new G4ITBrownianState();
127}
void SetInstantiateProcessState(G4bool flag)
virtual void StartTracking(G4Track *aTrack)
G4ProcessState * fpState

Member Data Documentation

◆ fNistWater

G4Material* G4DNABrownianTransportation::fNistWater
protected

Definition at line 92 of file G4DNABrownianTransportation.hh.

Referenced by G4DNABrownianTransportation().

◆ fpBrownianState

G4ITBrownianState* const& G4DNABrownianTransportation::fpBrownianState
protected

Definition at line 89 of file G4DNABrownianTransportation.hh.

◆ fpWaterDensity

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

◆ fUseMaximumTimeBeforeReachingBoundary

G4bool G4DNABrownianTransportation::fUseMaximumTimeBeforeReachingBoundary
protected

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