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

#include <G4FissionLibrary.hh>

Inheritance diagram for G4FissionLibrary:
G4NeutronHPFinalState

Public Member Functions

 G4FissionLibrary ()
 
 ~G4FissionLibrary ()
 
void Init (G4double A, G4double Z, G4int M, G4String &dirName, G4String &)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &theTrack)
 
G4NeutronHPFinalStateNew ()
 
- 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 76 of file G4FissionLibrary.hh.

Constructor & Destructor Documentation

G4FissionLibrary::G4FissionLibrary ( )

Definition at line 66 of file G4FissionLibrary.cc.

References G4NeutronHPFinalState::hasXsec.

Referenced by New().

67  : G4NeutronHPFinalState(), theIsotope(0), targetMass(0.0)
68 {
69  hasXsec = false;
70 }
G4FissionLibrary::~G4FissionLibrary ( )

Definition at line 72 of file G4FissionLibrary.cc.

73 {}

Member Function Documentation

G4HadFinalState * G4FissionLibrary::ApplyYourself ( const G4HadProjectile theTrack)
virtual

Reimplemented from G4NeutronHPFinalState.

Definition at line 138 of file G4FissionLibrary.cc.

References G4HadFinalState::AddSecondary(), G4HadFinalState::Clear(), G4Gamma::Gamma(), G4HadProjectile::Get4Momentum(), G4Nucleus::GetBiasedThermalNucleus(), G4HadProjectile::GetDefinition(), G4ReactionProduct::GetDefinition(), G4NeutronHPFissionERelease::GetFragmentKinetic(), G4HadProjectile::GetKineticEnergy(), G4ReactionProduct::GetKineticEnergy(), G4HadProjectile::GetMaterial(), G4ReactionProduct::GetMomentum(), getndircosu_(), getndircosv_(), getndircosw_(), getneng_(), G4ParticleDefinition::GetPDGMass(), getpdircosu_(), getpdircosv_(), getpdircosw_(), getpeng_(), G4Material::GetTemperature(), G4ReactionProduct::GetTotalMomentum(), G4DynamicParticle::GetTotalMomentum(), G4ReactionProduct::Lorentz(), python.hepunit::MeV, G4Neutron::Neutron(), G4ReactionProduct::SetDefinition(), G4DynamicParticle::SetDefinition(), G4ReactionProduct::SetKineticEnergy(), G4DynamicParticle::SetKineticEnergy(), G4HadFinalState::SetLocalEnergyDeposit(), G4ReactionProduct::SetMomentum(), G4DynamicParticle::SetMomentum(), G4NeutronHPAngular::SetNeutron(), G4HadFinalState::SetStatusChange(), G4NeutronHPAngular::SetTarget(), stopAndKill, G4NeutronHPFinalState::theResult, theTarget, and CLHEP::HepLorentzVector::vect().

