source: trunk/source/processes/hadronic/models/neutron_hp/src/G4NeutronHPElasticFS.cc@ 830

Last change on this file since 830 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 13.4 KB
Line 
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//
33#include "G4NeutronHPElasticFS.hh"
34#include "G4ReactionProduct.hh"
35#include "G4Nucleus.hh"
36#include "G4Proton.hh"
37#include "G4Deuteron.hh"
38#include "G4Triton.hh"
39#include "G4Alpha.hh"
40#include "G4ThreeVector.hh"
41#include "G4LorentzVector.hh"
42#include "G4ParticleTable.hh"
43#include "G4NeutronHPDataUsed.hh"
44
45 void G4NeutronHPElasticFS::Init (G4double A, G4double Z, G4String & dirName, G4String & )
46 {
47 G4String tString = "/FS/";
48 G4bool dbool;
49 G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), dirName, tString, dbool);
50 G4String filename = aFile.GetName();
51 theBaseA = aFile.GetA();
52 theBaseZ = aFile.GetZ();
53 if(!dbool)
54 {
55 hasAnyData = false;
56 hasFSData = false;
57 hasXsec = false;
58 return;
59 }
60 std::ifstream theData(filename, std::ios::in);
61 theData >> repFlag >> targetMass >> frameFlag;
62 if(repFlag==1)
63 {
64 G4int nEnergy;
65 theData >> nEnergy;
66 theCoefficients = new G4NeutronHPLegendreStore(nEnergy);
67 theCoefficients->InitInterpolation(theData);
68 G4double temp, energy;
69 G4int tempdep, nLegendre;
70 G4int i, ii;
71 for (i=0; i<nEnergy; i++)
72 {
73 theData >> temp >> energy >> tempdep >> nLegendre;
74 energy *=eV;
75 theCoefficients->Init(i, energy, nLegendre);
76 theCoefficients->SetTemperature(i, temp);
77 G4double coeff=0;
78 for(ii=0; ii<nLegendre; ii++)
79 {
80 // load legendre coefficients.
81 theData >> coeff;
82 theCoefficients->SetCoeff(i, ii+1, coeff); // @@@HPW@@@
83 }
84 }
85 }
86 else if (repFlag==2)
87 {
88 G4int nEnergy;
89 theData >> nEnergy;
90 theProbArray = new G4NeutronHPPartial(nEnergy, nEnergy);
91 theProbArray->InitInterpolation(theData);
92 G4double temp, energy;
93 G4int tempdep, nPoints;
94 for(G4int i=0; i<nEnergy; i++)
95 {
96 theData >> temp >> energy >> tempdep >> nPoints;
97 energy *= eV;
98 theProbArray->InitInterpolation(i, theData);
99 theProbArray->SetT(i, temp);
100 theProbArray->SetX(i, energy);
101 G4double prob, costh;
102 for(G4int ii=0; ii<nPoints; ii++)
103 {
104 // fill probability arrays.
105 theData >> costh >> prob;
106 theProbArray->SetX(i, ii, costh);
107 theProbArray->SetY(i, ii, prob);
108 }
109 }
110 }
111 else if ( repFlag==3 )
112 {
113 G4int nEnergy_Legendre;
114 theData >> nEnergy_Legendre;
115 theCoefficients = new G4NeutronHPLegendreStore( nEnergy_Legendre );
116 theCoefficients->InitInterpolation( theData );
117 G4double temp, energy;
118 G4int tempdep, nLegendre;
119 G4int i, ii;
120 for ( i = 0 ; i < nEnergy_Legendre ; i++ )
121 {
122 theData >> temp >> energy >> tempdep >> nLegendre;
123 energy *=eV;
124 theCoefficients->Init( i , energy , nLegendre );
125 theCoefficients->SetTemperature( i , temp );
126 G4double coeff = 0;
127 for ( ii = 0 ; ii < nLegendre ; ii++ )
128 {
129 // load legendre coefficients.
130 theData >> coeff;
131 theCoefficients->SetCoeff(i, ii+1, coeff); // @@@HPW@@@
132 }
133 }
134
135 tE_of_repFlag3 = energy;
136
137 G4int nEnergy_Prob;
138 theData >> nEnergy_Prob;
139 theProbArray = new G4NeutronHPPartial( nEnergy_Prob , nEnergy_Prob );
140 theProbArray->InitInterpolation( theData );
141 G4int nPoints;
142 for ( G4int i=0 ; i < nEnergy_Prob ; i++ )
143 {
144 theData >> temp >> energy >> tempdep >> nPoints;
145
146 energy *= eV;
147
148// consistensy check
149 if ( i == 0 )
150 if ( energy != tE_of_repFlag3 )
151 G4cout << "Warning Trangition Energy of repFlag3 is not consistent." << G4endl;
152
153 theProbArray->InitInterpolation( i , theData );
154 theProbArray->SetT( i , temp );
155 theProbArray->SetX( i , energy );
156 G4double prob, costh;
157 for( G4int ii = 0 ; ii < nPoints ; ii++ )
158 {
159 // fill probability arrays.
160 theData >> costh >> prob;
161 theProbArray->SetX( i , ii , costh );
162 theProbArray->SetY( i , ii , prob );
163 }
164 }
165 }
166 else if (repFlag==0)
167 {
168 theData >> frameFlag;
169 }
170 else
171 {
172 G4cout << "unusable number for repFlag: repFlag="<<repFlag<<G4endl;
173 throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPElasticFS::Init -- unusable number for repFlag");
174 }
175 theData.close();
176 }
177 G4HadFinalState * G4NeutronHPElasticFS::ApplyYourself(const G4HadProjectile & theTrack)
178 {
179// G4cout << "G4NeutronHPElasticFS::ApplyYourself+"<<G4endl;
180 theResult.Clear();
181 G4double eKinetic = theTrack.GetKineticEnergy();
182 const G4HadProjectile *incidentParticle = &theTrack;
183 G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
184 theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
185 theNeutron.SetKineticEnergy( eKinetic );
186// G4cout << "G4NeutronHPElasticFS::ApplyYourself++"<<eKinetic<<" "<<G4endl;
187// G4cout << "CMSVALUES 0 "<<theNeutron.GetTotalMomentum()<<G4endl;
188
189 G4ReactionProduct theTarget;
190 G4Nucleus aNucleus;
191 G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum();
192 theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
193// G4cout << "Nucleus-test"<<" "<<targetMass<<" ";
194// G4cout << theTarget.GetMomentum().x()<<" ";
195// G4cout << theTarget.GetMomentum().y()<<" ";
196// G4cout << theTarget.GetMomentum().z()<<G4endl;
197
198 // neutron and target defined as reaction products.
199
200// prepare lorentz-transformation to Lab.
201
202 G4ThreeVector the3Neutron = theNeutron.GetMomentum();
203 G4double nEnergy = theNeutron.GetTotalEnergy();
204 G4ThreeVector the3Target = theTarget.GetMomentum();
205// cout << "@@@" << the3Target<<G4endl;
206 G4double tEnergy = theTarget.GetTotalEnergy();
207 G4ReactionProduct theCMS;
208 G4double totE = nEnergy+tEnergy;
209 G4ThreeVector the3CMS = the3Target+the3Neutron;
210 theCMS.SetMomentum(the3CMS);
211 G4double cmsMom = std::sqrt(the3CMS*the3CMS);
212 G4double sqrts = std::sqrt((totE-cmsMom)*(totE+cmsMom));
213 theCMS.SetMass(sqrts);
214 theCMS.SetTotalEnergy(totE);
215
216 // data come as fcn of n-energy in nuclear rest frame
217 G4ReactionProduct boosted;
218 boosted.Lorentz(theNeutron, theTarget);
219 eKinetic = boosted.GetKineticEnergy(); // get kinetic energy for scattering
220 G4double cosTh = -2;
221 if(repFlag == 1)
222 {
223 cosTh = theCoefficients->SampleElastic(eKinetic);
224 }
225
226 else if (repFlag==2)
227 {
228 cosTh = theProbArray->Sample(eKinetic);
229 }
230 else if (repFlag==3)
231 {
232 if ( eKinetic <= tE_of_repFlag3 )
233 {
234 cosTh = theCoefficients->SampleElastic(eKinetic);
235 }
236 else
237 {
238 cosTh = theProbArray->Sample(eKinetic);
239 }
240 }
241 else if (repFlag==0)
242 {
243 cosTh = 2.*G4UniformRand()-1.;
244 }
245 else
246 {
247 G4cout << "unusable number for repFlag: repFlag="<<repFlag<<G4endl;
248 throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPElasticFS::Init -- unusable number for repFlag");
249 }
250 if(cosTh<-1.1) { return 0; }
251 G4double phi = twopi*G4UniformRand();
252 G4double theta = std::acos(cosTh);
253 G4double sinth = std::sin(theta);
254 if (frameFlag == 1) // final state data given in target rest frame.
255 {
256 // we have the scattering angle, now we need the energy, then do the
257 // boosting.
258 // relativistic elastic scattering energy angular correlation:
259 theNeutron.Lorentz(theNeutron, theTarget);
260 G4double e0 = theNeutron.GetTotalEnergy();
261 G4double p0 = theNeutron.GetTotalMomentum();
262 G4double mN = theNeutron.GetMass();
263 G4double mT = theTarget.GetMass();
264 G4double eE = e0+mT;
265 G4double ap = (mT+eE)*(mT-eE) + (p0+mN)*(p0-mN);
266 G4double a = 4*(eE+p0*cosTh)*(eE-p0*cosTh);
267 G4double b = 4*ap*p0*cosTh;
268 G4double c = (2.*eE*mN-ap)*(2.*eE*mN+ap);
269 G4double en = (-b+std::sqrt(b*b - 4*a*c) )/(2*a);
270 G4ThreeVector tempVector(en*sinth*std::cos(phi), en*sinth*std::sin(phi), en*std::cos(theta) );
271 theNeutron.SetMomentum(tempVector);
272 theNeutron.SetTotalEnergy(std::sqrt(en*en+theNeutron.GetMass()*theNeutron.GetMass()));
273 // first to lab
274 theNeutron.Lorentz(theNeutron, -1.*theTarget);
275 // now to CMS
276 theNeutron.Lorentz(theNeutron, theCMS);
277 theTarget.SetMomentum(-theNeutron.GetMomentum());
278 theTarget.SetTotalEnergy(theNeutron.GetTotalEnergy());
279 // and back to lab
280 theNeutron.Lorentz(theNeutron, -1.*theCMS);
281 theTarget.Lorentz(theTarget, -1.*theCMS);
282 }
283 else if (frameFlag == 2) // CMS
284 {
285 theNeutron.Lorentz(theNeutron, theCMS);
286 theTarget.Lorentz(theTarget, theCMS);
287 G4double en = theNeutron.GetTotalMomentum(); // already in CMS.
288 G4ThreeVector cmsMom=theNeutron.GetMomentum(); // for neutron direction in CMS
289 G4double cms_theta=cmsMom.theta();
290 G4double cms_phi=cmsMom.phi();
291 G4ThreeVector tempVector;
292 tempVector.setX(std::cos(theta)*std::sin(cms_theta)*std::cos(cms_phi)
293 +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::cos(cms_phi)
294 -std::sin(theta)*std::sin(phi)*std::sin(cms_phi) );
295 tempVector.setY(std::cos(theta)*std::sin(cms_theta)*std::sin(cms_phi)
296 +std::sin(theta)*std::cos(phi)*std::cos(cms_theta)*std::sin(cms_phi)
297 +std::sin(theta)*std::sin(phi)*std::cos(cms_phi) );
298 tempVector.setZ(std::cos(theta)*std::cos(cms_theta)
299 -std::sin(theta)*std::cos(phi)*std::sin(cms_theta) );
300 tempVector *= en;
301 theNeutron.SetMomentum(tempVector);
302 theTarget.SetMomentum(-tempVector);
303 G4double tP = theTarget.GetTotalMomentum();
304 G4double tM = theTarget.GetMass();
305 theTarget.SetTotalEnergy(std::sqrt((tP+tM)*(tP+tM)-2.*tP*tM));
306 theNeutron.Lorentz(theNeutron, -1.*theCMS);
307 theTarget.Lorentz(theTarget, -1.*theCMS);
308 }
309 else
310 {
311 G4cout <<"Value of frameFlag (1=LAB, 2=CMS): "<<frameFlag;
312 throw G4HadronicException(__FILE__, __LINE__, "G4NeutronHPElasticFS::ApplyYourSelf frameflag incorrect");
313 }
314 // now all in Lab
315// nun den recoil generieren...und energy change, momentum change angeben.
316 theResult.SetEnergyChange(theNeutron.GetKineticEnergy());
317 theResult.SetMomentumChange(theNeutron.GetMomentum().unit());
318 G4DynamicParticle* theRecoil = new G4DynamicParticle;
319 if(targetMass<4.5)
320 {
321 if(targetMass<1)
322 {
323 // proton
324 theRecoil->SetDefinition(G4Proton::Proton());
325 }
326 else if(targetMass<2 )
327 {
328 // deuteron
329 theRecoil->SetDefinition(G4Deuteron::Deuteron());
330 }
331 else if(targetMass<2.999 )
332 {
333 // 3He
334 theRecoil->SetDefinition(G4He3::He3());
335 }
336 else if(targetMass<3 )
337 {
338 // Triton
339 theRecoil->SetDefinition(G4Triton::Triton());
340 }
341 else
342 {
343 // alpha
344 theRecoil->SetDefinition(G4Alpha::Alpha());
345 }
346 }
347 else
348 {
349 theRecoil->SetDefinition(G4ParticleTable::GetParticleTable()
350 ->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA), 0, static_cast<G4int>(theBaseZ)));
351 }
352 theRecoil->SetMomentum(theTarget.GetMomentum());
353 theResult.AddSecondary(theRecoil);
354// G4cout << "G4NeutronHPElasticFS::ApplyYourself 10+"<<G4endl;
355 // postpone the tracking of the primary neutron
356 theResult.SetStatusChange(suspend);
357 return &theResult;
358 }
Note: See TracBrowser for help on using the repository browser.