Geant4-11
G4BOptnLeadingParticle.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//
28//---------------------------------------------------------------------
29//
30// G4BOptnLeadingParticle
31//
32// Class Description:
33// A G4VBiasingOperation that implements the so-called "Leading
34// particle biasing scheme". It is of interest in the shield problem
35// to estimate the flux leaking from the shield.
36// It works as follows:
37// - it is intented for hadronic inelastic interaction
38// - at each interaction, are kept:
39// - the most energetic particle (the leading particle)
40// - with unmodified weight
41// - randomly one particle of each species
42// - with this particle weight = n * primary_weight where
43// n is the number of particles of this species
44//---------------------------------------------------------------------
45// Initial version Nov. 2019 M. Verderi
46
47
48#ifndef G4BOptnLeadingParticle_hh
49#define G4BOptnLeadingParticle_hh 1
50
52#include "G4ParticleChange.hh"
53
55public:
56 // -- Constructor :
58 // -- destructor:
60
61public:
62 // -- Methods from G4VBiasingOperation interface:
63 // ----------------------------------------------
64 // -- Unused:
66 // -- Used:
67 virtual G4VParticleChange* ApplyFinalStateBiasing( const G4BiasingProcessInterface*, // -- Method used for this biasing. The related biasing operator
68 const G4Track*, // -- returns this biasing operation at the post step do it level
69 const G4Step*, // -- when the wrapped process has won the interaction length race.
70 G4bool& ); // -- The wrapped process final state is then trimmed.
71 // -- Unused:
74 G4ForceCondition*) {return 0;}
76 const G4Step* ) {return nullptr;}
77
78public:
79 // -- The possibility is given to further apply a Russian roulette on tracks that are accompagnying the leading particle
80 // -- after the classical leading particle biasing algorithm has been applied.
81 // -- This is of interest when applying the technique to e+ -> gamma gamma for example. Given one gamma is leading,
82 // -- the second one is alone in its category, hence selected. With the Russian roulette it is then possible to keep
83 // -- this one randomly. This is also of interest for pi0 decays, or for brem. e- -> e- gamma where the e- or gamma
84 // -- are alone in their category.
85 void SetFurtherKillingProbability( G4double p ) { fRussianRouletteKillingProbability = p; } // -- if p <= 0.0 the killing is ignored.
87
88private:
89 // -- Particle change used to return the trimmed final state:
92
93
94};
95
96#endif
G4ForceCondition
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
virtual G4double DistanceToApplyOperation(const G4Track *, G4double, G4ForceCondition *)
void SetFurtherKillingProbability(G4double p)
G4double GetFurtherKillingProbability() const
virtual const G4VBiasingInteractionLaw * ProvideOccurenceBiasingInteractionLaw(const G4BiasingProcessInterface *, G4ForceCondition &)
virtual G4VParticleChange * ApplyFinalStateBiasing(const G4BiasingProcessInterface *, const G4Track *, const G4Step *, G4bool &)
virtual G4VParticleChange * GenerateBiasingFinalState(const G4Track *, const G4Step *)
Definition: G4Step.hh:62
const char * name(G4int ptype)