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

#include <G4BOptnChangeCrossSection.hh>

+ Inheritance diagram for G4BOptnChangeCrossSection:

Public Member Functions

 G4BOptnChangeCrossSection (const G4String &name)
 
virtual ~G4BOptnChangeCrossSection ()
 
virtual const G4VBiasingInteractionLawProvideOccurenceBiasingInteractionLaw (const G4BiasingProcessInterface *, G4ForceCondition &proposeForceCondition)
 
virtual G4VParticleChangeApplyFinalStateBiasing (const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
 
virtual G4double DistanceToApplyOperation (const G4Track *, G4double, G4ForceCondition *)
 
virtual G4VParticleChangeGenerateBiasingFinalState (const G4Track *, const G4Step *)
 
G4InteractionLawPhysicalGetBiasedExponentialLaw ()
 
void SetBiasedCrossSection (G4double xst, G4bool updateInteractionLength=false)
 
G4double GetBiasedCrossSection () const
 
void Sample ()
 
void UpdateForStep (G4double stepLength)
 
G4bool GetInteractionOccured () const
 
void SetInteractionOccured ()
 
- Public Member Functions inherited from G4VBiasingOperation
 G4VBiasingOperation (const G4String &name)
 
virtual ~G4VBiasingOperation ()=default
 
virtual G4double ProposeAlongStepLimit (const G4BiasingProcessInterface *)
 
virtual G4GPILSelection ProposeGPILSelection (const G4GPILSelection wrappedProcessSelection)
 
virtual void AlongMoveBy (const G4BiasingProcessInterface *, const G4Step *, G4double)
 
const G4StringGetName () const
 
std::size_t GetUniqueID () const
 

Detailed Description

Definition at line 41 of file G4BOptnChangeCrossSection.hh.

Constructor & Destructor Documentation

◆ G4BOptnChangeCrossSection()

G4BOptnChangeCrossSection::G4BOptnChangeCrossSection ( const G4String & name)

Definition at line 32 of file G4BOptnChangeCrossSection.cc.

33 : G4VBiasingOperation( name )
34{
35 fBiasedExponentialLaw = new G4InteractionLawPhysical("LawForOperation"+name);
36}
G4VBiasingOperation(const G4String &name)

◆ ~G4BOptnChangeCrossSection()

G4BOptnChangeCrossSection::~G4BOptnChangeCrossSection ( )
virtual

Definition at line 38 of file G4BOptnChangeCrossSection.cc.

39{
40 delete fBiasedExponentialLaw;
41}

Member Function Documentation

◆ ApplyFinalStateBiasing()

virtual G4VParticleChange * G4BOptnChangeCrossSection::ApplyFinalStateBiasing ( const G4BiasingProcessInterface * ,
const G4Track * ,
const G4Step * ,
G4bool &  )
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 59 of file G4BOptnChangeCrossSection.hh.

60 { return nullptr; }

◆ DistanceToApplyOperation()

virtual G4double G4BOptnChangeCrossSection::DistanceToApplyOperation ( const G4Track * ,
G4double ,
G4ForceCondition *  )
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 61 of file G4BOptnChangeCrossSection.hh.

62 { return DBL_MAX; }
#define DBL_MAX
Definition templates.hh:62

◆ GenerateBiasingFinalState()

virtual G4VParticleChange * G4BOptnChangeCrossSection::GenerateBiasingFinalState ( const G4Track * ,
const G4Step *  )
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 63 of file G4BOptnChangeCrossSection.hh.

64 { return nullptr; }

◆ GetBiasedCrossSection()

G4double G4BOptnChangeCrossSection::GetBiasedCrossSection ( ) const

Definition at line 54 of file G4BOptnChangeCrossSection.cc.

55{
56 return fBiasedExponentialLaw->GetPhysicalCrossSection();
57}

◆ GetBiasedExponentialLaw()

G4InteractionLawPhysical * G4BOptnChangeCrossSection::GetBiasedExponentialLaw ( )
inline

Definition at line 69 of file G4BOptnChangeCrossSection.hh.

69{ return fBiasedExponentialLaw; }

◆ GetInteractionOccured()

G4bool G4BOptnChangeCrossSection::GetInteractionOccured ( ) const
inline

Definition at line 79 of file G4BOptnChangeCrossSection.hh.

79{ return fInteractionOccured; }

◆ ProvideOccurenceBiasingInteractionLaw()

const G4VBiasingInteractionLaw * G4BOptnChangeCrossSection::ProvideOccurenceBiasingInteractionLaw ( const G4BiasingProcessInterface * ,
G4ForceCondition & proposeForceCondition )
virtual

Implements G4VBiasingOperation.

Definition at line 43 of file G4BOptnChangeCrossSection.cc.

44{
45 return fBiasedExponentialLaw;
46}

◆ Sample()

void G4BOptnChangeCrossSection::Sample ( )

Definition at line 59 of file G4BOptnChangeCrossSection.cc.

60{
61 fInteractionOccured = false;
62 fBiasedExponentialLaw->Sample();
63}

◆ SetBiasedCrossSection()

void G4BOptnChangeCrossSection::SetBiasedCrossSection ( G4double xst,
G4bool updateInteractionLength = false )

Definition at line 48 of file G4BOptnChangeCrossSection.cc.

49{
50 fBiasedExponentialLaw->SetPhysicalCrossSection( xst );
51 if ( updateInteractionLength ) UpdateForStep( 0.0 );
52}
void UpdateForStep(G4double stepLength)

◆ SetInteractionOccured()

void G4BOptnChangeCrossSection::SetInteractionOccured ( )
inline

Definition at line 80 of file G4BOptnChangeCrossSection.hh.

80{ fInteractionOccured = true; }

◆ UpdateForStep()

void G4BOptnChangeCrossSection::UpdateForStep ( G4double stepLength)

Definition at line 65 of file G4BOptnChangeCrossSection.cc.

66{
67 fBiasedExponentialLaw->UpdateForStep( truePathLength );
68}

Referenced by SetBiasedCrossSection().


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