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