Geant4  10.01.p01
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;
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);
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);
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  theProbArray->DoneSetXY( i );
121  }
122  }
123  else if ( repFlag==3 )
124  {
125  G4int nEnergy_Legendre;
126  theData >> nEnergy_Legendre;
127  theCoefficients = new G4NeutronHPLegendreStore( nEnergy_Legendre );
129  G4double temp, energy;
130  G4int tempdep, nLegendre;
131  //G4int i, ii;
132  for ( G4int i = 0 ; i < nEnergy_Legendre ; i++ )
133  {
134  theData >> temp >> energy >> tempdep >> nLegendre;
135  energy *=eV;
136  theCoefficients->Init( i , energy , nLegendre );
137  theCoefficients->SetTemperature( i , temp );
138  G4double coeff = 0;
139  for (G4int ii = 0 ; ii < nLegendre ; ii++ )
140  {
141  // load legendre coefficients.
142  theData >> coeff;
143  theCoefficients->SetCoeff(i, ii+1, coeff); // @@@HPW@@@
144  }
145  tE_of_repFlag3 = energy; //the last one will be transition energy from Legendre to table
146  }
147 
148 
149  G4int nEnergy_Prob;
150  theData >> nEnergy_Prob;
151  theProbArray = new G4NeutronHPPartial( nEnergy_Prob , nEnergy_Prob );
152  theProbArray->InitInterpolation( theData );
153  G4int nPoints;
154  for ( G4int i=0 ; i < nEnergy_Prob ; i++ )
155  {
156  theData >> temp >> energy >> tempdep >> nPoints;
157 
158  energy *= eV;
159 
160 // consistency check
161  if ( i == 0 )
162  //if ( energy != tE_of_repFlag3 ) //110620TK This is too tight for 32bit machines
163  if ( std::abs( energy - tE_of_repFlag3 ) / tE_of_repFlag3 > 1.0e-15 )
164  G4cout << "Warning Transition Energy of repFlag3 is not consistent." << G4endl;
165 
166  theProbArray->InitInterpolation( i , theData );
167  theProbArray->SetT( i , temp );
168  theProbArray->SetX( i , energy );
169  G4double prob, costh;
170  for( G4int ii = 0 ; ii < nPoints ; ii++ )
171  {
172  // fill probability arrays.
173  theData >> costh >> prob;
174  theProbArray->SetX( i , ii , costh );
175  theProbArray->SetY( i , ii , prob );
176  }
177  theProbArray->DoneSetXY( i );
178  }
179  }
180  else if (repFlag==0)
181  {
182  theData >> frameFlag;
183  }
184  else
185  {
186  G4cout << "unusable number for repFlag: repFlag="<<repFlag<<G4endl;
187  throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPElasticFS::Init -- unusable number for repFlag");
188  }
189  //130205 For compressed data files(theData changed from ifstream to istringstream)
190  //theData.close();
191  }
193  {
194 // G4cout << "G4NeutronHPElasticFS::ApplyYourself+"<<G4endl;
195  if ( theResult.Get() == NULL ) theResult.Put( new G4HadFinalState );
196  theResult.Get()->Clear();
197  G4double eKinetic = theTrack.GetKineticEnergy();
198  const G4HadProjectile *incidentParticle = &theTrack;
199  G4ReactionProduct theNeutron( incidentParticle->GetDefinition() );
200  theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
201  theNeutron.SetKineticEnergy( eKinetic );
202 // G4cout << "G4NeutronHPElasticFS::ApplyYourself++"<<eKinetic<<" "<<G4endl;
203 // G4cout << "CMSVALUES 0 "<<theNeutron.GetTotalMomentum()<<G4endl;
204 
206  G4Nucleus aNucleus;
207  G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum();
208  theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
209 // G4cout << "Nucleus-test"<<" "<<targetMass<<" ";
210 // G4cout << theTarget.GetMomentum().x()<<" ";
211 // G4cout << theTarget.GetMomentum().y()<<" ";
212 // G4cout << theTarget.GetMomentum().z()<<G4endl;
213 
214  // neutron and target defined as reaction products.
215 
216 // prepare lorentz-transformation to Lab.
217 
218  G4ThreeVector the3Neutron = theNeutron.GetMomentum();
219  G4double nEnergy = theNeutron.GetTotalEnergy();
220  G4ThreeVector the3Target = theTarget.GetMomentum();
221 // cout << "@@@" << the3Target<<G4endl;
222  G4double tEnergy = theTarget.GetTotalEnergy();
223  G4ReactionProduct theCMS;
224  G4double totE = nEnergy+tEnergy;
225  G4ThreeVector the3CMS = the3Target+the3Neutron;
226  theCMS.SetMomentum(the3CMS);
227  G4double cmsMom = std::sqrt(the3CMS*the3CMS);
228  G4double sqrts = std::sqrt((totE-cmsMom)*(totE+cmsMom));
229  theCMS.SetMass(sqrts);
230  theCMS.SetTotalEnergy(totE);
231 
232  // data come as fcn of n-energy in nuclear rest frame
233  G4ReactionProduct boosted;
234  boosted.Lorentz(theNeutron, theTarget);
235  eKinetic = boosted.GetKineticEnergy(); // get kinetic energy for scattering
236  G4double cosTh = -2;
237  if(repFlag == 1)
238  {
239  cosTh = theCoefficients->SampleElastic(eKinetic);
240  }
241 
242  else if (repFlag==2)
243  {
244  cosTh = theProbArray->Sample(eKinetic);
245  }
246  else if (repFlag==3)
247  {
248  if ( eKinetic <= tE_of_repFlag3 )
249  {
250  cosTh = theCoefficients->SampleElastic(eKinetic);
251  }
252  else
253  {
254  cosTh = theProbArray->Sample(eKinetic);
255  }
256  }
257  else if (repFlag==0)
258  {
259  cosTh = 2.*G4UniformRand()-1.;
260  }
261  else
262  {
263  G4cout << "unusable number for repFlag: repFlag="<<repFlag<<G4endl;
264  throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPElasticFS::Init -- unusable number for repFlag");
265  }
266  if(cosTh<-1.1) { return 0; }
267  G4double phi = twopi*G4UniformRand();
268  G4double theta = std::acos(cosTh);
269  G4double sinth = std::sin(theta);
270  if (frameFlag == 1) // final state data given in target rest frame.
271  {
272  // we have the scattering angle, now we need the energy, then do the
273  // boosting.
274  // relativistic elastic scattering energy angular correlation:
275  theNeutron.Lorentz(theNeutron, theTarget);
276  G4double e0 = theNeutron.GetTotalEnergy();
277  G4double p0 = theNeutron.GetTotalMomentum();
278  G4double mN = theNeutron.GetMass();
279  G4double mT = theTarget.GetMass();
280  G4double eE = e0+mT;
281  G4double ap = (mT+eE)*(mT-eE) + (p0+mN)*(p0-mN);
282  G4double a = 4*(eE+p0*cosTh)*(eE-p0*cosTh);
283  G4double b = 4*ap*p0*cosTh;
284  G4double c = (2.*eE*mN-ap)*(2.*eE*mN+ap);
285  G4double en = (-b+std::sqrt(b*b - 4*a*c) )/(2*a);
286  G4ThreeVector tempVector(en*sinth*std::cos(phi), en*sinth*std::sin(phi), en*std::cos(theta) );
287  theNeutron.SetMomentum(tempVector);
288  theNeutron.SetTotalEnergy(std::sqrt(en*en+theNeutron.GetMass()*theNeutron.GetMass()));
289  // first to lab
290  theNeutron.Lorentz(theNeutron, -1.*theTarget);
291  // now to CMS
292  theNeutron.Lorentz(theNeutron, theCMS);
293  theTarget.SetMomentum(-theNeutron.GetMomentum());
294  theTarget.SetTotalEnergy(theNeutron.GetTotalEnergy());
295  // and back to lab
296  theNeutron.Lorentz(theNeutron, -1.*theCMS);
297  theTarget.Lorentz(theTarget, -1.*theCMS);
298 //111005 Protection for not producing 0 kinetic energy target
299  if ( theNeutron.GetKineticEnergy() <= 0 ) theNeutron.SetTotalEnergy ( theNeutron.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
300  if ( theTarget.GetKineticEnergy() <= 0 ) theTarget.SetTotalEnergy ( theTarget.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
301  }
302  else if (frameFlag == 2) // CMS
303  {
304  theNeutron.Lorentz(theNeutron, theCMS);
305  theTarget.Lorentz(theTarget, theCMS);
306  G4double en = theNeutron.GetTotalMomentum(); // already in CMS.
307  G4ThreeVector cmsMom_tmp=theNeutron.GetMomentum(); // for neutron direction in CMS
308  G4double cms_theta=cmsMom_tmp.theta();
309  G4double cms_phi=cmsMom_tmp.phi();
310  G4ThreeVector tempVector;
311  tempVector.setX(std::cos(theta)*std::sin(cms_theta)*std::cos(cms_phi)
312  +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::cos(cms_phi)
313  -std::sin(theta)*std::sin(phi)*std::sin(cms_phi) );
314  tempVector.setY(std::cos(theta)*std::sin(cms_theta)*std::sin(cms_phi)
315  +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::sin(cms_phi)
316  +std::sin(theta)*std::sin(phi)*std::cos(cms_phi) );
317  tempVector.setZ(std::cos(theta)*std::cos(cms_theta)
318  -std::sin(theta)*std::cos(phi)*std::sin(cms_theta) );
319  tempVector *= en;
320  theNeutron.SetMomentum(tempVector);
321  theTarget.SetMomentum(-tempVector);
322  G4double tP = theTarget.GetTotalMomentum();
323  G4double tM = theTarget.GetMass();
324  theTarget.SetTotalEnergy(std::sqrt((tP+tM)*(tP+tM)-2.*tP*tM));
325 
326 /*
327 For debug purpose.
328 Same transformation G4ReactionProduct.Lorentz() by 4vectors
329 {
330 G4LorentzVector n4p = G4LorentzVector ( theNeutron.GetMomentum() , theNeutron.GetKineticEnergy() + theNeutron.GetMass() );
331 G4cout << "before " << ( n4p.e() - n4p.m() ) / eV<< G4endl;
332 G4LorentzVector cm4p = G4LorentzVector ( theCMS.GetMomentum() , theCMS.GetKineticEnergy() + theCMS.GetMass() );
333 n4p.boost( cm4p.boostVector() );
334 G4cout << cm4p/eV << G4endl;
335 G4cout << "after " << ( n4p.e() - n4p.m() ) / eV<< G4endl;
336 }
337 */
338 
339  theNeutron.Lorentz(theNeutron, -1.*theCMS);
340 //080904 Add Protection for very low energy (1e-6eV) scattering
341  if ( theNeutron.GetKineticEnergy() <= 0 )
342  {
343  //theNeutron.SetMomentum( G4ThreeVector(0) );
344  //theNeutron.SetTotalEnergy ( theNeutron.GetMass() );
345 //110822 Protection for not producing 0 kinetic energy neutron
346  theNeutron.SetTotalEnergy ( theNeutron.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
347  }
348 
349  theTarget.Lorentz(theTarget, -1.*theCMS);
350 //080904 Add Protection for very low energy (1e-6eV) scattering
351  if ( theTarget.GetKineticEnergy() < 0 )
352  {
353  //theTarget.SetMomentum( G4ThreeVector(0) );
354  //theTarget.SetTotalEnergy ( theTarget.GetMass() );
355 //110822 Protection for not producing 0 kinetic energy target
356  theTarget.SetTotalEnergy ( theTarget.GetMass() * ( 1 + std::pow( 10 , -15.65 ) ) );
357  }
358  }
359  else
360  {
361  G4cout <<"Value of frameFlag (1=LAB, 2=CMS): "<<frameFlag;
362  throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPElasticFS::ApplyYourSelf frameflag incorrect");
363  }
364  // now all in Lab
365 // nun den recoil generieren...und energy change, momentum change angeben.
366  theResult.Get()->SetEnergyChange(theNeutron.GetKineticEnergy());
367  theResult.Get()->SetMomentumChange(theNeutron.GetMomentum().unit());
368  G4DynamicParticle* theRecoil = new G4DynamicParticle;
369  theRecoil->SetDefinition( G4IonTable::GetIonTable()->GetIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA), 0 ) );
370  theRecoil->SetMomentum(theTarget.GetMomentum());
371  theResult.Get()->AddSecondary(theRecoil);
372 // G4cout << "G4NeutronHPElasticFS::ApplyYourself 10+"<<G4endl;
373  // postpone the tracking of the primary neutron
375  return theResult.Get();
376  }
void SetMomentum(const G4ThreeVector &momentum)
G4double GetTotalMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4Cache< G4HadFinalState * > theResult
CLHEP::Hep3Vector G4ThreeVector
static G4NeutronHPManager * GetInstance()
void Init(G4int i, G4double e, G4int n)
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4NeutronHPLegendreStore * theCoefficients
G4NeutronHPPartial * theProbArray
value_type & Get() const
Definition: G4Cache.hh:253
void SetCoeff(G4int i, G4int l, G4double coeff)
G4double a
Definition: TRTMaterials.hh:39
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 SetStatusChange(G4HadFinalStateStatus aS)
void SetT(G4int i, G4double x)
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:108
void SetMass(const G4double mas)
#define G4UniformRand()
Definition: Randomize.hh:95
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
void InitInterpolation(std::istream &aDataFile)
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack)
static const G4double A[nN]
const G4LorentzVector & Get4Momentum() const
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
G4double GetKineticEnergy() const
static const double eV
Definition: G4SIunits.hh:194
void SetEnergyChange(G4double anEnergy)
G4double GetTotalEnergy() const
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType)
G4double GetPDGMass() const
G4double energy(const ThreeVector &p, const G4double m)
G4ThreeVector GetMomentum() const
G4double GetTemperature() const
Definition: G4Material.hh:182
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
void InitInterpolation(G4int i, std::istream &aDataFile)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetAZMs(G4double anA, G4double aZ, G4int aM, G4NeutronHPDataUsed used)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void Put(const value_type &val) const
Definition: G4Cache.hh:257
void SetMomentumChange(const G4ThreeVector &aV)
void SetTemperature(G4int i, G4double temp)
G4double GetMass() const
G4double Sample(G4double x)