00001 // 00002 // ******************************************************************** 00003 // * License and Disclaimer * 00004 // * * 00005 // * The Geant4 software is copyright of the Copyright Holders of * 00006 // * the Geant4 Collaboration. It is provided under the terms and * 00007 // * conditions of the Geant4 Software License, included in the file * 00008 // * LICENSE and available at http://cern.ch/geant4/license . These * 00009 // * include a list of copyright holders. * 00010 // * * 00011 // * Neither the authors of this software system, nor their employing * 00012 // * institutes,nor the agencies providing financial support for this * 00013 // * work make any representation or warranty, express or implied, * 00014 // * regarding this software system or assume any liability for its * 00015 // * use. Please see the license in the file LICENSE and URL above * 00016 // * for the full disclaimer and the limitation of liability. * 00017 // * * 00018 // * This code implementation is the result of the scientific and * 00019 // * technical work of the GEANT4 collaboration. * 00020 // * By using, copying, modifying or distributing the software (or * 00021 // * any work based on the software) you agree to acknowledge its * 00022 // * use in resulting scientific publications, and indicate your * 00023 // * acceptance of all terms of the Geant4 Software license. * 00024 // ******************************************************************** 00025 // 00026 // neutron_hp -- source file 00027 // J.P. Wellisch, Nov-1996 00028 // A prototype of the low energy neutron transport model. 00029 // 00030 #include "G4NeutronHPIsotropic.hh" 00031 #include "G4PhysicalConstants.hh" 00032 #include "G4SystemOfUnits.hh" 00033 #include "Randomize.hh" 00034 #include "G4Gamma.hh" 00035 #include "G4Electron.hh" 00036 #include "G4Positron.hh" 00037 #include "G4Neutron.hh" 00038 #include "G4Proton.hh" 00039 #include "G4Deuteron.hh" 00040 #include "G4Triton.hh" 00041 #include "G4He3.hh" 00042 #include "G4Alpha.hh" 00043 #include "G4ParticleTable.hh" 00044 00045 00046 void G4NeutronHPIsotropic::Init(std::ifstream & ) 00047 { 00048 } 00049 00050 G4ReactionProduct * G4NeutronHPIsotropic::Sample(G4double anEnergy, G4double massCode, G4double ) 00051 { 00052 G4ReactionProduct * result = new G4ReactionProduct; 00053 G4int Z = static_cast<G4int>(massCode/1000); 00054 G4int A = static_cast<G4int>(massCode-1000*Z); 00055 00056 if(massCode==0) 00057 { 00058 result->SetDefinition(G4Gamma::Gamma()); 00059 } 00060 else if(A==0) 00061 { 00062 result->SetDefinition(G4Electron::Electron()); 00063 if(Z==1) result->SetDefinition(G4Positron::Positron()); 00064 } 00065 else if(A==1) 00066 { 00067 result->SetDefinition(G4Neutron::Neutron()); 00068 if(Z==1) result->SetDefinition(G4Proton::Proton()); 00069 } 00070 else if(A==2) 00071 { 00072 result->SetDefinition(G4Deuteron::Deuteron()); 00073 } 00074 else if(A==3) 00075 { 00076 result->SetDefinition(G4Triton::Triton()); 00077 if(Z==2) result->SetDefinition(G4He3::He3()); 00078 } 00079 else if(A==4) 00080 { 00081 result->SetDefinition(G4Alpha::Alpha()); 00082 //110607 TK modified following parts for migration to G4NDL3.15 (ENDF VII.r0) 00083 //if(Z!=2) throw G4HadronicException(__FILE__, __LINE__, "Unknown ion case 1"); 00084 if(Z!=2) 00085 { 00086 result->SetDefinition( G4ParticleTable::GetParticleTable()->GetIon ( Z , A , 0.0 ) ); 00087 } 00088 } 00089 else 00090 { 00091 //110607 TK modified following parts for migration to G4NDL3.15 (ENDF VII.r0) 00092 result->SetDefinition( G4ParticleTable::GetParticleTable()->GetIon ( Z , A , 0.0 ) ); 00093 //throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPIsotropic: Unknown ion case 2"); 00094 } 00095 00096 G4double cosTh = G4UniformRand(); 00097 G4double phi = twopi*G4UniformRand(); 00098 G4double theta = std::acos(cosTh); 00099 G4double sinth = std::sin(theta); 00100 00101 // we need the the Q value of the reaction 00102 result->SetKineticEnergy(std::max(0.001*MeV, anEnergy+GetQValue())); 00103 G4double mtot = result->GetTotalMomentum(); 00104 G4ThreeVector tempVector(mtot*sinth*std::cos(phi), mtot*sinth*std::sin(phi), mtot*std::cos(theta) ); 00105 result->SetMomentum(tempVector); 00106 00107 return result; 00108 }