Geant4-11
Public Member Functions | Private Attributes
G4BOptnLeadingParticle Class Reference

#include <G4BOptnLeadingParticle.hh>

Inheritance diagram for G4BOptnLeadingParticle:
G4VBiasingOperation

Public Member Functions

virtual void AlongMoveBy (const G4BiasingProcessInterface *, const G4Step *, G4double)
 
virtual G4VParticleChangeApplyFinalStateBiasing (const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
 
virtual G4double DistanceToApplyOperation (const G4Track *, G4double, G4ForceCondition *)
 
 G4BOptnLeadingParticle (G4String name)
 
virtual G4VParticleChangeGenerateBiasingFinalState (const G4Track *, const G4Step *)
 
G4double GetFurtherKillingProbability () const
 
const G4StringGetName () const
 
std::size_t GetUniqueID () const
 
virtual G4double ProposeAlongStepLimit (const G4BiasingProcessInterface *)
 
virtual G4GPILSelection ProposeGPILSelection (const G4GPILSelection wrappedProcessSelection)
 
virtual const G4VBiasingInteractionLawProvideOccurenceBiasingInteractionLaw (const G4BiasingProcessInterface *, G4ForceCondition &)
 
void SetFurtherKillingProbability (G4double p)
 
virtual ~G4BOptnLeadingParticle ()
 

Private Attributes

const G4String fName
 
G4ParticleChange fParticleChange
 
G4double fRussianRouletteKillingProbability
 
std::size_t fUniqueID
 

Detailed Description

Definition at line 54 of file G4BOptnLeadingParticle.hh.

Constructor & Destructor Documentation

◆ G4BOptnLeadingParticle()

G4BOptnLeadingParticle::G4BOptnLeadingParticle ( G4String  name)

Definition at line 33 of file G4BOptnLeadingParticle.cc.

36{
37}
G4VBiasingOperation(G4String name)
const char * name(G4int ptype)

◆ ~G4BOptnLeadingParticle()

G4BOptnLeadingParticle::~G4BOptnLeadingParticle ( )
virtual

Definition at line 39 of file G4BOptnLeadingParticle.cc.

40{
41}

Member Function Documentation

◆ AlongMoveBy()

virtual void G4VBiasingOperation::AlongMoveBy ( const G4BiasingProcessInterface ,
const G4Step ,
G4double   
)
inlinevirtualinherited

Reimplemented in G4BOptnForceFreeFlight.

Definition at line 139 of file G4VBiasingOperation.hh.

141 {}

Referenced by G4BiasingProcessInterface::AlongStepDoIt().

◆ ApplyFinalStateBiasing()

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

G4cout << G4endl;

Implements G4VBiasingOperation.

Definition at line 43 of file G4BOptnLeadingParticle.cc.

47{
48 // -- collect wrapped process particle change:
49 auto wrappedProcessParticleChange = callingProcess->GetWrappedProcess()->PostStepDoIt(*track,*step);
50
51 // -- does nothing in case the primary stays alone or in weird situation where all are killed...
52 if ( wrappedProcessParticleChange->GetNumberOfSecondaries() == 0 ) return wrappedProcessParticleChange;
53 if ( wrappedProcessParticleChange->GetTrackStatus() == fKillTrackAndSecondaries ) return wrappedProcessParticleChange;
54 // -- ... else, collect the secondaries in a same vector (the primary is pushed in this vector, if surviving, later see [**]):
55 std::vector < G4Track* > secondariesAndPrimary;
56 for ( auto i = 0 ; i < wrappedProcessParticleChange->GetNumberOfSecondaries() ; i++ ) secondariesAndPrimary.push_back( wrappedProcessParticleChange->GetSecondary( i ) );
57
58
59 // -- If case the primary survives, need to collect its new state. In the general case of the base class G4VParticleChange
60 // -- this is tricky, as this class does not hold the primary changes (and we have to build a fake step and fake track
61 // -- caring about the touchables, etc.) So we limit here to the G4ParticleChange case, checking the reality of this
62 // -- class with a dynamic cast. If we don't have such an actual G4DynamicParticle object, we give up the biasing and return
63 // -- the untrimmed process final state.
64 // -- Note that this case does not happen often, as the technique is intended for inelastic processes. For case where several
65 // -- particles can be produced without killing the primary, we have for example the electron-/positron-nuclear
66 G4ParticleChange* castParticleChange ( nullptr );
67 G4Track* finalStatePrimary ( nullptr );
68 if ( ( wrappedProcessParticleChange->GetTrackStatus() != fStopAndKill ) )
69 {
70 // fFakePrimaryTrack->CopyTrackInfo( *track );
71 // fFakeStep ->InitializeStep( fFakePrimaryTrack );
72 // wrappedProcessParticleChange->UpdateStepForPostStep( fFakeStep );
73 // fFakeStep->UpdateTrack();
74 castParticleChange = dynamic_cast< G4ParticleChange* >( wrappedProcessParticleChange );
75 if ( castParticleChange == nullptr )
76 {
77 G4cout << " **** G4BOptnLeadingParticle::ApplyFinalStateBiasing(...) : can not bias for " << callingProcess->GetProcessName() << ", this is just a warning." << G4endl;
78 return wrappedProcessParticleChange;
79 }
80 finalStatePrimary = new G4Track( *track );
81 finalStatePrimary->SetKineticEnergy ( castParticleChange->GetEnergy() );
82 finalStatePrimary->SetWeight ( castParticleChange->GetWeight() );
83 finalStatePrimary->SetMomentumDirection( *(castParticleChange->GetMomentumDirection()) );
84 // -- [**] push the primary as the last track in the vector of tracks:
85 secondariesAndPrimary.push_back( finalStatePrimary );
86 }
87
88 // -- Ensure the secondaries all have the primary weight:
89 // ---- collect primary track weight, from updated by process if alive, or from original copy if died:
90 G4double primaryWeight;
91 if ( finalStatePrimary ) primaryWeight = finalStatePrimary->GetWeight();
92 else primaryWeight = track ->GetWeight();
93 // ---- now set this same weight to all secondaries:
94 for ( auto i = 0 ; i < wrappedProcessParticleChange->GetNumberOfSecondaries() ; i++ ) secondariesAndPrimary[ i ]->SetWeight( primaryWeight );
95
96
97 // -- finds the leading particle, initialize a map of surviving tracks, tag as surviving the leading track:
98 size_t leadingIDX = 0;
99 G4double leadingEnergy = -1;
100 std::map< G4Track*, G4bool > survivingMap;
101 for ( size_t idx = 0; idx < secondariesAndPrimary.size(); idx++ )
102 {
103 survivingMap[ secondariesAndPrimary[idx] ] = false;
104 if ( secondariesAndPrimary[idx]->GetKineticEnergy() > leadingEnergy )
105 {
106 leadingEnergy = secondariesAndPrimary[idx]->GetKineticEnergy();
107 leadingIDX = idx;
108 }
109 }
110 survivingMap[ secondariesAndPrimary[leadingIDX] ] = true; // -- tag as surviving the leading particle
111
112
113 // -- now make track vectors of given types ( choose type = abs(PDG) ), excluding the leading particle:
114 std::map < G4int, std::vector< G4Track* > > typesAndTracks;
115 for ( size_t idx = 0; idx < secondariesAndPrimary.size(); idx++ )
116 {
117 if ( idx == leadingIDX ) continue; // -- excludes the leading particle
118 auto currentTrack = secondariesAndPrimary[idx];
119 auto GROUP = std::abs( currentTrack->GetDefinition()->GetPDGEncoding() ); // -- merge particles and anti-particles in the same category -- §§ this might be proposed as an option in future
120 if ( currentTrack->GetDefinition()->GetBaryonNumber() >= 2 ) GROUP = -1000; // -- merge all baryons above proton/neutron in one same group -- §§ might be proposed as an option too
121
122 if ( typesAndTracks.find( GROUP ) == typesAndTracks.end() )
123 {
124 std::vector< G4Track* > v;
125 v.push_back( currentTrack );
126 typesAndTracks[ GROUP ] = v;
127 }
128 else
129 {
130 typesAndTracks[ GROUP ].push_back( currentTrack );
131 }
132 }
133 // -- and on these vectors, randomly select the surviving particles:
134 // ---- randomly select one surviving track per species
135 // ---- for this surviving track, further apply a Russian roulette
136 G4int nSecondaries = 0; // -- the number of secondaries to be returned
137 for ( auto typeAndTrack : typesAndTracks )
138 {
139 size_t nTracks = (typeAndTrack.second).size();
140 G4Track* keptTrack;
141 // -- select one track among ones in same species:
142 if ( nTracks > 1 )
143 {
144 auto keptTrackIDX = G4RandFlat::shootInt( nTracks );
145 keptTrack = (typeAndTrack.second)[keptTrackIDX];
146 keptTrack->SetWeight( keptTrack->GetWeight() * nTracks );
147 }
148 else
149 {
150 keptTrack = (typeAndTrack.second)[0];
151 }
152 // -- further apply a Russian Roulette on it:
153 G4bool keepTrack = false;
155 {
157 {
158 keptTrack->SetWeight( keptTrack->GetWeight() / (1. - fRussianRouletteKillingProbability) );
159 keepTrack = true;
160 }
161 }
162 else keepTrack = true;
163 if ( keepTrack )
164 {
165 survivingMap[ keptTrack ] = true;
166 if ( keptTrack != finalStatePrimary ) nSecondaries++;
167 }
168 }
169 // -- and if the leading is not the primary, we have to count it in nSecondaries:
170 if ( secondariesAndPrimary[leadingIDX] != finalStatePrimary ) nSecondaries++;
171
172 // -- verify if the primary is still alive or not after above selection:
173 G4bool primarySurvived = false;
174 if ( finalStatePrimary ) primarySurvived = survivingMap[ finalStatePrimary ];
175
176
177 // -- fill the trimmed particle change:
178 // ---- fill for the primary:
180 if ( primarySurvived )
181 {
182 fParticleChange.ProposeTrackStatus ( wrappedProcessParticleChange->GetTrackStatus() );
183 fParticleChange.ProposeParentWeight ( finalStatePrimary->GetWeight() ); // -- take weight from copy of primary, this one being updated in the random selection loop above
184 fParticleChange.ProposeEnergy ( finalStatePrimary->GetKineticEnergy() );
185 fParticleChange.ProposeMomentumDirection ( finalStatePrimary->GetMomentumDirection() );
186 }
187 else
188 {
192 }
193 // -- fill for surviving secondaries:
196 // ---- note we loop up to on the number of secondaries, which excludes the primary, last in secondariesAndPrimary vector:
198 for ( auto idx = 0 ; idx < wrappedProcessParticleChange->GetNumberOfSecondaries() ; idx++ )
199 {
200 G4Track* secondary = secondariesAndPrimary[idx];
201 // ********************
207 // ******************
208 if ( survivingMap[ secondary ] ) fParticleChange.AddSecondary( secondary );
209 else delete secondary;
210 }
212
213 // -- clean the wrapped process particle change:
214 wrappedProcessParticleChange->Clear();
215
216 if ( finalStatePrimary ) delete finalStatePrimary;
217
218 // -- finally, returns the trimmed particle change:
219 return &fParticleChange;
220
221}
@ fKillTrackAndSecondaries
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4VProcess * GetWrappedProcess() const
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
virtual void Initialize(const G4Track &)
G4double GetWeight() const
void SetWeight(G4double aValue)
void ProposeTrackStatus(G4TrackStatus status)
void SetSecondaryWeightByProcess(G4bool)
void SetNumberOfSecondaries(G4int totSecondaries)
void ProposeParentWeight(G4double finalWeight)
virtual G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &stepData)=0
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

