Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPCapture.cc
Go to the documentation of this file.
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// 070523 bug fix for G4FPE_DEBUG on by A. Howard ( and T. Koi)
31//
32// P. Arce, June-2014 Conversion neutron_hp to particle_hp
33//
36#include "G4SystemOfUnits.hh"
39#include "G4ParticleTable.hh"
40#include "G4IonTable.hh"
41#include "G4Threading.hh"
42
44 :G4HadronicInteraction("NeutronHPCapture")
45 ,theCapture(NULL)
46 ,numEle(0)
47 {
48 SetMinEnergy( 0.0 );
49 SetMaxEnergy( 20.*MeV );
50/*
51// G4cout << "Capture : start of construction!!!!!!!!"<<G4endl;
52 if(!G4FindDataDir("G4NEUTRONHPDATA"))
53 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
54 dirName = G4FindDataDir("G4NEUTRONHPDATA");
55 G4String tString = "/Capture";
56 dirName = dirName + tString;
57 numEle = G4Element::GetNumberOfElements();
58// G4cout << "+++++++++++++++++++++++++++++++++++++++++++++++++"<<G4endl;
59// G4cout <<"Disname="<<dirName<<" numEle="<<numEle<<G4endl;
60 //theCapture = new G4ParticleHPChannel[numEle];
61// G4cout <<"G4ParticleHPChannel constructed"<<G4endl;
62 G4ParticleHPCaptureFS * theFS = new G4ParticleHPCaptureFS;
63 //for (G4int i=0; i<numEle; i++)
64 //{
65// // G4cout << "initializing theCapture "<<i<<" "<< numEle<<G4endl;
66 // theCapture[i].Init((*(G4Element::GetElementTable()))[i], dirName);
67 // theCapture[i].Register(theFS);
68 //}
69 for ( G4int i = 0 ; i < numEle ; i++ )
70 {
71 theCapture.push_back( new G4ParticleHPChannel );
72 (*theCapture[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
73 (*theCapture[i]).Register(theFS);
74 }
75 delete theFS;
76// G4cout << "-------------------------------------------------"<<G4endl;
77// G4cout << "Leaving G4ParticleHPCapture::G4ParticleHPCapture"<<G4endl;
78*/
79 }
80
82 {
83 //delete [] theCapture;
84 //vector is shared, only master deletes
86 if ( theCapture != nullptr ) {
87 for ( std::vector<G4ParticleHPChannel*>::iterator
88 ite = theCapture->begin() ; ite != theCapture->end() ; ite++ ) {
89 delete *ite;
90 }
91 theCapture->clear();
92 }
93 }
94 }
95
98 {
99
100 //if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
101
103 const G4Material* theMaterial = aTrack.GetMaterial();
104 G4int n = (G4int)theMaterial->GetNumberOfElements();
105 std::size_t index = theMaterial->GetElement(0)->GetIndex();
106 if(n!=1)
107 {
108 G4double* xSec = new G4double[n];
109 G4double sum=0;
110 G4int i;
111 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
112 G4double rWeight;
113 G4ParticleHPThermalBoost aThermalE;
114 for (i=0; i<n; ++i)
115 {
116 index = theMaterial->GetElement(i)->GetIndex();
117 rWeight = NumAtomsPerVolume[i];
118 //xSec[i] = theCapture[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
119 xSec[i] = ((*theCapture)[index])->GetXsec(aThermalE.GetThermalEnergy(aTrack,
120 theMaterial->GetElement(i),
121 theMaterial->GetTemperature()));
122 xSec[i] *= rWeight;
123 sum+=xSec[i];
124 }
125 G4double random = G4UniformRand();
126 G4double running = 0;
127 for (i=0; i<n; ++i)
128 {
129 running += xSec[i];
130 index = theMaterial->GetElement(i)->GetIndex();
131 //if(random<=running/sum) break;
132 if( sum == 0 || random <= running/sum ) break;
133 }
134 if(i==n) i=std::max(0, n-1);
135 delete [] xSec;
136 }
137
138 //return theCapture[index].ApplyYourself(aTrack);
139 //G4HadFinalState* result = theCapture[index].ApplyYourself(aTrack);
140 G4HadFinalState* result = ((*theCapture)[index])->ApplyYourself(aTrack);
141
142 //Overwrite target parameters
144 const G4Element* target_element = (*G4Element::GetElementTable())[index];
145 const G4Isotope* target_isotope=nullptr;
146 G4int iele = (G4int)target_element->GetNumberOfIsotopes();
147 for ( G4int j = 0 ; j != iele ; ++j ) {
148 target_isotope=target_element->GetIsotope( j );
149 if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
150 }
151 //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
152 //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
153 //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
154 aNucleus.SetIsotope( target_isotope );
155
157 return result;
158 }
159
160const std::pair<G4double, G4double> G4ParticleHPCapture::GetFatalEnergyCheckLevels() const
161{
162 // max energy non-conservation is mass of heavy nucleus
163 return std::pair<G4double, G4double>(10.0*perCent, 350.0*CLHEP::GeV);
164}
165
166/*
167void G4ParticleHPCapture::addChannelForNewElement()
168{
169 G4ParticleHPCaptureFS* theFS = new G4ParticleHPCaptureFS;
170 for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
171 {
172 G4cout << "G4ParticleHPCapture Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
173 theCapture.push_back( new G4ParticleHPChannel );
174 (*theCapture[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
175 (*theCapture[i]).Register(theFS);
176 }
177 delete theFS;
178 numEle = (G4int)G4Element::GetNumberOfElements();
179}
180*/
181
183{
185}
187{
189}
190
192{
193
195
196 theCapture = hpmanager->GetCaptureFinalStates();
197
199
200 if ( theCapture == nullptr ) theCapture = new std::vector<G4ParticleHPChannel*>;
201
202 if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
203
204 if ( theCapture->size() == G4Element::GetNumberOfElements() ) {
206 return;
207 }
208
209 if ( !G4FindDataDir("G4NEUTRONHPDATA") )
210 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
211 dirName = G4FindDataDir("G4NEUTRONHPDATA");
212 G4String tString = "/Capture";
213 dirName = dirName + tString;
214
216 for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; ++i )
217 {
218 theCapture->push_back( new G4ParticleHPChannel );
219 ((*theCapture)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
220 ((*theCapture)[i])->Register(theFS);
221 }
222 delete theFS;
223 hpmanager->RegisterCaptureFinalStates( theCapture );
224 }
226}
227
228void G4ParticleHPCapture::ModelDescription(std::ostream& outFile) const
229{
230 outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for radiative capture reaction of neutrons below 20MeV\n";
231}
const char * G4FindDataDir(const char *)
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetIndex() const
Definition: G4Element.hh:182
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
const G4Material * GetMaterial() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
G4int GetN() const
Definition: G4Isotope.hh:93
G4double GetTemperature() const
Definition: G4Material.hh:177
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:197
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201
void SetParameters(const G4double A, const G4double Z, const G4int numberOfLambdas=0)
Definition: G4Nucleus.cc:307
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:114
void BuildPhysicsTable(const G4ParticleDefinition &)
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
virtual void ModelDescription(std::ostream &outFile) const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
std::vector< G4ParticleHPChannel * > * GetCaptureFinalStates()
static G4ParticleHPManager * GetInstance()
void RegisterCaptureFinalStates(std::vector< G4ParticleHPChannel * > *val)
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4bool IsWorkerThread()
Definition: G4Threading.cc:123
G4bool IsMasterThread()
Definition: G4Threading.cc:124