Geant4 Cross Reference |
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 /// \file GB05BOptnSplitAndKillByCrossSection.cc 28 /// \brief Implementation of the GB05BOptnSplitAndKillByCrossSection class 29 30 #include "GB05BOptnSplitAndKillByCrossSection.hh" 31 32 #include "Randomize.hh" 33 34 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 35 36 GB05BOptnSplitAndKillByCrossSection::GB05BOptnSplitAndKillByCrossSection(G4String name) 37 : G4VBiasingOperation(name), fParticleChange(), fInteractionLength(-1.0) 38 {} 39 40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 41 42 GB05BOptnSplitAndKillByCrossSection::~GB05BOptnSplitAndKillByCrossSection() {} 43 44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 45 46 G4double GB05BOptnSplitAndKillByCrossSection::DistanceToApplyOperation(const G4Track*, G4double, 47 G4ForceCondition* condition) 48 { 49 *condition = NotForced; 50 51 // -- Sample the exponential law using the total interaction length of processes 52 // -- to couterbalance for: 53 G4double proposedStepLength = -std::log(G4UniformRand()) * fInteractionLength; 54 55 return proposedStepLength; 56 } 57 58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 59 60 G4VParticleChange* 61 GB05BOptnSplitAndKillByCrossSection::GenerateBiasingFinalState(const G4Track* track, const G4Step*) 62 { 63 // -- This method is called if we have limited the step. 64 // -- We hence make the splitting or killing. 65 66 // Get track weight: 67 G4double initialWeight = track->GetWeight(); 68 69 // The "particle change" is the object to be used to communicate to 70 // the tracking the update of the primary state and/or creation 71 // secondary tracks. 72 fParticleChange.Initialize(*track); 73 74 // -- Splitting and killing factors. 75 // -- They are taken the same, but the killing factor can be make bigger. 76 G4double splittingFactor = 2.0; 77 G4double killingFactor = 2.0; 78 79 if (track->GetMomentumDirection().z() > 0) { 80 // -- We split if the track is moving forward: 81 82 // Define the tracks weight: 83 G4double weightOfTrack = initialWeight / splittingFactor; 84 85 // Ask currect track weight to be changed to new value: 86 fParticleChange.ProposeParentWeight(weightOfTrack); 87 // Now make clones of this track (this is the actual splitting): 88 // we will then have the primary and clone of it, hence the 89 // splitting by a factor 2: 90 G4Track* clone = new G4Track(*track); 91 clone->SetWeight(weightOfTrack); 92 fParticleChange.AddSecondary(clone); 93 // -- Below's call added for safety & illustration : inform particle change to not 94 // -- modify the clone (ie : daughter) weight to male it that of the 95 // -- primary. Here call not mandatory and both tracks have same weights. 96 fParticleChange.SetSecondaryWeightByProcess(true); 97 } 98 else { 99 // -- We apply Russian roulette if the track is moving backward: 100 101 // Shoot a random number (in ]0,1[ segment): 102 G4double random = G4UniformRand(); 103 G4double killingProbability = 1.0 - 1.0 / killingFactor; 104 if (random < killingProbability) { 105 // We ask for the the track to be killed: 106 fParticleChange.ProposeTrackStatus(fStopAndKill); 107 } 108 else { 109 // In this case, the track survives. We change its weight 110 // to conserve weight among killed and survival tracks: 111 fParticleChange.ProposeParentWeight(initialWeight * killingFactor); 112 } 113 } 114 115 return &fParticleChange; 116 } 117 118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 119