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