Geant4  10.02
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 //
34 #include "G4ParticleHPCapture.hh"
35 #include "G4ParticleHPManager.hh"
36 #include "G4SystemOfUnits.hh"
37 #include "G4ParticleHPCaptureFS.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(!getenv("G4NEUTRONHPDATA"))
53  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
54  dirName = getenv("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 // G4cout << "Leaving G4ParticleHPCapture::~G4ParticleHPCapture"<<G4endl;
85  for ( std::vector<G4ParticleHPChannel*>::iterator
86  ite = theCapture->begin() ; ite != theCapture->end() ; ite++ )
87  {
88  delete *ite;
89  }
90  theCapture->clear();
91  }
92 
95  {
96 
97  //if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
98 
100  if(getenv("NeutronHPCapture")) G4cout <<" ####### G4ParticleHPCapture called"<<G4endl;
101  const G4Material * theMaterial = aTrack.GetMaterial();
102  G4int n = theMaterial->GetNumberOfElements();
103  G4int index = theMaterial->GetElement(0)->GetIndex();
104  if(n!=1)
105  {
106  xSec = new G4double[n];
107  G4double sum=0;
108  G4int i;
109  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
110  G4double rWeight;
111  G4ParticleHPThermalBoost aThermalE;
112  for (i=0; i<n; i++)
113  {
114  index = theMaterial->GetElement(i)->GetIndex();
115  rWeight = NumAtomsPerVolume[i];
116  //xSec[i] = theCapture[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
117  xSec[i] = ((*theCapture)[index])->GetXsec(aThermalE.GetThermalEnergy(aTrack,
118  theMaterial->GetElement(i),
119  theMaterial->GetTemperature()));
120  xSec[i] *= rWeight;
121  sum+=xSec[i];
122  }
123  G4double random = G4UniformRand();
124  G4double running = 0;
125  for (i=0; i<n; i++)
126  {
127  running += xSec[i];
128  index = theMaterial->GetElement(i)->GetIndex();
129  //if(random<=running/sum) break;
130  if( sum == 0 || random <= running/sum ) break;
131  }
132  if(i==n) i=std::max(0, n-1);
133  delete [] xSec;
134  }
135 
136  //return theCapture[index].ApplyYourself(aTrack);
137  //G4HadFinalState* result = theCapture[index].ApplyYourself(aTrack);
138  G4HadFinalState* result = ((*theCapture)[index])->ApplyYourself(aTrack);
139 
140  //Overwrite target parameters
141  aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
142  const G4Element* target_element = (*G4Element::GetElementTable())[index];
143  const G4Isotope* target_isotope=NULL;
144  G4int iele = target_element->GetNumberOfIsotopes();
145  for ( G4int j = 0 ; j != iele ; j++ ) {
146  target_isotope=target_element->GetIsotope( j );
147  if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
148  }
149  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
150  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
151  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
152  aNucleus.SetIsotope( target_isotope );
153 
155  return result;
156  }
157 
158 const std::pair<G4double, G4double> G4ParticleHPCapture::GetFatalEnergyCheckLevels() const
159 {
160  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
161  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
162 }
163 
164 /*
165 void G4ParticleHPCapture::addChannelForNewElement()
166 {
167  G4ParticleHPCaptureFS* theFS = new G4ParticleHPCaptureFS;
168  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
169  {
170  G4cout << "G4ParticleHPCapture Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
171  theCapture.push_back( new G4ParticleHPChannel );
172  (*theCapture[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
173  (*theCapture[i]).Register(theFS);
174  }
175  delete theFS;
176  numEle = (G4int)G4Element::GetNumberOfElements();
177 }
178 */
179 
181 {
183 }
185 {
187 }
188 
190 {
191 
193 
194  theCapture = hpmanager->GetCaptureFinalStates();
195 
196  if ( G4Threading::IsMasterThread() ) {
197 
198  if ( theCapture == NULL ) theCapture = new std::vector<G4ParticleHPChannel*>;
199 
200  if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
201 
202  if ( theCapture->size() == G4Element::GetNumberOfElements() ) {
204  return;
205  }
206 
207  if ( !getenv("G4NEUTRONHPDATA") )
208  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
209  dirName = getenv("G4NEUTRONHPDATA");
210  G4String tString = "/Capture";
211  dirName = dirName + tString;
212 
214  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
215  {
216  theCapture->push_back( new G4ParticleHPChannel );
217  ((*theCapture)[i])->Init((*(G4Element::GetElementTable()))[i], dirName);
218  ((*theCapture)[i])->Register(theFS);
219  }
220  delete theFS;
222  }
224 }
225 
226 void G4ParticleHPCapture::ModelDescription(std::ostream& outFile) const
227 {
228  outFile << "High Precision model based on Evaluated Nuclear Data Files (ENDF) for radiative capture reaction of neutrons below 20MeV\n";
229 }
static G4ParticleHPManager * GetInstance()
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
std::vector< G4ParticleHPChannel * > * GetCaptureFinalStates()
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
void Init()
Definition: G4IonTable.cc:87
static const double MeV
Definition: G4SIunits.hh:211
std::vector< G4ParticleHPChannel * > * theCapture
G4int GetVerboseLevel() const
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
int G4int
Definition: G4Types.hh:78
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
void SetMinEnergy(G4double anEnergy)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
void Register(T *inst)
Definition: G4AutoDelete.hh:65
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
virtual void ModelDescription(std::ostream &outFile) const
size_t GetIndex() const
Definition: G4Element.hh:181
static const double perCent
Definition: G4SIunits.hh:329
void RegisterCaptureFinalStates(std::vector< G4ParticleHPChannel * > *val)
const G4int n
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetMaxEnergy(const G4double anEnergy)
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
G4double GetTemperature() const
Definition: G4Material.hh:182
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
G4bool IsMasterThread()
Definition: G4Threading.cc:130
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:211
void BuildPhysicsTable(const G4ParticleDefinition &)