Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4FissLib.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//
27// This software was developed by Lawrence Livermore National Laboratory.
28//
29// Redistribution and use in source and binary forms, with or without
30// modification, are permitted provided that the following conditions are met:
31//
32// 1. Redistributions of source code must retain the above copyright notice,
33// this list of conditions and the following disclaimer.
34// 2. Redistributions in binary form must reproduce the above copyright notice,
35// this list of conditions and the following disclaimer in the documentation
36// and/or other materials provided with the distribution.
37// 3. The name of the author may not be used to endorse or promote products
38// derived from this software without specific prior written permission.
39//
40// THIS SOFTWARE IS PROVIDED BY THE AUTHOR ``AS IS'' AND ANY EXPRESS OR IMPLIED
41// WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF
42// MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO
43// EVENT SHALL THE AUTHOR BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
44// SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
45// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
46// OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY,
47// WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR
48// OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF
49// ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
50//
51// Copyright (c) 2006 The Regents of the University of California.
52// All rights reserved.
53// UCRL-CODE-224807
54//
55//
56// $Id$
57//
58// neutron_hp -- source file
59// J.M. Verbeke, Jan-2007
60// A prototype of the low energy neutron transport model.
61//
62#include "G4FissLib.hh"
63#include "G4SystemOfUnits.hh"
64
66 :xSec(0)
67{
68 SetMinEnergy(0.0);
69 SetMaxEnergy(20.*MeV);
70 if(!getenv("G4NEUTRONHPDATA")) {
71 G4cout << "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files." << G4endl;
72 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
73 }
74 dirName = getenv("G4NEUTRONHPDATA");
75 G4String tString = "/Fission/";
76 dirName = dirName + tString;
78 theFission = new G4NeutronHPChannel[numEle];
79
80 for (G4int i=0; i<numEle; i++)
81 {
82// G4cout << "G4FissLib::G4FissLib(): element "<< i << " : " << (*(G4Element::GetElementTable()))[i]->GetZ()<< G4endl;
83 if((*(G4Element::GetElementTable()))[i]->GetZ()>89)
84 {
85 theFission[i].Init((*(G4Element::GetElementTable()))[i], dirName);
86 theFission[i].Register(&theLibrary);
87 }
88 }
89}
90
92{
93 delete [] theFission;
94}
95
98{
99 const G4Material* theMaterial = aTrack.GetMaterial();
100 G4int n = theMaterial->GetNumberOfElements();
101 G4int index = theMaterial->GetElement(0)->GetIndex();
102
103 if (n != 1) {
104 xSec = new G4double[n];
105 G4double sum = 0;
106 G4int i;
107 G4int imat;
108 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
109 G4double rWeight;
110 G4NeutronHPThermalBoost aThermalE;
111 for (i = 0; i < n; i++) {
112 imat = theMaterial->GetElement(i)->GetIndex();
113 rWeight = NumAtomsPerVolume[i];
114 xSec[i] = theFission[imat].GetXsec(aThermalE.GetThermalEnergy(aTrack,
115 theMaterial->GetElement(i),
116 theMaterial->GetTemperature()));
117 xSec[i] *= rWeight;
118 sum+=xSec[i];
119 }
120
121 G4double random = G4UniformRand();
122 G4double running = 0;
123 for (i = 0; i < n; i++) {
124 running += xSec[i];
125 index = theMaterial->GetElement(i)->GetIndex();
126 if(random<=running/sum) break;
127 }
128 delete [] xSec;
129 }
130
131 return theFission[index].ApplyYourself(aTrack);
132}
133
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
size_t GetIndex() const
Definition: G4Element.hh:182
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
Definition: G4FissLib.cc:97
~G4FissLib()
Definition: G4FissLib.cc:91
const G4Material * GetMaterial() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
G4double GetTemperature() const
Definition: G4Material.hh:181
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:201
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
G4bool Register(G4NeutronHPFinalState *theFS)
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)
G4double GetXsec(G4double energy)
void Init(G4Element *theElement, const G4String dirName)
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)