Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4NeutronHPFSFissionFS Class Reference

#include <G4NeutronHPFSFissionFS.hh>

Inheritance diagram for G4NeutronHPFSFissionFS:
G4NeutronHPFinalState

Public Member Functions

 G4NeutronHPFSFissionFS ()
 
 ~G4NeutronHPFSFissionFS ()
 
void Init (G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType)
 
G4DynamicParticleVectorApplyYourself (G4int Prompt, G4int delayed, G4double *decayconst)
 
G4NeutronHPFinalStateNew ()
 
G4double GetMass ()
 
void SampleNeutronMult (G4int &all, G4int &Prompt, G4int &delayed, G4double energy, G4int off)
 
void SetNeutron (const G4ReactionProduct &aNeutron)
 
void SetTarget (const G4ReactionProduct &aTarget)
 
G4DynamicParticleVectorGetPhotons ()
 
G4NeutronHPFissionEReleaseGetEnergyRelease ()
 
- Public Member Functions inherited from G4NeutronHPFinalState
 G4NeutronHPFinalState ()
 
virtual ~G4NeutronHPFinalState ()
 
void Init (G4double A, G4double Z, G4String &dirName, G4String &aFSType)
 
G4bool HasXsec ()
 
G4bool HasFSData ()
 
G4bool HasAnyData ()
 
virtual G4double GetXsec (G4double)
 
virtual G4NeutronHPVectorGetXsec ()
 
void SetA_Z (G4double anA, G4double aZ, G4int aM=0)
 
G4double GetZ ()
 
G4double GetN ()
 
G4int GetM ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4NeutronHPFinalState
void SetAZMs (G4double anA, G4double aZ, G4int aM, G4NeutronHPDataUsed used)
 
void adjust_final_state (G4LorentzVector)
 
G4bool DoNotAdjustFinalState ()
 
- Protected Attributes inherited from G4NeutronHPFinalState
G4bool hasXsec
 
G4bool hasFSData
 
G4bool hasAnyData
 
G4NeutronHPNames theNames
 
G4HadFinalState theResult
 
G4double theBaseA
 
G4double theBaseZ
 
G4int theBaseM
 
G4int theNDLDataZ
 
G4int theNDLDataA
 
G4int theNDLDataM
 

Detailed Description

Definition at line 43 of file G4NeutronHPFSFissionFS.hh.

Constructor & Destructor Documentation

G4NeutronHPFSFissionFS::G4NeutronHPFSFissionFS ( )
inline

Definition at line 47 of file G4NeutronHPFSFissionFS.hh.

References G4NeutronHPFinalState::hasXsec.

Referenced by New().

47 { hasXsec = true; }
G4NeutronHPFSFissionFS::~G4NeutronHPFSFissionFS ( )
inline

Definition at line 48 of file G4NeutronHPFSFissionFS.hh.

48 {}

Member Function Documentation

G4DynamicParticleVector * G4NeutronHPFSFissionFS::ApplyYourself ( G4int  Prompt,
G4int  delayed,
G4double decayconst 
)

Definition at line 101 of file G4NeutronHPFSFissionFS.cc.

References G4NeutronHPNeutronYield::GetDecayConstant(), G4ReactionProduct::GetKineticEnergy(), G4ReactionProduct::Lorentz(), G4Neutron::Neutron(), G4NeutronHPEnergyDistribution::Sample(), G4NeutronHPAngular::SampleAndUpdate(), G4ReactionProduct::SetDefinition(), G4DynamicParticle::SetDefinition(), G4ReactionProduct::SetKineticEnergy(), and G4DynamicParticle::SetMomentum().

Referenced by G4NeutronHPFissionFS::ApplyYourself().

103  {
104  G4int i;
106  G4ReactionProduct boosted;
107  boosted.Lorentz(theNeutron, theTarget);
108  G4double eKinetic = boosted.GetKineticEnergy();
109 
110 // Build neutrons
111  G4ReactionProduct * theNeutrons = new G4ReactionProduct[nPrompt+nDelayed];
112  for(i=0; i<nPrompt+nDelayed; i++)
113  {
114  theNeutrons[i].SetDefinition(G4Neutron::Neutron());
115  }
116 
117 // sample energies
118  G4int it, dummy;
119  G4double tempE;
120  for(i=0; i<nPrompt; i++)
121  {
122  tempE = thePromptNeutronEnDis.Sample(eKinetic, dummy); // energy distribution (file5) always in lab
123  theNeutrons[i].SetKineticEnergy(tempE);
124  }
125  for(i=nPrompt; i<nPrompt+nDelayed; i++)
126  {
127  theNeutrons[i].SetKineticEnergy(theDelayedNeutronEnDis.Sample(eKinetic, it)); // dito
128  if(it==0) theNeutrons[i].SetKineticEnergy(thePromptNeutronEnDis.Sample(eKinetic, dummy));
129  theDecayConst[i-nPrompt] = theFinalStateNeutrons.GetDecayConstant(it); // this is returned
130  }
131 
132 // sample neutron angular distribution
133  for(i=0; i<nPrompt+nDelayed; i++)
134  {
135  theNeutronAngularDis.SampleAndUpdate(theNeutrons[i]); // angular comes back in lab automatically
136  }
137 
138 // already in lab. Add neutrons to dynamic particle vector
139  for(i=0; i<nPrompt+nDelayed; i++)
140  {
142  dp->SetDefinition(theNeutrons[i].GetDefinition());
143  dp->SetMomentum(theNeutrons[i].GetMomentum());
144  aResult->push_back(dp);
145  }
146  delete [] theNeutrons;
147 // return the result
148  return aResult;
149  }
void SetMomentum(const G4ThreeVector &momentum)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
void SetKineticEnergy(const G4double en)
G4double Sample(G4double anEnergy, G4int &it)
int G4int
Definition: G4Types.hh:78
std::vector< G4DynamicParticle * > G4DynamicParticleVector
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double GetKineticEnergy() const
void SampleAndUpdate(G4ReactionProduct &aNeutron)
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4NeutronHPFissionERelease* G4NeutronHPFSFissionFS::GetEnergyRelease ( )
inline

Definition at line 82 of file G4NeutronHPFSFissionFS.hh.

Referenced by G4NeutronHPFissionFS::ApplyYourself().

83  {
84  return &theEnergyRelease;
85  }
G4double G4NeutronHPFSFissionFS::GetMass ( )
inline

Definition at line 60 of file G4NeutronHPFSFissionFS.hh.

Referenced by G4NeutronHPFissionFS::ApplyYourself().

60 { return targetMass; }
G4DynamicParticleVector * G4NeutronHPFSFissionFS::GetPhotons ( )

Definition at line 175 of file G4NeutronHPFSFissionFS.cc.

References G4ReactionProduct::GetKineticEnergy(), G4NeutronHPPhotonDist::GetPhotons(), G4ReactionProduct::Lorentz(), G4DynamicParticle::SetDefinition(), and G4DynamicParticle::SetMomentum().

Referenced by G4NeutronHPFissionFS::ApplyYourself().

176 {
177 // sample photons
179  G4ReactionProduct boosted;
180 // the photon distributions are in the Nucleus rest frame.
181  boosted.Lorentz(theNeutron, theTarget);
182  G4double anEnergy = boosted.GetKineticEnergy();
183  temp = theFinalStatePhotons.GetPhotons(anEnergy);
184  if(temp == 0) { return 0; }
185 
186 // lorentz transform, and add photons to final state
187  unsigned int i;
189  for(i=0; i<temp->size(); i++)
190  {
191  // back to lab
192  temp->operator[](i)->Lorentz(*(temp->operator[](i)), -1.*theTarget);
193  G4DynamicParticle * theOne = new G4DynamicParticle;
194  theOne->SetDefinition(temp->operator[](i)->GetDefinition());
195  theOne->SetMomentum(temp->operator[](i)->GetMomentum());
196  result->push_back(theOne);
197  delete temp->operator[](i);
198  }
199  delete temp;
200  return result;
201 }
G4ReactionProductVector * GetPhotons(G4double anEnergy)
void SetMomentum(const G4ThreeVector &momentum)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
std::vector< G4ReactionProduct * > G4ReactionProductVector
std::vector< G4DynamicParticle * > G4DynamicParticleVector
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void G4NeutronHPFSFissionFS::Init ( G4double  A,
G4double  Z,
G4int  M,
G4String dirName,
G4String aFSType 
)
virtual

Implements G4NeutronHPFinalState.

Definition at line 43 of file G4NeutronHPFSFissionFS.cc.

References G4cout, G4endl, G4NeutronHPManager::GetDataStream(), G4NeutronHPManager::GetInstance(), G4NeutronHPNames::GetName(), G4NeutronHPDataUsed::GetName(), G4NeutronHPNeutronYield::GetTargetMass(), G4NeutronHPFinalState::hasAnyData, G4NeutronHPFinalState::hasFSData, G4NeutronHPFinalState::hasXsec, G4NeutronHPFissionERelease::Init(), G4NeutronHPAngular::Init(), G4NeutronHPEnergyDistribution::Init(), G4NeutronHPPhotonDist::InitAngular(), G4NeutronHPNeutronYield::InitDelayed(), G4NeutronHPPhotonDist::InitEnergies(), G4NeutronHPNeutronYield::InitMean(), G4NeutronHPPhotonDist::InitMean(), G4NeutronHPNeutronYield::InitPrompt(), and G4NeutronHPFinalState::SetAZMs().

Referenced by G4NeutronHPFissionFS::Init().

