Geant4-11
Public Member Functions | Private Attributes
G4ParticleHPChannel Class Reference

#include <G4ParticleHPChannel.hh>

Public Member Functions

G4HadFinalStateApplyYourself (const G4HadProjectile &theTrack, G4int isoNumber=-1)
 
void DumpInfo ()
 
 G4ParticleHPChannel ()
 
 G4ParticleHPChannel (G4ParticleDefinition *projectile)
 
G4ParticleHPFinalState ** GetFinalStates () const
 
G4double GetFSCrossSection (G4double energy, G4int isoNumber)
 
G4String GetFSType () const
 
G4double GetM (G4int i)
 
G4double GetN (G4int i)
 
G4int GetNiso ()
 
G4double GetWeightedXsec (G4double energy, G4int isoNumber)
 
G4double GetXsec (G4double energy)
 
G4double GetZ (G4int i)
 
void Harmonise (G4ParticleHPVector *&theStore, G4ParticleHPVector *theNew)
 
G4bool HasAnyData (G4int isoNumber)
 
G4bool HasDataInAnyFinalState ()
 
G4bool HasFSData (G4int isoNumber)
 
void Init (G4Element *theElement, const G4String dirName)
 
void Init (G4Element *theElement, const G4String dirName, const G4String fsType)
 
G4bool IsActive (G4int isoNumber)
 
G4bool Register (G4ParticleHPFinalState *theFS)
 
void UpdateData (G4int A, G4int Z, G4int index, G4double abundance, G4ParticleDefinition *projectile)
 
void UpdateData (G4int A, G4int Z, G4int M, G4int index, G4double abundance, G4ParticleDefinition *projectile)
 
 ~G4ParticleHPChannel ()
 

Private Attributes

G4boolactive
 
G4int niso
 
G4int registerCount
 
G4ParticleHPVectortheBuffer
 
G4ParticleHPVectortheChannelData
 
G4String theDir
 
G4ElementtheElement
 
G4ParticleHPFinalState ** theFinalStates
 
G4String theFSType
 
G4ParticleHPIsoDatatheIsotopeWiseData
 
G4ParticleDefinitiontheProjectile
 
G4StableIsotopes theStableOnes
 
G4WendtFissionFragmentGenerator *const wendtFissionGenerator
 

Detailed Description

Definition at line 55 of file G4ParticleHPChannel.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPChannel() [1/2]

G4ParticleHPChannel::G4ParticleHPChannel ( G4ParticleDefinition projectile)
inline

Definition at line 59 of file G4ParticleHPChannel.hh.

59 :
60 wendtFissionGenerator( G4ParticleHPManager::GetInstance()->GetUseWendtFissionModel() ?
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);
69 theBuffer = 0;
72 active = 0;
73 registerCount = -1;
74 niso = -1;
75 theElement = NULL;
76 }
G4WendtFissionFragmentGenerator *const wendtFissionGenerator
G4ParticleHPFinalState ** theFinalStates
G4ParticleDefinition * theProjectile
G4ParticleHPVector * theChannelData
G4ParticleHPIsoData * theIsotopeWiseData
G4ParticleHPVector * theBuffer
void SetProduceFissionFragments(G4bool val)
static G4ParticleHPManager * GetInstance()
static G4WendtFissionFragmentGenerator * GetInstance()

References active, G4ParticleHPManager::GetInstance(), niso, registerCount, G4ParticleHPManager::SetProduceFissionFragments(), theBuffer, theChannelData, theElement, theFinalStates, theIsotopeWiseData, and theProjectile.

◆ G4ParticleHPChannel() [2/2]

G4ParticleHPChannel::G4ParticleHPChannel ( )
inline

Definition at line 78 of file G4ParticleHPChannel.hh.

78 : wendtFissionGenerator( G4ParticleHPManager::GetInstance()->GetUseWendtFissionModel() ?
80 {
81 if ( G4ParticleHPManager::GetInstance()->GetUseWendtFissionModel() ) {
82 // Make sure both fission fragment models are not active at same time
84 }
87 theBuffer = 0;
90 active = 0;
91 registerCount = -1;
92 niso = -1;
93 theElement = NULL;
94 }
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103

References active, G4ParticleHPManager::GetInstance(), G4Neutron::Neutron(), niso, registerCount, G4ParticleHPManager::SetProduceFissionFragments(), theBuffer, theChannelData, theElement, theFinalStates, theIsotopeWiseData, and theProjectile.

◆ ~G4ParticleHPChannel()

G4ParticleHPChannel::~G4ParticleHPChannel ( )
inline

Definition at line 96 of file G4ParticleHPChannel.hh.

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 }
int G4int
Definition: G4Types.hh:85

