Geant4 11.1.1
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)
 
virtual ~G4WeightWindowAlgorithm ()
 
virtual G4Nsplit_Weight Calculate (G4double init_w, G4double lowerWeightBound) const
 
- Public Member Functions inherited from G4VWeightWindowAlgorithm
 G4VWeightWindowAlgorithm ()
 
virtual ~G4VWeightWindowAlgorithm ()
 
virtual G4Nsplit_Weight Calculate (G4double init_w, G4double lowerWeightBound) const =0
 

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 ( )
virtual

Definition at line 43 of file G4WeightWindowAlgorithm.cc.

44{
45}

Member Function Documentation

◆ Calculate()

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

Implements G4VWeightWindowAlgorithm.

Definition at line 48 of file G4WeightWindowAlgorithm.cc.

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