44  {
45  G4String tString = "/FS/";
46  G4bool dbool;
47  G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, dirName, tString, dbool);
48  G4String filename = aFile.GetName();
49  SetAZMs( A, Z, M, aFile );
50  if(!dbool)
51  {
52  hasAnyData = false;
53  hasFSData = false;
54  hasXsec = false;
55  return;
56  }
57  //std::ifstream theData(filename, std::ios::in);
58  std::istringstream theData(std::ios::in);
60 
61  // here it comes
62  G4int infoType, dataType;
63  hasFSData = false;
64  while (theData >> infoType)
65  {
66  hasFSData = true;
67  theData >> dataType;
68  switch(infoType)
69  {
70  case 1:
71  if(dataType==4) theNeutronAngularDis.Init(theData);
72  if(dataType==5) thePromptNeutronEnDis.Init(theData);
73  if(dataType==12) theFinalStatePhotons.InitMean(theData);
74  if(dataType==14) theFinalStatePhotons.InitAngular(theData);
75  if(dataType==15) theFinalStatePhotons.InitEnergies(theData);
76  break;
77  case 2:
78  if(dataType==1) theFinalStateNeutrons.InitMean(theData);
79  break;
80  case 3:
81  if(dataType==1) theFinalStateNeutrons.InitDelayed(theData);
82  if(dataType==5) theDelayedNeutronEnDis.Init(theData);
83  break;
84  case 4:
85  if(dataType==1) theFinalStateNeutrons.InitPrompt(theData);
86  break;
87  case 5:
88  if(dataType==1) theEnergyRelease.Init(theData);
89  break;
90  default:
91  G4cout << "G4NeutronHPFSFissionFS::Init: unknown data type"<<dataType<<G4endl;
92  throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPFSFissionFS::Init: unknown data type");
93  break;
94  }
95  }
96  targetMass = theFinalStateNeutrons.GetTargetMass();
97  //theData.close();
98  }
void Init(std::istream &aDataFile)
static G4NeutronHPManager * GetInstance()
void InitDelayed(std::istream &aDataFile)
void GetDataStream(G4String, std::istringstream &iss)
G4NeutronHPDataUsed GetName(G4int A, G4int Z, G4String base, G4String rest, G4bool &active)
int G4int
Definition: G4Types.hh:78
void InitAngular(std::istream &aDataFile)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
void InitEnergies(std::istream &aDataFile)
void Init(std::istream &aDataFile)
G4bool InitMean(std::istream &aDataFile)
void InitPrompt(std::istream &aDataFile)
#define G4endl
Definition: G4ios.hh:61
void SetAZMs(G4double anA, G4double aZ, G4int aM, G4NeutronHPDataUsed used)
void InitMean(std::istream &aDataFile)
G4NeutronHPFinalState* G4NeutronHPFSFissionFS::New ( )
inlinevirtual

Implements G4NeutronHPFinalState.

Definition at line 54 of file G4NeutronHPFSFissionFS.hh.

References G4NeutronHPFSFissionFS().

void G4NeutronHPFSFissionFS::SampleNeutronMult ( G4int all,
G4int Prompt,
G4int delayed,
G4double  energy,
G4int  off 
)

Definition at line 151 of file G4NeutronHPFSFissionFS.cc.

References G4Poisson(), G4NeutronHPNeutronYield::GetDelayed(), G4NeutronHPNeutronYield::GetMean(), and G4NeutronHPNeutronYield::GetPrompt().

Referenced by G4NeutronHPFissionFS::ApplyYourself().

152 {
153  G4double promptNeutronMulti = 0;
154  promptNeutronMulti = theFinalStateNeutrons.GetPrompt(eKinetic);
155  G4double delayedNeutronMulti = 0;
156  delayedNeutronMulti = theFinalStateNeutrons.GetDelayed(eKinetic);
157 
158  if(delayedNeutronMulti==0&&promptNeutronMulti==0)
159  {
160  Prompt = 0;
161  delayed = 0;
162  G4double totalNeutronMulti = theFinalStateNeutrons.GetMean(eKinetic);
163  all = G4Poisson(totalNeutronMulti-off);
164  all += off;
165  }
166  else
167  {
168  Prompt = G4Poisson(promptNeutronMulti-off);
169  Prompt += off;
170  delayed = G4Poisson(delayedNeutronMulti);
171  all = Prompt+delayed;
172  }
173 }
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
G4double GetMean(G4double anEnergy)
G4double GetDelayed(G4double anEnergy)
double G4double
Definition: G4Types.hh:76
G4double GetPrompt(G4double anEnergy)
void G4NeutronHPFSFissionFS::SetNeutron ( const G4ReactionProduct aNeutron)
inline

Definition at line 68 of file G4NeutronHPFSFissionFS.hh.

References G4NeutronHPAngular::SetNeutron().

Referenced by G4NeutronHPFissionFS::ApplyYourself().

69  {
70  theNeutron = aNeutron;
71  theNeutronAngularDis.SetNeutron(aNeutron);
72  }
void SetNeutron(const G4ReactionProduct &aNeutron)
void G4NeutronHPFSFissionFS::SetTarget ( const G4ReactionProduct aTarget)
inline

Definition at line 74 of file G4NeutronHPFSFissionFS.hh.

References G4NeutronHPAngular::SetTarget().

Referenced by G4NeutronHPFissionFS::ApplyYourself().

75  {
76  theTarget = aTarget;
77  theNeutronAngularDis.SetTarget(aTarget);
78  }
void SetTarget(const G4ReactionProduct &aTarget)

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