| 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 | // neutron_hp -- source file
|
|---|
| 27 | // J.P. Wellisch, Nov-1996
|
|---|
| 28 | // A prototype of the low energy neutron transport model.
|
|---|
| 29 | //
|
|---|
| 30 | // 070618 fix memory leaking by T. Koi
|
|---|
| 31 | // 071002 enable cross section dump by T. Koi
|
|---|
| 32 | // 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
|
|---|
| 33 | // 081124 Protect invalid read which caused run time errors by T. Koi
|
|---|
| 34 |
|
|---|
| 35 | #include "G4NeutronHPFissionData.hh"
|
|---|
| 36 | #include "G4Neutron.hh"
|
|---|
| 37 | #include "G4ElementTable.hh"
|
|---|
| 38 | #include "G4NeutronHPData.hh"
|
|---|
| 39 |
|
|---|
| 40 | G4bool G4NeutronHPFissionData::IsApplicable(const G4DynamicParticle*aP, const G4Element*)
|
|---|
| 41 | {
|
|---|
| 42 | G4bool result = true;
|
|---|
| 43 | G4double eKin = aP->GetKineticEnergy();
|
|---|
| 44 | if(eKin>20*MeV||aP->GetDefinition()!=G4Neutron::Neutron()) result = false;
|
|---|
| 45 | return result;
|
|---|
| 46 | }
|
|---|
| 47 |
|
|---|
| 48 | G4NeutronHPFissionData::G4NeutronHPFissionData()
|
|---|
| 49 | {
|
|---|
| 50 | theCrossSections = 0;
|
|---|
| 51 | BuildPhysicsTable(*G4Neutron::Neutron());
|
|---|
| 52 | }
|
|---|
| 53 |
|
|---|
| 54 | G4NeutronHPFissionData::~G4NeutronHPFissionData()
|
|---|
| 55 | {
|
|---|
| 56 |
|
|---|
| 57 | // TKDB
|
|---|
| 58 | if ( theCrossSections != NULL )
|
|---|
| 59 | {
|
|---|
| 60 | theCrossSections->clearAndDestroy();
|
|---|
| 61 | delete theCrossSections;
|
|---|
| 62 | }
|
|---|
| 63 | }
|
|---|
| 64 |
|
|---|
| 65 | void G4NeutronHPFissionData::BuildPhysicsTable(const G4ParticleDefinition& aP)
|
|---|
| 66 | {
|
|---|
| 67 | if(&aP!=G4Neutron::Neutron())
|
|---|
| 68 | throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
|
|---|
| 69 | size_t numberOfElements = G4Element::GetNumberOfElements();
|
|---|
| 70 | //theCrossSections = new G4PhysicsTable( numberOfElements );
|
|---|
| 71 | // TKDB
|
|---|
| 72 | if ( theCrossSections == NULL ) theCrossSections = new G4PhysicsTable( numberOfElements );
|
|---|
| 73 |
|
|---|
| 74 | // make a PhysicsVector for each element
|
|---|
| 75 |
|
|---|
| 76 | static const G4ElementTable *theElementTable = G4Element::GetElementTable();
|
|---|
| 77 | for( size_t i=0; i<numberOfElements; ++i )
|
|---|
| 78 | {
|
|---|
| 79 | G4PhysicsVector* physVec = G4NeutronHPData::
|
|---|
| 80 | Instance()->MakePhysicsVector((*theElementTable)[i], this);
|
|---|
| 81 | theCrossSections->push_back(physVec);
|
|---|
| 82 | }
|
|---|
| 83 | }
|
|---|
| 84 |
|
|---|
| 85 | void G4NeutronHPFissionData::DumpPhysicsTable(const G4ParticleDefinition& aP)
|
|---|
| 86 | {
|
|---|
| 87 | if(&aP!=G4Neutron::Neutron())
|
|---|
| 88 | throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
|
|---|
| 89 |
|
|---|
| 90 | //
|
|---|
| 91 | // Dump element based cross section
|
|---|
| 92 | // range 10e-5 eV to 20 MeV
|
|---|
| 93 | // 10 point per decade
|
|---|
| 94 | // in barn
|
|---|
| 95 | //
|
|---|
| 96 |
|
|---|
| 97 | G4cout << G4endl;
|
|---|
| 98 | G4cout << G4endl;
|
|---|
| 99 | G4cout << "Fission Cross Section of Neutron HP"<< G4endl;
|
|---|
| 100 | G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
|
|---|
| 101 | G4cout << G4endl;
|
|---|
| 102 | G4cout << "Name of Element" << G4endl;
|
|---|
| 103 | G4cout << "Energy[eV] XS[barn]" << G4endl;
|
|---|
| 104 | G4cout << G4endl;
|
|---|
| 105 |
|
|---|
| 106 | size_t numberOfElements = G4Element::GetNumberOfElements();
|
|---|
| 107 | static const G4ElementTable *theElementTable = G4Element::GetElementTable();
|
|---|
| 108 |
|
|---|
| 109 | for ( size_t i = 0 ; i < numberOfElements ; ++i )
|
|---|
| 110 | {
|
|---|
| 111 |
|
|---|
| 112 | G4cout << (*theElementTable)[i]->GetName() << G4endl;
|
|---|
| 113 |
|
|---|
| 114 | if ( (*((*theCrossSections)(i))).GetVectorLength() == 0 )
|
|---|
| 115 | {
|
|---|
| 116 | G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
|
|---|
| 117 | G4cout << G4endl;
|
|---|
| 118 | continue;
|
|---|
| 119 | }
|
|---|
| 120 |
|
|---|
| 121 | G4int ie = 0;
|
|---|
| 122 |
|
|---|
| 123 | for ( ie = 0 ; ie < 130 ; ie++ )
|
|---|
| 124 | {
|
|---|
| 125 | G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
|
|---|
| 126 | G4bool outOfRange = false;
|
|---|
| 127 |
|
|---|
| 128 | if ( eKinetic < 20*MeV )
|
|---|
| 129 | {
|
|---|
| 130 | G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
|
|---|
| 131 | }
|
|---|
| 132 |
|
|---|
| 133 | }
|
|---|
| 134 |
|
|---|
| 135 | G4cout << G4endl;
|
|---|
| 136 | }
|
|---|
| 137 |
|
|---|
| 138 | //G4cout << "G4NeutronHPFissionData::DumpPhysicsTable still to be implemented"<<G4endl;
|
|---|
| 139 | }
|
|---|
| 140 |
|
|---|
| 141 | #include "G4NucleiProperties.hh"
|
|---|
| 142 |
|
|---|
| 143 | G4double G4NeutronHPFissionData::
|
|---|
| 144 | GetCrossSection(const G4DynamicParticle* aP, const G4Element*anE, G4double aT)
|
|---|
| 145 | {
|
|---|
| 146 | G4double result = 0;
|
|---|
| 147 | if(anE->GetZ()<90) return result;
|
|---|
| 148 | G4bool outOfRange;
|
|---|
| 149 | G4int index = anE->GetIndex();
|
|---|
| 150 |
|
|---|
| 151 | // prepare neutron
|
|---|
| 152 | G4double eKinetic = aP->GetKineticEnergy();
|
|---|
| 153 | G4ReactionProduct theNeutron( aP->GetDefinition() );
|
|---|
| 154 | theNeutron.SetMomentum( aP->GetMomentum() );
|
|---|
| 155 | theNeutron.SetKineticEnergy( eKinetic );
|
|---|
| 156 |
|
|---|
| 157 | // prepare thermal nucleus
|
|---|
| 158 | G4Nucleus aNuc;
|
|---|
| 159 | G4double eps = 0.0001;
|
|---|
| 160 | G4double theA = anE->GetN();
|
|---|
| 161 | G4double theZ = anE->GetZ();
|
|---|
| 162 | G4double eleMass;
|
|---|
| 163 | eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
|
|---|
| 164 | ) / G4Neutron::Neutron()->GetPDGMass();
|
|---|
| 165 |
|
|---|
| 166 | G4ReactionProduct boosted;
|
|---|
| 167 | G4double aXsection;
|
|---|
| 168 |
|
|---|
| 169 | // MC integration loop
|
|---|
| 170 | G4int counter = 0;
|
|---|
| 171 | G4double buffer = 0;
|
|---|
| 172 | G4int size = G4int(std::max(10., aT/60*kelvin));
|
|---|
| 173 | G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
|
|---|
| 174 | G4double neutronVMag = neutronVelocity.mag();
|
|---|
| 175 |
|
|---|
| 176 | while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
|
|---|
| 177 | {
|
|---|
| 178 | if(counter) buffer = result/counter;
|
|---|
| 179 | while (counter<size)
|
|---|
| 180 | {
|
|---|
| 181 | counter ++;
|
|---|
| 182 | G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
|
|---|
| 183 | boosted.Lorentz(theNeutron, aThermalNuc);
|
|---|
| 184 | G4double theEkin = boosted.GetKineticEnergy();
|
|---|
| 185 | aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
|
|---|
| 186 | // velocity correction.
|
|---|
| 187 | G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
|
|---|
| 188 | aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
|
|---|
| 189 | result += aXsection;
|
|---|
| 190 | }
|
|---|
| 191 | size += size;
|
|---|
| 192 | }
|
|---|
| 193 | result /= counter;
|
|---|
| 194 | return result;
|
|---|
| 195 | }
|
|---|