Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParticleHPChannel.hh
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// Hadronic Process: Very Low Energy Neutron X-Sections
27// original by H.P. Wellisch, TRIUMF, 14-Feb-97
28// Builds and has the Cross-section data for one element and channel.
29//
30// Bug fixes and workarounds in the destructor, F.W.Jones 06-Jul-1999
31// 070612 Fix memory leaking by T. Koi
32//
33// 080520 Delete unnecessary dependencies by T. Koi
34
35// P. Arce, June-2014 Conversion neutron_hp to particle_hp
36//
37#ifndef G4ParticleHPChannel_h
38#define G4ParticleHPChannel_h 1
39
40#include "globals.hh"
42#include "G4ParticleHPVector.hh"
43#include "G4Material.hh"
44#include "G4HadProjectile.hh"
45#include "G4StableIsotopes.hh"
48#include "G4Element.hh"
51
53
54
56{
57public:
58
60 wendtFissionGenerator( G4ParticleHPManager::GetInstance()->GetUseWendtFissionModel() ?
61 G4WendtFissionFragmentGenerator::GetInstance() : nullptr )
62 {
63 if ( G4ParticleHPManager::GetInstance()->GetUseWendtFissionModel() ) {
64 // Make sure both fission fragment models are not active at same time
66 }
67 theProjectile = const_cast<G4ParticleDefinition*>(projectile);
68 theChannelData = new G4ParticleHPVector;
69 theBuffer = 0;
70 theIsotopeWiseData = 0;
71 theFinalStates = 0;
72 active = 0;
73 registerCount = -1;
74 niso = -1;
75 theElement = NULL;
76 }
77
78 G4ParticleHPChannel() : wendtFissionGenerator( G4ParticleHPManager::GetInstance()->GetUseWendtFissionModel() ?
79 G4WendtFissionFragmentGenerator::GetInstance() : nullptr )
80 {
81 if ( G4ParticleHPManager::GetInstance()->GetUseWendtFissionModel() ) {
82 // Make sure both fission fragment models are not active at same time
84 }
85 theProjectile = G4Neutron::Neutron();
86 theChannelData = new G4ParticleHPVector;
87 theBuffer = 0;
88 theIsotopeWiseData = 0;
89 theFinalStates = 0;
90 active = 0;
91 registerCount = -1;
92 niso = -1;
93 theElement = NULL;
94 }
95
97 {
98 delete theChannelData;
99 // Following statement disabled to avoid SEGV
100 // theBuffer is also deleted as "theChannelData" in
101 // ~G4ParticleHPIsoData. FWJ 06-Jul-1999
102 //if(theBuffer != 0) delete theBuffer;
103 if(theIsotopeWiseData != 0) delete [] theIsotopeWiseData;
104 // Deletion of FinalStates disabled to avoid endless looping
105 // in the destructor heirarchy. FWJ 06-Jul-1999
106 //if(theFinalStates != 0)
107 //{
108 // for(i=0; i<niso; i++)
109 // {
110 // delete theFinalStates[i];
111 // }
112 // delete [] theFinalStates;
113 //}
114 // FWJ experiment
115 //if(active!=0) delete [] active;
116 // T.K.
117 if ( theFinalStates != 0 )
118 {
119 for ( G4int i = 0 ; i < niso ; i++ )
120 {
121 delete theFinalStates[i];
122 }
123 delete [] theFinalStates;
124 }
125 if ( active != 0 ) delete [] active;
126 }
127
128 G4double GetXsec(G4double energy);
129
130 G4double GetWeightedXsec(G4double energy, G4int isoNumber);
131
132 G4double GetFSCrossSection(G4double energy, G4int isoNumber);
133
134 inline G4bool IsActive(G4int isoNumber) { return active[isoNumber]; }
135
136 inline G4bool HasFSData(G4int isoNumber) { return theFinalStates[isoNumber]->HasFSData(); }
137
138 inline G4bool HasAnyData(G4int isoNumber) { return theFinalStates[isoNumber]->HasAnyData(); }
139
141
142 void Init(G4Element * theElement, const G4String dirName);
143
144 void Init(G4Element * theElement, const G4String dirName, const G4String fsType);
145
146 //void UpdateData(G4int A, G4int Z, G4int index, G4double abundance);
147 void UpdateData(G4int A, G4int Z, G4int index, G4double abundance, G4ParticleDefinition* projectile) { G4int M = 0; UpdateData( A, Z, M, index, abundance, projectile); };
148 void UpdateData(G4int A, G4int Z, G4int M, G4int index, G4double abundance, G4ParticleDefinition* projectile);
149
150 void Harmonise(G4ParticleHPVector *& theStore, G4ParticleHPVector * theNew);
151
152 G4HadFinalState * ApplyYourself(const G4HadProjectile & theTrack, G4int isoNumber=-1);
153
154 inline G4int GetNiso() {return niso;}
155
156 inline G4double GetN(G4int i) {return theFinalStates[i]->GetN();}
157 inline G4double GetZ(G4int i) {return theFinalStates[i]->GetZ();}
158 inline G4double GetM(G4int i) {return theFinalStates[i]->GetM();}
159
161 {
162 G4bool result = false;
163 G4int i;
164 for(i=0; i<niso; i++)
165 {
166 if(theFinalStates[i]->HasAnyData()) result = true;
167 }
168 return result;
169 }
170
171 void DumpInfo();
172
174 return theFSType;
175 }
176
178 return theFinalStates;
179 }
180
181private:
182 G4ParticleDefinition* theProjectile;
183
184 G4ParticleHPVector * theChannelData; // total (element) cross-section for this channel
185 G4ParticleHPVector * theBuffer;
186
187 G4ParticleHPIsoData * theIsotopeWiseData; // these are the isotope-wise cross-sections for each final state.
188 G4ParticleHPFinalState ** theFinalStates; // also these are isotope-wise pionters, parallel to the above.
189 G4bool * active;
190 G4int niso;
191
192 G4StableIsotopes theStableOnes;
193
194 G4String theDir;
195 G4String theFSType;
196 G4Element * theElement;
197
198 G4int registerCount;
199
200 G4WendtFissionFragmentGenerator* const wendtFissionGenerator;
201
202};
203
204#endif
#define M(row, col)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4bool HasAnyData(G4int isoNumber)
void UpdateData(G4int A, G4int Z, G4int index, G4double abundance, G4ParticleDefinition *projectile)
G4bool HasFSData(G4int isoNumber)
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)
G4bool IsActive(G4int isoNumber)
G4String GetFSType() const
void Init(G4Element *theElement, const G4String dirName)
G4ParticleHPChannel(G4ParticleDefinition *projectile)
void Harmonise(G4ParticleHPVector *&theStore, G4ParticleHPVector *theNew)
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
G4double GetM(G4int i)
G4double GetN(G4int i)
G4ParticleHPFinalState ** GetFinalStates() const
G4bool Register(G4ParticleHPFinalState *theFS)
G4double GetXsec(G4double energy)
G4double GetZ(G4int i)
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
void SetProduceFissionFragments(G4bool val)
static G4ParticleHPManager * GetInstance()