Geant4-11
G4EmBiasingManager.hh
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class header file
30//
31//
32// File name: G4EmBiasingManager
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 28.07.2011
37//
38// Modifications:
39//
40// Class Description:
41//
42// It is a class providing step limit for forced process biasing
43
44// -------------------------------------------------------------------
45//
46
47#ifndef G4EmBiasingManager_h
48#define G4EmBiasingManager_h 1
49
50#include "globals.hh"
52#include "G4DynamicParticle.hh"
53#include "Randomize.hh"
54#include <vector>
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57
58class G4Region;
59class G4Track;
61class G4VEmModel;
65
67{
68public:
69
71
73
74 void Initialise(const G4ParticleDefinition& part,
75 const G4String& procName, G4int verbose);
76
77 // default parameters are possible
78 void ActivateForcedInteraction(G4double length = 0.0,
79 const G4String& r = "");
80
81 // no default parameters
82 void ActivateSecondaryBiasing(const G4String& region, G4double factor,
83 G4double energyLimit);
84
85 // return forced step limit
86 G4double GetStepLimit(G4int coupleIdx, G4double previousStep);
87
88 // return weight of splitting or Russian roulette
89 // G4DynamicParticle may be deleted
90 // two functions are required because of the different ParticleChange
91 // ApplySecondaryBiasing() are wrappers
92
93 // for G4VEmProcess
94 G4double ApplySecondaryBiasing(std::vector<G4DynamicParticle*>&,
95 const G4Track& track,
96 G4VEmModel* currentModel,
97 G4ParticleChangeForGamma* pParticleChange,
98 G4double& eloss,
99 G4int coupleIdx,
100 G4double tcut,
101 G4double safety = 0.0);
102
103 // for G4VEnergyLossProcess
104 G4double ApplySecondaryBiasing(std::vector<G4DynamicParticle*>&,
105 const G4Track& track,
106 G4VEmModel* currentModel,
107 G4ParticleChangeForLoss* pParticleChange,
108 G4double& eloss,
109 G4int coupleIdx,
110 G4double tcut,
111 G4double safety = 0.0);
112
113 // for G4VEnergyLossProcess
114 G4double ApplySecondaryBiasing(std::vector<G4Track*>&,
115 G4int coupleIdx);
116
117 inline G4bool SecondaryBiasingRegion(G4int coupleIdx);
118
119 inline G4bool ForcedInteractionRegion(G4int coupleIdx);
120
121 inline void ResetForcedInteraction();
122
124
127
133
134 // hide copy constructor and assignment operator
137
138private:
139
140 void ApplyRangeCut(std::vector<G4DynamicParticle*>& vd,
141 const G4Track& track,
142 G4double& eloss,
143 G4double safety);
144
145 G4double ApplySplitting(std::vector<G4DynamicParticle*>& vd,
146 const G4Track& track,
147 G4VEmModel* currentModel,
148 G4int index,
149 G4double tcut);
150
151 G4double ApplyDirectionalSplitting(std::vector<G4DynamicParticle*>& vd,
152 const G4Track& track,
153 G4VEmModel* currentModel,
154 G4int index,
155 G4double tcut,
156 G4ParticleChangeForGamma* partChange);
157
158 G4double ApplyDirectionalSplitting(std::vector<G4DynamicParticle*>& vd,
159 const G4Track& track,
160 G4VEmModel* currentModel,
161 G4int index,
162 G4double tcut);
163
164 inline G4double ApplyRussianRoulette(std::vector<G4DynamicParticle*>& vd,
165 G4int index);
166
168
171
175
178
181
183 std::vector<G4double> fDirectionalSplittingWeights;
184
185 std::vector<G4double> lengthForRegion;
186 std::vector<G4double> secBiasedWeight;
187 std::vector<G4double> secBiasedEnegryLimit;
188
189 std::vector<const G4Region*> forcedRegions;
190 std::vector<const G4Region*> secBiasedRegions;
191
192 std::vector<G4int> nBremSplitting;
193 std::vector<G4int> idxForcedCouple;
194 std::vector<G4int> idxSecBiasedCouple;
195
196 std::vector<G4DynamicParticle*> tmpSecondaries;
197};
198
199inline G4bool
201{
202 G4bool res = false;
203 if(nSecBiasedRegions > 0) {
204 if(idxSecBiasedCouple[coupleIdx] >= 0) { res = true; }
205 }
206 return res;
207}
208
210{
211 G4bool res = false;
212 if(nForcedRegions > 0) {
213 if(idxForcedCouple[coupleIdx] >= 0) { res = true; }
214 }
215 return res;
216}
217
219{
220 startTracking = true;
221}
222
223//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
224
225inline G4double
226G4EmBiasingManager::ApplyRussianRoulette(std::vector<G4DynamicParticle*>& vd,
227 G4int index)
228{
229 size_t n = vd.size();
230 G4double weight = secBiasedWeight[index];
231 for(size_t k=0; k<n; ++k) {
232 if(G4UniformRand()*weight > 1.0) {
233 const G4DynamicParticle* dp = vd[k];
234 delete dp;
235 vd[k] = nullptr;
236 }
237 }
238 return weight;
239}
240
241//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
242
243#endif
static const G4double pos
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
G4double ApplySplitting(std::vector< G4DynamicParticle * > &vd, const G4Track &track, G4VEmModel *currentModel, G4int index, G4double tcut)
G4EmBiasingManager(G4EmBiasingManager &)=delete
std::vector< G4DynamicParticle * > tmpSecondaries
std::vector< const G4Region * > forcedRegions
std::vector< G4double > secBiasedEnegryLimit
void SetDirectionalSplittingRadius(G4double r)
void SetDirectionalSplitting(G4bool v)
G4bool ForcedInteractionRegion(G4int coupleIdx)
std::vector< G4double > secBiasedWeight
G4double ApplySecondaryBiasing(std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForGamma *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
G4VEnergyLossProcess * eIonisation
void ApplyRangeCut(std::vector< G4DynamicParticle * > &vd, const G4Track &track, G4double &eloss, G4double safety)
const G4ParticleDefinition * theGamma
std::vector< G4int > idxSecBiasedCouple
void Initialise(const G4ParticleDefinition &part, const G4String &procName, G4int verbose)
G4double GetWeight(G4int i)
std::vector< const G4Region * > secBiasedRegions
std::vector< G4int > nBremSplitting
std::vector< G4double > lengthForRegion
std::vector< G4double > fDirectionalSplittingWeights
G4ThreeVector fDirectionalSplittingTarget
const G4ParticleDefinition * theElectron
G4bool CheckDirection(G4ThreeVector pos, G4ThreeVector momdir) const
G4EmBiasingManager & operator=(const G4EmBiasingManager &right)=delete
G4double fDirectionalSplittingRadius
G4double ApplyDirectionalSplitting(std::vector< G4DynamicParticle * > &vd, const G4Track &track, G4VEmModel *currentModel, G4int index, G4double tcut, G4ParticleChangeForGamma *partChange)
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="")
G4double ApplyRussianRoulette(std::vector< G4DynamicParticle * > &vd, G4int index)
void SetDirectionalSplittingTarget(G4ThreeVector v)
G4bool SecondaryBiasingRegion(G4int coupleIdx)
G4double GetStepLimit(G4int coupleIdx, G4double previousStep)
std::vector< G4int > idxForcedCouple