Geant4  10.01.p02
G4NeutronHPCapture.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 #include "G4NeutronHPCapture.hh"
33 #include "G4NeutronHPManager.hh"
34 #include "G4SystemOfUnits.hh"
35 #include "G4NeutronHPCaptureFS.hh"
36 #include "G4NeutronHPDeExGammas.hh"
37 #include "G4ParticleTable.hh"
38 #include "G4IonTable.hh"
39 #include "G4Threading.hh"
40 
42  :G4HadronicInteraction("NeutronHPCapture")
43  ,theCapture(NULL)
44  ,numEle(0)
45  {
46  SetMinEnergy( 0.0 );
47  SetMaxEnergy( 20.*MeV );
48 /*
49 // G4cout << "Capture : start of construction!!!!!!!!"<<G4endl;
50  if(!getenv("G4NEUTRONHPDATA"))
51  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
52  dirName = getenv("G4NEUTRONHPDATA");
53  G4String tString = "/Capture";
54  dirName = dirName + tString;
55  numEle = G4Element::GetNumberOfElements();
56 // G4cout << "+++++++++++++++++++++++++++++++++++++++++++++++++"<<G4endl;
57 // G4cout <<"Disname="<<dirName<<" numEle="<<numEle<<G4endl;
58  //theCapture = new G4NeutronHPChannel[numEle];
59 // G4cout <<"G4NeutronHPChannel constructed"<<G4endl;
60  G4NeutronHPCaptureFS * theFS = new G4NeutronHPCaptureFS;
61  //for (G4int i=0; i<numEle; i++)
62  //{
63 // // G4cout << "initializing theCapture "<<i<<" "<< numEle<<G4endl;
64  // theCapture[i].Init((*(G4Element::GetElementTable()))[i], dirName);
65  // theCapture[i].Register(theFS);
66  //}
67  for ( G4int i = 0 ; i < numEle ; i++ )
68  {
69  theCapture.push_back( new G4NeutronHPChannel );
70  (*theCapture[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
71  (*theCapture[i]).Register(theFS);
72  }
73  delete theFS;
74 // G4cout << "-------------------------------------------------"<<G4endl;
75 // G4cout << "Leaving G4NeutronHPCapture::G4NeutronHPCapture"<<G4endl;
76 */
77  }
78 
80  {
81  //delete [] theCapture;
82 // G4cout << "Leaving G4NeutronHPCapture::~G4NeutronHPCapture"<<G4endl;
83  for ( std::vector<G4NeutronHPChannel*>::iterator
84  ite = (*theCapture).begin() ; ite != (*theCapture).end() ; ite++ )
85  {
86  delete *ite;
87  }
88  (*theCapture).clear();
89  }
90 
93  {
94 
95  //if ( numEle < (G4int)G4Element::GetNumberOfElements() ) addChannelForNewElement();
96 
98  if(getenv("NeutronHPCapture")) G4cout <<" ####### G4NeutronHPCapture called"<<G4endl;
99  const G4Material * theMaterial = aTrack.GetMaterial();
100  G4int n = theMaterial->GetNumberOfElements();
101  G4int index = theMaterial->GetElement(0)->GetIndex();
102  if(n!=1)
103  {
104  xSec = new G4double[n];
105  G4double sum=0;
106  G4int i;
107  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
108  G4double rWeight;
109  G4NeutronHPThermalBoost aThermalE;
110  for (i=0; i<n; i++)
111  {
112  index = theMaterial->GetElement(i)->GetIndex();
113  rWeight = NumAtomsPerVolume[i];
114  //xSec[i] = theCapture[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
115  xSec[i] = (*(*theCapture)[index]).GetXsec(aThermalE.GetThermalEnergy(aTrack,
116  theMaterial->GetElement(i),
117  theMaterial->GetTemperature()));
118  xSec[i] *= rWeight;
119  sum+=xSec[i];
120  }
121  G4double random = G4UniformRand();
122  G4double running = 0;
123  for (i=0; i<n; i++)
124  {
125  running += xSec[i];
126  index = theMaterial->GetElement(i)->GetIndex();
127  //if(random<=running/sum) break;
128  if( sum == 0 || random <= running/sum ) break;
129  }
130  if(i==n) i=std::max(0, n-1);
131  delete [] xSec;
132  }
133 
134  //return theCapture[index].ApplyYourself(aTrack);
135  //G4HadFinalState* result = theCapture[index].ApplyYourself(aTrack);
136  G4HadFinalState* result = (*(*theCapture)[index]).ApplyYourself(aTrack);
137 
138  //Overwrite target parameters
139  aNucleus.SetParameters(G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
140  const G4Element* target_element = (*G4Element::GetElementTable())[index];
141  const G4Isotope* target_isotope=NULL;
142  G4int iele = target_element->GetNumberOfIsotopes();
143  for ( G4int j = 0 ; j != iele ; j++ ) {
144  target_isotope=target_element->GetIsotope( j );
145  if ( target_isotope->GetN() == G4NeutronHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
146  }
147  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
148  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
149  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
150  aNucleus.SetIsotope( target_isotope );
151 
153  return result;
154  }
155 
156 const std::pair<G4double, G4double> G4NeutronHPCapture::GetFatalEnergyCheckLevels() const
157 {
158  //return std::pair<G4double, G4double>(10*perCent,10*GeV);
159  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
160 }
161 
162 /*
163 void G4NeutronHPCapture::addChannelForNewElement()
164 {
165  G4NeutronHPCaptureFS* theFS = new G4NeutronHPCaptureFS;
166  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
167  {
168  G4cout << "G4NeutronHPCapture Prepairing Data for the new element of " << (*(G4Element::GetElementTable()))[i]->GetName() << G4endl;
169  theCapture.push_back( new G4NeutronHPChannel );
170  (*theCapture[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
171  (*theCapture[i]).Register(theFS);
172  }
173  delete theFS;
174  numEle = (G4int)G4Element::GetNumberOfElements();
175 }
176 */
177 
179 {
181 }
183 {
185 }
187 {
189 
190  theCapture = hpmanager->GetCaptureFinalStates();
191 
192  if ( !G4Threading::IsWorkerThread() ) {
193 
194  if ( theCapture == NULL ) theCapture = new std::vector<G4NeutronHPChannel*>;
195 
196  if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
197 
198  if ( theCapture->size() == G4Element::GetNumberOfElements() ) {
200  return;
201  }
202 
203  if ( !getenv("G4NEUTRONHPDATA") )
204  throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
205  dirName = getenv("G4NEUTRONHPDATA");
206  G4String tString = "/Capture";
207  dirName = dirName + tString;
208 
210  for ( G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements() ; i++ )
211  {
212  (*theCapture).push_back( new G4NeutronHPChannel );
213  (*(*theCapture)[i]).Init((*(G4Element::GetElementTable()))[i], dirName);
214  (*(*theCapture)[i]).Register(theFS);
215  }
216  delete theFS;
218 
219  }
221 
222 }
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
void SetIsotope(const G4Isotope *iso)
Definition: G4Nucleus.hh:122
void Init()
Definition: G4IonTable.cc:87
void BuildPhysicsTable(const G4ParticleDefinition &)
static const double MeV
Definition: G4SIunits.hh:193
std::vector< G4NeutronHPChannel * > * GetCaptureFinalStates()
static G4NeutronHPManager * GetInstance()
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
virtual const std::pair< G4double, G4double > GetFatalEnergyCheckLevels() const
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
int G4int
Definition: G4Types.hh:78
G4int GetVerboseLevel() const
void SetMinEnergy(G4double anEnergy)
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
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
size_t GetIndex() const
Definition: G4Element.hh:181
static const double perCent
Definition: G4SIunits.hh:296
const G4int n
G4bool IsWorkerThread()
Definition: G4Threading.cc:128
void SetVerboseLevel(G4int i)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
std::vector< G4NeutronHPChannel * > * theCapture
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
#define DBL_MAX
Definition: templates.hh:83
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:198
void RegisterCaptureFinalStates(std::vector< G4NeutronHPChannel * > *val)