Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 ()
 
G4bool CheckDirection (G4ThreeVector pos, G4ThreeVector momdir) const
 
G4bool GetDirectionalSplitting ()
 
void SetDirectionalSplitting (G4bool v)
 
void SetDirectionalSplittingTarget (G4ThreeVector v)
 
void SetDirectionalSplittingRadius (G4double r)
 
G4double GetWeight (G4int i)
 

Detailed Description

Definition at line 66 of file G4EmBiasingManager.hh.

Constructor & Destructor Documentation

◆ G4EmBiasingManager()

G4EmBiasingManager::G4EmBiasingManager ( )

Definition at line 67 of file G4EmBiasingManager.cc.

68 : nForcedRegions(0),nSecBiasedRegions(0),eIonisation(nullptr),
69 currentStepLimit(0.0),startTracking(true)
70{
71 fSafetyMin = 1.e-6*mm;
72 theElectron = G4Electron::Electron();
73 theGamma = G4Gamma::Gamma();
74
75 fDirectionalSplitting = false;
76 fDirectionalSplittingRadius = 0.;
77 fDirectionalSplittingTarget = G4ThreeVector(0.,0.,0.);
78 fDirectionalSplittingWeights.clear();
79}
CLHEP::Hep3Vector G4ThreeVector
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85

◆ ~G4EmBiasingManager()

G4EmBiasingManager::~G4EmBiasingManager ( )

Definition at line 83 of file G4EmBiasingManager.cc.

84{}

Member Function Documentation

◆ ActivateForcedInteraction()

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

Definition at line 169 of file G4EmBiasingManager.cc.

171{
173 G4String name = rname;
174 if(name == "" || name == "world" || name == "World") {
175 name = "DefaultRegionForTheWorld";
176 }
177 const G4Region* reg = regionStore->GetRegion(name, false);
178 if(!reg) {
179 G4cout << "### G4EmBiasingManager::ForcedInteraction WARNING: "
180 << " G4Region <"
181 << rname << "> is unknown" << G4endl;
182 return;
183 }
184
185 // the region is in the list
186 if (0 < nForcedRegions) {
187 for (G4int i=0; i<nForcedRegions; ++i) {
188 if (reg == forcedRegions[i]) {
189 lengthForRegion[i] = val;
190 return;
191 }
192 }
193 }
194 if(val < 0.0) {
195 G4cout << "### G4EmBiasingManager::ForcedInteraction WARNING: "
196 << val << " < 0.0, so no activation for the G4Region <"
197 << rname << ">" << G4endl;
198 return;
199 }
200
201 // new region
202 forcedRegions.push_back(reg);
203 lengthForRegion.push_back(val);
204 ++nForcedRegions;
205}
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
const char * name(G4int ptype)

Referenced by G4VEnergyLossProcess::ActivateForcedInteraction(), and G4VEmProcess::ActivateForcedInteraction().

◆ ActivateSecondaryBiasing()

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

Definition at line 210 of file G4EmBiasingManager.cc.

213{
214 //G4cout << "G4EmBiasingManager::ActivateSecondaryBiasing: "
215 // << rname << " F= " << factor << " E(MeV)= " << energyLimit/MeV
216 // << G4endl;
218 G4String name = rname;
219 if(name == "" || name == "world" || name == "World") {
220 name = "DefaultRegionForTheWorld";
221 }
222 const G4Region* reg = regionStore->GetRegion(name, false);
223 if(!reg) {
224 G4cout << "### G4EmBiasingManager::ActivateBremsstrahlungSplitting "
225 << "WARNING: G4Region <"
226 << rname << "> is unknown" << G4endl;
227 return;
228 }
229
230 // Range cut
231 G4int nsplit = 0;
232 G4double w = factor;
233
234 // splitting
235 if(factor >= 1.0) {
236 nsplit = G4lrint(factor);
237 w = 1.0/G4double(nsplit);
238
239 // Russian roulette
240 } else if(0.0 < factor) {
241 nsplit = 1;
242 w = 1.0/factor;
243 }
244
245 // the region is in the list - overwrite parameters
246 if (0 < nSecBiasedRegions) {
247 for (G4int i=0; i<nSecBiasedRegions; ++i) {
248 if (reg == secBiasedRegions[i]) {
249 secBiasedWeight[i] = w;
250 nBremSplitting[i] = nsplit;
251 secBiasedEnegryLimit[i] = energyLimit;
252 return;
253 }
254 }
255 }
256 /*
257 G4cout << "### G4EmBiasingManager::ActivateSecondaryBiasing: "
258 << " nsplit= " << nsplit << " for the G4Region <"
259 << rname << ">" << G4endl;
260 */
261
262 // new region
263 secBiasedRegions.push_back(reg);
264 secBiasedWeight.push_back(w);
265 nBremSplitting.push_back(nsplit);
266 secBiasedEnegryLimit.push_back(energyLimit);
267 ++nSecBiasedRegions;
268 //G4cout << "nSecBiasedRegions= " << nSecBiasedRegions << G4endl;
269}
double G4double
Definition: G4Types.hh:83
int G4lrint(double ad)
Definition: templates.hh:134

