Geant4  10.01.p03
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 
43  :G4HadronicInteraction("NeutronHPCapture")
44  {
45  SetMinEnergy( 0.0 );
46  SetMaxEnergy( 20.*MeV );
47 // G4cout << "Capture : start of construction!!!!!!!!"<<G4endl;
48  if(!getenv("G4NEUTRONHPDATA"))
49  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
50  dirName = getenv("G4NEUTRONHPDATA");
51  G4String tString = "/Capture";
52  dirName = dirName + tString;
54 // G4cout << "+++++++++++++++++++++++++++++++++++++++++++++++++"<<G4endl;
55 // G4cout <<"Disname="<<dirName<<" numEle="<<numEle<<G4endl;
56  //theCapture = new G4ParticleHPChannel[numEle];
57 // G4cout <<"G4ParticleHPChannel constructed"<<G4endl;
59  //for (G4int i=0; i<numEle; i++)
60  //{
61 // // G4cout << "initializing theCapture "<<i<<" "<< numEle<<G4endl;
62  // theCapture[i].Init((*(G4Element::GetElementTable()))[i], dirName);
63  // theCapture[i].Register(theFS);
64  //}
65  for ( G4int i = 0 ; i < numEle ; i++ )
66  {
67  theCapture.push_back( new G4ParticleHPChannel );
69  (*theCapture[i]).Register(theFS);
70  }
71  delete theFS;
72 // G4cout << "-------------------------------------------------"<<G4endl;
73 // G4cout << "Leaving G4ParticleHPCapture::G4ParticleHPCapture"<<G4endl;
74  }
75 
77  {
78  //delete [] theCapture;
79 // G4cout << "Leaving G4ParticleHPCapture::~G4ParticleHPCapture"<<G4endl;
80  for ( std::vector<G4ParticleHPChannel*>::iterator
81  ite = theCapture.begin() ; ite != theCapture.end() ; ite++ )
82  {
83  delete *ite;
84  }
85  theCapture.clear();
86  }
87 
90  {
91 
93 
95  if(getenv("NeutronHPCapture")) G4cout <<" ####### G4ParticleHPCapture called"<<G4endl;
96  const G4Material * theMaterial = aTrack.GetMaterial();
97  G4int n = theMaterial->GetNumberOfElements();
98  G4int index = theMaterial->GetElement(0)->GetIndex();
99  if(n!=1)
100  {
101  xSec = new G4double[n];
102  G4double sum=0;
103  G4int i;
104  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
105  G4double rWeight;
106  G4ParticleHPThermalBoost aThermalE;
107  for (i=0; i<n; i++)
108  {
109  index = theMaterial->GetElement(i)->GetIndex();
110  rWeight = NumAtomsPerVolume[i];
111  //xSec[i] = theCapture[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
112  xSec[i] = (*theCapture[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
113  theMaterial->GetElement(i),
114  theMaterial->GetTemperature()));
115  xSec[i] *= rWeight;
116  sum+=xSec[i];
117  }
118  G4double random = G4UniformRand();
119  G4double running = 0;
120  for (i=0; i<n; i++)
121  {
122  running += xSec[i];
123  index = theMaterial->GetElement(i)->GetIndex();
124  //if(random<=running/sum) break;
125  if( sum == 0 || random <= running/sum ) break;
126  }
127  if(i==n) i=std::max(0, n-1);
128  delete [] xSec;
129  }
130 
131  //return theCapture[index].ApplyYourself(aTrack);
132  //G4HadFinalState* result = theCapture[index].ApplyYourself(aTrack);
133  G4HadFinalState* result = (*theCapture[index]).ApplyYourself(aTrack);
134 
135  //Overwrite target parameters
136  aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
137  const G4Element* target_element = (*G4Element::GetElementTable())[index];
138  const G4Isotope* target_isotope=NULL;
139  G4int iele = target_element->GetNumberOfIsotopes();
140  for ( G4int j = 0 ; j != iele ; j++ ) {
141  target_isotope=target_element->GetIsotope( j );
142  if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
143  }
144  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
145  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
146  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
147  aNucleus.SetIsotope( target_isotope );
148 
150  return result;
151  }
152 
153 const std::pair<G4double, G4double> G4ParticleHPCapture::GetFatalEnergyCheckLevels() const
154 {
155  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
156  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
157 }
158 
160 {
162  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
163  {
164  G4cout << "G4ParticleHPCapture Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
165  theCapture.push_back( new G4ParticleHPChannel );
167  (*theCapture[i]).Register(theFS);
168  }
169  delete theFS;
171 }
172 
174 {
176 }
178 {
180 }
static G4ParticleHPManager * GetInstance()
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
void Init()
Definition: G4IonTable.cc:87
static const double MeV
Definition: G4SIunits.hh:193
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:93
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
size_t GetIndex() const
Definition: G4Element.hh:181
static const double perCent
Definition: G4SIunits.hh:296
const G4int n
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
std::vector< G4ParticleHPChannel * > theCapture
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
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198