References active, niso, theChannelData, theFinalStates, and theIsotopeWiseData.

Member Function Documentation

◆ ApplyYourself()

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

Definition at line 228 of file G4ParticleHPChannel.cc.

230 {
231// G4cout << "G4ParticleHPChannel::ApplyYourself+"<<niso<<G4endl;
232 if ( anIsotope != -1 && anIsotope != -2 )
233 {
234 //Inelastic Case
235 //G4cout << "G4ParticleHPChannel Inelastic Case"
236 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
239 return theFinalStates[anIsotope]->ApplyYourself(theTrack);
240 }
241 G4double sum=0;
242 G4int it=0;
243 G4double * xsec = new G4double[niso];
244 G4ParticleHPThermalBoost aThermalE;
245 for (G4int i=0; i<niso; i++)
246 {
247 if(theFinalStates[i]->HasAnyData())
248 {
249 xsec[i] = theIsotopeWiseData[i].GetXsec(aThermalE.GetThermalEnergy(theTrack,
250 theFinalStates[i]->GetN(),
251 theFinalStates[i]->GetZ(),
252 theTrack.GetMaterial()->GetTemperature()));
253 sum += xsec[i];
254 }
255 else
256 {
257 xsec[i]=0;
258 }
259 }
260 if(sum == 0)
261 {
262// G4cout << "G4ParticleHPChannel::ApplyYourself theFinalState->Initialize+"<<G4endl;
263// G4cout << "G4ParticleHPChannel::ApplyYourself theFinalState->Initialize-"<<G4endl;
264 it = static_cast<G4int>(niso*G4UniformRand());
265 }
266 else
267 {
268// G4cout << "Are we still here? "<<sum<<G4endl;
269// G4cout << "TESTHP 23 NISO="<<niso<<G4endl;
270 G4double random = G4UniformRand();
271 G4double running=0;
272// G4cout << "G4ParticleHPChannel::ApplyYourself Done the sum"<<niso<<G4endl;
273// G4cout << "TESTHP 24 NISO="<<niso<<G4endl;
274 for (G4int ix=0; ix<niso; ix++)
275 {
276 running += xsec[ix];
277 //if(random<=running/sum)
278 if( sum == 0 || random <= running/sum )
279 {
280 it = ix;
281 break;
282 }
283 }
284 if(it==niso) it--;
285 }
286 delete [] xsec;
287 G4HadFinalState * theFinalState=0;
288 const G4int A = (G4int)this->GetN(it);
289 const G4int Z = (G4int)this->GetZ(it);
290 const G4int M = (G4int)this->GetM(it);
291
292 //-2:Marker for Fission
293 if(wendtFissionGenerator&&anIsotope==-2)
294 {
295 theFinalState = wendtFissionGenerator->ApplyYourself(theTrack, Z, A);
296 }
297
298 // Use the standard procedure if the G4FissionFragmentGenerator model fails
299 if (!theFinalState)
300 {
301
302 G4int icounter=0;
303 G4int icounter_max=1024;
304 while(theFinalState==0) // Loop checking, 11.05.2015, T. Koi
305 {
306 icounter++;
307 if ( icounter > icounter_max ) {
308 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
309 break;
310 }
311// G4cout << "TESTHP 24 it="<<it<<G4endl;
312 theFinalState = theFinalStates[it]->ApplyYourself(theTrack);
313 }
314 }
315
316 //G4cout <<"THE IMPORTANT RETURN"<<G4endl;
317 //G4cout << "TK G4ParticleHPChannel Elastic, Capture and Fission Cases "
318 //<< " Z= " << this->GetZ(it) << " A = " << this->GetN(it) << G4endl;
322
323 return theFinalState;
324 }
#define M(row, col)
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
const G4Material * GetMaterial() const
G4double GetTemperature() const
Definition: G4Material.hh:178
G4bool HasAnyData(G4int isoNumber)
G4double GetM(G4int i)
G4double GetN(G4int i)
G4double GetZ(G4int i)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &)
G4double GetXsec(G4double energy)
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4HadFinalState * ApplyYourself(const G4HadProjectile &projectile, G4int Z, G4int A)

References A, G4ParticleHPFinalState::ApplyYourself(), G4WendtFissionFragmentGenerator::ApplyYourself(), G4cout, G4endl, G4UniformRand, G4ParticleHPManager::GetInstance(), GetM(), G4HadProjectile::GetMaterial(), G4ParticleHPFinalState::GetN(), GetN(), G4ParticleHPManager::GetReactionWhiteBoard(), G4Material::GetTemperature(), G4ParticleHPThermalBoost::GetThermalEnergy(), G4ParticleHPIsoData::GetXsec(), G4ParticleHPFinalState::GetZ(), GetZ(), HasAnyData(), M, niso, G4ParticleHPReactionWhiteBoard::SetTargA(), G4ParticleHPReactionWhiteBoard::SetTargM(), G4ParticleHPReactionWhiteBoard::SetTargZ(), theFinalStates, theIsotopeWiseData, wendtFissionGenerator, and Z.

Referenced by G4ParticleHPChannelList::ApplyYourself(), and G4FissLib::ApplyYourself().

◆ DumpInfo()

void G4ParticleHPChannel::DumpInfo ( )

Definition at line 327 of file G4ParticleHPChannel.cc.

327 {
328
329 G4cout<<" Element: "<<theElement->GetName()<<G4endl;
330 G4cout<<" Directory name: "<<theDir<<G4endl;
331 G4cout<<" FS name: "<<theFSType<<G4endl;
332 G4cout<<" Number of Isotopes: "<<niso<<G4endl;
333 G4cout<<" Have cross sections: "<<G4endl;
334 for(int i=0;i<niso;i++){
335 G4cout<<theFinalStates[i]->HasXsec()<<" ";
336 }
337 G4cout<<G4endl;
338 if(theChannelData){
339 G4cout<<" Cross Section (total for this channel):"<<G4endl;
341 G4cout<<np<<G4endl;
342 for(int i=0;i<np;i++){
344 }
345 }
346
347}
static constexpr double eV
Definition: G4SIunits.hh:201
const G4String & GetName() const
Definition: G4Element.hh:127
G4double GetXsec(G4int i)
G4double GetEnergy(G4int i) const
G4int GetVectorLength() const

References eV, G4cout, G4endl, G4ParticleHPVector::GetEnergy(), G4Element::GetName(), G4ParticleHPVector::GetVectorLength(), G4ParticleHPVector::GetXsec(), G4ParticleHPFinalState::HasXsec(), niso, theChannelData, theDir, theElement, theFinalStates, and theFSType.

Referenced by G4ParticleHPChannelList::DumpInfo().

◆ GetFinalStates()

G4ParticleHPFinalState ** G4ParticleHPChannel::GetFinalStates ( ) const
inline

Definition at line 177 of file G4ParticleHPChannel.hh.

177 {
178 return theFinalStates;
179 }

References theFinalStates.

◆ GetFSCrossSection()

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

Definition at line 61 of file G4ParticleHPChannel.cc.

62{
63 return theFinalStates[isoNumber]->GetXsec(energy);
64}
virtual G4double GetXsec(G4double)
G4double energy(const ThreeVector &p, const G4double m)

References G4INCL::KinematicsUtils::energy(), G4ParticleHPFinalState::GetXsec(), and theFinalStates.

Referenced by G4ParticleHPChannelList::ApplyYourself().

◆ GetFSType()

G4String G4ParticleHPChannel::GetFSType ( ) const
inline

Definition at line 173 of file G4ParticleHPChannel.hh.

173 {
174 return theFSType;
175 }

References theFSType.

◆ GetM()

G4double G4ParticleHPChannel::GetM ( G4int  i)
inline

Definition at line 158 of file G4ParticleHPChannel.hh.

References G4ParticleHPFinalState::GetM(), and theFinalStates.

Referenced by ApplyYourself().

◆ GetN()

G4double G4ParticleHPChannel::GetN ( G4int  i)
inline

◆ GetNiso()

G4int G4ParticleHPChannel::GetNiso ( )
inline

Definition at line 154 of file G4ParticleHPChannel.hh.

154{return niso;}

References niso.

Referenced by G4ParticleHPChannelList::ApplyYourself().

◆ GetWeightedXsec()

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

◆ GetXsec()

G4double G4ParticleHPChannel::GetXsec ( G4double  energy)

