Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4EmBiasingManager Class Reference

#include <G4EmBiasingManager.hh>

Public Member Functions

 G4EmBiasingManager ()
 
 ~G4EmBiasingManager ()
 
void Initialise (const G4ParticleDefinition &part, const G4String &procName, G4int verbose)
 
void ActivateForcedInteraction (G4double length=0.0, const G4String &r="")
 
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
 
G4double GetStepLimit (G4int coupleIdx, G4double previousStep)
 
G4double ApplySecondaryBiasing (std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForGamma *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
 
G4double ApplySecondaryBiasing (std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForLoss *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
 
G4double ApplySecondaryBiasing (std::vector< G4Track * > &, G4int coupleIdx)
 
G4bool SecondaryBiasingRegion (G4int coupleIdx)
 
G4bool ForcedInteractionRegion (G4int coupleIdx)
 
void ResetForcedInteraction ()
 

Detailed Description

Definition at line 67 of file G4EmBiasingManager.hh.

Constructor & Destructor Documentation

G4EmBiasingManager::G4EmBiasingManager ( )

Definition at line 65 of file G4EmBiasingManager.cc.

66  : nForcedRegions(0),nSecBiasedRegions(0),eIonisation(nullptr),
67  currentStepLimit(0.0),startTracking(true)
68 {
69  fSafetyMin = 1.e-6*mm;
70  theElectron = G4Electron::Electron();
71 }
static constexpr double mm
Definition: G4SIunits.hh:115
static G4Electron * Electron()
Definition: G4Electron.cc:94

Here is the call graph for this function:

G4EmBiasingManager::~G4EmBiasingManager ( )

Definition at line 75 of file G4EmBiasingManager.cc.

76 {}

Member Function Documentation

void G4EmBiasingManager::ActivateForcedInteraction ( G4double  length = 0.0,
const G4String r = "" 
)

Definition at line 146 of file G4EmBiasingManager.cc.

148 {
149  G4RegionStore* regionStore = G4RegionStore::GetInstance();
150  G4String name = rname;
151  if(name == "" || name == "world" || name == "World") {
152  name = "DefaultRegionForTheWorld";
153  }
154  const G4Region* reg = regionStore->GetRegion(name, false);
155  if(!reg) {
156  G4cout << "### G4EmBiasingManager::ForcedInteraction WARNING: "
157  << " G4Region <"
158  << rname << "> is unknown" << G4endl;
159  return;
160  }
161 
162  // the region is in the list
163  if (0 < nForcedRegions) {
164  for (G4int i=0; i<nForcedRegions; ++i) {
165  if (reg == forcedRegions[i]) {
166  lengthForRegion[i] = val;
167  return;
168  }
169  }
170  }
171  if(val < 0.0) {
172  G4cout << "### G4EmBiasingManager::ForcedInteraction WARNING: "
173  << val << " < 0.0, so no activation for the G4Region <"
174  << rname << ">" << G4endl;
175  return;
176  }
177 
178  // new region
179  forcedRegions.push_back(reg);
180  lengthForRegion.push_back(val);
181  ++nForcedRegions;
182 }
const XML_Char * name
Definition: expat.h:151
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
int G4int
Definition: G4Types.hh:78
static const G4double reg
static G4RegionStore * GetInstance()
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmBiasingManager::ActivateSecondaryBiasing ( const G4String region,
G4double  factor,
G4double  energyLimit 
)

Definition at line 187 of file G4EmBiasingManager.cc.

190 {
191  //G4cout << "G4EmBiasingManager::ActivateSecondaryBiasing: "
192  // << rname << " F= " << factor << " E(MeV)= " << energyLimit/MeV
193  // << G4endl;
194  G4RegionStore* regionStore = G4RegionStore::GetInstance();
195  G4String name = rname;
196  if(name == "" || name == "world" || name == "World") {
197  name = "DefaultRegionForTheWorld";
198  }
199  const G4Region* reg = regionStore->GetRegion(name, false);
200  if(!reg) {
201  G4cout << "### G4EmBiasingManager::ActivateBremsstrahlungSplitting "
202  << "WARNING: G4Region <"
203  << rname << "> is unknown" << G4endl;
204  return;
205  }
206 
207  // Range cut
208  G4int nsplit = 0;
209  G4double w = factor;
210 
211  // splitting
212  if(factor >= 1.0) {
213  nsplit = G4lrint(factor);
214  w = 1.0/G4double(nsplit);
215 
216  // Russian roulette
217  } else if(0.0 < factor) {
218  nsplit = 1;
219  w = 1.0/factor;
220  }
221 
222  // the region is in the list - overwrite parameters
223  if (0 < nSecBiasedRegions) {
224  for (G4int i=0; i<nSecBiasedRegions; ++i) {
225  if (reg == secBiasedRegions[i]) {
226  secBiasedWeight[i] = w;
227  nBremSplitting[i] = nsplit;
228  secBiasedEnegryLimit[i] = energyLimit;
229  return;
230  }
231  }
232  }
233  /*
234  G4cout << "### G4EmBiasingManager::ActivateSecondaryBiasing: "
235  << " nsplit= " << nsplit << " for the G4Region <"
236  << rname << ">" << G4endl;
237  */
238 
239  // new region
240  secBiasedRegions.push_back(reg);
241  secBiasedWeight.push_back(w);
242  nBremSplitting.push_back(nsplit);
243  secBiasedEnegryLimit.push_back(energyLimit);
244  ++nSecBiasedRegions;
245  //G4cout << "nSecBiasedRegions= " << nSecBiasedRegions << G4endl;
246 }
const XML_Char * name
Definition: expat.h:151
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
int G4int
Definition: G4Types.hh:78
static const G4double reg
static G4RegionStore * GetInstance()
G4GLOB_DLL std::ostream G4cout
int G4lrint(double ad)
Definition: templates.hh:163
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmBiasingManager::ApplySecondaryBiasing ( std::vector< G4DynamicParticle * > &  vd,
const G4Track track,
G4VEmModel currentModel,
G4ParticleChangeForGamma pParticleChange,
G4double eloss,
G4int  coupleIdx,
G4double  tcut,
G4double  safety = 0.0 
)

Definition at line 320 of file G4EmBiasingManager.cc.

329 {
330  G4int index = idxSecBiasedCouple[coupleIdx];
331  G4double weight = 1.0;
332  if(0 <= index) {
333  size_t n = vd.size();
334 
335  // the check cannot be applied per secondary particle
336  // because weight correction is common, so the first
337  // secondary is checked
338  if(0 < n && vd[0]->GetKineticEnergy() < secBiasedEnegryLimit[index]) {
339 
340  G4int nsplit = nBremSplitting[index];
341 
342  // Range cut
343  if(0 == nsplit) {
344  if(safety > fSafetyMin) { ApplyRangeCut(vd, track, eloss, safety); }
345 
346  // Russian Roulette
347  } else if(1 == nsplit) {
348  weight = ApplyRussianRoulette(vd, index);
349 
350  // Splitting
351  } else {
352  G4double tmpEnergy = pPartChange->GetProposedKineticEnergy();
353  G4ThreeVector tmpMomDir = pPartChange->GetProposedMomentumDirection();
354 
355  weight = ApplySplitting(vd, track, currentModel, index, tcut);
356 
357  pPartChange->SetProposedKineticEnergy(tmpEnergy);
358  pPartChange->ProposeMomentumDirection(tmpMomDir);
359  }
360  }
361  }
362  return weight;
363 }
int G4int
Definition: G4Types.hh:78
const G4int n
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4EmBiasingManager::ApplySecondaryBiasing ( std::vector< G4DynamicParticle * > &  vd,
const G4Track track,
G4VEmModel currentModel,
G4ParticleChangeForLoss pParticleChange,
G4double eloss,
G4int  coupleIdx,
G4double  tcut,
G4double  safety = 0.0 
)

Definition at line 272 of file G4EmBiasingManager.cc.

281 {
282  G4int index = idxSecBiasedCouple[coupleIdx];
283  G4double weight = 1.0;
284  if(0 <= index) {
285  size_t n = vd.size();
286 
287  // the check cannot be applied per secondary particle
288  // because weight correction is common, so the first
289  // secondary is checked
290  if(0 < n && vd[0]->GetKineticEnergy() < secBiasedEnegryLimit[index]) {
291 
292  G4int nsplit = nBremSplitting[index];
293 
294  // Range cut
295  if(0 == nsplit) {
296  if(safety > fSafetyMin) { ApplyRangeCut(vd, track, eloss, safety); }
297 
298  // Russian Roulette
299  } else if(1 == nsplit) {
300  weight = ApplyRussianRoulette(vd, index);
301 
302  // Splitting
303  } else {
304  G4double tmpEnergy = pPartChange->GetProposedKineticEnergy();
305  G4ThreeVector tmpMomDir = pPartChange->GetProposedMomentumDirection();
306 
307  weight = ApplySplitting(vd, track, currentModel, index, tcut);
308 
309  pPartChange->SetProposedKineticEnergy(tmpEnergy);
310  pPartChange->ProposeMomentumDirection(tmpMomDir);
311  }
312  }
313  }
314  return weight;
315 }
int G4int
Definition: G4Types.hh:78
const G4int n
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4EmBiasingManager::ApplySecondaryBiasing ( std::vector< G4Track * > &  track,
G4int  coupleIdx 
)

Definition at line 368 of file G4EmBiasingManager.cc.

370 {
371  G4int index = idxSecBiasedCouple[coupleIdx];
372  G4double weight = 1.0;
373  if(0 <= index) {
374  size_t n = track.size();
375 
376  // the check cannot be applied per secondary particle
377  // because weight correction is common, so the first
378  // secondary is checked
379  if(0 < n && track[0]->GetKineticEnergy() < secBiasedEnegryLimit[index]) {
380 
381  G4int nsplit = nBremSplitting[index];
382 
383  // Russian Roulette only
384  if(1 == nsplit) {
385  weight = secBiasedWeight[index];
386  for(size_t k=0; k<n; ++k) {
387  if(G4UniformRand()*weight > 1.0) {
388  const G4Track* t = track[k];
389  delete t;
390  track[k] = 0;
391  }
392  }
393  }
394  }
395  }
396  return weight;
397 }
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
const G4int n
double G4double
Definition: G4Types.hh:76
G4bool G4EmBiasingManager::ForcedInteractionRegion ( G4int  coupleIdx)
inline

Definition at line 177 of file G4EmBiasingManager.hh.

178 {
179  G4bool res = false;
180  if(nForcedRegions > 0) {
181  if(idxForcedCouple[coupleIdx] >= 0) { res = true; }
182  }
183  return res;
184 }
bool G4bool
Definition: G4Types.hh:79

Here is the caller graph for this function:

G4double G4EmBiasingManager::GetStepLimit ( G4int  coupleIdx,
G4double  previousStep 
)

Definition at line 250 of file G4EmBiasingManager.cc.

252 {
253  if(startTracking) {
254  startTracking = false;
255  G4int i = idxForcedCouple[coupleIdx];
256  if(i < 0) {
257  currentStepLimit = DBL_MAX;
258  } else {
259  currentStepLimit = lengthForRegion[i];
260  if(currentStepLimit > 0.0) { currentStepLimit *= G4UniformRand(); }
261  }
262  } else {
263  currentStepLimit -= previousStep;
264  }
265  if(currentStepLimit < 0.0) { currentStepLimit = 0.0; }
266  return currentStepLimit;
267 }
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:97
#define DBL_MAX
Definition: templates.hh:83

Here is the caller graph for this function:

void G4EmBiasingManager::Initialise ( const G4ParticleDefinition part,
const G4String procName,
G4int  verbose 
)

Definition at line 80 of file G4EmBiasingManager.cc.

82 {
83  //G4cout << "G4EmBiasingManager::Initialise for "
84  // << part.GetParticleName()
85  // << " and " << procName << G4endl;
86  const G4ProductionCutsTable* theCoupleTable=
88  size_t numOfCouples = theCoupleTable->GetTableSize();
89 
90  if(0 < nForcedRegions) { idxForcedCouple.resize(numOfCouples, -1); }
91  if(0 < nSecBiasedRegions) { idxSecBiasedCouple.resize(numOfCouples, -1); }
92 
93  // Deexcitation
94  for (size_t j=0; j<numOfCouples; ++j) {
95  const G4MaterialCutsCouple* couple =
96  theCoupleTable->GetMaterialCutsCouple(j);
97  const G4ProductionCuts* pcuts = couple->GetProductionCuts();
98  if(0 < nForcedRegions) {
99  for(G4int i=0; i<nForcedRegions; ++i) {
100  if(forcedRegions[i]) {
101  if(pcuts == forcedRegions[i]->GetProductionCuts()) {
102  idxForcedCouple[j] = i;
103  break;
104  }
105  }
106  }
107  }
108  if(0 < nSecBiasedRegions) {
109  for(G4int i=0; i<nSecBiasedRegions; ++i) {
110  if(secBiasedRegions[i]) {
111  if(pcuts == secBiasedRegions[i]->GetProductionCuts()) {
112  idxSecBiasedCouple[j] = i;
113  break;
114  }
115  }
116  }
117  }
118  }
119  if (nForcedRegions > 0 && 0 < verbose) {
120  G4cout << " Forced Interaction is activated for "
121  << part.GetParticleName() << " and "
122  << procName
123  << " inside G4Regions: " << G4endl;
124  for (G4int i=0; i<nForcedRegions; ++i) {
125  const G4Region* r = forcedRegions[i];
126  if(r) { G4cout << " " << r->GetName() << G4endl; }
127  }
128  }
129  if (nSecBiasedRegions > 0 && 0 < verbose) {
130  G4cout << " Secondary biasing is activated for "
131  << part.GetParticleName() << " and "
132  << procName
133  << " inside G4Regions: " << G4endl;
134  for (G4int i=0; i<nSecBiasedRegions; ++i) {
135  const G4Region* r = secBiasedRegions[i];
136  if(r) {
137  G4cout << " " << r->GetName()
138  << " BiasingWeight= " << secBiasedWeight[i] << G4endl;
139  }
140  }
141  }
142 }
const G4String & GetName() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
#define G4endl
Definition: G4ios.hh:61
G4ProductionCuts * GetProductionCuts() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4EmBiasingManager::ResetForcedInteraction ( )
inline

Definition at line 186 of file G4EmBiasingManager.hh.

187 {
188  startTracking = true;
189 }

Here is the caller graph for this function:

G4bool G4EmBiasingManager::SecondaryBiasingRegion ( G4int  coupleIdx)
inline

Definition at line 168 of file G4EmBiasingManager.hh.

169 {
170  G4bool res = false;
171  if(nSecBiasedRegions > 0) {
172  if(idxSecBiasedCouple[coupleIdx] >= 0) { res = true; }
173  }
174  return res;
175 }
bool G4bool
Definition: G4Types.hh:79

Here is the caller graph for this function:


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