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 GB03BOptnSplitOrKillOnBoundary.cc 28 /// \brief Implementation of the GB03BOptnSplitOrKillOnBoundary class 29 30 #include "GB03BOptnSplitOrKillOnBoundary.hh" 31 32 #include "Randomize.hh" 33 34 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 35 36 GB03BOptnSplitOrKillOnBoundary::GB03BOptnSplitOrKillOnBoundary(G4String name) 37 : G4VBiasingOperation(name), fParticleChange(), fParticleChangeForNothing() 38 {} 39 40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 41 42 GB03BOptnSplitOrKillOnBoundary::~GB03BOptnSplitOrKillOnBoundary() {} 43 44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 45 46 G4double GB03BOptnSplitOrKillOnBoundary::DistanceToApplyOperation(const G4Track*, G4double, 47 G4ForceCondition* condition) 48 { 49 // -- return "infinite" distance for interaction, but asks for GenerateBiasingFinalState 50 // -- being called anyway at the end of the step, by returning the "Forced" condition 51 // -- flag. 52 *condition = Forced; 53 return DBL_MAX; 54 } 55 56 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 57 58 G4VParticleChange* GB03BOptnSplitOrKillOnBoundary::GenerateBiasingFinalState(const G4Track* track, 59 const G4Step* step) 60 { 61 // Check if step is limited by the geometry: as we attached the biasing operator 62 // to the absorber layer, this volume boundary is the one of the absorber. 63 // (check of current step # of track is inelegant, but is to fix a "feature" 64 // that a cloned track can wrongly be seen in the wrong volume, because of numerical 65 // precision issue. In this case it makes a tiny step, which should be disregarded). 66 if ((step->GetPostStepPoint()->GetStepStatus() == fGeomBoundary) 67 && (track->GetCurrentStepNumber() != 1)) 68 { 69 // -- Before deciding for killing or splitting, we make decision on applying 70 // -- the technique or not: 71 G4double trial = G4UniformRand(); // -- Note: G4UniformRand() is thread-safe 72 // -- engine for random numbers 73 if (trial <= fApplyProbability) { 74 // -- Get z component of track, to see if it moves forward or backward: 75 G4double pz = track->GetMomentum().z(); 76 77 if (pz > 0.0) { 78 // ------------------------------------------------- 79 // Here, we are moving "forward". We do "splitting": 80 // ------------------------------------------------- 81 82 // Get track weight: 83 G4double initialWeight = track->GetWeight(); 84 // Define the tracks weight: 85 G4double weightOfTrack = initialWeight / fSplittingFactor; 86 87 // The "particle change" is the object to be used to communicate to 88 // the tracking the update of the primary state and/or creation 89 // secondary tracks. 90 fParticleChange.Initialize(*track); 91 92 // ask currect track weight to be changed to new value: 93 fParticleChange.ProposeParentWeight(weightOfTrack); 94 95 // Now make clones of this track (this is the actual splitting): 96 // we will then have the primary and N-1 clones of it, hence the 97 // splitting by a factor N: 98 fParticleChange.SetNumberOfSecondaries(fSplittingFactor - 1); 99 for (G4int iSplit = 1; iSplit < fSplittingFactor; iSplit++) { 100 G4Track* clone = new G4Track(*track); 101 clone->SetWeight(weightOfTrack); 102 fParticleChange.AddSecondary(clone); 103 } 104 fParticleChange.SetSecondaryWeightByProcess(true); // -- tricky 105 // -- take it as is ;) [though not necessary here, put for safety] 106 107 // this new final state is returned to the tracking; 108 return &fParticleChange; 109 } 110 111 else { 112 // -------------------------------------------------------------- 113 // Here, we are moving backward. We do killing, playing a russian 114 // roulette, killing 1/fSplittingFactor of the tracks in average: 115 // -------------------------------------------------------------- 116 117 // Get track weight: 118 G4double initialWeight = track->GetWeight(); 119 120 // The "particle change" is the object to be used to communicate to 121 // the tracking the update of the primary state and/or creation 122 // secondary tracks. 123 fParticleChange.Initialize(*track); 124 125 // Shoot a random number (in ]0,1[ segment): 126 G4double random = G4UniformRand(); 127 128 // Decide to kill, keeping 1/fSplittingFactor of tracks: 129 G4double survivingProbability = 1.0 / fSplittingFactor; 130 if (random > survivingProbability) { 131 // We ask for the the track to be killed: 132 fParticleChange.ProposeTrackStatus(fStopAndKill); 133 } 134 else { 135 // In this case, the track survives. We change its weight 136 // to conserve weight among killed and survival tracks: 137 fParticleChange.ProposeParentWeight(initialWeight * fSplittingFactor); 138 } 139 140 // this new final state is returned to the tracking; 141 return &fParticleChange; 142 } 143 } // -- end of : if ( trial > probaForApplying ) 144 } // -- end of : if ( ( step->GetPostStepPoint()->GetStepStatus() == 145 // fGeomBoundary ) ... 146 147 // Here, the step was not limited by the geometry (but certainly by a physics 148 // process). We do nothing: ie we make no change to the current track. 149 fParticleChangeForNothing.Initialize(*track); 150 return &fParticleChangeForNothing; 151 } 152 153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 154