Geant4-11
G4WeightWindowAlgorithm.cc
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// G4WeightWindowAlgorithm implementation
27//
28// Author: Michael Dressel (CERN), 2003
29// ----------------------------------------------------------------------
31#include "Randomize.hh"
32
34G4WeightWindowAlgorithm(G4double upperLimitFactor,
35 G4double survivalFactor,
36 G4int maxNumberOfSplits)
37 : fUpperLimitFactor(upperLimitFactor),
38 fSurvivalFactor(survivalFactor),
39 fMaxNumberOfSplits(maxNumberOfSplits)
40{
41}
42
44{
45}
46
49 G4double lowerWeightBound) const
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;
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}
112
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
virtual G4Nsplit_Weight Calculate(G4double init_w, G4double lowerWeightBound) const
G4WeightWindowAlgorithm(G4double upperLimitFactor=5, G4double survivalFactor=3, G4int maxNumberOfSplits=5)
T max(const T t1, const T t2)
brief Return the largest of the two arguments