Referenced by G4VEmProcess::ActivateSecondaryBiasing(), and G4VEnergyLossProcess::ActivateSecondaryBiasing().

◆ ApplySecondaryBiasing() [1/3]

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 348 of file G4EmBiasingManager.cc.

357{
358 G4int index = idxSecBiasedCouple[coupleIdx];
359 G4double weight = 1.;
360 if(0 <= index) {
361 size_t n = vd.size();
362
363 // the check cannot be applied per secondary particle
364 // because weight correction is common, so the first
365 // secondary is checked
366 if((0 < n && vd[0]->GetKineticEnergy() < secBiasedEnegryLimit[index])
367 || fDirectionalSplitting) {
368
369 G4int nsplit = nBremSplitting[index];
370
371 // Range cut
372 if(0 == nsplit) {
373 if(safety > fSafetyMin) { ApplyRangeCut(vd, track, eloss, safety); }
374
375 // Russian Roulette
376 } else if(1 == nsplit) {
377 weight = ApplyRussianRoulette(vd, index);
378
379 // Splitting
380 } else {
381 if (fDirectionalSplitting) {
382 weight = ApplyDirectionalSplitting(vd, track, currentModel,
383 index, tcut, pPartChange);
384 } else {
385 G4double tmpEnergy = pPartChange->GetProposedKineticEnergy();
386 G4ThreeVector tmpMomDir = pPartChange->GetProposedMomentumDirection();
387
388 weight = ApplySplitting(vd, track, currentModel, index, tcut);
389
390 pPartChange->SetProposedKineticEnergy(tmpEnergy);
391 pPartChange->ProposeMomentumDirection(tmpMomDir);
392 }
393 }
394 }
395 }
396 return weight;
397}

Referenced by G4eplusAnnihilation::AtRestDoIt(), G4VEmProcess::PostStepDoIt(), and G4VEnergyLossProcess::PostStepDoIt().

◆ ApplySecondaryBiasing() [2/3]

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 295 of file G4EmBiasingManager.cc.

304{
305 G4int index = idxSecBiasedCouple[coupleIdx];
306 G4double weight = 1.;
307 if(0 <= index) {
308 size_t n = vd.size();
309
310 // the check cannot be applied per secondary particle
311 // because weight correction is common, so the first
312 // secondary is checked
313 if((0 < n && vd[0]->GetKineticEnergy() < secBiasedEnegryLimit[index])
314 || fDirectionalSplitting) {
315
316 G4int nsplit = nBremSplitting[index];
317
318 // Range cut
319 if(0 == nsplit) {
320 if(safety > fSafetyMin) { ApplyRangeCut(vd, track, eloss, safety); }
321
322 // Russian Roulette
323 } else if(1 == nsplit) {
324 weight = ApplyRussianRoulette(vd, index);
325
326 // Splitting
327 } else {
328 if (fDirectionalSplitting) {
329 weight = ApplyDirectionalSplitting(vd, track, currentModel, index, tcut);
330 } else {
331 G4double tmpEnergy = pPartChange->GetProposedKineticEnergy();
332 G4ThreeVector tmpMomDir = pPartChange->GetProposedMomentumDirection();
333
334 weight = ApplySplitting(vd, track, currentModel, index, tcut);
335
336 pPartChange->SetProposedKineticEnergy(tmpEnergy);
337 pPartChange->ProposeMomentumDirection(tmpMomDir);
338 }
339 }
340 }
341 }
342 return weight;
343}

◆ ApplySecondaryBiasing() [3/3]

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

Definition at line 402 of file G4EmBiasingManager.cc.