Definition at line 51 of file G4ParticleHPChannel.cc.

52{
54}
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References G4INCL::KinematicsUtils::energy(), G4ParticleHPVector::GetXsec(), G4INCL::Math::max(), and theChannelData.

Referenced by G4FissLib::ApplyYourself().

◆ GetZ()

G4double G4ParticleHPChannel::GetZ ( G4int  i)
inline

◆ Harmonise()

void G4ParticleHPChannel::Harmonise ( G4ParticleHPVector *&  theStore,
G4ParticleHPVector theNew 
)

Definition at line 184 of file G4ParticleHPChannel.cc.

185 {
186 G4int s_tmp = 0, n=0, m_tmp=0;
187 G4ParticleHPVector * theMerge = new G4ParticleHPVector;
188 G4ParticleHPVector * anActive = theStore;
189 G4ParticleHPVector * aPassive = theNew;
190 G4ParticleHPVector * tmp;
191 G4int a = s_tmp, p = n, t;
192 while (a<anActive->GetVectorLength()&&p<aPassive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
193 {
194 if(anActive->GetEnergy(a) <= aPassive->GetEnergy(p))
195 {
196 G4double xa = anActive->GetEnergy(a);
197 theMerge->SetData(m_tmp, xa, anActive->GetXsec(a)+std::max(0., aPassive->GetXsec(xa)) );
198 m_tmp++;
199 a++;
200 G4double xp = aPassive->GetEnergy(p);
201 if( std::abs(std::abs(xp-xa)/xa)<0.001 )
202 {
203 p++;
204 }
205 } else {
206 tmp = anActive; t=a;
207 anActive = aPassive; a=p;
208 aPassive = tmp; p=t;
209 }
210 }
211 while (a!=anActive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
212 {
213 theMerge->SetData(m_tmp++, anActive->GetEnergy(a), anActive->GetXsec(a));
214 a++;
215 }
216 while (p!=aPassive->GetVectorLength()) // Loop checking, 11.05.2015, T. Koi
217 {
218 if(std::abs(theMerge->GetEnergy(std::max(0,m_tmp-1))-aPassive->GetEnergy(p))/aPassive->GetEnergy(p)>0.001)
219 theMerge->SetData(m_tmp++, aPassive->GetEnergy(p), aPassive->GetXsec(p));
220 p++;
221 }
222 delete theStore;
223 theStore = theMerge;
224 }
void SetData(G4int i, G4double x, G4double y)

References G4ParticleHPVector::GetEnergy(), G4ParticleHPVector::GetVectorLength(), G4ParticleHPVector::GetXsec(), G4INCL::Math::max(), CLHEP::detail::n, and G4ParticleHPVector::SetData().

Referenced by UpdateData().

◆ HasAnyData()

G4bool G4ParticleHPChannel::HasAnyData ( G4int  isoNumber)
inline

◆ HasDataInAnyFinalState()

G4bool G4ParticleHPChannel::HasDataInAnyFinalState ( )
inline

Definition at line 160 of file G4ParticleHPChannel.hh.

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 }
bool G4bool
Definition: G4Types.hh:86

References HasAnyData(), niso, and theFinalStates.

Referenced by Register().

◆ HasFSData()

G4bool G4ParticleHPChannel::HasFSData ( G4int  isoNumber)
inline

Definition at line 136 of file G4ParticleHPChannel.hh.

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

References G4ParticleHPFinalState::HasFSData(), and theFinalStates.

◆ Init() [1/2]

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

Definition at line 73 of file G4ParticleHPChannel.cc.

74 {
75 theDir = dirName;
76 theElement = anElement;
77 }

References theDir, and theElement.

Referenced by G4FissLib::G4FissLib(), Init(), and G4ParticleHPChannelList::Register().

◆ Init() [2/2]

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

Definition at line 66 of file G4ParticleHPChannel.cc.

68 {
69 theFSType = aFSType;
70 Init(anElement, dirName);
71 }
void Init(G4Element *theElement, const G4String dirName)

References Init(), and theFSType.

◆ IsActive()

G4bool G4ParticleHPChannel::IsActive ( G4int  isoNumber)
inline

Definition at line 134 of file G4ParticleHPChannel.hh.

134{ return active[isoNumber]; }

References active.

◆ Register()

G4bool G4ParticleHPChannel::Register ( G4ParticleHPFinalState theFS)

Definition at line 79 of file G4ParticleHPChannel.cc.

80 {
83
85 if ( registerCount > 5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material"); // for Elastic, Capture, Fission case
86 if ( Z < 1 ) return false;
87/*
88 if(registerCount<5)
89 {
90 Z = Z-registerCount;
91 }
92*/
93 //if(Z=theElement->GetZ()-5) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
94 // Bug fix by TK on behalf of AH
95 //if ( Z <=theElement->GetZ()-5 ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this material");
96 G4int count = 0;
98 if(count == 0||registerCount!=0) count +=
100 niso = count;
101 delete [] theIsotopeWiseData;
103 delete [] active;
104 active = new G4bool[niso];
105
106 delete [] theFinalStates;
108 delete theChannelData;
110 for(G4int i=0; i<niso; i++)
111 {
112 theFinalStates[i] = theFS->New();
114 }
115 count = 0;
116 G4int nIsos = niso;
118 {
119 for (G4int i1=0; i1<nIsos; i1++)
120 {
121 // G4cout <<" Init: normal case"<<G4endl;
122 G4int A = theElement->GetIsotope(i1)->GetN();
123 G4int M = theElement->GetIsotope(i1)->Getm();
125 //theFinalStates[i1]->SetA_Z(A, Z);
126 //UpdateData(A, Z, count++, frac);
127 theFinalStates[i1]->SetA_Z(A, Z, M);
128 UpdateData(A, Z, M, count++, frac, theProjectile);
129 }
130 } else {
131 //G4cout <<" Init: mean case: "
132 // <<theStableOnes.GetNumberOfIsotopes(Z)<<" "
133 // <<Z<<" "<<theElement
134 // << G4endl;
136 for(G4int i1=0;
138 i1++)
139 {
141 G4double frac = theStableOnes.GetAbundance(first+i1);
142 theFinalStates[i1]->SetA_Z(A, Z);
143 UpdateData(A, Z, count++, frac, theProjectile);
144 }
145 }
147
148 //To avoid issuing hash by worker threads
149 if ( result ) theChannelData->Hash();
150
151 return result;
152 }
static constexpr double perCent
Definition: G4SIunits.hh:325
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
G4double GetZ() const
Definition: G4Element.hh:131
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4int Getm() const
Definition: G4Isotope.hh:99
G4int GetN() const
Definition: G4Isotope.hh:93
void UpdateData(G4int A, G4int Z, G4int index, G4double abundance, G4ParticleDefinition *projectile)
G4StableIsotopes theStableOnes
void SetA_Z(G4double anA, G4double aZ, G4int aM=0)
virtual G4ParticleHPFinalState * New()=0
void SetProjectile(G4ParticleDefinition *projectile)
G4double GetAbundance(G4int number)
G4int GetFirstIsotope(G4int Z)
G4int GetNumberOfIsotopes(G4int Z)
G4int GetIsotopeNucleonCount(G4int number)
int G4lrint(double ad)
Definition: templates.hh:134

References A, active, G4lrint(), G4StableIsotopes::GetAbundance(), G4StableIsotopes::GetFirstIsotope(), G4Element::GetIsotope(), G4StableIsotopes::GetIsotopeNucleonCount(), G4Isotope::Getm(), G4Isotope::GetN(), G4Element::GetNumberOfIsotopes(), G4StableIsotopes::GetNumberOfIsotopes(), G4Element::GetRelativeAbundanceVector(), G4Element::GetZ(), HasDataInAnyFinalState(), G4ParticleHPVector::Hash(), M, G4ParticleHPFinalState::New(), niso, perCent, registerCount, G4ParticleHPFinalState::SetA_Z(), G4ParticleHPFinalState::SetProjectile(), theChannelData, theElement, theFinalStates, theIsotopeWiseData, theProjectile, theStableOnes, UpdateData(), and Z.

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

◆ UpdateData() [1/2]

void G4ParticleHPChannel::UpdateData ( G4int  A,
G4int  Z,
G4int  index,
G4double  abundance,
G4ParticleDefinition projectile 
)
inline

Definition at line 147 of file G4ParticleHPChannel.hh.

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

References A, M, UpdateData(), and Z.

Referenced by Register(), and UpdateData().

◆ UpdateData() [2/2]

void G4ParticleHPChannel::UpdateData ( G4int  A,
G4int  Z,
G4int  M,
G4int  index,
G4double  abundance,
G4ParticleDefinition projectile 
)

Definition at line 155 of file G4ParticleHPChannel.cc.

156 {
157 // Initialze the G4FissionFragment generator for this isomer if needed
159 {
161 }
162
163 theFinalStates[index]->Init(A, Z, M, theDir, theFSType, projectile);
164 if(!theFinalStates[index]->HasAnyData()) return; // nothing there for exactly this isotope.
165
166 // the above has put the X-sec into the FS
167 theBuffer = 0;
168 if(theFinalStates[index]->HasXsec())
169 {
170 theBuffer = theFinalStates[index]->GetXsec();
171 theBuffer->Times(abundance/100.);
173 }
174 else // get data from CrossSection directory
175 {
176 G4String tString = "/CrossSection";
177 //active[index] = theIsotopeWiseData[index].Init(A, Z, abundance, theDir, tString);
178 active[index] = theIsotopeWiseData[index].Init(A, Z, M, abundance, theDir, tString);
180 }
182 }
void Harmonise(G4ParticleHPVector *&theStore, G4ParticleHPVector *theNew)
void Init(G4double A, G4double Z, G4String &dirName, G4String &aFSType, G4ParticleDefinition *projectile)
G4ParticleHPVector * MakeChannelData()
G4bool Init(G4int A, G4int Z, G4double abun, G4String dirName, G4String aFSType)
void FillChannelData(G4ParticleHPVector *aBuffer)
void Times(G4double factor)
void InitializeANucleus(const G4int A, const G4int Z, const G4int M, const G4String &dataDirectory)

References A, active, G4ParticleHPIsoData::FillChannelData(), G4ParticleHPFinalState::GetXsec(), Harmonise(), HasAnyData(), G4ParticleHPFinalState::Init(), G4ParticleHPIsoData::Init(), G4WendtFissionFragmentGenerator::InitializeANucleus(), M, G4ParticleHPIsoData::MakeChannelData(), theBuffer, theChannelData, theDir, theFinalStates, theFSType, theIsotopeWiseData, G4ParticleHPVector::Times(), wendtFissionGenerator, and Z.

Field Documentation

◆ active

G4bool* G4ParticleHPChannel::active
private

◆ niso

G4int G4ParticleHPChannel::niso
private

◆ registerCount

G4int G4ParticleHPChannel::registerCount
private

Definition at line 198 of file G4ParticleHPChannel.hh.

Referenced by G4ParticleHPChannel(), and Register().

◆ theBuffer

G4ParticleHPVector* G4ParticleHPChannel::theBuffer
private

Definition at line 185 of file G4ParticleHPChannel.hh.

Referenced by G4ParticleHPChannel(), and UpdateData().

◆ theChannelData

G4ParticleHPVector* G4ParticleHPChannel::theChannelData
private

◆ theDir

G4String G4ParticleHPChannel::theDir
private

Definition at line 194 of file G4ParticleHPChannel.hh.

Referenced by DumpInfo(), Init(), and UpdateData().

◆ theElement

G4Element* G4ParticleHPChannel::theElement
private

Definition at line 196 of file G4ParticleHPChannel.hh.

Referenced by DumpInfo(), G4ParticleHPChannel(), Init(), and Register().

◆ theFinalStates

G4ParticleHPFinalState** G4ParticleHPChannel::theFinalStates
private

◆ theFSType

G4String G4ParticleHPChannel::theFSType
private

Definition at line 195 of file G4ParticleHPChannel.hh.

Referenced by DumpInfo(), GetFSType(), Init(), and UpdateData().

◆ theIsotopeWiseData

G4ParticleHPIsoData* G4ParticleHPChannel::theIsotopeWiseData
private

◆ theProjectile

G4ParticleDefinition* G4ParticleHPChannel::theProjectile
private

Definition at line 182 of file G4ParticleHPChannel.hh.

Referenced by G4ParticleHPChannel(), and Register().

◆ theStableOnes

G4StableIsotopes G4ParticleHPChannel::theStableOnes
private

Definition at line 192 of file G4ParticleHPChannel.hh.

Referenced by Register().

◆ wendtFissionGenerator

G4WendtFissionFragmentGenerator* const G4ParticleHPChannel::wendtFissionGenerator
private

Definition at line 200 of file G4ParticleHPChannel.hh.

Referenced by ApplyYourself(), and UpdateData().


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