source: trunk/source/processes/hadronic/models/low_energy/src/G4LEProtonInelastic.cc@ 1196

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

import all except CVS

File size: 16.5 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//
27//
28// Hadronic Process: Low Energy Proton Inelastic Process
29// J.L. Chuma, TRIUMF, 19-Nov-1996
30
31#include "G4LEProtonInelastic.hh"
32#include "Randomize.hh"
33
34 G4HadFinalState *
35 G4LEProtonInelastic::ApplyYourself( const G4HadProjectile &aTrack,
36 G4Nucleus &targetNucleus )
37 {
38 theParticleChange.Clear();
39 const G4HadProjectile *originalIncident = &aTrack;
40 if (originalIncident->GetKineticEnergy()<= 0.1*MeV)
41 {
42 theParticleChange.SetStatusChange(isAlive);
43 theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
44 theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
45 return &theParticleChange;
46 }
47 //
48 // create the target particle
49 //
50 G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
51 if( verboseLevel > 1 )
52 {
53 const G4Material *targetMaterial = aTrack.GetMaterial();
54 G4cout << "G4LEProtonInelastic::ApplyYourself called" << G4endl;
55 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
56 G4cout << "target material = " << targetMaterial->GetName() << ", ";
57 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
58 << G4endl;
59 }
60 if( originalIncident->GetKineticEnergy()/GeV < 0.01+2.*G4UniformRand()/9. )
61 {
62 SlowProton( originalIncident, targetNucleus );
63 delete originalTarget;
64 return &theParticleChange;
65 }
66 //
67 // Fermi motion and evaporation
68 // As of Geant3, the Fermi energy calculation had not been Done
69 //
70 G4double ek = originalIncident->GetKineticEnergy()/MeV;
71 G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
72 G4ReactionProduct modifiedOriginal;
73 modifiedOriginal = *originalIncident;
74
75 G4double tkin = targetNucleus.Cinema( ek );
76 ek += tkin;
77 modifiedOriginal.SetKineticEnergy( ek*MeV );
78 G4double et = ek + amas;
79 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
80 G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
81 if( pp > 0.0 )
82 {
83 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
84 modifiedOriginal.SetMomentum( momentum * (p/pp) );
85 }
86 //
87 // calculate black track energies
88 //
89 tkin = targetNucleus.EvaporationEffects( ek );
90 ek -= tkin;
91 modifiedOriginal.SetKineticEnergy( ek*MeV );
92 et = ek + amas;
93 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
94 pp = modifiedOriginal.GetMomentum().mag()/MeV;
95 if( pp > 0.0 )
96 {
97 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
98 modifiedOriginal.SetMomentum( momentum * (p/pp) );
99 }
100 const G4double cutOff = 0.1;
101 if( modifiedOriginal.GetKineticEnergy()/MeV <= cutOff )
102 {
103 SlowProton( originalIncident, targetNucleus );
104 delete originalTarget;
105 return &theParticleChange;
106 }
107 G4ReactionProduct currentParticle = modifiedOriginal;
108 G4ReactionProduct targetParticle;
109 targetParticle = *originalTarget;
110 currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
111 targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
112 G4bool incidentHasChanged = false;
113 G4bool targetHasChanged = false;
114 G4bool quasiElastic = false;
115 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the sec. particles
116 G4int vecLen = 0;
117 vec.Initialize( 0 );
118
119 Cascade( vec, vecLen,
120 originalIncident, currentParticle, targetParticle,
121 incidentHasChanged, targetHasChanged, quasiElastic );
122
123 CalculateMomenta( vec, vecLen,
124 originalIncident, originalTarget, modifiedOriginal,
125 targetNucleus, currentParticle, targetParticle,
126 incidentHasChanged, targetHasChanged, quasiElastic );
127
128 SetUpChange( vec, vecLen,
129 currentParticle, targetParticle,
130 incidentHasChanged );
131
132 delete originalTarget;
133 return &theParticleChange;
134 }
135
136 void
137 G4LEProtonInelastic::SlowProton(
138 const G4HadProjectile *originalIncident,
139 G4Nucleus &targetNucleus )
140 {
141 const G4double A = targetNucleus.GetN(); // atomic weight
142 const G4double Z = targetNucleus.GetZ(); // atomic number
143// G4double currentKinetic = originalIncident->GetKineticEnergy()/MeV;
144 //
145 // calculate Q-value of reactions
146 //
147 G4double theAtomicMass = targetNucleus.AtomicMass( A, Z );
148 G4double massVec[9];
149 massVec[0] = targetNucleus.AtomicMass( A+1.0, Z+1.0 );
150 massVec[1] = 0.;
151 if (A > Z+1.0)
152 massVec[1] = targetNucleus.AtomicMass( A , Z+1.0 );
153 massVec[2] = theAtomicMass;
154 massVec[3] = 0.;
155 if (A > 1.0 && A-1.0 > Z)
156 massVec[3] = targetNucleus.AtomicMass( A-1.0, Z );
157 massVec[4] = 0.;
158 if (A > 2.0 && A-2.0 > Z)
159 massVec[4] = targetNucleus.AtomicMass( A-2.0, Z );
160 massVec[5] = 0.;
161 if (A > 3.0 && Z > 1.0 && A-3.0 > Z-1.0)
162 massVec[5] = targetNucleus.AtomicMass( A-3.0, Z-1.0 );
163 massVec[6] = 0.;
164 if (A > 1.0 && A-1.0 > Z+1.0)
165 massVec[6] = targetNucleus.AtomicMass( A-1.0, Z+1.0 );
166 massVec[7] = massVec[3];
167 massVec[8] = 0.;
168 if (A > 1.0 && Z > 1.0)
169 massVec[8] = targetNucleus.AtomicMass( A-1.0, Z-1.0 );
170
171 G4FastVector<G4ReactionProduct,4> vec; // vec will contain the secondary particles
172 G4int vecLen = 0;
173 vec.Initialize( 0 );
174
175 theReactionDynamics.NuclearReaction( vec, vecLen, originalIncident,
176 targetNucleus, theAtomicMass, massVec );
177
178 theParticleChange.SetStatusChange( stopAndKill );
179 theParticleChange.SetEnergyChange( 0.0 );
180
181 G4DynamicParticle *pd;
182 for( G4int i=0; i<vecLen; ++i )
183 {
184 pd = new G4DynamicParticle();
185 pd->SetDefinition( vec[i]->GetDefinition() );
186 pd->SetMomentum( vec[i]->GetMomentum() );
187 theParticleChange.AddSecondary( pd );
188 delete vec[i];
189 }
190 }
191
192 void
193 G4LEProtonInelastic::Cascade(
194 G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
195 G4int &vecLen,
196 const G4HadProjectile *originalIncident,
197 G4ReactionProduct &currentParticle,
198 G4ReactionProduct &targetParticle,
199 G4bool &incidentHasChanged,
200 G4bool &targetHasChanged,
201 G4bool &quasiElastic )
202 {
203 // derived from original FORTRAN code CASP by H. Fesefeldt (13-Sep-1987)
204 //
205 // Proton undergoes interaction with nucleon within a nucleus. Check if it is
206 // energetically possible to produce pions/kaons. In not, assume nuclear excitation
207 // occurs and input particle is degraded in energy. No other particles are produced.
208 // If reaction is possible, find the correct number of pions/protons/neutrons
209 // produced using an interpolation to multiplicity data. Replace some pions or
210 // protons/neutrons by kaons or strange baryons according to the average
211 // multiplicity per Inelastic reaction.
212 //
213 // the center of mass energy is based on the initial energy, before
214 // Fermi motion and evaporation effects are taken into account
215 //
216 const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
217 const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
218 const G4double targetMass = targetParticle.GetMass()/MeV;
219 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
220 targetMass*targetMass +
221 2.0*targetMass*etOriginal );
222 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
223 if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
224 { // not energetically possible to produce pion(s)
225 quasiElastic = true;
226 return;
227 }
228 static G4bool first = true;
229 const G4int numMul = 1200;
230 const G4int numSec = 60;
231 static G4double protmul[numMul], protnorm[numSec]; // proton constants
232 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
233 // np = number of pi+, nm = number of pi-, nz = number of pi0
234 G4int counter, nt=0, np=0, nm=0, nz=0;
235 const G4double c = 1.25;
236 const G4double b[] = { 0.70, 0.35 };
237 if( first ) // compute normalization constants, this will only be Done once
238 {
239 first = false;
240 G4int i;
241 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
242 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
243 counter = -1;
244 for( np=0; np<(numSec/3); ++np )
245 {
246 for( nm=std::max(0,np-2); nm<=np; ++nm )
247 {
248 for( nz=0; nz<numSec/3; ++nz )
249 {
250 if( ++counter < numMul )
251 {
252 nt = np+nm+nz;
253 if( nt>0 && nt<=numSec )
254 {
255 protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c) /
256 ( theReactionDynamics.Factorial(2-np+nm)*
257 theReactionDynamics.Factorial(np-nm) );
258 protnorm[nt-1] += protmul[counter];
259 }
260 }
261 }
262 }
263 }
264 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
265 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
266 counter = -1;
267 for( np=0; np<numSec/3; ++np )
268 {
269 for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
270 {
271 for( nz=0; nz<numSec/3; ++nz )
272 {
273 if( ++counter < numMul )
274 {
275 nt = np+nm+nz;
276 if( nt>0 && nt<=numSec )
277 {
278 neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c) /
279 ( theReactionDynamics.Factorial(1-np+nm)*
280 theReactionDynamics.Factorial(1+np-nm) );
281 neutnorm[nt-1] += neutmul[counter];
282 }
283 }
284 }
285 }
286 }
287 for( i=0; i<numSec; ++i )
288 {
289 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
290 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
291 }
292 } // end of initialization
293
294 const G4double expxu = 82.; // upper bound for arg. of exp
295 const G4double expxl = -expxu; // lower bound for arg. of exp
296 G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
297 G4ParticleDefinition *aProton = G4Proton::Proton();
298 G4int ieab = static_cast<G4int>(availableEnergy*5.0/GeV);
299 const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
300 G4double test, w0, wp, wt, wm;
301 if( (availableEnergy < 2.0*GeV) && (G4UniformRand() >= supp[ieab]) )
302 {
303 // suppress high multiplicity events at low momentum
304 // only one pion will be produced
305
306 np = nm = nz = 0;
307 if( targetParticle.GetDefinition() == aProton )
308 {
309 test = std::exp( std::min( expxu, std::max(
310 expxl, -(1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
311 w0 = test/2.0;
312 wp = test;
313 if( G4UniformRand() < w0/(w0+wp) )
314 nz = 1;
315 else
316 np = 1;
317 }
318 else // target is a neutron
319 {
320 test = std::exp( std::min( expxu, std::max(
321 expxl, -(1.0+b[1])*(1.0+b[1])/(2.0*c*c) ) ) );
322 w0 = test;
323 wp = test/2.0;
324 test = std::exp( std::min( expxu, std::max(
325 expxl, -(-1.0+b[1])*(-1.0+b[1])/(2.0*c*c) ) ) );
326 wm = test/2.0;
327 wt = w0+wp+wm;
328 wp += w0;
329 G4double ran = G4UniformRand();
330 if( ran < w0/wt )
331 nz = 1;
332 else if( ran < wp/wt )
333 np = 1;
334 else
335 nm = 1;
336 }
337 }
338 else // (availableEnergy >= 2.0*GeV) || (random number < supp[ieab])
339 {
340 G4double n, anpn;
341 GetNormalizationConstant( availableEnergy, n, anpn );
342 G4double ran = G4UniformRand();
343 G4double dum, excs = 0.0;
344 if( targetParticle.GetDefinition() == aProton )
345 {
346 counter = -1;
347 for( np=0; np<numSec/3 && ran>=excs; ++np )
348 {
349 for( nm=std::max(0,np-2); nm<=np && ran>=excs; ++nm )
350 {
351 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
352 {
353 if( ++counter < numMul )
354 {
355 nt = np+nm+nz;
356 if( nt>0 && nt<=numSec )
357 {
358 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
359 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
360 if( std::fabs(dum) < 1.0 )
361 {
362 if( test >= 1.0e-10 )excs += dum*test;
363 } else {
364 excs += dum*test;
365 }
366 }
367 }
368 }
369 }
370 }
371 if( ran >= excs ) // 3 previous loops continued to the end
372 {
373 quasiElastic = true;
374 return;
375 }
376 }
377 else // target must be a neutron
378 {
379 counter = -1;
380 for( np=0; np<numSec/3 && ran>=excs; ++np )
381 {
382 for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm )
383 {
384 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
385 {
386 if( ++counter < numMul )
387 {
388 nt = np+nm+nz;
389 if( nt>0 && nt<=numSec )
390 {
391 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
392 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
393 if( std::fabs(dum) < 1.0 )
394 {
395 if( test >= 1.0e-10 )excs += dum*test;
396 } else {
397 excs += dum*test;
398 }
399 }
400 }
401 }
402 }
403 }
404 if( ran >= excs ) // 3 previous loops continued to the end
405 {
406 quasiElastic = true;
407 return;
408 }
409 }
410 np--; nm--; nz--;
411 }
412 if( targetParticle.GetDefinition() == aProton )
413 {
414 switch( np-nm )
415 {
416 case 1:
417 if( G4UniformRand() < 0.5 )
418 {
419 targetParticle.SetDefinitionAndUpdateE( aNeutron );
420 targetHasChanged = true;
421 } else {
422 currentParticle.SetDefinitionAndUpdateE( aNeutron );
423 incidentHasChanged = true;
424 }
425 break;
426 case 2:
427 currentParticle.SetDefinitionAndUpdateE( aNeutron );
428 targetParticle.SetDefinitionAndUpdateE( aNeutron );
429 incidentHasChanged = true;
430 targetHasChanged = true;
431 break;
432 default:
433 break;
434 }
435 }
436 else // target is a neutron
437 {
438 switch( np-nm )
439 {
440 case 0:
441 if( G4UniformRand() < 0.333333 )
442 {
443 currentParticle.SetDefinitionAndUpdateE( aNeutron );
444 targetParticle.SetDefinitionAndUpdateE( aProton );
445 incidentHasChanged = true;
446 targetHasChanged = true;
447 }
448 break;
449 case 1:
450 currentParticle.SetDefinitionAndUpdateE( aNeutron );
451 incidentHasChanged = true;
452 break;
453 default:
454 targetParticle.SetDefinitionAndUpdateE( aProton );
455 targetHasChanged = true;
456 break;
457 }
458 }
459 SetUpPions( np, nm, nz, vec, vecLen );
460 return;
461 }
462
463 /* end of file */
464
Note: See TracBrowser for help on using the repository browser.