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 // G4ILawTruncatedExp 27 // -------------------------------------------------------------------- 28 29 #include "G4ILawTruncatedExp.hh" 30 #include "Randomize.hh" 31 #include "G4Track.hh" 32 33 34 G4ILawTruncatedExp::G4ILawTruncatedExp(const G4String& name) 35 : G4VBiasingInteractionLaw(name) 36 {} 37 38 G4ILawTruncatedExp::~G4ILawTruncatedExp() 39 {} 40 41 void G4ILawTruncatedExp::SetForceCrossSection(G4double crossSection) 42 { 43 if (crossSection < 0.0) 44 { 45 G4Exception("G4ILawTruncatedExp::SetForceCrossSection(..)", 46 "BIAS.GEN.09", JustWarning, 47 "Cross-section value passed is negative. It is set to zero !"); 48 fIsSingular = true; 49 crossSection = 0.0; 50 } 51 fIsSingular = false; 52 fCrossSectionDefined = true; 53 fCrossSection = crossSection; 54 } 55 56 G4double G4ILawTruncatedExp:: 57 ComputeEffectiveCrossSectionAt(G4double distance) const 58 { 59 if ( !fCrossSectionDefined ) 60 { 61 G4Exception("G4ILawTruncatedExp::ComputeEffectiveCrossSection(..)", 62 "BIAS.GEN.10", JustWarning, 63 "Cross-section value requested, but has not been defined yet. Assumes 0 !"); 64 // -- zero cross-section, returns the limit form of the effective cross-section: 65 return 1.0 / ( fMaximumDistance - distance ); 66 } 67 68 G4double denum = 1.0 - std::exp(-fCrossSection *(fMaximumDistance-distance)); 69 return fCrossSection/denum; 70 } 71 72 G4double G4ILawTruncatedExp:: 73 ComputeNonInteractionProbabilityAt(G4double distance) const 74 { 75 if (!fCrossSectionDefined) 76 { 77 G4Exception("G4ILawTruncatedExp::ComputeNonInteractionProbability(..)", 78 "BIAS.GEN.11", JustWarning, 79 "Non interaction probability value requested, but cross section has not been defined yet. Assumes it to be 0 !"); 80 // -- return limit case of null cross-section: 81 return 1.0 - distance / fMaximumDistance; 82 } 83 G4double num = 1.0 - std::exp( -fCrossSection*distance ); 84 G4double denum = 1.0 - std::exp( -fCrossSection*fMaximumDistance ); 85 return 1.0 - num/denum; 86 } 87 88 G4double G4ILawTruncatedExp::SampleInteractionLength() 89 { 90 if ( !fCrossSectionDefined ) 91 { 92 G4Exception("G4ILawTruncatedExp::Sample(..)", 93 "BIAS.GEN.12", JustWarning, 94 "Trying to sample while cross-section is not defined, assuming 0 !"); 95 fInteractionDistance = G4UniformRand() * fMaximumDistance; 96 return fInteractionDistance; 97 } 98 fInteractionDistance = -std::log(1.0-G4UniformRand()*(1.0-std::exp(-fCrossSection*fMaximumDistance)))/fCrossSection; 99 return fInteractionDistance; 100 } 101 102 G4double G4ILawTruncatedExp:: 103 UpdateInteractionLengthForStep(G4double truePathLength) 104 { 105 fInteractionDistance -= truePathLength; 106 fMaximumDistance -= truePathLength; 107 108 if ( fInteractionDistance < 0 ) 109 { 110 G4ExceptionDescription ed; 111 ed << " Negative number of interaction length for `" << GetName() 112 << "' " << fInteractionDistance << ", set it to zero !" << G4endl; 113 G4Exception("G4ILawTruncatedExp::UpdateInteractionLengthForStep(...)", 114 "BIAS.GEN.13", JustWarning, 115 "Trying to sample while cross-section is not defined, assuming 0 !"); 116 fInteractionDistance = 0.; 117 } 118 119 return fInteractionDistance; 120 } 121