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

#include <G4WeightWindowAlgorithm.hh>

+ Inheritance diagram for G4WeightWindowAlgorithm:

Public Member Functions

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

Detailed Description

Definition at line 53 of file G4WeightWindowAlgorithm.hh.

Constructor & Destructor Documentation

◆ G4WeightWindowAlgorithm()

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

Definition at line 33 of file G4WeightWindowAlgorithm.cc.

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

◆ ~G4WeightWindowAlgorithm()

G4WeightWindowAlgorithm::~G4WeightWindowAlgorithm ( )
overridedefault

Member Function Documentation

◆ Calculate()

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

Implements G4VWeightWindowAlgorithm.

Definition at line 46 of file G4WeightWindowAlgorithm.cc.

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