References G4ParticleChange::AddSecondary(), fKillTrackAndSecondaries, fParticleChange, fRussianRouletteKillingProbability, fStopAndKill, G4cout, G4endl, G4UniformRand, G4ParticleChange::GetEnergy(), G4Track::GetKineticEnergy(), G4ParticleChange::GetMomentumDirection(), G4Track::GetMomentumDirection(), G4VProcess::GetProcessName(), G4Track::GetWeight(), G4VParticleChange::GetWeight(), G4BiasingProcessInterface::GetWrappedProcess(), G4ParticleChange::Initialize(), G4VProcess::PostStepDoIt(), G4ParticleChange::ProposeEnergy(), G4ParticleChange::ProposeMomentumDirection(), G4VParticleChange::ProposeParentWeight(), G4VParticleChange::ProposeTrackStatus(), G4Track::SetKineticEnergy(), G4Track::SetMomentumDirection(), G4VParticleChange::SetNumberOfSecondaries(), G4VParticleChange::SetSecondaryWeightByProcess(), and G4Track::SetWeight().

◆ DistanceToApplyOperation()

virtual G4double G4BOptnLeadingParticle::DistanceToApplyOperation ( const G4Track ,
G4double  ,
G4ForceCondition  
)
inlinevirtual

Implements G4VBiasingOperation.

