![]() |
|
|||
File indexing completed on 2025-02-23 09:20:45
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 /// \file GB03BOptnSplitOrKillOnBoundary.cc 0028 /// \brief Implementation of the GB03BOptnSplitOrKillOnBoundary class 0029 0030 #include "GB03BOptnSplitOrKillOnBoundary.hh" 0031 0032 #include "Randomize.hh" 0033 0034 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0035 0036 GB03BOptnSplitOrKillOnBoundary::GB03BOptnSplitOrKillOnBoundary(G4String name) 0037 : G4VBiasingOperation(name), fParticleChange(), fParticleChangeForNothing() 0038 {} 0039 0040 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0041 0042 GB03BOptnSplitOrKillOnBoundary::~GB03BOptnSplitOrKillOnBoundary() {} 0043 0044 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0045 0046 G4double GB03BOptnSplitOrKillOnBoundary::DistanceToApplyOperation(const G4Track*, G4double, 0047 G4ForceCondition* condition) 0048 { 0049 // -- return "infinite" distance for interaction, but asks for GenerateBiasingFinalState 0050 // -- being called anyway at the end of the step, by returning the "Forced" condition 0051 // -- flag. 0052 *condition = Forced; 0053 return DBL_MAX; 0054 } 0055 0056 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 0057 0058 G4VParticleChange* GB03BOptnSplitOrKillOnBoundary::GenerateBiasingFinalState(const G4Track* track, 0059 const G4Step* step) 0060 { 0061 // Check if step is limited by the geometry: as we attached the biasing operator 0062 // to the absorber layer, this volume boundary is the one of the absorber. 0063 // (check of current step # of track is inelegant, but is to fix a "feature" 0064 // that a cloned track can wrongly be seen in the wrong volume, because of numerical 0065 // precision issue. In this case it makes a tiny step, which should be disregarded). 0066 if ((step->GetPostStepPoint()->GetStepStatus() == fGeomBoundary) 0067 && (track->GetCurrentStepNumber() != 1)) 0068 { 0069 // -- Before deciding for killing or splitting, we make decision on applying 0070 // -- the technique or not: 0071 G4double trial = G4UniformRand(); // -- Note: G4UniformRand() is thread-safe 0072 // -- engine for random numbers 0073 if (trial <= fApplyProbability) { 0074 // -- Get z component of track, to see if it moves forward or backward: 0075 G4double pz = track->GetMomentum().z(); 0076 0077 if (pz > 0.0) { 0078 // ------------------------------------------------- 0079 // Here, we are moving "forward". We do "splitting": 0080 // ------------------------------------------------- 0081 0082 // Get track weight: 0083 G4double initialWeight = track->GetWeight(); 0084 // Define the tracks weight: 0085 G4double weightOfTrack = initialWeight / fSplittingFactor; 0086 0087 // The "particle change" is the object to be used to communicate to 0088 // the tracking the update of the primary state and/or creation 0089 // secondary tracks. 0090 fParticleChange.Initialize(*track); 0091 0092 // ask currect track weight to be changed to new value: 0093 fParticleChange.ProposeParentWeight(weightOfTrack); 0094 0095 // Now make clones of this track (this is the actual splitting): 0096 // we will then have the primary and N-1 clones of it, hence the 0097 // splitting by a factor N: 0098 fParticleChange.SetNumberOfSecondaries(fSplittingFactor - 1); 0099 for (G4int iSplit = 1; iSplit < fSplittingFactor; iSplit++) { 0100 G4Track* clone = new G4Track(*track); 0101 clone->SetWeight(weightOfTrack); 0102 fParticleChange.AddSecondary(clone); 0103 } 0104 fParticleChange.SetSecondaryWeightByProcess(true); // -- tricky 0105 // -- take it as is ;) [though not necessary here, put for safety] 0106 0107 // this new final state is returned to the tracking; 0108 return &fParticleChange; 0109 } 0110 0111 else { 0112 // -------------------------------------------------------------- 0113 // Here, we are moving backward. We do killing, playing a russian 0114 // roulette, killing 1/fSplittingFactor of the tracks in average: 0115 // -------------------------------------------------------------- 0116 0117 // Get track weight: 0118 G4double initialWeight = track->GetWeight(); 0119 0120 // The "particle change" is the object to be used to communicate to 0121 // the tracking the update of the primary state and/or creation 0122 // secondary tracks. 0123 fParticleChange.Initialize(*track); 0124 0125 // Shoot a random number (in ]0,1[ segment): 0126 G4double random = G4UniformRand(); 0127 0128 // Decide to kill, keeping 1/fSplittingFactor of tracks: 0129 G4double survivingProbability = 1.0 / fSplittingFactor; 0130 if (random > survivingProbability) { 0131 // We ask for the the track to be killed: 0132 fParticleChange.ProposeTrackStatus(fStopAndKill); 0133 } 0134 else { 0135 // In this case, the track survives. We change its weight 0136 // to conserve weight among killed and survival tracks: 0137 fParticleChange.ProposeParentWeight(initialWeight * fSplittingFactor); 0138 } 0139 0140 // this new final state is returned to the tracking; 0141 return &fParticleChange; 0142 } 0143 } // -- end of : if ( trial > probaForApplying ) 0144 } // -- end of : if ( ( step->GetPostStepPoint()->GetStepStatus() == 0145 // fGeomBoundary ) ... 0146 0147 // Here, the step was not limited by the geometry (but certainly by a physics 0148 // process). We do nothing: ie we make no change to the current track. 0149 fParticleChangeForNothing.Initialize(*track); 0150 return &fParticleChangeForNothing; 0151 } 0152 0153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |