Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronHPChannel Class Reference

#include <G4NeutronHPChannel.hh>

Public Member Functions

 G4NeutronHPChannel ()
 
 ~G4NeutronHPChannel ()
 
G4double GetXsec (G4double energy)
 
G4double GetWeightedXsec (G4double energy, G4int isoNumber)
 
G4double GetFSCrossSection (G4double energy, G4int isoNumber)
 
G4bool IsActive (G4int isoNumber)
 
G4bool HasFSData (G4int isoNumber)
 
G4bool HasAnyData (G4int isoNumber)
 
G4bool Register (G4NeutronHPFinalState *theFS)
 
void Init (G4Element *theElement, const G4String dirName)
 
void Init (G4Element *theElement, const G4String dirName, const G4String fsType)
 
void UpdateData (G4int A, G4int Z, G4int index, G4double abundance)
 
void UpdateData (G4int A, G4int Z, G4int M, G4int index, G4double abundance)
 
void Harmonise (G4NeutronHPVector *&theStore, G4NeutronHPVector *theNew)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &theTrack, G4int isoNumber=-1)
 
G4int GetNiso ()
 
G4double GetN (G4int i)
 
G4double GetZ (G4int i)
 
G4bool HasDataInAnyFinalState ()
 

Detailed Description

Definition at line 54 of file G4NeutronHPChannel.hh.

Constructor & Destructor Documentation

◆ G4NeutronHPChannel()

G4NeutronHPChannel::G4NeutronHPChannel ( )
inline

Definition at line 58 of file G4NeutronHPChannel.hh.

59 {
60 theChannelData = new G4NeutronHPVector;
61 theBuffer = 0;
62 theIsotopeWiseData = 0;
63 theFinalStates = 0;
64 active = 0;
65 registerCount = -1;
66 }

◆ ~G4NeutronHPChannel()

G4NeutronHPChannel::~G4NeutronHPChannel ( )
inline

Definition at line 68 of file G4NeutronHPChannel.hh.

69 {
70 delete theChannelData;
71 // Following statement disabled to avoid SEGV
72 // theBuffer is also deleted as "theChannelData" in
73 // ~G4NeutronHPIsoData. FWJ 06-Jul-1999
74 //if(theBuffer != 0) delete theBuffer;
75 if(theIsotopeWiseData != 0) delete [] theIsotopeWiseData;
76 // Deletion of FinalStates disabled to avoid endless looping
77 // in the destructor heirarchy. FWJ 06-Jul-1999
78 //if(theFinalStates != 0)
79 //{
80 // for(i=0; i<niso; i++)
81 // {
82 // delete theFinalStates[i];
83 // }
84 // delete [] theFinalStates;
85 //}
86 // FWJ experiment
87 //if(active!=0) delete [] active;
88// T.K.
89 if ( theFinalStates != 0 )
90 {
91 for ( G4int i = 0 ; i < niso ; i++ )
92 {
93 delete theFinalStates[i];
94 }
95 delete [] theFinalStates;
96 }
97 if ( active != 0 ) delete [] active;
98
99 }
int G4int
Definition: G4Types.hh:66

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4NeutronHPChannel::ApplyYourself ( const G4HadProjectile theTrack,
G4int  isoNumber = -1 
)

Definition at line 210 of file G4NeutronHPChannel.cc.

212 {
213// G4cout << "G4NeutronHPChannel::ApplyYourself+"<<niso<<G4endl;
214 if ( anIsotope != -1 )
215 {
216 //Inelastic Case
217 //G4cout << "G4NeutronHPChannel Inelastic Case"
218 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
221 return theFinalStates[anIsotope]->ApplyYourself(theTrack);
222 }
223 G4double sum=0;
224 G4int it=0;
225 G4double * xsec = new G4double[niso];
226 G4NeutronHPThermalBoost aThermalE;
227 for (G4int i=0; i<niso; i++)
228 {
229 if(theFinalStates[i]->HasAnyData())
230 {
231 xsec[i] = theIsotopeWiseData[i].GetXsec(aThermalE.GetThermalEnergy(theTrack,
232 theFinalStates[i]->GetN(),
233 theFinalStates[i]->GetZ(),
234 theTrack.GetMaterial()->GetTemperature()));
235 sum += xsec[i];
236 }
237 else
238 {
239 xsec[i]=0;
240 }
241 }
242 if(sum == 0)
243 {
244// G4cout << "G4NeutronHPChannel::ApplyYourself theFinalState->Initialize+"<<G4endl;
245// G4cout << "G4NeutronHPChannel::ApplyYourself theFinalState->Initialize-"<<G4endl;
246 it = static_cast<G4int>(niso*G4UniformRand());
247 }
248 else
249 {
250// G4cout << "Are we still here? "<<sum<<G4endl;
251// G4cout << "TESTHP 23 NISO="<<niso<<G4endl;
252 G4double random = G4UniformRand();
253 G4double running=0;
254// G4cout << "G4NeutronHPChannel::ApplyYourself Done the sum"<<niso<<G4endl;
255// G4cout << "TESTHP 24 NISO="<<niso<<G4endl;
256 for (G4int ix=0; ix<niso; ix++)
257 {
258 running += xsec[ix];
259 //if(random<=running/sum)
260 if( sum == 0 || random <= running/sum )
261 {
262 it = ix;
263 break;
264 }
265 }
266 if(it==niso) it--;
267 }
268 delete [] xsec;
269 G4HadFinalState * theFinalState=0;
270 while(theFinalState==0)
271 {
272// G4cout << "TESTHP 24 it="<<it<<G4endl;
273 theFinalState = theFinalStates[it]->ApplyYourself(theTrack);
274 }
275// G4cout <<"THE IMPORTANT RETURN"<<G4endl;
276 //G4cout << "TK G4NeutronHPChannel Elastic, Capture and Fission Cases "
277 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
280 return theFinalState;
281 }
double G4double
Definition: G4Types.hh:64
#define G4UniformRand()
Definition: Randomize.hh:53
const G4Material * GetMaterial() const
G4double GetTemperature() const
Definition: G4Material.hh:181
G4double GetN(G4int i)
G4bool HasAnyData(G4int isoNumber)
G4double GetZ(G4int i)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &)
G4double GetXsec(G4double energy)
static G4NeutronHPManager * GetInstance()
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)

Referenced by G4NeutronHPChannelList::ApplyYourself(), G4FissLib::ApplyYourself(), G4NeutronHPorLCapture::ApplyYourself(), G4NeutronHPorLElastic::ApplyYourself(), and G4NeutronHPorLFission::ApplyYourself().

◆ GetFSCrossSection()

G4double G4NeutronHPChannel::GetFSCrossSection ( G4double  energy,
G4int  isoNumber 
)

Definition at line 53 of file G4NeutronHPChannel.cc.

54 {
55 return theFinalStates[isoNumber]->GetXsec(energy);
56 }
virtual G4double GetXsec(G4double)

Referenced by G4NeutronHPChannelList::ApplyYourself().

◆ GetN()

G4double G4NeutronHPChannel::GetN ( G4int  i)
inline

Definition at line 129 of file G4NeutronHPChannel.hh.

129{return theFinalStates[i]->GetN();}

Referenced by G4NeutronHPChannelList::ApplyYourself(), and ApplyYourself().

◆ GetNiso()

G4int G4NeutronHPChannel::GetNiso ( )
inline

◆ GetWeightedXsec()

G4double G4NeutronHPChannel::GetWeightedXsec ( G4double  energy,
G4int  isoNumber 
)

Definition at line 48 of file G4NeutronHPChannel.cc.

49 {
50 return theIsotopeWiseData[isoNumber].GetXsec(energy);
51 }

Referenced by G4NeutronHPChannelList::ApplyYourself().

◆ GetXsec()

◆ GetZ()

G4double G4NeutronHPChannel::GetZ ( G4int  i)
inline

Definition at line 130 of file G4NeutronHPChannel.hh.

130{return theFinalStates[i]->GetZ();}

Referenced by G4NeutronHPChannelList::ApplyYourself(), ApplyYourself(), and Register().

◆ Harmonise()

void G4NeutronHPChannel::Harmonise ( G4NeutronHPVector *&  theStore,
G4NeutronHPVector theNew 
)

Definition at line 166 of file G4NeutronHPChannel.cc.

167 {
168 G4int s_tmp = 0, n=0, m_tmp=0;
169 G4NeutronHPVector * theMerge = new G4NeutronHPVector;
170 G4NeutronHPVector * anActive = theStore;
171 G4NeutronHPVector * aPassive = theNew;
172 G4NeutronHPVector * tmp;
173 G4int a = s_tmp, p = n, t;
174 while (a<anActive->GetVectorLength()&&p<aPassive->GetVectorLength())
175 {
176 if(anActive->GetEnergy(a) <= aPassive->GetEnergy(p))
177 {
178 G4double xa = anActive->GetEnergy(a);
179 theMerge->SetData(m_tmp, xa, anActive->GetXsec(a)+std::max(0., aPassive->GetXsec(xa)) );
180 m_tmp++;
181 a++;
182 G4double xp = aPassive->GetEnergy(p);
183 if( std::abs(std::abs(xp-xa)/xa)<0.001 )
184 {
185 p++;
186 }
187 } else {
188 tmp = anActive; t=a;
189 anActive = aPassive; a=p;
190 aPassive = tmp; p=t;
191 }
192 }
193 while (a!=anActive->GetVectorLength())
194 {
195 theMerge->SetData(m_tmp++, anActive->GetEnergy(a), anActive->GetXsec(a));
196 a++;
197 }
198 while (p!=aPassive->GetVectorLength())
199 {
200 if(std::abs(theMerge->GetEnergy(std::max(0,m_tmp-1))-aPassive->GetEnergy(p))/aPassive->GetEnergy(p)>0.001)
201 theMerge->SetData(m_tmp++, aPassive->GetEnergy(p), aPassive->GetXsec(p));
202 p++;
203 }
204 delete theStore;
205 theStore = theMerge;
206 }
G4int GetVectorLength() const
G4double GetEnergy(G4int i) const
void SetData(G4int i, G4double x, G4double y)

Referenced by UpdateData().

◆ HasAnyData()

G4bool G4NeutronHPChannel::HasAnyData ( G4int  isoNumber)
inline

Definition at line 111 of file G4NeutronHPChannel.hh.

111{ return theFinalStates[isoNumber]->HasAnyData(); }

Referenced by ApplyYourself(), HasDataInAnyFinalState(), and UpdateData().

◆ HasDataInAnyFinalState()

G4bool G4NeutronHPChannel::HasDataInAnyFinalState ( )
inline

Definition at line 132 of file G4NeutronHPChannel.hh.

133 {
134 G4bool result = false;
135 G4int i;
136 for(i=0; i<niso; i++)
137 {
138 if(theFinalStates[i]->HasAnyData()) result = true;
139 }
140 return result;
141 }
bool G4bool
Definition: G4Types.hh:67

Referenced by Register().

◆ HasFSData()

G4bool G4NeutronHPChannel::HasFSData ( G4int  isoNumber)
inline

Definition at line 109 of file G4NeutronHPChannel.hh.

109{ return theFinalStates[isoNumber]->HasFSData(); }

◆ Init() [1/2]

void G4NeutronHPChannel::Init ( G4Element theElement,
const G4String  dirName 
)

◆ Init() [2/2]

void G4NeutronHPChannel::Init ( G4Element theElement,
const G4String  dirName,
const G4String  fsType 
)

Definition at line 58 of file G4NeutronHPChannel.cc.

60 {
61 theFSType = aFSType;
62 Init(anElement, dirName);
63 }
void Init(G4Element *theElement, const G4String dirName)

◆ IsActive()

G4bool G4NeutronHPChannel::IsActive ( G4int  isoNumber)
inline

Definition at line 107 of file G4NeutronHPChannel.hh.

107{ return active[isoNumber]; }

◆ Register()

G4bool G4NeutronHPChannel::Register ( G4NeutronHPFinalState theFS)

Definition at line 71 of file G4NeutronHPChannel.cc.

