Geant4 9.6.0
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 upperLimitFaktor=5, G4double survivalFaktor=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 58 of file G4WeightWindowAlgorithm.hh.

Constructor & Destructor Documentation

◆ G4WeightWindowAlgorithm()

G4WeightWindowAlgorithm::G4WeightWindowAlgorithm ( G4double  upperLimitFaktor = 5,
G4double  survivalFaktor = 3,
G4int  maxNumberOfSplits = 5 
)

Definition at line 39 of file G4WeightWindowAlgorithm.cc.

42 :
43 fUpperLimitFaktor(upperLimitFaktor),
44 fSurvivalFaktor(survivalFaktor),
45 fMaxNumberOfSplits(maxNumberOfSplits)
46{}

◆ ~G4WeightWindowAlgorithm()

G4WeightWindowAlgorithm::~G4WeightWindowAlgorithm ( )
virtual

Definition at line 48 of file G4WeightWindowAlgorithm.cc.

49{}

Member Function Documentation

◆ Calculate()

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

Implements G4VWeightWindowAlgorithm.

Definition at line 54 of file G4WeightWindowAlgorithm.cc.

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

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