Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4BOptnChangeCrossSection Class Reference

#include <G4BOptnChangeCrossSection.hh>

Inheritance diagram for G4BOptnChangeCrossSection:
G4VBiasingOperation

Public Member Functions

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

Detailed Description

Definition at line 47 of file G4BOptnChangeCrossSection.hh.

Constructor & Destructor Documentation

G4BOptnChangeCrossSection::G4BOptnChangeCrossSection ( G4String  name)

Definition at line 31 of file G4BOptnChangeCrossSection.cc.

32  : G4VBiasingOperation(name)
33 {
34  fBiasedExponentialLaw = new G4InteractionLawPhysical("LawForOperation"+name);
35 }
G4VBiasingOperation(G4String name)
G4BOptnChangeCrossSection::~G4BOptnChangeCrossSection ( )
virtual

Definition at line 37 of file G4BOptnChangeCrossSection.cc.

38 {}

Member Function Documentation

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

Implements G4VBiasingOperation.

Definition at line 61 of file G4BOptnChangeCrossSection.hh.

63  {return 0;}
virtual G4double G4BOptnChangeCrossSection::DistanceToApplyOperation ( const G4Track ,
G4double  ,
G4ForceCondition  
)
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 64 of file G4BOptnChangeCrossSection.hh.

References DBL_MAX.

66  {return DBL_MAX;}
#define DBL_MAX
Definition: templates.hh:83
virtual G4VParticleChange* G4BOptnChangeCrossSection::GenerateBiasingFinalState ( const G4Track ,
const G4Step  
)
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 67 of file G4BOptnChangeCrossSection.hh.

68  {return 0;}
G4double G4BOptnChangeCrossSection::GetBiasedCrossSection ( ) const

Definition at line 50 of file G4BOptnChangeCrossSection.cc.

References G4InteractionLawPhysical::GetPhysicalCrossSection().

51 {
52  return fBiasedExponentialLaw->GetPhysicalCrossSection();
53 }
G4double GetPhysicalCrossSection() const
G4InteractionLawPhysical* G4BOptnChangeCrossSection::GetBiasedExponentialLaw ( )
inline

Definition at line 74 of file G4BOptnChangeCrossSection.hh.

74 {return fBiasedExponentialLaw;}
G4bool G4BOptnChangeCrossSection::GetInteractionOccured ( ) const
inline

Definition at line 84 of file G4BOptnChangeCrossSection.hh.

84 { return fInteractionOccured; }
const G4VBiasingInteractionLaw * G4BOptnChangeCrossSection::ProvideOccurenceBiasingInteractionLaw ( const G4BiasingProcessInterface )
virtual

Implements G4VBiasingOperation.

Definition at line 40 of file G4BOptnChangeCrossSection.cc.

41 {
42  return fBiasedExponentialLaw;
43 }
void G4BOptnChangeCrossSection::Sample ( )

Definition at line 55 of file G4BOptnChangeCrossSection.cc.

References G4VBiasingInteractionLaw::Sample().

56 {
57  fInteractionOccured = false;
58  fBiasedExponentialLaw->Sample();
59 }
void G4BOptnChangeCrossSection::SetBiasedCrossSection ( G4double  xst)

Definition at line 45 of file G4BOptnChangeCrossSection.cc.

References G4InteractionLawPhysical::SetPhysicalCrossSection().

46 {
47  fBiasedExponentialLaw->SetPhysicalCrossSection( xst );
48 }
void SetPhysicalCrossSection(G4double crossSection)
void G4BOptnChangeCrossSection::SetInteractionOccured ( )
inline

Definition at line 85 of file G4BOptnChangeCrossSection.hh.

85 { fInteractionOccured = true; }
void G4BOptnChangeCrossSection::UpdateForStep ( G4double  stepLength)

Definition at line 61 of file G4BOptnChangeCrossSection.cc.

References G4VBiasingInteractionLaw::UpdateForStep().

62 {
63  fBiasedExponentialLaw->UpdateForStep( truePathLength );
64 }
G4double UpdateForStep(G4double truePathLength)

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