Definition at line 72 of file G4BOptnLeadingParticle.hh.

74 {return 0;}

◆ GenerateBiasingFinalState()

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

Implements G4VBiasingOperation.

Definition at line 75 of file G4BOptnLeadingParticle.hh.

76 {return nullptr;}

◆ GetFurtherKillingProbability()

G4double G4BOptnLeadingParticle::GetFurtherKillingProbability ( ) const
inline

◆ GetName()

const G4String & G4VBiasingOperation::GetName ( ) const
inlineinherited

Definition at line 198 of file G4VBiasingOperation.hh.

198{return fName;}

References G4VBiasingOperation::fName.

◆ GetUniqueID()

std::size_t G4VBiasingOperation::GetUniqueID ( ) const
inlineinherited

Definition at line 199 of file G4VBiasingOperation.hh.

199{return fUniqueID;}

References G4VBiasingOperation::fUniqueID.

◆ ProposeAlongStepLimit()

virtual G4double G4VBiasingOperation::ProposeAlongStepLimit ( const G4BiasingProcessInterface )
inlinevirtualinherited

Reimplemented in G4BOptnForceCommonTruncatedExp.

Definition at line 127 of file G4VBiasingOperation.hh.

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

References DBL_MAX.

Referenced by G4BiasingProcessInterface::AlongStepGetPhysicalInteractionLength().