72 {
73 registerCount++;
74 G4int Z = G4lrint(theElement->GetZ());
75
76 Z = Z-registerCount;
77 if ( registerCount > 5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material"); // for Elastic, Capture, Fission case
78 if ( Z < 1 ) return false;
79/*
80 if(registerCount<5)
81 {
82 Z = Z-registerCount;
83 }
84*/
85 //if(Z=theElement->GetZ()-5) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
86 // Bug fix by TK on behalf of AH
87 if ( Z <=theElement->GetZ()-5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
88 G4int count = 0;
89 if(registerCount==0) count = theElement->GetNumberOfIsotopes();
90 if(count == 0||registerCount!=0) count +=
91 theStableOnes.GetNumberOfIsotopes(Z);
92 niso = count;
93 delete [] theIsotopeWiseData;
94 theIsotopeWiseData = new G4NeutronHPIsoData [niso];
95 delete [] active;
96 active = new G4bool[niso];
97
98 delete [] theFinalStates;
99 theFinalStates = new G4NeutronHPFinalState * [niso];
100 delete theChannelData;
101 theChannelData = new G4NeutronHPVector;
102 for(G4int i=0; i<niso; i++)
103 {
104 theFinalStates[i] = theFS->New();
105 }
106 count = 0;
107 G4int nIsos = niso;
108 if(theElement->GetNumberOfIsotopes()!=0&&registerCount==0)
109 {
110 for (G4int i1=0; i1<nIsos; i1++)
111 {
112 // G4cout <<" Init: normal case"<<G4endl;
113 G4int A = theElement->GetIsotope(i1)->GetN();
114 G4int M = theElement->GetIsotope(i1)->Getm();
115 G4double frac = theElement->GetRelativeAbundanceVector()[i1]/perCent;
116 //theFinalStates[i1]->SetA_Z(A, Z);
117 //UpdateData(A, Z, count++, frac);
118 theFinalStates[i1]->SetA_Z(A, Z, M);
119 UpdateData(A, Z, M, count++, frac);
120 }
121 } else {
122 //G4cout <<" Init: mean case: "
123 // <<theStableOnes.GetNumberOfIsotopes(Z)<<" "
124 // <<Z<<" "<<theElement
125 // << G4endl;
126 G4int first = theStableOnes.GetFirstIsotope(Z);
127 for(G4int i1=0;
128 i1<theStableOnes.GetNumberOfIsotopes(Z);
129 i1++)
130 {
131 G4int A = theStableOnes.GetIsotopeNucleonCount(first+i1);
132 G4double frac = theStableOnes.GetAbundance(first+i1);
133 theFinalStates[i1]->SetA_Z(A, Z);
134 UpdateData(A, Z, count++, frac);
135 }
136 }
138 return result;
139 }
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
G4double GetZ() const
Definition: G4Element.hh:131
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4int Getm() const
Definition: G4Isotope.hh:100
G4int GetN() const
Definition: G4Isotope.hh:94
void UpdateData(G4int A, G4int Z, G4int index, G4double abundance)
virtual G4NeutronHPFinalState * New()=0
void SetA_Z(G4double anA, G4double aZ, G4int aM=0)
G4double GetAbundance(G4int number)
G4int GetFirstIsotope(G4int Z)
G4int GetNumberOfIsotopes(G4int Z)
G4int GetIsotopeNucleonCount(G4int number)
int G4lrint(double ad)
Definition: templates.hh:163

Referenced by G4FissLib::G4FissLib(), and G4NeutronHPChannelList::Register().

◆ UpdateData() [1/2]

void G4NeutronHPChannel::UpdateData ( G4int  A,
G4int  Z,
G4int  index,
G4double  abundance 
)
inline

Definition at line 120 of file G4NeutronHPChannel.hh.

120{ G4int M = 0; UpdateData( A, Z, M, index, abundance); };

Referenced by Register(), and UpdateData().

◆ UpdateData() [2/2]

void G4NeutronHPChannel::UpdateData ( G4int  A,
G4int  Z,
G4int  M,
G4int  index,
G4double  abundance 
)

Definition at line 142 of file G4NeutronHPChannel.cc.

143 {
144 //theFinalStates[index]->Init(A, Z, theDir, theFSType);
145 theFinalStates[index]->Init(A, Z, M, theDir, theFSType);
146 if(!theFinalStates[index]->HasAnyData()) return; // nothing there for exactly this isotope.
147
148 // the above has put the X-sec into the FS
149 theBuffer = 0;
150 if(theFinalStates[index]->HasXsec())
151 {
152 theBuffer = theFinalStates[index]->GetXsec();
153 theBuffer->Times(abundance/100.);
154 theIsotopeWiseData[index].FillChannelData(theBuffer);
155 }
156 else // get data from CrossSection directory
157 {
158 G4String tString = "/CrossSection";
159 //active[index] = theIsotopeWiseData[index].Init(A, Z, abundance, theDir, tString);
160 active[index] = theIsotopeWiseData[index].Init(A, Z, M, abundance, theDir, tString);
161 if(active[index]) theBuffer = theIsotopeWiseData[index].MakeChannelData();
162 }
163 if(theBuffer != 0) Harmonise(theChannelData, theBuffer);
164 }
void Harmonise(G4NeutronHPVector *&theStore, G4NeutronHPVector *theNew)
void Init(G4double A, G4double Z, G4String &dirName, G4String &aFSType)
G4NeutronHPVector * MakeChannelData()
G4bool Init(G4int A, G4int Z, G4double abun, G4String dirName, G4String aFSType)
void FillChannelData(G4NeutronHPVector *aBuffer)
void Times(G4double factor)

The documentation for this class was generated from the following files: