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 // 27 // ------------------------------------------------------------------- 28 // 29 // Geant4 source file 30 // 31 // File name: G4ParticleHPInelasticURR.cc 32 // 33 // Authors: Marek Zmeskal (CTU, Czech Technical University in Prague, Czech Republic) 34 // Loic Thulliez (CEA France) 35 // 36 // Creation date: 4 June 2024 37 // 38 // Description: Class to handle URR range, can be omitted once 39 // the proper isotope cross-section is stored in 40 // ParticleHP. 41 // 42 // Modifications: 43 // 44 // ------------------------------------------------------------------- 45 // 46 // 47 48 #include "G4ParticleHPInelasticURR.hh" 49 #include "G4ParticleHPManager.hh" 50 #include "G4HadronicParameters.hh" 51 #include "G4ParticleHPChannel.hh" 52 #include "G4ParticleHPInelastic.hh" 53 #include "G4ParticleHPProbabilityTablesStore.hh" 54 #include "G4SystemOfUnits.hh" 55 #include "G4Threading.hh" 56 57 58 G4ParticleHPInelasticURR::G4ParticleHPInelasticURR() : G4HadronicInteraction( "NeutronHPInelasticURR" ) { 59 SetMinEnergy( 0.0 * CLHEP::eV ); 60 SetMaxEnergy( 20.0 * CLHEP::MeV ); 61 particleHPinelastic = new G4ParticleHPInelastic( G4Neutron::Neutron(), "NeutronHPInelastic" ); 62 } 63 64 65 G4ParticleHPInelasticURR::~G4ParticleHPInelasticURR() {} 66 67 68 G4HadFinalState* G4ParticleHPInelasticURR::ApplyYourself( const G4HadProjectile& aTrack, G4Nucleus& aNucleus ) { 69 if ( doNOTusePTforInelastic ) { 70 return particleHPinelastic->ApplyYourself( aTrack, aNucleus ); 71 } 72 const G4Material* theMaterial = aTrack.GetMaterial(); 73 G4double kineticEnergy = aTrack.GetKineticEnergy(); 74 G4HadFinalState* theFinalState = nullptr; 75 if ( kineticEnergy < (*URRlimits).back().first || kineticEnergy > (*URRlimits).back().second ) { 76 return particleHPinelastic->ApplyYourself( aTrack, aNucleus ); 77 } 78 G4int elementI = -1; 79 G4int isotopeJ = -1; 80 G4int A = aNucleus.GetA_asInt(); 81 G4int Z = aNucleus.GetZ_asInt(); 82 G4ParticleHPManager::GetInstance()->OpenReactionWhiteBoard(); 83 // finds the element and isotope of the selected target aNucleus 84 for ( G4int i = 0; i < (G4int)theMaterial->GetNumberOfElements(); ++i ) { 85 if ( Z == theMaterial->GetElement(i)->GetZasInt() ) { 86 for ( G4int j = 0; j < (G4int)theMaterial->GetElement(i)->GetNumberOfIsotopes(); ++j ) { 87 if ( A == theMaterial->GetElement(i)->GetIsotope(j)->GetN() ) { 88 isotopeJ = j; 89 break; 90 } 91 } 92 // the loop cannot be ended here because the material can have two elements with same Z but different isotopic composition 93 if ( isotopeJ != -1 ) { 94 // isotope was found and for loop is ended 95 elementI = (G4int)theMaterial->GetElement(i)->GetIndex(); 96 break; 97 } 98 } // end if find element 99 } // end element loop 100 // Check whether the energy is out of the URR limits for the given element 101 if ( kineticEnergy < (*URRlimits).at(elementI).first || kineticEnergy > (*URRlimits).at(elementI).second ) { 102 // Call inelastic final state in G4ParicleHPChannel and SELECT ISOTOPE (to be improved in the future) 103 const G4Element* target_element = (*G4Element::GetElementTable())[elementI]; 104 theFinalState = (*G4ParticleHPManager::GetInstance()->GetInelasticFinalStates( aTrack.GetDefinition() ))[elementI] 105 ->ApplyYourself( target_element, aTrack ); 106 // Update target nucleus information according to the selected isotope 107 G4int selectedIsotope_A = G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(); 108 aNucleus.SetParameters( selectedIsotope_A, Z ); 109 const G4Isotope* target_isotope = nullptr; 110 // Find the selected isotope among in the element 111 for ( G4int j = 0; j < (G4int)target_element->GetNumberOfIsotopes(); ++j ) { 112 target_isotope = target_element->GetIsotope(j); 113 if ( target_isotope->GetN() == selectedIsotope_A ) break; 114 } 115 aNucleus.SetIsotope( target_isotope ); 116 } else { 117 // the energy is inside the limits of the URR and the isotope has to be found, calls the final state for the found element and isotope 118 theFinalState = (*G4ParticleHPManager::GetInstance()->GetInelasticFinalStates( aTrack.GetDefinition() ))[elementI] 119 ->ApplyYourself( isotopeJ, Z, A, aTrack ); 120 } 121 G4ParticleHPManager::GetInstance()->CloseReactionWhiteBoard(); 122 return theFinalState; 123 } 124 125 126 void G4ParticleHPInelasticURR::BuildPhysicsTable( const G4ParticleDefinition& ) { 127 particleHPinelastic->BuildPhysicsTable( *(G4Neutron::Neutron()) ); 128 if ( G4HadronicParameters::Instance()->GetTypeTablePT() == "njoy" ) { 129 doNOTusePTforInelastic = true; 130 } else if ( G4HadronicParameters::Instance()->GetTypeTablePT() == "calendf" ) { 131 doNOTusePTforInelastic = false; 132 // in the case of calendf probability tables, it sets the limits of the URR 133 URRlimits = G4ParticleHPManager::GetInstance()->GetURRlimits(); 134 if ( URRlimits == nullptr ) { 135 G4ParticleHPProbabilityTablesStore::GetInstance()->InitURRlimits(); 136 URRlimits = G4ParticleHPProbabilityTablesStore::GetInstance()->GetURRlimits(); 137 G4ParticleHPManager::GetInstance()->RegisterURRlimits( URRlimits ); 138 } 139 } 140 } 141 142 143 const std::pair< G4double, G4double > G4ParticleHPInelasticURR::GetFatalEnergyCheckLevels() const { 144 // max energy non-conservation is mass of heavy nucleus 145 return std::pair< G4double, G4double >( 10.0 * perCent, 350.0 * CLHEP::GeV ); 146 } 147 148 149 G4int G4ParticleHPInelasticURR::GetVerboseLevel() const { 150 return G4ParticleHPManager::GetInstance()->GetVerboseLevel(); 151 } 152 153 154 void G4ParticleHPInelasticURR::SetVerboseLevel( G4int newValue ) { 155 G4ParticleHPManager::GetInstance()->SetVerboseLevel( newValue ); 156 } 157 158 159 void G4ParticleHPInelasticURR::ModelDescription( std::ostream& outFile ) const { 160 outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for Inelastic reaction of neutrons in the unresolved resonance region."; 161 } 162