◆ ProposeGPILSelection()

virtual G4GPILSelection G4VBiasingOperation::ProposeGPILSelection ( const G4GPILSelection  wrappedProcessSelection)
inlinevirtualinherited

Reimplemented in G4BOptnForceCommonTruncatedExp.

Definition at line 131 of file G4VBiasingOperation.hh.

132 {return wrappedProcessSelection;}

Referenced by G4BiasingProcessInterface::AlongStepGetPhysicalInteractionLength().

◆ ProvideOccurenceBiasingInteractionLaw()

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

Implements G4VBiasingOperation.

Definition at line 65 of file G4BOptnLeadingParticle.hh.

65{return nullptr;}

◆ SetFurtherKillingProbability()

void G4BOptnLeadingParticle::SetFurtherKillingProbability ( G4double  p)
inline

Definition at line 85 of file G4BOptnLeadingParticle.hh.

85{ fRussianRouletteKillingProbability = p; } // -- if p <= 0.0 the killing is ignored.

References fRussianRouletteKillingProbability.

Field Documentation

◆ fName

const G4String G4VBiasingOperation::fName
privateinherited

Definition at line 203 of file G4VBiasingOperation.hh.

Referenced by G4VBiasingOperation::GetName().

◆ fParticleChange

G4ParticleChange G4BOptnLeadingParticle::fParticleChange
private

Definition at line 90 of file G4BOptnLeadingParticle.hh.

Referenced by ApplyFinalStateBiasing().

◆ fRussianRouletteKillingProbability

G4double G4BOptnLeadingParticle::fRussianRouletteKillingProbability
private

◆ fUniqueID

std::size_t G4VBiasingOperation::fUniqueID
privateinherited

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