Geant4 11.1.1
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)
 
 G4EmBiasingManager (G4EmBiasingManager &)=delete
 
G4EmBiasingManageroperator= (const G4EmBiasingManager &right)=delete
 

Detailed Description

Definition at line 66 of file G4EmBiasingManager.hh.

Constructor & Destructor Documentation

◆ G4EmBiasingManager() [1/2]

G4EmBiasingManager::G4EmBiasingManager ( )

Definition at line 67 of file G4EmBiasingManager.cc.

68 : fDirectionalSplittingTarget(0.0,0.0,0.0)
69{
70 fSafetyMin = 1.e-6*mm;
71 theElectron = G4Electron::Electron();
72 theGamma = G4Gamma::Gamma();
73}
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85

◆ ~G4EmBiasingManager()

G4EmBiasingManager::~G4EmBiasingManager ( )
default

◆ G4EmBiasingManager() [2/2]

G4EmBiasingManager::G4EmBiasingManager ( G4EmBiasingManager )
delete

Member Function Documentation

◆ ActivateForcedInteraction()

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

Definition at line 162 of file G4EmBiasingManager.cc.

164{
166 G4String name = rname;
167 if(name == "" || name == "world" || name == "World") {
168 name = "DefaultRegionForTheWorld";
169 }
170 const G4Region* reg = regionStore->GetRegion(name, false);
171 if(!reg) {
172 G4cout << "### G4EmBiasingManager::ForcedInteraction WARNING: "
173 << " G4Region <"
174 << rname << "> is unknown" << G4endl;
175 return;
176 }
177
178 // the region is in the list
179 if (0 < nForcedRegions) {
180 for (G4int i=0; i<nForcedRegions; ++i) {
181 if (reg == forcedRegions[i]) {
182 lengthForRegion[i] = val;
183 return;
184 }
185 }
186 }
187 if(val < 0.0) {
188 G4cout << "### G4EmBiasingManager::ForcedInteraction WARNING: "
189 << val << " < 0.0, so no activation for the G4Region <"
190 << rname << ">" << G4endl;
191 return;
192 }
193
194 // new region
195 forcedRegions.push_back(reg);
196 lengthForRegion.push_back(val);
197 ++nForcedRegions;
198}
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 203 of file G4EmBiasingManager.cc.

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

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

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

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

◆ ApplySecondaryBiasing() [3/3]

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

Definition at line 395 of file G4EmBiasingManager.cc.

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

◆ CheckDirection()

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

Definition at line 455 of file G4EmBiasingManager.cc.

457{
458 G4ThreeVector delta = fDirectionalSplittingTarget - pos;
459 G4double angle = momdir.angle(delta);
460 G4double dist = delta.cross(momdir).mag();
461 if (dist <= fDirectionalSplittingRadius && angle < halfpi) {
462 return true;
463 }
464 return false;
465}
Hep3Vector cross(const Hep3Vector &) const
double angle(const Hep3Vector &) const
double mag() const

◆ ForcedInteractionRegion()

G4bool G4EmBiasingManager::ForcedInteractionRegion ( G4int  coupleIdx)
inline

Definition at line 209 of file G4EmBiasingManager.hh.

210{
211 G4bool res = false;
212 if(nForcedRegions > 0) {
213 if(idxForcedCouple[coupleIdx] >= 0) { res = true; }
214 }
215 return res;
216}
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 266 of file G4EmBiasingManager.cc.

268{
269 if(startTracking) {
270 startTracking = false;
271 G4int i = idxForcedCouple[coupleIdx];
272 if(i < 0) {
273 currentStepLimit = DBL_MAX;
274 } else {
275 currentStepLimit = lengthForRegion[i];
276 if(currentStepLimit > 0.0) { currentStepLimit *= G4UniformRand(); }
277 }
278 } else {
279 currentStepLimit -= previousStep;
280 }
281 if(currentStepLimit < 0.0) { currentStepLimit = 0.0; }
282 return currentStepLimit;
283}
#define DBL_MAX
Definition: templates.hh:62

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

◆ GetWeight()

G4double G4EmBiasingManager::GetWeight ( G4int  i)

