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

G4WeightWindowAlgorithm is a concrete implementation of a weight window algorithm. More...

#include <G4WeightWindowAlgorithm.hh>

Inheritance diagram for G4WeightWindowAlgorithm:

Public Member Functions

 G4WeightWindowAlgorithm (G4double upperLimitFactor=5, G4double survivalFactor=3, G4int maxNumberOfSplits=5)
 ~G4WeightWindowAlgorithm () override=default
G4Nsplit_Weight Calculate (G4double init_w, G4double lowerWeightBound) const override
Public Member Functions inherited from G4VWeightWindowAlgorithm
 G4VWeightWindowAlgorithm ()=default
virtual ~G4VWeightWindowAlgorithm ()=default

Detailed Description

G4WeightWindowAlgorithm is a concrete implementation of a weight window algorithm.

Definition at line 58 of file G4WeightWindowAlgorithm.hh.

Constructor & Destructor Documentation

◆ G4WeightWindowAlgorithm()

G4WeightWindowAlgorithm::G4WeightWindowAlgorithm ( G4double upperLimitFactor = 5,
G4double survivalFactor = 3,
G4int maxNumberOfSplits = 5 )

Constructor. The arguments configure the algorithm. In case of upperLimitFactor=survivalFactor=1 the algorithm becomes the expected weight algorithm of the importance sampling technique.

Parameters
[in]upperLimitFactorThe factor defining the upper weight limit W_u = upperLimitFactor * W_l (W_l lower weight bound).
[in]survivalFactoruUsed in calculating the survival weight W_s = survivalFactor * W_l.
[in]maxNumberOfSplitsThe maximal number of splits allowed to be created in one go, and the reciprocal of the minimal survival probability in case of Russian roulette.

Definition at line 33 of file G4WeightWindowAlgorithm.cc.

37 : fUpperLimitFactor(upperLimitFactor),
38 fSurvivalFactor(survivalFactor),
39 fMaxNumberOfSplits(maxNumberOfSplits)
40{
41}

◆ ~G4WeightWindowAlgorithm()

G4WeightWindowAlgorithm::~G4WeightWindowAlgorithm ( )
overridedefault

Default Destructor.

Member Function Documentation

◆ Calculate()

G4Nsplit_Weight G4WeightWindowAlgorithm::Calculate ( G4double init_w,
G4double lowerWeightBound ) const
overridevirtual

Calculates the number of tracks and their weight according to the initial track weight and the lower energy bound.

Implements G4VWeightWindowAlgorithm.

Definition at line 44 of file G4WeightWindowAlgorithm.cc.

46{
47 G4double survivalWeight = lowerWeightBound * fSurvivalFactor;
48 G4double upperWeight = lowerWeightBound * fUpperLimitFactor;
49
50 // initialize return value for case weight in window
51 //
52 G4Nsplit_Weight nw;
53 nw.fN = 1;
54 nw.fW = init_w;
55
56 if (init_w > upperWeight)
57 {
58 // splitting
59 //
60 G4double temp_wi_ws = init_w/upperWeight;
61 auto split_i = static_cast<G4int>(temp_wi_ws);
62 if(split_i != temp_wi_ws) { ++split_i; }
63 G4double wi_ws = init_w/split_i;
64
65 // values in case integer mode or in csae of double
66 // mode and the lower number of splits will be diced
67 //
68 nw.fN = split_i;
69 nw.fW = wi_ws;
70
71//TB if (wi_ws <= fMaxNumberOfSplits) {
72//TB if (wi_ws > int_wi_ws) {
73//TB // double mode
74//TB G4double p2 = wi_ws - int_wi_ws;
75//TB G4double r = G4UniformRand();
76//TB if (r<p2) {
77//TB nw.fN = int_wi_ws + 1;
78//TB }
79//TB }
80//TB }
81//TB else {
82//TB // fMaxNumberOfSplits < wi_ws
83//TB nw.fN = fMaxNumberOfSplits;
84//TB nw.fW = init_w/fMaxNumberOfSplits;
85//TB }
86
87 }
88 else if (init_w < lowerWeightBound) // Russian roulette
89 {
90 G4double wi_ws = init_w/survivalWeight;
91 G4double p = std::max(wi_ws,1./fMaxNumberOfSplits);
93 if (r<p)
94 {
95 nw.fW = init_w/p;
96 nw.fN = 1;
97 }
98 else
99 {
100 nw.fW = 0;
101 nw.fN = 0;
102 }
103 }
104 // else do nothing
105
106 return nw;
107}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52

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