Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronHPFissionFS.cc
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 // neutron_hp -- source file
27 // J.P. Wellisch, Nov-1996
28 // A prototype of the low energy neutron transport model.
29 //
30 // 12-Apr-06 fix in delayed neutron and photon emission without FS data by T. Koi
31 // 07-Sep-11 M. Kelsey -- Follow change to G4HadFinalState interface
32 
33 #include "G4NeutronHPFissionFS.hh"
34 #include "G4PhysicalConstants.hh"
35 #include "G4Nucleus.hh"
38 #include "G4IonTable.hh"
39 
40  void G4NeutronHPFissionFS::Init (G4double A, G4double Z, G4int M, G4String & dirName, G4String & aFSType)
41  {
42  //G4cout << "G4NeutronHPFissionFS::Init " << A << " " << Z << " " << M << G4endl;
43  theFS.Init(A, Z, M, dirName, aFSType);
44  theFC.Init(A, Z, M, dirName, aFSType);
45  theSC.Init(A, Z, M, dirName, aFSType);
46  theTC.Init(A, Z, M, dirName, aFSType);
47  theLC.Init(A, Z, M, dirName, aFSType);
48 
49  theFF.Init(A, Z, M, dirName, aFSType);
50  if ( getenv("G4NEUTRONHP_PRODUCE_FISSION_FRAGMENTS") && theFF.HasFSData() )
51  {
52  G4cout << "Activate Fission Fragments Production for the target isotope of "
53  << "Z = " << (G4int)Z
54  << ", A = " << (G4int)A
55  //<< "M = " << M
56  << G4endl;
57  G4cout << "As the result, delayed neutrons are omitted and they should be taken care by RadioaActiveDecay."
58  << G4endl;
59  produceFissionFragments = true;
60  }
61  }
63  {
64  //G4cout << "G4NeutronHPFissionFS::ApplyYourself " << G4endl;
65 // prepare neutron
66  theResult.Clear();
67  G4double eKinetic = theTrack.GetKineticEnergy();
68  const G4HadProjectile *incidentParticle = &theTrack;
69  G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
70  theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
71  theNeutron.SetKineticEnergy( eKinetic );
72 
73 // prepare target
74  G4Nucleus aNucleus;
76  G4double targetMass = theFS.GetMass();
77  G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum();
78  theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
79 
80 // set neutron and target in the FS classes
81  theFS.SetNeutron(theNeutron);
82  theFS.SetTarget(theTarget);
83  theFC.SetNeutron(theNeutron);
84  theFC.SetTarget(theTarget);
85  theSC.SetNeutron(theNeutron);
86  theSC.SetTarget(theTarget);
87  theTC.SetNeutron(theNeutron);
88  theTC.SetTarget(theTarget);
89  theLC.SetNeutron(theNeutron);
90  theLC.SetTarget(theTarget);
91 
92 
93  theFF.SetNeutron(theNeutron);
94  theFF.SetTarget(theTarget);
95 
96 //TKWORK 120531
97 //G4cout << theTarget.GetDefinition() << G4endl; this should be NULL
98 //G4cout << "Z = " << theBaseZ << ", A = " << theBaseA << ", M = " << theBaseM << G4endl;
99 // theNDLDataZ,A,M should be filled in each FS (theFS, theFC, theSC, theTC, theLC and theFF)
100 ////G4cout << "Z = " << theNDLDataZ << ", A = " << theNDLDataA << ", M = " << theNDLDataM << G4endl;
101 
102 // boost to target rest system and decide on channel.
103  theNeutron.Lorentz(theNeutron, -1*theTarget);
104 
105 // dice the photons
106 
107  G4DynamicParticleVector * thePhotons;
108  thePhotons = theFS.GetPhotons();
109 
110 // select the FS in charge
111 
112  eKinetic = theNeutron.GetKineticEnergy();
113  G4double xSec[4];
114  xSec[0] = theFC.GetXsec(eKinetic);
115  xSec[1] = xSec[0]+theSC.GetXsec(eKinetic);
116  xSec[2] = xSec[1]+theTC.GetXsec(eKinetic);
117  xSec[3] = xSec[2]+theLC.GetXsec(eKinetic);
118  G4int it;
119  unsigned int i=0;
120  G4double random = G4UniformRand();
121  if(xSec[3]==0)
122  {
123  it=-1;
124  }
125  else
126  {
127  for(i=0; i<4; i++)
128  {
129  it =i;
130  if(random<xSec[i]/xSec[3]) break;
131  }
132  }
133 
134 // dice neutron multiplicities, energies and momenta in Lab. @@
135 // no energy conservation on an event-to-event basis. we rely on the data to be ok. @@
136 // also for mean, we rely on the consistancy of the data. @@
137 
138  G4int Prompt=0, delayed=0, all=0;
139  G4DynamicParticleVector * theNeutrons = 0;
140  switch(it) // check logic, and ask, if partials can be assumed to correspond to individual particles @@@
141  {
142  case 0:
143  theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 0);
144  if(Prompt==0&&delayed==0) Prompt=all;
145  theNeutrons = theFC.ApplyYourself(Prompt); // delayed always in FS
146  // take 'U' into account explicitely (see 5.4) in the sampling of energy @@@@
147  break;
148  case 1:
149  theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 1);
150  if(Prompt==0&&delayed==0) Prompt=all;
151  theNeutrons = theSC.ApplyYourself(Prompt); // delayed always in FS, off done in FSFissionFS
152  break;
153  case 2:
154  theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 2);
155  if(Prompt==0&&delayed==0) Prompt=all;
156  theNeutrons = theTC.ApplyYourself(Prompt); // delayed always in FS
157  break;
158  case 3:
159  theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 3);
160  if(Prompt==0&&delayed==0) Prompt=all;
161  theNeutrons = theLC.ApplyYourself(Prompt); // delayed always in FS
162  break;
163  default:
164  break;
165  }
166 
167 // dice delayed neutrons and photons, and fallback
168 // for Prompt in case channel had no FS data; add all paricles to FS.
169 
170  //TKWORK120531
171  if ( produceFissionFragments ) delayed=0;
172 
173  G4double * theDecayConstants;
174 
175  if( theNeutrons != 0)
176  {
177  theDecayConstants = new G4double[delayed];
178  //
179  //110527TKDB Unused codes, Detected by gcc4.6 compiler
180  //G4int nPhotons = 0;
181  //if(thePhotons!=0) nPhotons = thePhotons->size();
182  for(i=0; i<theNeutrons->size(); i++)
183  {
184  theResult.AddSecondary(theNeutrons->operator[](i));
185  }
186  delete theNeutrons;
187 
188  G4DynamicParticleVector * theDelayed = 0;
189 // G4cout << "delayed" << G4endl;
190  theDelayed = theFS.ApplyYourself(0, delayed, theDecayConstants);
191  for(i=0; i<theDelayed->size(); i++)
192  {
193  G4double time = -std::log(G4UniformRand())/theDecayConstants[i];
194  time += theTrack.GetGlobalTime();
195  theResult.AddSecondary(theDelayed->operator[](i));
197  }
198  delete theDelayed;
199  }
200  else
201  {
202 // cout << " all = "<<all<<G4endl;
203  theFS.SampleNeutronMult(all, Prompt, delayed, eKinetic, 0);
204  theDecayConstants = new G4double[delayed];
205  if(Prompt==0&&delayed==0) Prompt=all;
206  theNeutrons = theFS.ApplyYourself(Prompt, delayed, theDecayConstants);
207  //110527TKDB Unused codes, Detected by gcc4.6 compiler
208  //G4int nPhotons = 0;
209  //if(thePhotons!=0) nPhotons = thePhotons->size();
210  G4int i0;
211  for(i0=0; i0<Prompt; i0++)
212  {
213  theResult.AddSecondary(theNeutrons->operator[](i0));
214  }
215 
216 //G4cout << "delayed" << G4endl;
217  for(i0=Prompt; i0<Prompt+delayed; i0++)
218  {
219  G4double time = -std::log(G4UniformRand())/theDecayConstants[i0-Prompt];
220  time += theTrack.GetGlobalTime();
221  theResult.AddSecondary(theNeutrons->operator[](i0));
223  }
224  delete theNeutrons;
225  }
226  delete [] theDecayConstants;
227 // cout << "all delayed "<<delayed<<G4endl;
228  unsigned int nPhotons = 0;
229  if(thePhotons!=0)
230  {
231  nPhotons = thePhotons->size();
232  for(i=0; i<nPhotons; i++)
233  {
234  theResult.AddSecondary(thePhotons->operator[](i));
235  }
236  delete thePhotons;
237  }
238 
239 // finally deal with local energy depositions.
240 // G4cout <<"Number of secondaries = "<<theResult.GetNumberOfSecondaries()<< G4endl;
241 // G4cout <<"Number of photons = "<<nPhotons<<G4endl;
242 // G4cout <<"Number of Prompt = "<<Prompt<<G4endl;
243 // G4cout <<"Number of delayed = "<<delayed<<G4endl;
244 
245  G4NeutronHPFissionERelease * theERelease = theFS.GetEnergyRelease();
246  G4double eDepByFragments = theERelease->GetFragmentKinetic();
247  //theResult.SetLocalEnergyDeposit(eDepByFragments);
248  if ( !produceFissionFragments ) theResult.SetLocalEnergyDeposit(eDepByFragments);
249 // cout << "local energy deposit" << eDepByFragments<<G4endl;
250 // clean up the primary neutron
252  //G4cout << "Prompt = " << Prompt << ", Delayed = " << delayed << ", All= " << all << G4endl;
253  //G4cout << "local energy deposit " << eDepByFragments/MeV << "MeV " << G4endl;
254 
255  //TKWORK120531
256  if ( produceFissionFragments )
257  {
258  G4int fragA_Z=0;
259  G4int fragA_A=0;
260  G4int fragA_M=0;
261  // System is traget rest!
262  theFF.GetAFissionFragment(eKinetic,fragA_Z,fragA_A,fragA_M);
263  G4int fragB_Z=(G4int)theBaseZ-fragA_Z;
264  G4int fragB_A=(G4int)theBaseA-fragA_A-Prompt;
265  //fragA_M ignored
266  //G4int fragB_M=theBaseM-fragA_M;
267  //G4cout << fragA_Z << " " << fragA_A << " " << fragA_M << G4endl;
268  //G4cout << fragB_Z << " " << fragB_A << G4endl;
269 
271  //Excitation energy is not taken into account
272  G4ParticleDefinition* pdA = pt->GetIon( fragA_Z , fragA_A , 0.0 );
273  G4ParticleDefinition* pdB = pt->GetIon( fragB_Z , fragB_A , 0.0 );
274 
275  //Isotropic Distribution
276  G4double phi = twopi*G4UniformRand();
277  G4double theta = pi*G4UniformRand();
278  G4double sinth = std::sin(theta);
279  G4ThreeVector direction (sinth*std::cos(phi) , sinth*std::sin(phi), std::cos(theta) );
280 
281  // Just use ENDF value for this
282  G4double ER = eDepByFragments;
283  G4double ma = pdA->GetPDGMass();
284  G4double mb = pdB->GetPDGMass();
285  G4double EA = ER / ( 1 + ma/mb);
286  G4double EB = ER - EA;
287  G4DynamicParticle* dpA = new G4DynamicParticle( pdA , direction , EA);
288  G4DynamicParticle* dpB = new G4DynamicParticle( pdB , -direction , EB);
289  theResult.AddSecondary(dpA);
290  theResult.AddSecondary(dpB);
291  }
292  //TKWORK 120531 END
293 
294  return &theResult;
295  }
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType)
void SetNeutron(const G4ReactionProduct &aNeutron)
virtual G4double GetXsec(G4double anEnergy)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4HadSecondary * GetSecondary(size_t i)
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType)
void SetKineticEnergy(const G4double en)
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
G4NeutronHPFissionERelease * GetEnergyRelease()
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType)
void SetTarget(const G4ReactionProduct &aTarget)
int G4int
Definition: G4Types.hh:78
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType)
void SetStatusChange(G4HadFinalStateStatus aS)
void SetTarget(const G4ReactionProduct &aTarget)
G4DynamicParticleVector * ApplyYourself(G4int NNeutrons)
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:108
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack)
const G4ParticleDefinition * GetDefinition() const
G4DynamicParticleVector * ApplyYourself(G4int Prompt, G4int delayed, G4double *decayconst)
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType)
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType)
G4double GetKineticEnergy() const
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType)
G4double GetGlobalTime() const
G4ErrorTarget * theTarget
Definition: errprop.cc:59
G4DynamicParticleVector * ApplyYourself(G4int NNeutrons)
G4DynamicParticleVector * GetPhotons()
std::vector< G4DynamicParticle * > G4DynamicParticleVector
void GetAFissionFragment(G4double, G4int &, G4int &, G4int &)
const G4LorentzVector & Get4Momentum() const
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:80
G4double GetKineticEnergy() const
G4double GetPDGMass() const
G4DynamicParticleVector * ApplyYourself(G4int nNeutrons)
void SampleNeutronMult(G4int &all, G4int &Prompt, G4int &delayed, G4double energy, G4int off)
G4ThreeVector GetMomentum() const
G4double GetTemperature() const
Definition: G4Material.hh:180
void SetLocalEnergyDeposit(G4double aE)
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
void SetNeutron(const G4ReactionProduct &aNeutron)
double G4double
Definition: G4Types.hh:76
G4int GetNumberOfSecondaries() const
void AddSecondary(G4DynamicParticle *aP)
G4DynamicParticleVector * ApplyYourself(G4int NNeutrons)