Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
GB04BOptnBremSplitting Class Reference

#include <Doxymodules_biasing.h>

Inheritance diagram for GB04BOptnBremSplitting:
G4VBiasingOperation

Public Member Functions

 GB04BOptnBremSplitting (G4String name)
 
virtual ~GB04BOptnBremSplitting ()
 
virtual const G4VBiasingInteractionLawProvideOccurenceBiasingInteractionLaw (const G4BiasingProcessInterface *, G4ForceCondition &)
 
virtual G4VParticleChangeApplyFinalStateBiasing (const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
 
virtual G4double DistanceToApplyOperation (const G4Track *, G4double, G4ForceCondition *)
 
virtual G4VParticleChangeGenerateBiasingFinalState (const G4Track *, const G4Step *)
 
void SetSplittingFactor (G4int splittingFactor)
 
G4int GetSplittingFactor () const
 

Private Attributes

G4int fSplittingFactor
 
G4ParticleChange fParticleChange
 

Detailed Description

Definition at line 109 of file Doxymodules_biasing.h.

Constructor & Destructor Documentation

◆ GB04BOptnBremSplitting()

GB04BOptnBremSplitting::GB04BOptnBremSplitting ( G4String  name)

◆ ~GB04BOptnBremSplitting()

GB04BOptnBremSplitting::~GB04BOptnBremSplitting ( )
virtual

Definition at line 46 of file GB04BOptnBremSplitting.cc.

47{
48}

Member Function Documentation

◆ ProvideOccurenceBiasingInteractionLaw()

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

Definition at line 60 of file GB04BOptnBremSplitting.hh.

62 { return 0; }

◆ ApplyFinalStateBiasing()

G4VParticleChange * GB04BOptnBremSplitting::ApplyFinalStateBiasing ( const G4BiasingProcessInterface callingProcess,
const G4Track track,
const G4Step step,
G4bool &   
)
virtual

Definition at line 51 of file GB04BOptnBremSplitting.cc.

56{
57
58 // -- Collect brem. process (wrapped process) final state:
59 G4VParticleChange* processFinalState =
60 callingProcess->GetWrappedProcess()->PostStepDoIt(*track, *step);
61
62 // -- if no splitting requested, let the brem. process to return directly its
63 // -- generated final state:
64 if ( fSplittingFactor == 1 ) return processFinalState;
65
66 // -- a special case here: the brem. process corrects for cross-section change
67 // -- over the step due to energy loss by sometimes "abandoning" the interaction,
68 // -- returning an unchanged incoming electron/positron.
69 // -- We respect this correction, and if no secondary is produced, its means this
70 // -- case is happening:
71 if ( processFinalState->GetNumberOfSecondaries() == 0 ) return processFinalState;
72
73 // -- Now start the biasing:
74 // -- - the electron state will be taken as the first one produced by the brem.
75 // -- process, hence the one stored in above processFinalState particle change.
76 // -- This state will be stored in our fParticleChange object.
77 // -- - the photon accompagnying the electron will be stored also this way.
78 // -- - we will then do fSplittingFactor - 1 call to the brem. process to collect
79 // -- fSplittingFactor - 1 additionnal gammas. All these will be stored in our
80 // -- fParticleChange object.
81
82 // -- We called the brem. process above. Its concrete particle change is indeed
83 // -- a "G4ParticleChangeForLoss" object. We cast this particle change to access
84 // -- methods of the concrete G4ParticleChangeForLoss type:
85 G4ParticleChangeForLoss* actualParticleChange =
86 ( G4ParticleChangeForLoss* ) processFinalState ;
87
88 fParticleChange.Initialize(*track);
89
90 // -- Store electron final state:
92 ProposeTrackStatus ( actualParticleChange->GetTrackStatus() );
94 ProposeEnergy ( actualParticleChange->GetProposedKineticEnergy() );
96 ProposeMomentumDirection( actualParticleChange->GetProposedMomentumDirection() );
97
98 // -- Now deal with the gamma's:
99 // -- their common weight:
100 G4double gammaWeight = track->GetWeight() / fSplittingFactor;
101
102 // -- inform we will have fSplittingFactor gamma's:
103 fParticleChange.SetNumberOfSecondaries( fSplittingFactor );
104
105 // -- inform we take care of secondaries weight (otherwise these
106 // -- secondaries are by default given the primary weight).
107 fParticleChange.SetSecondaryWeightByProcess(true);
108
109 // -- Store first gamma:
110 G4Track* gammaTrack = actualParticleChange->GetSecondary(0);
111 gammaTrack->SetWeight( gammaWeight );
112 fParticleChange.AddSecondary( gammaTrack );
113 // -- and clean-up the brem. process particle change:
114 actualParticleChange->Clear();
115
116 // -- now start the fSplittingFactor-1 calls to the brem. process to store each
117 // -- related gamma:
118 G4int nCalls = 1;
119 while ( nCalls < fSplittingFactor )
120 {
121 // ( note: we don't need to cast to actual type here, as methods for accessing
122 // secondary particles are from base class G4VParticleChange )
123 processFinalState = callingProcess->GetWrappedProcess()->PostStepDoIt(*track, *step);
124 if ( processFinalState->GetNumberOfSecondaries() == 1 )
125 {
126 gammaTrack = processFinalState->GetSecondary(0);
127 gammaTrack->SetWeight( gammaWeight );
128 fParticleChange.AddSecondary( gammaTrack );
129 nCalls++;
130 }
131 // -- very rare special case: we ignore for now.
132 else if ( processFinalState->GetNumberOfSecondaries() > 1 )
133 {
134 for ( G4int i = 0 ; i < processFinalState->GetNumberOfSecondaries() ; i++)
135 delete processFinalState->GetSecondary(i);
136 }
137 processFinalState->Clear();
138 }
139
140 // -- we are done:
141 return &fParticleChange;
142}

◆ DistanceToApplyOperation()

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

Definition at line 71 of file GB04BOptnBremSplitting.hh.

74 {return DBL_MAX;}

◆ GenerateBiasingFinalState()

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

Definition at line 75 of file GB04BOptnBremSplitting.hh.

77 {return 0;}

◆ SetSplittingFactor()

void GB04BOptnBremSplitting::SetSplittingFactor ( G4int  splittingFactor)
inline

Definition at line 85 of file GB04BOptnBremSplitting.hh.

86 { fSplittingFactor = splittingFactor; }

◆ GetSplittingFactor()

G4int GB04BOptnBremSplitting::GetSplittingFactor ( ) const
inline

Definition at line 87 of file GB04BOptnBremSplitting.hh.

87{ return fSplittingFactor; }

Member Data Documentation

◆ fSplittingFactor

G4int GB04BOptnBremSplitting::fSplittingFactor
private

Definition at line 90 of file GB04BOptnBremSplitting.hh.

◆ fParticleChange

G4ParticleChange GB04BOptnBremSplitting::fParticleChange
private

Definition at line 91 of file GB04BOptnBremSplitting.hh.


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

Applications | User Support | Publications | Collaboration