139 {
140  theResult.Clear();
141 
142 // prepare neutron
143  G4double eKinetic = theTrack.GetKineticEnergy();
144  const G4HadProjectile *incidentParticle = &theTrack;
145  G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
146  theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
147  theNeutron.SetKineticEnergy( eKinetic );
148 
149 // prepare target
150  G4Nucleus aNucleus;
152  G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum();
153  theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
154 
155 // set neutron and target in the FS classes
156  theNeutronAngularDis.SetNeutron(theNeutron);
157  theNeutronAngularDis.SetTarget(theTarget);
158 
159 // boost to target rest system
160  theNeutron.Lorentz(theNeutron, -1*theTarget);
161 
162  eKinetic = theNeutron.GetKineticEnergy();
163 
164 // dice neutron and gamma multiplicities, energies and momenta in Lab. @@
165 // no energy conservation on an event-to-event basis. we rely on the data to be ok. @@
166 // also for mean, we rely on the consistency of the data. @@
167 
168  G4int nPrompt=0, gPrompt=0;
169  SampleMult(theTrack, &nPrompt, &gPrompt, eKinetic);
170 
171 // Build neutrons and add them to dynamic particle vector
172  G4double momentum;
173  for(G4int i=0; i<nPrompt; i++)
174  {
177  it->SetKineticEnergy(getneng_(&i)*MeV);
178  momentum = it->GetTotalMomentum();
179  G4ThreeVector temp(momentum*getndircosu_(&i),
180  momentum*getndircosv_(&i),
181  momentum*getndircosw_(&i));
182  it->SetMomentum( temp );
183 // it->SetGlobalTime(getnage_(&i)*second);
185 // G4cout <<"G4FissionLibrary::ApplyYourself: energy of prompt neutron " << i << " = " << it->GetKineticEnergy()<<G4endl;
186  }
187 
188 // Build gammas, lorentz transform them, and add them to dynamic particle vector
189  for(G4int i=0; i<gPrompt; i++)
190  {
191  G4ReactionProduct * thePhoton = new G4ReactionProduct;
192  thePhoton->SetDefinition(G4Gamma::Gamma());
193  thePhoton->SetKineticEnergy(getpeng_(&i)*MeV);
194  momentum = thePhoton->GetTotalMomentum();
195  G4ThreeVector temp(momentum*getpdircosu_(&i),
196  momentum*getpdircosv_(&i),
197  momentum*getpdircosw_(&i));
198  thePhoton->SetMomentum( temp );
199  thePhoton->Lorentz(*thePhoton, -1.*theTarget);
200 
202  it->SetDefinition(thePhoton->GetDefinition());
203  it->SetMomentum(thePhoton->GetMomentum());
204 // it->SetGlobalTime(getpage_(&i)*second);
205 // G4cout <<"G4FissionLibrary::ApplyYourself: energy of prompt photon " << i << " = " << it->GetKineticEnergy()<<G4endl;
207  delete thePhoton;
208  }
209 // G4cout <<"G4FissionLibrary::ApplyYourself: Number of secondaries = "<<theResult.GetNumberOfSecondaries()<< G4endl;
210 // G4cout <<"G4FissionLibrary::ApplyYourself: Number of induced prompt neutron = "<<nPrompt<<G4endl;
211 // G4cout <<"G4FissionLibrary::ApplyYourself: Number of induced prompt photons = "<<gPrompt<<G4endl;
212 
213 // finally deal with local energy depositions.
214  G4double eDepByFragments = theEnergyRelease.GetFragmentKinetic();
215  theResult.SetLocalEnergyDeposit(eDepByFragments);
216 // G4cout << "G4FissionLibrary::local energy deposit" << eDepByFragments<<G4endl;
217 // clean up the primary neutron
219  return &theResult;
220 }
void SetNeutron(const G4ReactionProduct &aNeutron)
G4double getndircosu_(G4int *index)
void SetMomentum(const G4ThreeVector &momentum)
G4double GetTotalMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double getneng_(G4int *index)
void SetKineticEnergy(const G4double en)
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTarget(const G4ReactionProduct &aTarget)
int G4int
Definition: G4Types.hh:78
G4ParticleDefinition * GetDefinition() const
G4double GetTotalMomentum() const
void SetStatusChange(G4HadFinalStateStatus aS)
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:108
Hep3Vector vect() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ErrorTarget * theTarget
Definition: errprop.cc:59
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4LorentzVector & Get4Momentum() const
void SetKineticEnergy(G4double aEnergy)
G4double getpdircosw_(G4int *index)
G4double getpeng_(G4int *index)
G4double getpdircosu_(G4int *index)
G4double GetPDGMass() const
G4double getpdircosv_(G4int *index)
void SetDefinition(G4ParticleDefinition *aParticleDefinition)
G4double getndircosv_(G4int *index)
G4ThreeVector GetMomentum() const
G4double GetTemperature() const
Definition: G4Material.hh:180
void SetLocalEnergyDeposit(G4double aE)
const G4Material * GetMaterial() const
double G4double
Definition: G4Types.hh:76
G4double getndircosw_(G4int *index)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void AddSecondary(G4DynamicParticle *aP)
void G4FissionLibrary::Init ( G4double  A,
G4double  Z,
G4int  M,
G4String dirName,
G4String  
)
virtual

Implements G4NeutronHPFinalState.

Definition at line 82 of file G4FissionLibrary.cc.

References G4cout, G4endl, 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::theNames.

83 {
84  G4String tString = "/FS/";
85  G4bool dbool;
86  theIsotope = static_cast<G4int>(1000*Z+A);
87  G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, dirName, tString, dbool);
88  G4String filename = aFile.GetName();
89 
90  if(!dbool)
91  {
92  hasAnyData = false;
93  hasFSData = false;
94  hasXsec = false;
95  return;
96  }
97  std::ifstream theData(filename, std::ios::in);
98 
99  // here it comes
100  G4int infoType, dataType;
101  hasFSData = false;
102  while (theData >> infoType)
103  {
104  hasFSData = true;
105  theData >> dataType;
106  switch(infoType)
107  {
108  case 1:
109  if(dataType==4) theNeutronAngularDis.Init(theData);
110  if(dataType==5) thePromptNeutronEnDis.Init(theData);
111  if(dataType==12) theFinalStatePhotons.InitMean(theData);
112  if(dataType==14) theFinalStatePhotons.InitAngular(theData);
113  if(dataType==15) theFinalStatePhotons.InitEnergies(theData);
114  break;
115  case 2:
116  if(dataType==1) theFinalStateNeutrons.InitMean(theData);
117  break;
118  case 3:
119  if(dataType==1) theFinalStateNeutrons.InitDelayed(theData);
120  if(dataType==5) theDelayedNeutronEnDis.Init(theData);
121  break;
122  case 4:
123  if(dataType==1) theFinalStateNeutrons.InitPrompt(theData);
124  break;
125  case 5:
126  if(dataType==1) theEnergyRelease.Init(theData);
127  break;
128  default:
129  G4cout << "G4FissionLibrary::Init: unknown data type"<<dataType<<G4endl;
130  throw G4HadronicException(__FILE__, __LINE__, "G4FissionLibrary::Init: unknown data type");
131  break;
132  }
133  }
134  targetMass = theFinalStateNeutrons.GetTargetMass();
135  theData.close();
136 }
void Init(std::istream &aDataFile)
void InitDelayed(std::istream &aDataFile)
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 InitMean(std::istream &aDataFile)
G4NeutronHPFinalState * G4FissionLibrary::New ( )
virtual

Implements G4NeutronHPFinalState.

Definition at line 75 of file G4FissionLibrary.cc.

References G4FissionLibrary().

76 {
77  G4FissionLibrary * theNew = new G4FissionLibrary;
78  return theNew;
79 }

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