404{
405 G4int index = idxSecBiasedCouple[coupleIdx];
406 G4double weight = 1.;
407 if(0 <= index) {
408 size_t n = track.size();
409
410 // the check cannot be applied per secondary particle
411 // because weight correction is common, so the first
412 // secondary is checked
413 if(0 < n && track[0]->GetKineticEnergy() < secBiasedEnegryLimit[index]) {
414
415 G4int nsplit = nBremSplitting[index];
416
417 // Russian Roulette only
418 if(1 == nsplit) {
419 weight = secBiasedWeight[index];
420 for(size_t k=0; k<n; ++k) {
421 if(G4UniformRand()*weight > 1.0) {
422 const G4Track* t = track[k];
423 delete t;
424 track[k] = 0;
425 }
426 }
427 }
428 }
429 }
430 return weight;
431}
#define G4UniformRand()
Definition: Randomize.hh:52

◆ CheckDirection()

G4bool G4EmBiasingManager::CheckDirection ( G4ThreeVector  pos,
G4ThreeVector  momdir 
) const

Definition at line 463 of file G4EmBiasingManager.cc.

465{
466 G4ThreeVector delta = fDirectionalSplittingTarget - pos;
467 G4double angle = momdir.angle(delta);
468 G4double dist = delta.cross(momdir).mag();
469 if (dist <= fDirectionalSplittingRadius && angle < halfpi) {
470 return true;
471 }
472 return false;
473}
Hep3Vector cross(const Hep3Vector &) const
double angle(const Hep3Vector &) const
double mag() const

◆ ForcedInteractionRegion()

G4bool G4EmBiasingManager::ForcedInteractionRegion ( G4int  coupleIdx)
inline

Definition at line 206 of file G4EmBiasingManager.hh.

207{
208 G4bool res = false;
209 if(nForcedRegions > 0) {
210 if(idxForcedCouple[coupleIdx] >= 0) { res = true; }
211 }
212 return res;
213}
bool G4bool
Definition: G4Types.hh:86

Referenced by G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), and G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength().

◆ GetDirectionalSplitting()

G4bool G4EmBiasingManager::GetDirectionalSplitting ( )
inline

Definition at line 125 of file G4EmBiasingManager.hh.

125{ return fDirectionalSplitting; }

◆ GetStepLimit()

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

Definition at line 273 of file G4EmBiasingManager.cc.

275{
276 if(startTracking) {
277 startTracking = false;
278 G4int i = idxForcedCouple[coupleIdx];
279 if(i < 0) {
280 currentStepLimit = DBL_MAX;
281 } else {
282 currentStepLimit = lengthForRegion[i];
283 if(currentStepLimit > 0.0) { currentStepLimit *= G4UniformRand(); }
284 }
285 } else {
286 currentStepLimit -= previousStep;
287 }
288 if(currentStepLimit < 0.0) { currentStepLimit = 0.0; }
289 return currentStepLimit;
290}
#define DBL_MAX
Definition: templates.hh:62

Referenced by G4VEmProcess::PostStepGetPhysicalInteractionLength(), and G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength().

◆ GetWeight()

G4double G4EmBiasingManager::GetWeight ( G4int  i)

Definition at line 636 of file G4EmBiasingManager.cc.

637{
638 // normally return 1. If a directionally split particle survives RR,
639 // return 1./(splitting factor)
640 if (fDirectionalSplittingWeights.size() >= (unsigned int)(i+1) ) {
641 G4double w = fDirectionalSplittingWeights[i];
642 fDirectionalSplittingWeights[i] = 1.; // ensure it's not used again
643 return w;
644 } else {
645 return 1.;
646 }
647}

Referenced by G4eplusAnnihilation::AtRestDoIt(), G4VEmProcess::PostStepDoIt(), and G4VEnergyLossProcess::PostStepDoIt().

◆ Initialise()

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

Definition at line 88 of file G4EmBiasingManager.cc.

