Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
45 G4double lowerWeightBound) const
46{
47 G4double survivalWeight = lowerWeightBound * fSurvivalFactor;
48 G4double upperWeight = lowerWeightBound * fUpperLimitFactor;
49
50 // initialize return value for case weight in window
51 //
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}
108
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
G4Nsplit_Weight is a class (struct) used by importance sampling. It contains the number of tracks a m...
G4Nsplit_Weight Calculate(G4double init_w, G4double lowerWeightBound) const override
G4WeightWindowAlgorithm(G4double upperLimitFactor=5, G4double survivalFactor=3, G4int maxNumberOfSplits=5)