// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // neutron_hp -- source file // J.P. Wellisch, Nov-1996 // A prototype of the low energy neutron transport model. // // 09-May-06 fix in Sample by T. Koi // 080318 Fix Compilation warnings - gcc-4.3.0 by T. Koi // (This fix has a real effect to the code.) // 080409 Fix div0 error with G4FPE by T. Koi // #include "G4NeutronHPContAngularPar.hh" #include "G4NeutronHPLegendreStore.hh" #include "G4Gamma.hh" #include "G4Electron.hh" #include "G4Positron.hh" #include "G4Neutron.hh" #include "G4Proton.hh" #include "G4Deuteron.hh" #include "G4Triton.hh" #include "G4He3.hh" #include "G4Alpha.hh" #include "G4NeutronHPVector.hh" #include "G4NucleiPropertiesTable.hh" #include "G4NeutronHPKallbachMannSyst.hh" #include "G4ParticleTable.hh" void G4NeutronHPContAngularPar::Init(std::ifstream & aDataFile) { aDataFile >> theEnergy >> nEnergies >> nDiscreteEnergies >> nAngularParameters; theEnergy *= eV; theAngular = new G4NeutronHPList [nEnergies]; for(G4int i=0; i> sEnergy; sEnergy*=eV; theAngular[i].SetLabel(sEnergy); theAngular[i].Init(aDataFile, nAngularParameters, 1.); } } G4ReactionProduct * G4NeutronHPContAngularPar::Sample(G4double anEnergy, G4double massCode, G4double /*targetMass*/, G4int angularRep, G4int /*interpolE*/) { G4ReactionProduct * result = new G4ReactionProduct; G4int Z = static_cast(massCode/1000); G4int A = static_cast(massCode-1000*Z); if(massCode==0) { result->SetDefinition(G4Gamma::Gamma()); } else if(A==0) { result->SetDefinition(G4Electron::Electron()); if(Z==1) result->SetDefinition(G4Positron::Positron()); } else if(A==1) { result->SetDefinition(G4Neutron::Neutron()); if(Z==1) result->SetDefinition(G4Proton::Proton()); } else if(A==2) { result->SetDefinition(G4Deuteron::Deuteron()); } else if(A==3) { result->SetDefinition(G4Triton::Triton()); if(Z==2) result->SetDefinition(G4He3::He3()); } else if(A==4) { result->SetDefinition(G4Alpha::Alpha()); if(Z!=2) throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPContAngularPar: Unknown ion case 1"); } else { result->SetDefinition(G4ParticleTable::GetParticleTable()->FindIon(Z,A,0,Z)); } G4int i(0); G4int it(0); G4double fsEnergy(0); G4double cosTh(0); if(angularRep==1) { G4double random = G4UniformRand(); G4double * running = new G4double[nEnergies]; running[0]=0; G4double weighted = 0; for(i=1; iGetPDGMass(); G4double productEnergy = fsEnergy; G4double productMass = result->GetMass(); G4int targetZ = G4int(theTargetCode/1000); G4int targetA = G4int(theTargetCode-1000*targetZ); // To correspond to natural composition (-nat-) data files. if ( targetA == 0 ) targetA = int ( theTarget->GetMass()/amu_c2 + 0.5 ); G4double targetMass = theTarget->GetMass(); G4int residualA = targetA+1-A; G4int residualZ = targetZ-Z; G4double residualMass = residualZ*G4Proton::Proton()->GetPDGMass(); residualMass +=(residualA-residualZ)*G4Neutron::Neutron()->GetPDGMass(); residualMass -= G4NucleiPropertiesTable::GetBindingEnergy(residualZ, residualA); G4NeutronHPKallbachMannSyst theKallbach(compoundFraction, incidentEnergy, incidentMass, productEnergy, productMass, residualMass, residualA, residualZ, targetMass, targetA, targetZ); cosTh = theKallbach.Sample(anEnergy); } else if(angularRep>10&&angularRep<16) { G4double random = G4UniformRand(); G4double * running = new G4double[nEnergies]; running[0]=0; G4double weighted = 0; for(i=1; iSetKineticEnergy(fsEnergy); G4double phi = twopi*G4UniformRand(); G4double theta = std::acos(cosTh); G4double sinth = std::sin(theta); G4double mtot = result->GetTotalMomentum(); G4ThreeVector tempVector(mtot*sinth*std::cos(phi), mtot*sinth*std::sin(phi), mtot*std::cos(theta) ); result->SetMomentum(tempVector); // return the result. return result; }