[819] | 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 | // 05-11-21 NeutronHP or Low Energy Parameterization Models |
---|
| 28 | // Implemented by T. Koi (SLAC/SCCS) |
---|
| 29 | // If NeutronHP data do not available for an element, then Low Energy |
---|
| 30 | // Parameterization models handle the interactions of the element. |
---|
[962] | 31 | // 081024 G4NucleiPropertiesTable:: to G4NucleiProperties:: |
---|
[819] | 32 | // |
---|
| 33 | |
---|
| 34 | #include "G4NeutronHPorLElasticData.hh" |
---|
| 35 | #include "G4Neutron.hh" |
---|
| 36 | #include "G4ElementTable.hh" |
---|
| 37 | #include "G4NeutronHPData.hh" |
---|
| 38 | |
---|
| 39 | #include "G4PhysicsVector.hh" |
---|
| 40 | |
---|
| 41 | |
---|
| 42 | |
---|
| 43 | G4NeutronHPorLElasticData::G4NeutronHPorLElasticData( G4NeutronHPChannel* pChannel , std::set< G4String >* pSet ) |
---|
| 44 | { |
---|
| 45 | theElasticChannel = pChannel; |
---|
| 46 | unavailable_elements = pSet; |
---|
| 47 | } |
---|
| 48 | |
---|
| 49 | |
---|
| 50 | |
---|
| 51 | G4bool G4NeutronHPorLElasticData::IsApplicable(const G4DynamicParticle*aP, const G4Element* anElement) |
---|
| 52 | { |
---|
| 53 | G4bool result = true; |
---|
| 54 | G4double eKin = aP->GetKineticEnergy(); |
---|
| 55 | if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false; |
---|
| 56 | if ( unavailable_elements->find( anElement->GetName() ) != unavailable_elements->end() ) result = false; |
---|
| 57 | return result; |
---|
| 58 | } |
---|
| 59 | |
---|
| 60 | G4NeutronHPorLElasticData::G4NeutronHPorLElasticData() |
---|
| 61 | { |
---|
| 62 | // BuildPhysicsTable(*G4Neutron::Neutron()); |
---|
| 63 | } |
---|
| 64 | |
---|
| 65 | |
---|
| 66 | |
---|
| 67 | G4NeutronHPorLElasticData::~G4NeutronHPorLElasticData() |
---|
| 68 | { |
---|
| 69 | // delete theCrossSections; |
---|
| 70 | } |
---|
| 71 | |
---|
| 72 | |
---|
| 73 | |
---|
| 74 | void G4NeutronHPorLElasticData::BuildPhysicsTable( const G4ParticleDefinition& aP ) |
---|
| 75 | { |
---|
| 76 | if( &aP!=G4Neutron::Neutron() ) |
---|
| 77 | throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!"); |
---|
| 78 | } |
---|
| 79 | |
---|
| 80 | |
---|
| 81 | |
---|
| 82 | void G4NeutronHPorLElasticData::DumpPhysicsTable(const G4ParticleDefinition& aP) |
---|
| 83 | { |
---|
| 84 | if(&aP!=G4Neutron::Neutron()) |
---|
| 85 | throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!"); |
---|
| 86 | // G4cout << "G4NeutronHPorLElasticData::DumpPhysicsTable still to be implemented"<<G4endl; |
---|
| 87 | } |
---|
| 88 | |
---|
| 89 | |
---|
| 90 | |
---|
| 91 | #include "G4Nucleus.hh" |
---|
[962] | 92 | #include "G4NucleiProperties.hh" |
---|
[819] | 93 | #include "G4Neutron.hh" |
---|
| 94 | #include "G4Electron.hh" |
---|
| 95 | |
---|
| 96 | G4double G4NeutronHPorLElasticData:: |
---|
| 97 | GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT) |
---|
| 98 | { |
---|
| 99 | |
---|
| 100 | //G4cout << "Choice G4NeutronHPorLElasticData for element " << anE->GetName() << G4endl; |
---|
| 101 | G4double result = 0; |
---|
| 102 | // G4bool outOfRange; |
---|
| 103 | G4int index = anE->GetIndex(); |
---|
| 104 | |
---|
| 105 | // prepare neutron |
---|
| 106 | G4double eKinetic = aP->GetKineticEnergy(); |
---|
| 107 | G4ReactionProduct theNeutron( aP->GetDefinition() ); |
---|
| 108 | theNeutron.SetMomentum( aP->GetMomentum() ); |
---|
| 109 | theNeutron.SetKineticEnergy( eKinetic ); |
---|
| 110 | |
---|
| 111 | // prepare thermal nucleus |
---|
| 112 | G4Nucleus aNuc; |
---|
| 113 | G4double eps = 0.0001; |
---|
| 114 | G4double theA = anE->GetN(); |
---|
| 115 | G4double theZ = anE->GetZ(); |
---|
| 116 | G4double eleMass; |
---|
[962] | 117 | eleMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps)) |
---|
[819] | 118 | ) / G4Neutron::Neutron()->GetPDGMass(); |
---|
| 119 | |
---|
| 120 | G4ReactionProduct boosted; |
---|
| 121 | G4double aXsection; |
---|
| 122 | |
---|
| 123 | // MC integration loop |
---|
| 124 | G4int counter = 0; |
---|
| 125 | G4double buffer = 0; |
---|
| 126 | G4int size = G4int(std::max(10., aT/60*kelvin)); |
---|
| 127 | G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum(); |
---|
| 128 | G4double neutronVMag = neutronVelocity.mag(); |
---|
| 129 | while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer) |
---|
| 130 | { |
---|
| 131 | if(counter) buffer = result/counter; |
---|
| 132 | while (counter<size) |
---|
| 133 | { |
---|
| 134 | counter ++; |
---|
| 135 | G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT); |
---|
| 136 | boosted.Lorentz(theNeutron, aThermalNuc); |
---|
| 137 | G4double theEkin = boosted.GetKineticEnergy(); |
---|
| 138 | //aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange); |
---|
| 139 | aXsection = theElasticChannel[index].GetXsec( theEkin ); |
---|
| 140 | // velocity correction. |
---|
| 141 | G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum(); |
---|
| 142 | aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag; |
---|
| 143 | result += aXsection; |
---|
| 144 | } |
---|
| 145 | size += size; |
---|
| 146 | } |
---|
| 147 | result /= counter; |
---|
| 148 | return result; |
---|
| 149 | } |
---|