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

#include <G4BOptnCloning.hh>

+ Inheritance diagram for G4BOptnCloning:

Public Member Functions

 G4BOptnCloning (const G4String &name)
 
virtual ~G4BOptnCloning ()
 
virtual const G4VBiasingInteractionLawProvideOccurenceBiasingInteractionLaw (const G4BiasingProcessInterface *, G4ForceCondition &)
 
virtual G4VParticleChangeApplyFinalStateBiasing (const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
 
virtual G4double DistanceToApplyOperation (const G4Track *, G4double, G4ForceCondition *condition)
 
virtual G4VParticleChangeGenerateBiasingFinalState (const G4Track *, const G4Step *)
 
void SetCloneWeights (G4double clone1Weight, G4double clone2Weight)
 
G4TrackGetCloneTrack () const
 
- 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 G4BOptnCloning.hh.

Constructor & Destructor Documentation

◆ G4BOptnCloning()

G4BOptnCloning::G4BOptnCloning ( const G4String & name)

Definition at line 31 of file G4BOptnCloning.cc.

32 : G4VBiasingOperation( name ),
33 fParticleChange()
34{}
G4VBiasingOperation(const G4String &name)

◆ ~G4BOptnCloning()

G4BOptnCloning::~G4BOptnCloning ( )
virtual

Definition at line 36 of file G4BOptnCloning.cc.

37{}

Member Function Documentation

◆ ApplyFinalStateBiasing()

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

Implements G4VBiasingOperation.

Definition at line 57 of file G4BOptnCloning.hh.

58 { return nullptr; }

◆ DistanceToApplyOperation()

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

Implements G4VBiasingOperation.

Definition at line 61 of file G4BOptnCloning.hh.

63 {
64 *condition = NotForced; return 0.; // -- acts immediately
65 }
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced

◆ GenerateBiasingFinalState()

G4VParticleChange * G4BOptnCloning::GenerateBiasingFinalState ( const G4Track * track,
const G4Step *  )
virtual

Implements G4VBiasingOperation.

Definition at line 40 of file G4BOptnCloning.cc.

42{
43 fParticleChange.Initialize(*track);
44 fParticleChange.ProposeParentWeight( fClone1W );
45 fParticleChange.SetSecondaryWeightByProcess(true);
46 fParticleChange.SetNumberOfSecondaries(1);
47 fCloneTrack = new G4Track( *track );
48 fCloneTrack->SetWeight( fClone2W );
49 fParticleChange.AddSecondary( fCloneTrack );
50 return &fParticleChange;
51}

◆ GetCloneTrack()

G4Track * G4BOptnCloning::GetCloneTrack ( ) const
inline

Definition at line 77 of file G4BOptnCloning.hh.

77{ return fCloneTrack; }

◆ ProvideOccurenceBiasingInteractionLaw()

virtual const G4VBiasingInteractionLaw * G4BOptnCloning::ProvideOccurenceBiasingInteractionLaw ( const G4BiasingProcessInterface * ,
G4ForceCondition &  )
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 54 of file G4BOptnCloning.hh.

55 { return nullptr; }

◆ SetCloneWeights()

void G4BOptnCloning::SetCloneWeights ( G4double clone1Weight,
G4double clone2Weight )
inline

Definition at line 71 of file G4BOptnCloning.hh.

72 {
73 fClone1W = clone1Weight;
74 fClone2W = clone2Weight;
75 }

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