68 : fDirectionalSplittingTarget(0.0,0.0,0.0)
70 fSafetyMin = 1.e-6*mm;
91 if(0 < nForcedRegions) { idxForcedCouple.resize(numOfCouples, -1); }
92 if(0 < nSecBiasedRegions) { idxSecBiasedCouple.resize(numOfCouples, -1); }
95 for (
G4int j=0; j<numOfCouples; ++j) {
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;
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;
123 if (fDirectionalSplitting) {
128 if (nForcedRegions > 0 && 0 < verbose) {
129 G4cout <<
" Forced Interaction is activated for "
132 <<
" inside G4Regions: " <<
G4endl;
133 for (
G4int i=0; i<nForcedRegions; ++i) {
134 const G4Region* r = forcedRegions[i];
138 if (nSecBiasedRegions > 0 && 0 < verbose) {
139 G4cout <<
" Secondary biasing is activated for "
142 <<
" inside G4Regions: " <<
G4endl;
143 for (
G4int i=0; i<nSecBiasedRegions; ++i) {
144 const G4Region* r = secBiasedRegions[i];
147 <<
" BiasingWeight= " << secBiasedWeight[i] <<
G4endl;
150 if (fDirectionalSplitting) {
151 G4cout <<
" Directional splitting activated, with target position: "
152 << fDirectionalSplittingTarget/cm
154 << fDirectionalSplittingRadius/cm
167 if(name ==
"" || name ==
"world" || name ==
"World") {
168 name =
"DefaultRegionForTheWorld";
172 G4cout <<
"### G4EmBiasingManager::ForcedInteraction WARNING: "
174 << rname <<
"> is unknown" <<
G4endl;
179 if (0 < nForcedRegions) {
180 for (
G4int i=0; i<nForcedRegions; ++i) {
181 if (reg == forcedRegions[i]) {
182 lengthForRegion[i] = val;
188 G4cout <<
"### G4EmBiasingManager::ForcedInteraction WARNING: "
189 << val <<
" < 0.0, so no activation for the G4Region <"
190 << rname <<
">" <<
G4endl;
195 forcedRegions.push_back(reg);
196 lengthForRegion.push_back(val);
212 if(name ==
"" || name ==
"world" || name ==
"World") {
213 name =
"DefaultRegionForTheWorld";
217 G4cout <<
"### G4EmBiasingManager::ActivateBremsstrahlungSplitting "
218 <<
"WARNING: G4Region <"
219 << rname <<
"> is unknown" <<
G4endl;
233 }
else if(0.0 < factor) {
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;
256 secBiasedRegions.push_back(reg);
257 secBiasedWeight.push_back(w);
258 nBremSplitting.push_back(nsplit);
259 secBiasedEnegryLimit.push_back(energyLimit);
270 startTracking =
false;
271 G4int i = idxForcedCouple[coupleIdx];
275 currentStepLimit = lengthForRegion[i];
276 if(currentStepLimit > 0.0) { currentStepLimit *=
G4UniformRand(); }
279 currentStepLimit -= previousStep;
281 if(currentStepLimit < 0.0) { currentStepLimit = 0.0; }
282 return currentStepLimit;
289 std::vector<G4DynamicParticle*>& vd,
298 G4int index = idxSecBiasedCouple[coupleIdx];
301 std::size_t n = vd.size();
306 if((0 < n && vd[0]->GetKineticEnergy() < secBiasedEnegryLimit[index])
307 || fDirectionalSplitting) {
309 G4int nsplit = nBremSplitting[index];
313 if(safety > fSafetyMin) { ApplyRangeCut(vd, track, eloss, safety); }
316 }
else if(1 == nsplit) {
317 weight = ApplyRussianRoulette(vd, index);
321 if (fDirectionalSplitting) {
322 weight = ApplyDirectionalSplitting(vd, track, currentModel, index, tcut);
327 weight = ApplySplitting(vd, track, currentModel, index, tcut);
342 std::vector<G4DynamicParticle*>& vd,
351 G4int index = idxSecBiasedCouple[coupleIdx];
354 std::size_t n = vd.size();
359 if((0 < n && vd[0]->GetKineticEnergy() < secBiasedEnegryLimit[index])
360 || fDirectionalSplitting) {
362 G4int nsplit = nBremSplitting[index];
366 if(safety > fSafetyMin) { ApplyRangeCut(vd, track, eloss, safety); }
369 }
else if(1 == nsplit) {
370 weight = ApplyRussianRoulette(vd, index);
374 if (fDirectionalSplitting) {
375 weight = ApplyDirectionalSplitting(vd, track, currentModel,
376 index, tcut, pPartChange);
381 weight = ApplySplitting(vd, track, currentModel, index, tcut);
398 G4int index = idxSecBiasedCouple[coupleIdx];
401 std::size_t n = track.size();
406 if(0 < n && track[0]->GetKineticEnergy() < secBiasedEnegryLimit[index]) {
408 G4int nsplit = nBremSplitting[index];
412 weight = secBiasedWeight[index];
413 for(std::size_t k=0; k<n; ++k) {
429G4EmBiasingManager::ApplyRangeCut(std::vector<G4DynamicParticle*>& vd,
433 std::size_t n = vd.size();
439 for(std::size_t k=0; k<
n; ++k) {
461 if (dist <= fDirectionalSplittingRadius && angle < halfpi) {
470G4EmBiasingManager::ApplySplitting(std::vector<G4DynamicParticle*>& vd,
479 std::size_t n = vd.size();
480 G4double w = secBiasedWeight[index];
482 if(1 != n || 1.0 <= w) {
return weight; }
487 G4int nsplit = nBremSplitting[index];
490 if(1 < nsplit && trackWeight>w) {
493 if(nsplit > (
G4int)tmpSecondaries.size()) {
494 tmpSecondaries.reserve(nsplit);
498 for(
G4int k=1; k<nsplit; ++k) {
499 tmpSecondaries.clear();
502 for (std::size_t kk=0; kk<tmpSecondaries.size(); ++kk) {
503 vd.push_back(tmpSecondaries[kk]);
513G4EmBiasingManager::ApplyDirectionalSplitting(
514 std::vector<G4DynamicParticle*>& vd,
525 G4double w = secBiasedWeight[index];
527 fDirectionalSplittingWeights.clear();
529 fDirectionalSplittingWeights.push_back(weight);
534 G4int nsplit = nBremSplitting[index];
537 if(1 < nsplit && trackWeight>w) {
542 G4bool foundPrimaryParticle =
false;
545 G4double primaryWeight = trackWeight;
550 for (
G4int k=0; k<nsplit; ++k) {
552 tmpSecondaries.clear();
558 for (std::size_t kk=0; kk<tmpSecondaries.size(); ++kk) {
559 if (tmpSecondaries[kk]->GetParticleDefinition() == theGamma) {
560 if (
CheckDirection(pos, tmpSecondaries[kk]->GetMomentumDirection())){
561 vd.push_back(tmpSecondaries[kk]);
562 fDirectionalSplittingWeights.push_back(1.);
564 vd.push_back(tmpSecondaries[kk]);
565 fDirectionalSplittingWeights.push_back(1./weight);
567 delete tmpSecondaries[kk];
568 tmpSecondaries[kk] =
nullptr;
571 vd.push_back(tmpSecondaries[kk]);
572 fDirectionalSplittingWeights.push_back(1./weight);
574 delete tmpSecondaries[kk];
575 tmpSecondaries[kk] =
nullptr;
585 if (!foundPrimaryParticle) {
587 primaryMomdir = momdir;
588 foundPrimaryParticle =
true;
589 primaryWeight = weight;
595 fDirectionalSplittingWeights.push_back(1.);
598 if (!foundPrimaryParticle) {
599 foundPrimaryParticle =
true;
601 primaryMomdir = momdir;
608 fDirectionalSplittingWeights.push_back(1./weight);
618 for (std::size_t i = 0; i < vd.size(); ++i) {
619 fDirectionalSplittingWeights.push_back(1.);
632 if (fDirectionalSplittingWeights.size() >= (
unsigned int)(i+1) ) {
633 G4double w = fDirectionalSplittingWeights[i];
634 fDirectionalSplittingWeights[i] = 1.;
642G4EmBiasingManager::ApplyDirectionalSplitting(
643 std::vector<G4DynamicParticle*>& vd,
652 G4double w = secBiasedWeight[index];
654 fDirectionalSplittingWeights.clear();
656 fDirectionalSplittingWeights.push_back(weight);
661 G4int nsplit = nBremSplitting[index];
664 if(1 < nsplit && trackWeight>w) {
672 for (
G4int k=0; k<nsplit; ++k) {
674 tmpSecondaries.clear();
680 for (std::size_t kk=0; kk < tmpSecondaries.size(); ++kk) {
681 if (
CheckDirection(pos, tmpSecondaries[kk]->GetMomentumDirection())) {
682 vd.push_back(tmpSecondaries[kk]);
683 fDirectionalSplittingWeights.push_back(1.);
685 vd.push_back(tmpSecondaries[kk]);
686 fDirectionalSplittingWeights.push_back(1./weight);
688 delete tmpSecondaries[kk];
689 tmpSecondaries[kk] =
nullptr;
694 for (std::size_t i = 0; i < vd.size(); ++i) {
695 fDirectionalSplittingWeights.push_back(1.0);
G4GLOB_DLL std::ostream G4cout
Hep3Vector cross(const Hep3Vector &) const
double angle(const Hep3Vector &) const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
void SetDirectionalSplittingRadius(G4double r)
void SetDirectionalSplitting(G4bool v)
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 ®ion, G4double factor, G4double energyLimit)
void Initialise(const G4ParticleDefinition &part, const G4String &procName, G4int verbose)
G4double GetWeight(G4int i)
G4bool CheckDirection(G4ThreeVector pos, G4ThreeVector momdir) const
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="")
void SetDirectionalSplittingTarget(G4ThreeVector v)
G4double GetStepLimit(G4int coupleIdx, G4double previousStep)
static G4EmParameters * Instance()
G4double GetDirectionalSplittingRadius()
G4bool GetDirectionalSplitting() const
G4ThreeVector GetDirectionalSplittingTarget() const
static G4LossTableManager * Instance()
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
G4ProductionCuts * GetProductionCuts() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4double GetProposedKineticEnergy() const
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
const G4ThreeVector & GetProposedMomentumDirection() const
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4double GetProposedKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4ThreeVector & GetProposedMomentumDirection() const
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)