Definition at line 628 of file G4EmBiasingManager.cc.

629{
630 // normally return 1. If a directionally split particle survives RR,
631 // return 1./(splitting factor)
632 if (fDirectionalSplittingWeights.size() >= (unsigned int)(i+1) ) {
633 G4double w = fDirectionalSplittingWeights[i];
634 fDirectionalSplittingWeights[i] = 1.; // ensure it's not used again
635 return w;
636 } else {
637 return 1.;
638 }
639}

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

◆ Initialise()

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

Definition at line 81 of file G4EmBiasingManager.cc.

83{
84 //G4cout << "G4EmBiasingManager::Initialise for "
85 // << part.GetParticleName()
86 // << " and " << procName << G4endl;
87 const G4ProductionCutsTable* theCoupleTable=
89 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
90
91 if(0 < nForcedRegions) { idxForcedCouple.resize(numOfCouples, -1); }
92 if(0 < nSecBiasedRegions) { idxSecBiasedCouple.resize(numOfCouples, -1); }
93
94 // Deexcitation
95 for (G4int j=0; j<numOfCouples; ++j) {
96 const G4MaterialCutsCouple* couple =
97 theCoupleTable->GetMaterialCutsCouple(j);
98 const G4ProductionCuts* pcuts = couple->GetProductionCuts();
99 if(0 < nForcedRegions) {
100 for(G4int i=0; i<nForcedRegions; ++i) {
101 if(forcedRegions[i]) {
102 if(pcuts == forcedRegions[i]->GetProductionCuts()) {
103 idxForcedCouple[j] = i;
104 break;
105 }
106 }
107 }
108 }
109 if(0 < nSecBiasedRegions) {
110 for(G4int i=0; i<nSecBiasedRegions; ++i) {
111 if(secBiasedRegions[i]) {
112 if(pcuts == secBiasedRegions[i]->GetProductionCuts()) {
113 idxSecBiasedCouple[j] = i;
114 break;
115 }
116 }
117 }
118 }
119 }
120
123 if (fDirectionalSplitting) {
126 }
127
128 if (nForcedRegions > 0 && 0 < verbose) {
129 G4cout << " Forced Interaction is activated for "
130 << part.GetParticleName() << " and "
131 << procName
132 << " inside G4Regions: " << G4endl;
133 for (G4int i=0; i<nForcedRegions; ++i) {
134 const G4Region* r = forcedRegions[i];
135 if(r) { G4cout << " " << r->GetName() << G4endl; }
136 }
137 }
138 if (nSecBiasedRegions > 0 && 0 < verbose) {
139 G4cout << " Secondary biasing is activated for "
140 << part.GetParticleName() << " and "
141 << procName
142 << " inside G4Regions: " << G4endl;
143 for (G4int i=0; i<nSecBiasedRegions; ++i) {
144 const G4Region* r = secBiasedRegions[i];
145 if(r) {
146 G4cout << " " << r->GetName()
147 << " BiasingWeight= " << secBiasedWeight[i] << G4endl;
148 }
149 }
150 if (fDirectionalSplitting) {
151 G4cout << " Directional splitting activated, with target position: "
152 << fDirectionalSplittingTarget/cm
153 << " cm; radius: "
154 << fDirectionalSplittingRadius/cm
155 << "cm." << G4endl;
156 }
157 }
158}
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().

◆ operator=()

G4EmBiasingManager & G4EmBiasingManager::operator= ( const G4EmBiasingManager right)
delete

◆ ResetForcedInteraction()

void G4EmBiasingManager::ResetForcedInteraction ( )
inline

Definition at line 218 of file G4EmBiasingManager.hh.

219{
220 startTracking = true;
221}

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

◆ SecondaryBiasingRegion()

G4bool G4EmBiasingManager::SecondaryBiasingRegion ( G4int  coupleIdx)
inline

Definition at line 200 of file G4EmBiasingManager.hh.

201{
202 G4bool res = false;
203 if(nSecBiasedRegions > 0) {
204 if(idxSecBiasedCouple[coupleIdx] >= 0) { res = true; }
205 }
206 return res;
207}

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: