Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmBiasingManager.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4EmBiasingManager
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 28.07.2011
37//
38// Modifications:
39//
40// 31-05-12 D. Sawkey put back in high energy limit for brem, russian roulette
41// 30-05-12 D. Sawkey brem split gammas are unique; do weight tests for
42// brem, russian roulette
43// -------------------------------------------------------------------
44//
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
47
48#include "G4EmBiasingManager.hh"
49#include "G4SystemOfUnits.hh"
53#include "G4ProductionCuts.hh"
54#include "G4Region.hh"
55#include "G4RegionStore.hh"
56#include "G4Track.hh"
57#include "G4Electron.hh"
58#include "G4Gamma.hh"
59#include "G4VEmModel.hh"
60#include "G4LossTableManager.hh"
63#include "G4EmParameters.hh"
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
68 : fDirectionalSplittingTarget(0.0,0.0,0.0)
69{
70 fSafetyMin = 1.e-6*mm;
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
78{}
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81
83 const G4String& procName, G4int verbose)
84{
85 //G4cout << "G4EmBiasingManager::Initialise for "
86 // << part.GetParticleName()
87 // << " and " << procName << G4endl;
88 const G4ProductionCutsTable* theCoupleTable=
90 size_t numOfCouples = theCoupleTable->GetTableSize();
91
92 if(0 < nForcedRegions) { idxForcedCouple.resize(numOfCouples, -1); }
93 if(0 < nSecBiasedRegions) { idxSecBiasedCouple.resize(numOfCouples, -1); }
94
95 // Deexcitation
96 for (size_t j=0; j<numOfCouples; ++j) {
97 const G4MaterialCutsCouple* couple =
98 theCoupleTable->GetMaterialCutsCouple(j);
99 const G4ProductionCuts* pcuts = couple->GetProductionCuts();
100 if(0 < nForcedRegions) {
101 for(G4int i=0; i<nForcedRegions; ++i) {
102 if(forcedRegions[i]) {
103 if(pcuts == forcedRegions[i]->GetProductionCuts()) {
104 idxForcedCouple[j] = i;
105 break;
106 }
107 }
108 }
109 }
110 if(0 < nSecBiasedRegions) {
111 for(G4int i=0; i<nSecBiasedRegions; ++i) {
112 if(secBiasedRegions[i]) {
113 if(pcuts == secBiasedRegions[i]->GetProductionCuts()) {
114 idxSecBiasedCouple[j] = i;
115 break;
116 }
117 }
118 }
119 }
120 }
121
127 }
128
129 if (nForcedRegions > 0 && 0 < verbose) {
130 G4cout << " Forced Interaction is activated for "
131 << part.GetParticleName() << " and "
132 << procName
133 << " inside G4Regions: " << G4endl;
134 for (G4int i=0; i<nForcedRegions; ++i) {
135 const G4Region* r = forcedRegions[i];
136 if(r) { G4cout << " " << r->GetName() << G4endl; }
137 }
138 }
139 if (nSecBiasedRegions > 0 && 0 < verbose) {
140 G4cout << " Secondary biasing is activated for "
141 << part.GetParticleName() << " and "
142 << procName
143 << " inside G4Regions: " << G4endl;
144 for (G4int i=0; i<nSecBiasedRegions; ++i) {
145 const G4Region* r = secBiasedRegions[i];
146 if(r) {
147 G4cout << " " << r->GetName()
148 << " BiasingWeight= " << secBiasedWeight[i] << G4endl;
149 }
150 }
152 G4cout << " Directional splitting activated, with target position: "
154 << " cm; radius: "
156 << "cm." << G4endl;
157 }
158 }
159}
160
161//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
162
164 const G4String& rname)
165{
167 G4String name = rname;
168 if(name == "" || name == "world" || name == "World") {
169 name = "DefaultRegionForTheWorld";
170 }
171 const G4Region* reg = regionStore->GetRegion(name, false);
172 if(!reg) {
173 G4cout << "### G4EmBiasingManager::ForcedInteraction WARNING: "
174 << " G4Region <"
175 << rname << "> is unknown" << G4endl;
176 return;
177 }
178
179 // the region is in the list
180 if (0 < nForcedRegions) {
181 for (G4int i=0; i<nForcedRegions; ++i) {
182 if (reg == forcedRegions[i]) {
183 lengthForRegion[i] = val;
184 return;
185 }
186 }
187 }
188 if(val < 0.0) {
189 G4cout << "### G4EmBiasingManager::ForcedInteraction WARNING: "
190 << val << " < 0.0, so no activation for the G4Region <"
191 << rname << ">" << G4endl;
192 return;
193 }
194
195 // new region
196 forcedRegions.push_back(reg);
197 lengthForRegion.push_back(val);
199}
200
201//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
202
203void
205 G4double factor,
206 G4double energyLimit)
207{
208 //G4cout << "G4EmBiasingManager::ActivateSecondaryBiasing: "
209 // << rname << " F= " << factor << " E(MeV)= " << energyLimit/MeV
210 // << G4endl;
212 G4String name = rname;
213 if(name == "" || name == "world" || name == "World") {
214 name = "DefaultRegionForTheWorld";
215 }
216 const G4Region* reg = regionStore->GetRegion(name, false);
217 if(!reg) {
218 G4cout << "### G4EmBiasingManager::ActivateBremsstrahlungSplitting "
219 << "WARNING: G4Region <"
220 << rname << "> is unknown" << G4endl;
221 return;
222 }
223
224 // Range cut
225 G4int nsplit = 0;
226 G4double w = factor;
227
228 // splitting
229 if(factor >= 1.0) {
230 nsplit = G4lrint(factor);
231 w = 1.0/G4double(nsplit);
232
233 // Russian roulette
234 } else if(0.0 < factor) {
235 nsplit = 1;
236 w = 1.0/factor;
237 }
238
239 // the region is in the list - overwrite parameters
240 if (0 < nSecBiasedRegions) {
241 for (G4int i=0; i<nSecBiasedRegions; ++i) {
242 if (reg == secBiasedRegions[i]) {
243 secBiasedWeight[i] = w;
244 nBremSplitting[i] = nsplit;
245 secBiasedEnegryLimit[i] = energyLimit;
246 return;
247 }
248 }
249 }
250 /*
251 G4cout << "### G4EmBiasingManager::ActivateSecondaryBiasing: "
252 << " nsplit= " << nsplit << " for the G4Region <"
253 << rname << ">" << G4endl;
254 */
255
256 // new region
257 secBiasedRegions.push_back(reg);
258 secBiasedWeight.push_back(w);
259 nBremSplitting.push_back(nsplit);
260 secBiasedEnegryLimit.push_back(energyLimit);
262 //G4cout << "nSecBiasedRegions= " << nSecBiasedRegions << G4endl;
263}
264
265//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
266
268 G4double previousStep)
269{
270 if(startTracking) {
271 startTracking = false;
272 G4int i = idxForcedCouple[coupleIdx];
273 if(i < 0) {
275 } else {
278 }
279 } else {
280 currentStepLimit -= previousStep;
281 }
282 if(currentStepLimit < 0.0) { currentStepLimit = 0.0; }
283 return currentStepLimit;
284}
285
286//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
287
290 std::vector<G4DynamicParticle*>& vd,
291 const G4Track& track,
292 G4VEmModel* currentModel,
293 G4ParticleChangeForLoss* pPartChange,
294 G4double& eloss,
295 G4int coupleIdx,
296 G4double tcut,
297 G4double safety)
298{
299 G4int index = idxSecBiasedCouple[coupleIdx];
300 G4double weight = 1.;
301 if(0 <= index) {
302 size_t n = vd.size();
303
304 // the check cannot be applied per secondary particle
305 // because weight correction is common, so the first
306 // secondary is checked
307 if((0 < n && vd[0]->GetKineticEnergy() < secBiasedEnegryLimit[index])
309
310 G4int nsplit = nBremSplitting[index];
311
312 // Range cut
313 if(0 == nsplit) {
314 if(safety > fSafetyMin) { ApplyRangeCut(vd, track, eloss, safety); }
315
316 // Russian Roulette
317 } else if(1 == nsplit) {
318 weight = ApplyRussianRoulette(vd, index);
319
320 // Splitting
321 } else {
323 weight = ApplyDirectionalSplitting(vd, track, currentModel, index, tcut);
324 } else {
325 G4double tmpEnergy = pPartChange->GetProposedKineticEnergy();
326 G4ThreeVector tmpMomDir = pPartChange->GetProposedMomentumDirection();
327
328 weight = ApplySplitting(vd, track, currentModel, index, tcut);
329
330 pPartChange->SetProposedKineticEnergy(tmpEnergy);
331 pPartChange->ProposeMomentumDirection(tmpMomDir);
332 }
333 }
334 }
335 }
336 return weight;
337}
338
339//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
340
343 std::vector<G4DynamicParticle*>& vd,
344 const G4Track& track,
345 G4VEmModel* currentModel,
346 G4ParticleChangeForGamma* pPartChange,
347 G4double& eloss,
348 G4int coupleIdx,
349 G4double tcut,
350 G4double safety)
351{
352 G4int index = idxSecBiasedCouple[coupleIdx];
353 G4double weight = 1.;
354 if(0 <= index) {
355 size_t n = vd.size();
356
357 // the check cannot be applied per secondary particle
358 // because weight correction is common, so the first
359 // secondary is checked
360 if((0 < n && vd[0]->GetKineticEnergy() < secBiasedEnegryLimit[index])
362
363 G4int nsplit = nBremSplitting[index];
364
365 // Range cut
366 if(0 == nsplit) {
367 if(safety > fSafetyMin) { ApplyRangeCut(vd, track, eloss, safety); }
368
369 // Russian Roulette
370 } else if(1 == nsplit) {
371 weight = ApplyRussianRoulette(vd, index);
372
373 // Splitting
374 } else {
376 weight = ApplyDirectionalSplitting(vd, track, currentModel,
377 index, tcut, pPartChange);
378 } else {
379 G4double tmpEnergy = pPartChange->GetProposedKineticEnergy();
380 G4ThreeVector tmpMomDir = pPartChange->GetProposedMomentumDirection();
381
382 weight = ApplySplitting(vd, track, currentModel, index, tcut);
383
384 pPartChange->SetProposedKineticEnergy(tmpEnergy);
385 pPartChange->ProposeMomentumDirection(tmpMomDir);
386 }
387 }
388 }
389 }
390 return weight;
391}
392
393//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
394
396G4EmBiasingManager::ApplySecondaryBiasing(std::vector<G4Track*>& track,
397 G4int coupleIdx)
398{
399 G4int index = idxSecBiasedCouple[coupleIdx];
400 G4double weight = 1.;
401 if(0 <= index) {
402 size_t n = track.size();
403
404 // the check cannot be applied per secondary particle
405 // because weight correction is common, so the first
406 // secondary is checked
407 if(0 < n && track[0]->GetKineticEnergy() < secBiasedEnegryLimit[index]) {
408
409 G4int nsplit = nBremSplitting[index];
410
411 // Russian Roulette only
412 if(1 == nsplit) {
413 weight = secBiasedWeight[index];
414 for(size_t k=0; k<n; ++k) {
415 if(G4UniformRand()*weight > 1.0) {
416 const G4Track* t = track[k];
417 delete t;
418 track[k] = 0;
419 }
420 }
421 }
422 }
423 }
424 return weight;
425}
426
427//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
428
429void
430G4EmBiasingManager::ApplyRangeCut(std::vector<G4DynamicParticle*>& vd,
431 const G4Track& track,
432 G4double& eloss, G4double safety)
433{
434 size_t n = vd.size();
435 if(!eIonisation) {
436 eIonisation =
438 }
439 if(eIonisation) {
440 for(size_t k=0; k<n; ++k) {
441 const G4DynamicParticle* dp = vd[k];
442 if(dp->GetDefinition() == theElectron) {
443 G4double e = dp->GetKineticEnergy();
444 if(eIonisation->GetRange(e, track.GetMaterialCutsCouple()) < safety) {
445 eloss += e;
446 delete dp;
447 vd[k] = 0;
448 }
449 }
450 }
451 }
452}
453
454//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
455
457 G4ThreeVector momdir) const
458{
460 G4double angle = momdir.angle(delta);
461 G4double dist = delta.cross(momdir).mag();
462 if (dist <= fDirectionalSplittingRadius && angle < halfpi) {
463 return true;
464 }
465 return false;
466}
467
468//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
469
471G4EmBiasingManager::ApplySplitting(std::vector<G4DynamicParticle*>& vd,
472 const G4Track& track,
473 G4VEmModel* currentModel,
474 G4int index,
475 G4double tcut)
476{
477 // method is applied only if 1 secondary created PostStep
478 // in the case of many secondaries there is a contradiction
479 G4double weight = 1.;
480 size_t n = vd.size();
481 G4double w = secBiasedWeight[index];
482
483 if(1 != n || 1.0 <= w) { return weight; }
484
485 G4double trackWeight = track.GetWeight();
486 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
487
488 G4int nsplit = nBremSplitting[index];
489
490 // double splitting is suppressed
491 if(1 < nsplit && trackWeight>w) {
492
493 weight = w;
494 if(nsplit > (G4int)tmpSecondaries.size()) {
495 tmpSecondaries.reserve(nsplit);
496 }
497 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
498 // start from 1, because already one secondary created
499 for(G4int k=1; k<nsplit; ++k) {
500 tmpSecondaries.clear();
501 currentModel->SampleSecondaries(&tmpSecondaries, couple, dynParticle,
502 tcut);
503 for (size_t kk=0; kk<tmpSecondaries.size(); ++kk) {
504 vd.push_back(tmpSecondaries[kk]);
505 }
506 }
507 }
508 return weight;
509}
510
511//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
512
515 std::vector<G4DynamicParticle*>& vd,
516 const G4Track& track,
517 G4VEmModel* currentModel,
518 G4int index,
519 G4double tcut,
520 G4ParticleChangeForGamma* partChange)
521{
522 // primary is gamma. do splitting/RR as appropriate
523 // method applied for any number of secondaries
524
525 G4double weight = 1.0;
526 G4double w = secBiasedWeight[index];
527
529 if(1.0 <= w) {
530 fDirectionalSplittingWeights.push_back(weight);
531 return weight;
532 }
533
534 G4double trackWeight = track.GetWeight();
535 G4int nsplit = nBremSplitting[index];
536
537 // double splitting is suppressed
538 if(1 < nsplit && trackWeight>w) {
539
540 weight = w;
541 const G4ThreeVector pos = track.GetPosition();
542
543 G4bool foundPrimaryParticle = false;
544 G4double primaryEnergy = 0.;
545 G4ThreeVector primaryMomdir(0.,0.,0.);
546 G4double primaryWeight = trackWeight;
547
548 tmpSecondaries = vd;
549 vd.clear();
550 vd.reserve(nsplit);
551 for (G4int k=0; k<nsplit; ++k) {
552 if (k>0) { // for k==0, SampleSecondaries has already been called
553 tmpSecondaries.clear();
554 // SampleSecondaries modifies primary info stored in partChange
555 currentModel->SampleSecondaries(&tmpSecondaries,
556 track.GetMaterialCutsCouple(),
557 track.GetDynamicParticle(), tcut);
558 }
559 for (size_t kk=0; kk<tmpSecondaries.size(); ++kk) {
560 if (tmpSecondaries[kk]->GetParticleDefinition() == theGamma) {
561 if (CheckDirection(pos, tmpSecondaries[kk]->GetMomentumDirection())){
562 vd.push_back(tmpSecondaries[kk]);
563 fDirectionalSplittingWeights.push_back(1.);
564 } else if (G4UniformRand() < w) {
565 vd.push_back(tmpSecondaries[kk]);
566 fDirectionalSplittingWeights.push_back(1./weight);
567 } else {
568 delete tmpSecondaries[kk];
569 tmpSecondaries[kk] = nullptr;
570 }
571 } else if (k==0) { // keep charged 2ry from first splitting
572 vd.push_back(tmpSecondaries[kk]);
573 fDirectionalSplittingWeights.push_back(1./weight);
574 } else {
575 delete tmpSecondaries[kk];
576 tmpSecondaries[kk] = nullptr;
577 }
578 }
579
580 // primary
581 G4double en = partChange->GetProposedKineticEnergy();
582 if (en>0.) { // don't add if kinetic energy = 0
583 G4ThreeVector momdir = partChange->GetProposedMomentumDirection();
584 if (CheckDirection(pos,momdir)) {
585 // keep only one primary; others are secondaries
586 if (!foundPrimaryParticle) {
587 primaryEnergy = en;
588 primaryMomdir = momdir;
589 foundPrimaryParticle = true;
590 primaryWeight = weight;
591 } else {
593 partChange->GetProposedMomentumDirection(),
594 partChange->GetProposedKineticEnergy());
595 vd.push_back(dp);
596 fDirectionalSplittingWeights.push_back(1.);
597 }
598 } else if (G4UniformRand()<w) { // not going to target. play RR.
599 if (!foundPrimaryParticle) {
600 foundPrimaryParticle = true;
601 primaryEnergy = en;
602 primaryMomdir = momdir;
603 primaryWeight = 1.;
604 } else {
606 partChange->GetProposedMomentumDirection(),
607 partChange->GetProposedKineticEnergy());
608 vd.push_back(dp);
609 fDirectionalSplittingWeights.push_back(1./weight);
610 }
611 }
612 }
613 } // end of loop over nsplit
614
615 partChange->ProposeWeight(primaryWeight);
616 partChange->SetProposedKineticEnergy(primaryEnergy);
617 partChange->ProposeMomentumDirection(primaryMomdir);
618 } else {
619 for (size_t i = 0; i < vd.size(); ++i) {
620 fDirectionalSplittingWeights.push_back(1.);
621 }
622 }
623
624 return weight;
625}
626
627//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
628
630{
631 // normally return 1. If a directionally split particle survives RR,
632 // return 1./(splitting factor)
633 if (fDirectionalSplittingWeights.size() >= (unsigned int)(i+1) ) {
635 fDirectionalSplittingWeights[i] = 1.; // ensure it's not used again
636 return w;
637 } else {
638 return 1.;
639 }
640}
641
644 std::vector<G4DynamicParticle*>& vd,
645 const G4Track& track,
646 G4VEmModel* currentModel,
647 G4int index,
648 G4double tcut)
649{
650 // primary is not a gamma. Do nothing with primary
651
652 G4double weight = 1.0;
653 G4double w = secBiasedWeight[index];
654
656 if(1.0 <= w) {
657 fDirectionalSplittingWeights.push_back(weight);
658 return weight;
659 }
660
661 G4double trackWeight = track.GetWeight();
662 G4int nsplit = nBremSplitting[index];
663
664 // double splitting is suppressed
665 if(1 < nsplit && trackWeight>w) {
666
667 weight = w;
668 const G4ThreeVector pos = track.GetPosition();
669
670 tmpSecondaries = vd;
671 vd.clear();
672 vd.reserve(nsplit);
673 for (G4int k=0; k<nsplit; ++k) {
674 if (k>0) {
675 tmpSecondaries.clear();
676 currentModel->SampleSecondaries(&tmpSecondaries,
677 track.GetMaterialCutsCouple(),
678 track.GetDynamicParticle(), tcut);
679 }
680 //for (auto sec : tmpSecondaries) {
681 for (size_t kk=0; kk < tmpSecondaries.size(); ++kk) {
682 if (CheckDirection(pos, tmpSecondaries[kk]->GetMomentumDirection())) {
683 vd.push_back(tmpSecondaries[kk]);
684 fDirectionalSplittingWeights.push_back(1.);
685 } else if (G4UniformRand()<w) {
686 vd.push_back(tmpSecondaries[kk]);
687 fDirectionalSplittingWeights.push_back(1./weight);
688 } else {
689 delete tmpSecondaries[kk];
690 tmpSecondaries[kk] = nullptr;
691 }
692 }
693 } // end of loop over nsplit
694 } else { // no splitting was done; still need weights
695 for (size_t i = 0; i < vd.size(); ++i) {
696 fDirectionalSplittingWeights.push_back(1.0);
697 }
698 }
699 return weight;
700}
static const G4double pos
static const G4double reg
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double halfpi
Definition: G4SIunits.hh:57
static constexpr double cm
Definition: G4SIunits.hh:99
static const G4double angle[DIMMOTT]
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector cross(const Hep3Vector &) const
double angle(const Hep3Vector &) const
double mag() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4double ApplySplitting(std::vector< G4DynamicParticle * > &vd, const G4Track &track, G4VEmModel *currentModel, G4int index, G4double tcut)
std::vector< G4DynamicParticle * > tmpSecondaries
std::vector< const G4Region * > forcedRegions
std::vector< G4double > secBiasedEnegryLimit
void SetDirectionalSplittingRadius(G4double r)
void SetDirectionalSplitting(G4bool v)
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
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)
G4double GetStepLimit(G4int coupleIdx, G4double previousStep)
std::vector< G4int > idxForcedCouple
static G4EmParameters * Instance()
G4double GetDirectionalSplittingRadius()
G4bool GetDirectionalSplitting() const
G4ThreeVector GetDirectionalSplittingTarget() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4LossTableManager * Instance()
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
G4ProductionCuts * GetProductionCuts() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4ThreeVector & GetProposedMomentumDirection() const
G4double GetProposedKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4ThreeVector & GetProposedMomentumDirection() const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
const G4String & GetName() const
G4double GetWeight() const
const G4ThreeVector & GetPosition() const
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
G4double GetRange(G4double kineticEnergy, const G4MaterialCutsCouple *)
void ProposeWeight(G4double finalWeight)
const char * name(G4int ptype)
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62