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: G4ParticleHPFissionURR.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 the 39 // proper isotope cross-section is stored in ParticleHP. 40 // 41 // Modifications: 42 // 43 // ------------------------------------------------------------------- 44 // 45 // 46 47 #include "G4ParticleHPFissionURR.hh" 48 #include "G4ParticleHPManager.hh" 49 #include "G4ParticleHPChannel.hh" 50 #include "G4ParticleHPFission.hh" 51 #include "G4WendtFissionFragmentGenerator.hh" 52 #include "G4ParticleHPProbabilityTablesStore.hh" 53 #include "G4SystemOfUnits.hh" 54 #include "G4Threading.hh" 55 56 57 G4ParticleHPFissionURR::G4ParticleHPFissionURR() : G4HadronicInteraction( "NeutronHPFissionURR" ) { 58 SetMinEnergy( 0.0 * CLHEP::eV ); 59 SetMaxEnergy( 20.0 * CLHEP::MeV ); 60 particleHPfission = new G4ParticleHPFission; 61 } 62 63 64 G4ParticleHPFissionURR::~G4ParticleHPFissionURR() {} 65 66 67 G4HadFinalState* G4ParticleHPFissionURR::ApplyYourself( const G4HadProjectile& aTrack, G4Nucleus& aNucleus ) { 68 const G4Material* theMaterial = aTrack.GetMaterial(); 69 G4double kineticEnergy = aTrack.GetKineticEnergy(); 70 G4HadFinalState* theFinalState = nullptr; 71 if ( kineticEnergy < (*URRlimits).back().first || kineticEnergy > (*URRlimits).back().second ) { 72 return particleHPfission->ApplyYourself( aTrack, aNucleus ); 73 } 74 G4int elementI = -1; 75 G4int isotopeJ = -1; 76 G4int A = aNucleus.GetA_asInt(); 77 G4int Z = aNucleus.GetZ_asInt(); 78 // finds the element and isotope of the selected target aNucleus 79 for ( G4int i = 0; i < (G4int)theMaterial->GetNumberOfElements(); ++i ) { 80 if ( Z == theMaterial->GetElement(i)->GetZasInt() ) { 81 for ( G4int j = 0; j < (G4int)theMaterial->GetElement(i)->GetNumberOfIsotopes(); ++j ) { 82 if ( A == theMaterial->GetElement(i)->GetIsotope(j)->GetN() ) { 83 isotopeJ = j; 84 break; 85 } 86 } 87 // the loop cannot be ended here because the material can have two elements with same Z but different isotopic composition 88 if ( isotopeJ != -1 ) { 89 // isotope was found and for loop is ended 90 elementI = (G4int)theMaterial->GetElement(i)->GetIndex(); 91 break; 92 } 93 } // end if find element 94 } // end element loop 95 // Check whether the energy is out of the URR limits for the given element 96 if ( kineticEnergy < (*URRlimits).at(elementI).first || kineticEnergy > (*URRlimits).at(elementI).second ) { 97 // Call fission final state in G4ParicleHPChannel and SELECT ISOTOPE (to be improved in the future) 98 G4ParticleHPManager::GetInstance()->OpenReactionWhiteBoard(); 99 theFinalState = (*G4ParticleHPManager::GetInstance()->GetFissionFinalStates())[elementI]->ApplyYourself( aTrack, -2 ); 100 // Update target nucleus information according to the selected isotope 101 G4int selectedIsotope_A = G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(); 102 aNucleus.SetParameters( selectedIsotope_A, Z ); 103 const G4Element* target_element = (*G4Element::GetElementTable())[elementI]; 104 const G4Isotope* target_isotope = nullptr; 105 // Find the selected isotope among in the element 106 for ( G4int j = 0; j < (G4int)target_element->GetNumberOfIsotopes(); ++j ) { 107 target_isotope = target_element->GetIsotope(j); 108 if ( target_isotope->GetN() == selectedIsotope_A ) break; 109 } 110 aNucleus.SetIsotope( target_isotope ); 111 G4ParticleHPManager::GetInstance()->CloseReactionWhiteBoard(); 112 } else { 113 // the energy is inside the limits of the URR, part copied from G4ParticleHPChannel::ApplyYourself 114 if ( G4ParticleHPManager::GetInstance()->GetUseWendtFissionModel() ) { 115 if ( (*G4ParticleHPManager::GetInstance()->GetFissionFinalStates())[elementI]->GetWendtFissionGenerator() ) { 116 theFinalState = (*G4ParticleHPManager::GetInstance()->GetFissionFinalStates())[elementI]->GetWendtFissionGenerator()->ApplyYourself( aTrack, Z, A ); 117 } 118 } 119 if ( ! theFinalState ) { 120 G4int icounter = 0; 121 G4int icounter_max = 1024; 122 while ( theFinalState == nullptr ) { 123 icounter++; 124 if ( icounter > icounter_max ) { 125 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl; 126 break; 127 } 128 // calls the final state for the found element and isotope 129 theFinalState = ((*G4ParticleHPManager::GetInstance()->GetFissionFinalStates())[elementI]->GetFinalStates())[isotopeJ]->ApplyYourself( aTrack ); 130 } 131 } 132 } 133 return theFinalState; 134 } 135 136 137 void G4ParticleHPFissionURR::BuildPhysicsTable( const G4ParticleDefinition& ) { 138 particleHPfission->BuildPhysicsTable( *(G4Neutron::Neutron()) ); 139 URRlimits = G4ParticleHPManager::GetInstance()->GetURRlimits(); 140 if ( URRlimits == nullptr ) { 141 G4ParticleHPProbabilityTablesStore::GetInstance()->InitURRlimits(); 142 URRlimits = G4ParticleHPProbabilityTablesStore::GetInstance()->GetURRlimits(); 143 G4ParticleHPManager::GetInstance()->RegisterURRlimits( URRlimits ); 144 } 145 } 146 147 148 const std::pair< G4double, G4double > G4ParticleHPFissionURR::GetFatalEnergyCheckLevels() const { 149 // max energy non-conservation is mass of heavy nucleus 150 return std::pair< G4double, G4double >( 10.0 * perCent, 350.0 * CLHEP::GeV ); 151 } 152 153 154 G4int G4ParticleHPFissionURR::GetVerboseLevel() const { 155 return G4ParticleHPManager::GetInstance()->GetVerboseLevel(); 156 } 157 158 159 void G4ParticleHPFissionURR::SetVerboseLevel( G4int newValue ) { 160 G4ParticleHPManager::GetInstance()->SetVerboseLevel( newValue ); 161 } 162 163 164 void G4ParticleHPFissionURR::ModelDescription( std::ostream& outFile ) const { 165 outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for fission reaction of neutrons in the unresolved resonance region."; 166 } 167