Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4NeutronHPElasticFS.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 // 25-08-06 New Final State type (refFlag==3 , Legendre (Low Energy) + Probability (High Energy) )
31 // is added by T. KOI
32 // 080904 Add Protection for negative energy results in very low energy ( 1E-6 eV ) scattering by T. Koi
33 //
34 #include "G4NeutronHPElasticFS.hh"
35 #include "G4NeutronHPManager.hh"
36 
37 #include "G4PhysicalConstants.hh"
38 #include "G4SystemOfUnits.hh"
39 #include "G4ReactionProduct.hh"
40 #include "G4Nucleus.hh"
41 #include "G4Proton.hh"
42 #include "G4Deuteron.hh"
43 #include "G4Triton.hh"
44 #include "G4Alpha.hh"
45 #include "G4ThreeVector.hh"
46 #include "G4LorentzVector.hh"
47 #include "G4IonTable.hh"
48 #include "G4NeutronHPDataUsed.hh"
49 
50 #include "zlib.h"
51 
53  {
54  G4String tString = "/FS";
55  G4bool dbool;
56  G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, dirName, tString, dbool);
57  G4String filename = aFile.GetName();
58  SetAZMs( A, Z, M, aFile );
59  //theBaseA = aFile.GetA();
60  //theBaseZ = aFile.GetZ();
61  if(!dbool)
62  {
63  hasAnyData = false;
64  hasFSData = false;
65  hasXsec = false;
66  return;
67  }
68  //130205 For compressed data files
69  std::istringstream theData(std::ios::in);
71  //130205 END
72  theData >> repFlag >> targetMass >> frameFlag;
73  if(repFlag==1)
74  {
75  G4int nEnergy;
76  theData >> nEnergy;
77  theCoefficients = new G4NeutronHPLegendreStore(nEnergy);
78  theCoefficients->InitInterpolation(theData);
79  G4double temp, energy;
80  G4int tempdep, nLegendre;
81  G4int i, ii;
82  for (i=0; i<nEnergy; i++)
83  {
84  theData >> temp >> energy >> tempdep >> nLegendre;
85  energy *=eV;
86  theCoefficients->Init(i, energy, nLegendre);
87  theCoefficients->SetTemperature(i, temp);
88  G4double coeff=0;
89  for(ii=0; ii<nLegendre; ii++)
90  {
91  // load legendre coefficients.
92  theData >> coeff;
93  theCoefficients->SetCoeff(i, ii+1, coeff); // @@@HPW@@@
94  }
95  }
96  }
97  else if (repFlag==2)
98  {
99  G4int nEnergy;
100  theData >> nEnergy;
101  theProbArray = new G4NeutronHPPartial(nEnergy, nEnergy);
102  theProbArray->InitInterpolation(theData);
103  G4double temp, energy;
104  G4int tempdep, nPoints;
105  for(G4int i=0; i<nEnergy; i++)
106  {
107  theData >> temp >> energy >> tempdep >> nPoints;
108  energy *= eV;
109  theProbArray->InitInterpolation(i, theData);
110  theProbArray->SetT(i, temp);
111  theProbArray->SetX(i, energy);
112  G4double prob, costh;
113  for(G4int ii=0; ii<nPoints; ii++)
114  {
115  // fill probability arrays.
116  theData >> costh >> prob;
117  theProbArray->SetX(i, ii, costh);
118  theProbArray->SetY(i, ii, prob);
119  }
120  }
121  }
122  else if ( repFlag==3 )
123  {
124  G4int nEnergy_Legendre;
125  theData >> nEnergy_Legendre;
126  theCoefficients = new G4NeutronHPLegendreStore( nEnergy_Legendre );
127  theCoefficients->InitInterpolation( theData );
128  G4double temp, energy;
129  G4int tempdep, nLegendre;
130  //G4int i, ii;
131  for ( G4int i = 0 ; i < nEnergy_Legendre ; i++ )
132  {
133  theData >> temp >> energy >> tempdep >> nLegendre;
134  energy *=eV;
135  theCoefficients->Init( i , energy , nLegendre );
136  theCoefficients->SetTemperature( i , temp );
137  G4double coeff = 0;
138  for (G4int ii = 0 ; ii < nLegendre ; ii++ )
139  {
140  // load legendre coefficients.
141  theData >> coeff;
142  theCoefficients->SetCoeff(i, ii+1, coeff); // @@@HPW@@@
143  }
144  }
145 
146  tE_of_repFlag3 = energy;
147 
148  G4int nEnergy_Prob;
149  theData >> nEnergy_Prob;
150  theProbArray = new G4NeutronHPPartial( nEnergy_Prob , nEnergy_Prob );
151  theProbArray->InitInterpolation( theData );
152  G4int nPoints;
153  for ( G4int i=0 ; i < nEnergy_Prob ; i++ )
154  {
155  theData >> temp >> energy >> tempdep >> nPoints;
156 
157  energy *= eV;
158 
159 // consistency check
160  if ( i == 0 )
161  //if ( energy != tE_of_repFlag3 ) //110620TK This is too tight for 32bit machines
162  if ( std::abs( energy - tE_of_repFlag3 ) / tE_of_repFlag3 > 1.0e-15 )
163  G4cout << "Warning Transition Energy of repFlag3 is not consistent." << G4endl;
164 
165  theProbArray->InitInterpolation( i , theData );
166  theProbArray->SetT( i , temp );
167  theProbArray->SetX( i , energy );
168  G4double prob, costh;
169  for( G4int ii = 0 ; ii < nPoints ; ii++ )
170  {
171  // fill probability arrays.
172  theData >> costh >> prob;
173  theProbArray->SetX( i , ii , costh );
174  theProbArray->SetY( i , ii , prob );
175  }
176  }
177  }
178  else if (repFlag==0)
179  {
180  theData >> frameFlag;
181  }
182  else
183  {
184  G4cout << "unusable number for repFlag: repFlag="<<repFlag<<G4endl;
185  throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPElasticFS::Init -- unusable number for repFlag");
186  }
187  //130205 For compressed data files(theData changed from ifstream to istringstream)
188  //theData.close();
189  }
191  {
192 // G4cout << "G4NeutronHPElasticFS::ApplyYourself+"<<G4endl;
193  theResult.Clear();
194  G4double eKinetic = theTrack.GetKineticEnergy();
195  const G4HadProjectile *incidentParticle = &theTrack;
196  G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
197  theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
198  theNeutron.SetKineticEnergy( eKinetic );
199 // G4cout << "G4NeutronHPElasticFS::ApplyYourself++"<<eKinetic<<" "<<G4endl;
200 // G4cout << "CMSVALUES 0 "<<theNeutron.GetTotalMomentum()<<G4endl;
201 
203  G4Nucleus aNucleus;
204  G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum();
205  theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
206 // G4cout << "Nucleus-test"<<" "<<targetMass<<" ";
207 // G4cout << theTarget.GetMomentum().x()<<" ";
208 // G4cout << theTarget.GetMomentum().y()<<" ";
209 // G4cout << theTarget.GetMomentum().z()<<G4endl;
210 
211  // neutron and target defined as reaction products.
212 
213 // prepare lorentz-transformation to Lab.
214 
215  G4ThreeVector the3Neutron = theNeutron.GetMomentum();
216  G4double nEnergy = theNeutron.GetTotalEnergy();
217  G4ThreeVector the3Target = theTarget.GetMomentum();
218 // cout << "@@@" << the3Target<<G4endl;
219  G4double tEnergy = theTarget.GetTotalEnergy();
220  G4ReactionProduct theCMS;
221  G4double totE = nEnergy+tEnergy;
222  G4ThreeVector the3CMS = the3Target+the3Neutron;
223  theCMS.SetMomentum(the3CMS);
224  G4double cmsMom = std::sqrt(the3CMS*the3CMS);
225  G4double sqrts = std::sqrt((totE-cmsMom)*(totE+cmsMom));
226  theCMS.SetMass(sqrts);
227  theCMS.SetTotalEnergy(totE);
228 
229  // data come as fcn of n-energy in nuclear rest frame
230  G4ReactionProduct boosted;
231  boosted.Lorentz(theNeutron, theTarget);
232  eKinetic = boosted.GetKineticEnergy(); // get kinetic energy for scattering
233  G4double cosTh = -2;
234  if(repFlag == 1)
235  {
236  cosTh = theCoefficients->SampleElastic(eKinetic);
237  }
238 
239  else if (repFlag==2)
240  {
241  cosTh = theProbArray->Sample(eKinetic);
242  }
243  else if (repFlag==3)
244  {
245  if ( eKinetic <= tE_of_repFlag3 )
246  {
247  cosTh = theCoefficients->SampleElastic(eKinetic);
248  }
249  else
250  {
251  cosTh = theProbArray->Sample(eKinetic);
252  }
253  }
254  else if (repFlag==0)
255  {
256  cosTh = 2.*G4UniformRand()-1.;
257  }
258  else
259  {
260  G4cout << "unusable number for repFlag: repFlag="<<repFlag<<G4endl;
261  throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPElasticFS::Init -- unusable number for repFlag");
262  }
263  if(cosTh<-1.1) { return 0; }
264  G4double phi = twopi*G4UniformRand();
265  G4double theta = std::acos(cosTh);
266  G4double sinth = std::sin(theta);
267  if (frameFlag == 1) // final state data given in target rest frame.
268  {
269  // we have the scattering angle, now we need the energy, then do the
270  // boosting.
271  // relativistic elastic scattering energy angular correlation:
272  theNeutron.Lorentz(theNeutron, theTarget);
273  G4double e0 = theNeutron.GetTotalEnergy();
274  G4double p0 = theNeutron.GetTotalMomentum();
275  G4double mN = theNeutron.GetMass();
276  G4double mT = theTarget.GetMass();
277  G4double eE = e0+mT;
278  G4double ap = (mT+eE)*(mT-eE) + (p0+mN)*(p0-mN);
279  G4double a = 4*(eE+p0*cosTh)*(eE-p0*cosTh);
280  G4double b = 4*ap*p0*cosTh;
281  G4double c = (2.*eE*mN-ap)*(2.*eE*mN+ap);
282  G4double en = (-b+std::sqrt(b*b - 4*a*c) )/(2*a);
283  G4ThreeVector tempVector(en*sinth*std::cos(phi), en*sinth*std::sin(phi), en*std::cos(theta) );
284  theNeutron.SetMomentum(tempVector);
285  theNeutron.SetTotalEnergy(std::sqrt(en*en+theNeutron.GetMass()*theNeutron.GetMass()));
286  // first to lab
287  theNeutron.Lorentz(theNeutron, -1.*theTarget);
288  // now to CMS
289  theNeutron.Lorentz(theNeutron, theCMS);
290  theTarget.SetMomentum(-theNeutron.GetMomentum());
291  theTarget.SetTotalEnergy(theNeutron.GetTotalEnergy());
292  // and back to lab
293  theNeutron.Lorentz(theNeutron, -1.*theCMS);
294  theTarget.Lorentz(theTarget, -1.*theCMS);
295 //111005 Protection for not producing 0 kinetic energy target
296  if ( theNeutron.GetKineticEnergy() <= 0 ) theNeutron.SetTotalEnergy ( theNeutron.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
297  if ( theTarget.GetKineticEnergy() <= 0 ) theTarget.SetTotalEnergy ( theTarget.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
298  }
299  else if (frameFlag == 2) // CMS
300  {
301  theNeutron.Lorentz(theNeutron, theCMS);
302  theTarget.Lorentz(theTarget, theCMS);
303  G4double en = theNeutron.GetTotalMomentum(); // already in CMS.
304  G4ThreeVector cmsMom_tmp=theNeutron.GetMomentum(); // for neutron direction in CMS
305  G4double cms_theta=cmsMom_tmp.theta();
306  G4double cms_phi=cmsMom_tmp.phi();
307  G4ThreeVector tempVector;
308  tempVector.setX(std::cos(theta)*std::sin(cms_theta)*std::cos(cms_phi)
309  +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::cos(cms_phi)
310  -std::sin(theta)*std::sin(phi)*std::sin(cms_phi) );
311  tempVector.setY(std::cos(theta)*std::sin(cms_theta)*std::sin(cms_phi)
312  +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::sin(cms_phi)
313  +std::sin(theta)*std::sin(phi)*std::cos(cms_phi) );
314  tempVector.setZ(std::cos(theta)*std::cos(cms_theta)
315  -std::sin(theta)*std::cos(phi)*std::sin(cms_theta) );
316  tempVector *= en;
317  theNeutron.SetMomentum(tempVector);
318  theTarget.SetMomentum(-tempVector);
319  G4double tP = theTarget.GetTotalMomentum();
320  G4double tM = theTarget.GetMass();
321  theTarget.SetTotalEnergy(std::sqrt((tP+tM)*(tP+tM)-2.*tP*tM));
322 
323 /*
324 For debug purpose.
325 Same transformation G4ReactionProduct.Lorentz() by 4vectors
326 {
327 G4LorentzVector n4p = G4LorentzVector ( theNeutron.GetMomentum() , theNeutron.GetKineticEnergy() + theNeutron.GetMass() );
328 G4cout << "before " << ( n4p.e() - n4p.m() ) / eV<< G4endl;
329 G4LorentzVector cm4p = G4LorentzVector ( theCMS.GetMomentum() , theCMS.GetKineticEnergy() + theCMS.GetMass() );
330 n4p.boost( cm4p.boostVector() );
331 G4cout << cm4p/eV << G4endl;
332 G4cout << "after " << ( n4p.e() - n4p.m() ) / eV<< G4endl;
333 }
334 */
335 
336  theNeutron.Lorentz(theNeutron, -1.*theCMS);
337 //080904 Add Protection for very low energy (1e-6eV) scattering
338  if ( theNeutron.GetKineticEnergy() <= 0 )
339  {
340  //theNeutron.SetMomentum( G4ThreeVector(0) );
341  //theNeutron.SetTotalEnergy ( theNeutron.GetMass() );
342 //110822 Protection for not producing 0 kinetic energy neutron
343  theNeutron.SetTotalEnergy ( theNeutron.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
344  }
345 
346  theTarget.Lorentz(theTarget, -1.*theCMS);
347 //080904 Add Protection for very low energy (1e-6eV) scattering
348  if ( theTarget.GetKineticEnergy() < 0 )
349  {
350  //theTarget.SetMomentum( G4ThreeVector(0) );
351  //theTarget.SetTotalEnergy ( theTarget.GetMass() );
352 //110822 Protection for not producing 0 kinetic energy target
353  theTarget.SetTotalEnergy ( theTarget.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
354  }
355  }
356  else
357  {
358  G4cout <<"Value of frameFlag (1=LAB, 2=CMS): "<<frameFlag;
359  throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPElasticFS::ApplyYourSelf frameflag incorrect");
360  }
361  // now all in Lab
362 // nun den recoil generieren...und energy change, momentum change angeben.
365  G4DynamicParticle* theRecoil = new G4DynamicParticle;
366  if(targetMass<4.5)
367  {
368  if(targetMass<1)
369  {
370  // proton
371  theRecoil->SetDefinition(G4Proton::Proton());
372  }
373  else if(targetMass<2 )
374  {
375  // deuteron
376  theRecoil->SetDefinition(G4Deuteron::Deuteron());
377  }
378  else if(targetMass<2.999 )
379  {
380  // 3He
381  theRecoil->SetDefinition(G4He3::He3());
382  }
383  else if(targetMass<3 )
384  {
385  // Triton
386  theRecoil->SetDefinition(G4Triton::Triton());
387  }
388  else
389  {
390  // alpha
391  theRecoil->SetDefinition(G4Alpha::Alpha());
392  }
393  }
394  else
395  {
396  //theRecoil->SetDefinition(G4ParticleTable::GetParticleTable()
397  // ->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA), 0, static_cast<G4int>(theBaseZ)));
399  ->GetIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA), 0 ));
400  }
401  theRecoil->SetMomentum(theTarget.GetMomentum());
402  theResult.AddSecondary(theRecoil);
403 // G4cout << "G4NeutronHPElasticFS::ApplyYourself 10+"<<G4endl;
404  // postpone the tracking of the primary neutron
406  return &theResult;
407  }
void SetMomentum(const G4ThreeVector &momentum)
G4double GetTotalMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
void SetKineticEnergy(const G4double en)
static G4NeutronHPManager * GetInstance()
void Init(G4int i, G4double e, G4int n)
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetCoeff(G4int i, G4int l, G4double coeff)
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 setY(double)
void setZ(double)
void setX(double)
void SetStatusChange(G4HadFinalStateStatus aS)
void SetT(G4int i, G4double x)
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:108
void SetMass(const G4double mas)
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4double SampleElastic(G4double anEnergy)
const G4ParticleDefinition * GetDefinition() const
void SetY(G4int i, G4int j, G4double y)
bool G4bool
Definition: G4Types.hh:79
void SetX(G4int i, G4double x)
G4double GetKineticEnergy() const
void SetTotalEnergy(const G4double en)
G4ErrorTarget * theTarget
Definition: errprop.cc:59
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
void InitInterpolation(std::istream &aDataFile)
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack)
double phi() const
const G4LorentzVector & Get4Momentum() const
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:80
double theta() const
G4double GetKineticEnergy() const
void SetEnergyChange(G4double anEnergy)
G4double GetTotalEnergy() const
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType)
G4double GetPDGMass() const
Hep3Vector unit() const
G4ThreeVector GetMomentum() const
G4double GetTemperature() const
Definition: G4Material.hh:180
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
void InitInterpolation(G4int i, std::istream &aDataFile)
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
void SetAZMs(G4double anA, G4double aZ, G4int aM, G4NeutronHPDataUsed used)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
static G4He3 * He3()
Definition: G4He3.cc:94
void SetMomentumChange(const G4ThreeVector &aV)
void SetTemperature(G4int i, G4double temp)
G4double GetMass() const
void AddSecondary(G4DynamicParticle *aP)
G4double Sample(G4double x)