90{
91 //G4cout << "G4EmBiasingManager::Initialise for "
92 // << part.GetParticleName()
93 // << " and " << procName << G4endl;
94 const G4ProductionCutsTable* theCoupleTable=
96 size_t numOfCouples = theCoupleTable->GetTableSize();
97
98 if(0 < nForcedRegions) { idxForcedCouple.resize(numOfCouples, -1); }
99 if(0 < nSecBiasedRegions) { idxSecBiasedCouple.resize(numOfCouples, -1); }
100
101 // Deexcitation
102 for (size_t j=0; j<numOfCouples; ++j) {
103 const G4MaterialCutsCouple* couple =
104 theCoupleTable->GetMaterialCutsCouple(j);
105 const G4ProductionCuts* pcuts = couple->GetProductionCuts();
106 if(0 < nForcedRegions) {
107 for(G4int i=0; i<nForcedRegions; ++i) {
108 if(forcedRegions[i]) {
109 if(pcuts == forcedRegions[i]->GetProductionCuts()) {
110 idxForcedCouple[j] = i;
111 break;
112 }
113 }
114 }
115 }
116 if(0 < nSecBiasedRegions) {
117 for(G4int i=0; i<nSecBiasedRegions; ++i) {
118 if(secBiasedRegions[i]) {
119 if(pcuts == secBiasedRegions[i]->GetProductionCuts()) {
120 idxSecBiasedCouple[j] = i;
121 break;
122 }
123 }
124 }
125 }
126 }
127
130 if (fDirectionalSplitting) {
133 }
134
135 if (nForcedRegions > 0 && 0 < verbose) {
136 G4cout << " Forced Interaction is activated for "
137 << part.GetParticleName() << " and "
138 << procName
139 << " inside G4Regions: " << G4endl;
140 for (G4int i=0; i<nForcedRegions; ++i) {
141 const G4Region* r = forcedRegions[i];
142 if(r) { G4cout << " " << r->GetName() << G4endl; }
143 }
144 }
145 if (nSecBiasedRegions > 0 && 0 < verbose) {
146 G4cout << " Secondary biasing is activated for "
147 << part.GetParticleName() << " and "
148 << procName
149 << " inside G4Regions: " << G4endl;
150 for (G4int i=0; i<nSecBiasedRegions; ++i) {
151 const G4Region* r = secBiasedRegions[i];
152 if(r) {
153 G4cout << " " << r->GetName()
154 << " BiasingWeight= " << secBiasedWeight[i] << G4endl;
155 }
156 }
157 if (fDirectionalSplitting) {
158 G4cout << " Directional splitting activated, with target position: "
159 << fDirectionalSplittingTarget/cm
160 << " cm; radius: "
161 << fDirectionalSplittingRadius/cm
162 << "cm." << G4endl;
163 }
164 }
165}
void SetDirectionalSplittingRadius(G4double r)
void SetDirectionalSplitting(G4bool v)
void SetDirectionalSplittingTarget(G4ThreeVector v)
static G4EmParameters * Instance()
G4double GetDirectionalSplittingRadius()
G4bool GetDirectionalSplitting() const
G4ThreeVector GetDirectionalSplittingTarget() const
G4ProductionCuts * GetProductionCuts() const
const G4String & GetParticleName() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
const G4String & GetName() const

Referenced by G4VEmProcess::PreparePhysicsTable(), and G4VEnergyLossProcess::PreparePhysicsTable().

◆ ResetForcedInteraction()

void G4EmBiasingManager::ResetForcedInteraction ( )
inline

Definition at line 215 of file G4EmBiasingManager.hh.

216{
217 startTracking = true;
218}

Referenced by G4VEmProcess::StartTracking(), and G4VEnergyLossProcess::StartTracking().

◆ SecondaryBiasingRegion()

G4bool G4EmBiasingManager::SecondaryBiasingRegion ( G4int  coupleIdx)
inline

Definition at line 197 of file G4EmBiasingManager.hh.

198{
199 G4bool res = false;
200 if(nSecBiasedRegions > 0) {
201 if(idxSecBiasedCouple[coupleIdx] >= 0) { res = true; }
202 }
203 return res;
204}

Referenced by G4eplusAnnihilation::AtRestDoIt(), G4VEmProcess::PostStepDoIt(), and G4VEnergyLossProcess::PostStepDoIt().

◆ SetDirectionalSplitting()

void G4EmBiasingManager::SetDirectionalSplitting ( G4bool  v)
inline

Definition at line 126 of file G4EmBiasingManager.hh.

126{ fDirectionalSplitting = v; }

Referenced by Initialise().

◆ SetDirectionalSplittingRadius()

void G4EmBiasingManager::SetDirectionalSplittingRadius ( G4double  r)
inline

Definition at line 130 of file G4EmBiasingManager.hh.

131 { fDirectionalSplittingRadius = r; }

Referenced by Initialise().

◆ SetDirectionalSplittingTarget()

void G4EmBiasingManager::SetDirectionalSplittingTarget ( G4ThreeVector  v)
inline

Definition at line 128 of file G4EmBiasingManager.hh.

129 { fDirectionalSplittingTarget = v; }

Referenced by Initialise().


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