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

#include <G4BOptnLeadingParticle.hh>

+ Inheritance diagram for G4BOptnLeadingParticle:

Public Member Functions

 G4BOptnLeadingParticle (const G4String &name)
 
virtual ~G4BOptnLeadingParticle ()
 
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 SetFurtherKillingProbability (G4double p)
 
G4double GetFurtherKillingProbability () 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 51 of file G4BOptnLeadingParticle.hh.

Constructor & Destructor Documentation

◆ G4BOptnLeadingParticle()

G4BOptnLeadingParticle::G4BOptnLeadingParticle ( const G4String & name)

Definition at line 36 of file G4BOptnLeadingParticle.cc.

37 : G4VBiasingOperation( name )
38{
39}
G4VBiasingOperation(const G4String &name)

◆ ~G4BOptnLeadingParticle()

G4BOptnLeadingParticle::~G4BOptnLeadingParticle ( )
virtual

Definition at line 41 of file G4BOptnLeadingParticle.cc.

42{
43}

Member Function Documentation

◆ ApplyFinalStateBiasing()

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

G4cout << G4endl;

Implements G4VBiasingOperation.

Definition at line 45 of file G4BOptnLeadingParticle.cc.

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

◆ DistanceToApplyOperation()

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

Implements G4VBiasingOperation.

Definition at line 74 of file G4BOptnLeadingParticle.hh.

74{ return 0.0; }

◆ GenerateBiasingFinalState()

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

Implements G4VBiasingOperation.

Definition at line 76 of file G4BOptnLeadingParticle.hh.

76{ return nullptr; }

◆ GetFurtherKillingProbability()

G4double G4BOptnLeadingParticle::GetFurtherKillingProbability ( ) const
inline

Definition at line 88 of file G4BOptnLeadingParticle.hh.

89 {
90 return fRussianRouletteKillingProbability;
91 }

◆ ProvideOccurenceBiasingInteractionLaw()

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

Implements G4VBiasingOperation.

Definition at line 64 of file G4BOptnLeadingParticle.hh.

65 { return nullptr; }

◆ SetFurtherKillingProbability()

void G4BOptnLeadingParticle::SetFurtherKillingProbability ( G4double p)
inline

Definition at line 84 of file G4BOptnLeadingParticle.hh.

85 {
86 fRussianRouletteKillingProbability = p;
87 }

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