source: trunk/source/processes/hadronic/models/rpg/src/G4RPGAntiNeutronInelastic.cc@ 1199

Last change on this file since 1199 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 18.5 KB
RevLine 
[819]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// $Id: G4RPGAntiNeutronInelastic.cc,v 1.1 2007/07/18 21:04:20 dennis Exp $
[1196]27// GEANT4 tag $Name: geant4-09-03-cand-01 $
[819]28//
29
30#include "G4RPGAntiNeutronInelastic.hh"
31#include "Randomize.hh"
32
33G4HadFinalState*
34G4RPGAntiNeutronInelastic::ApplyYourself( const G4HadProjectile &aTrack,
35 G4Nucleus &targetNucleus )
36{
37 const G4HadProjectile *originalIncident = &aTrack;
38
39 // create the target particle
40
41 G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
42
43 if( verboseLevel > 1 )
44 {
45 const G4Material *targetMaterial = aTrack.GetMaterial();
46 G4cout << "G4RPGAntiNeutronInelastic::ApplyYourself called" << G4endl;
47 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
48 G4cout << "target material = " << targetMaterial->GetName() << ", ";
49 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
50 << G4endl;
51 }
52 //
53 // Fermi motion and evaporation
54 // As of Geant3, the Fermi energy calculation had not been Done
55 //
56 G4double ek = originalIncident->GetKineticEnergy()/MeV;
57 G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
58 G4ReactionProduct modifiedOriginal;
59 modifiedOriginal = *originalIncident;
60
61 G4double tkin = targetNucleus.Cinema( ek );
62 ek += tkin;
63 modifiedOriginal.SetKineticEnergy( ek*MeV );
64 G4double et = ek + amas;
65 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
66 G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
67 if( pp > 0.0 )
68 {
69 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
70 modifiedOriginal.SetMomentum( momentum * (p/pp) );
71 }
72 //
73 // calculate black track energies
74 //
75 tkin = targetNucleus.EvaporationEffects( ek );
76 ek -= tkin;
77 modifiedOriginal.SetKineticEnergy( ek*MeV );
78 et = ek + amas;
79 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
80 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 G4ReactionProduct currentParticle = modifiedOriginal;
88 G4ReactionProduct targetParticle;
89 targetParticle = *originalTarget;
90 currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
91 targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
92 G4bool incidentHasChanged = false;
93 G4bool targetHasChanged = false;
94 G4bool quasiElastic = false;
95 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
96 G4int vecLen = 0;
97 vec.Initialize( 0 );
98
99 const G4double cutOff = 0.1*MeV;
100 const G4double anni = std::min( 1.3*currentParticle.GetTotalMomentum()/GeV, 0.4 );
101
102 if( (currentParticle.GetKineticEnergy()/MeV > cutOff) ||
103 (G4UniformRand() > anni) )
104 Cascade( vec, vecLen,
105 originalIncident, currentParticle, targetParticle,
106 incidentHasChanged, targetHasChanged, quasiElastic );
107 else
108 quasiElastic = true;
109
110 CalculateMomenta( vec, vecLen,
111 originalIncident, originalTarget, modifiedOriginal,
112 targetNucleus, currentParticle, targetParticle,
113 incidentHasChanged, targetHasChanged, quasiElastic );
114
115 SetUpChange( vec, vecLen,
116 currentParticle, targetParticle,
117 incidentHasChanged );
118
119 delete originalTarget;
120 return &theParticleChange;
121}
122
123
124void G4RPGAntiNeutronInelastic::Cascade(
125 G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
126 G4int& vecLen,
127 const G4HadProjectile *originalIncident,
128 G4ReactionProduct &currentParticle,
129 G4ReactionProduct &targetParticle,
130 G4bool &incidentHasChanged,
131 G4bool &targetHasChanged,
132 G4bool &quasiElastic )
133{
134 // Derived from H. Fesefeldt's original FORTRAN code CASNB
135 // AntiNeutron undergoes interaction with nucleon within a nucleus. Check if it is
136 // energetically possible to produce pions/kaons. In not, assume nuclear excitation
137 // occurs and input particle is degraded in energy. No other particles are produced.
138 // If reaction is possible, find the correct number of pions/protons/neutrons
139 // produced using an interpolation to multiplicity data. Replace some pions or
140 // protons/neutrons by kaons or strange baryons according to the average
141 // multiplicity per Inelastic reaction.
142
143 const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
144 const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
145 const G4double pOriginal = originalIncident->GetTotalMomentum()/MeV;
146 const G4double targetMass = targetParticle.GetMass()/MeV;
147 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
148 targetMass*targetMass +
149 2.0*targetMass*etOriginal );
150 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
151
152 static G4bool first = true;
153 const G4int numMul = 1200;
154 const G4int numMulA = 400;
155 const G4int numSec = 60;
156 static G4double protmul[numMul], protnorm[numSec]; // proton constants
157 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
158 static G4double protmulA[numMulA], protnormA[numSec]; // proton constants
159 static G4double neutmulA[numMulA], neutnormA[numSec]; // neutron constants
160 // np = number of pi+, nm = number of pi-, nz = number of pi0
161 G4int counter, nt=0, np=0, nm=0, nz=0;
162 G4double test;
163 const G4double c = 1.25;
164 const G4double b[] = { 0.70, 0.70 };
165 if( first ) // compute normalization constants, this will only be Done once
166 {
167 first = false;
168 G4int i;
169 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
170 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
171 counter = -1;
172 for( np=0; np<(numSec/3); ++np )
173 {
174 for( nm=std::max(0,np-2); nm<=np; ++nm )
175 {
176 for( nz=0; nz<numSec/3; ++nz )
177 {
178 if( ++counter < numMul )
179 {
180 nt = np+nm+nz;
181 if( nt>0 && nt<=numSec )
182 {
183 protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
184 protnorm[nt-1] += protmul[counter];
185 }
186 }
187 }
188 }
189 }
190 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
191 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
192 counter = -1;
193 for( np=0; np<numSec/3; ++np )
194 {
195 for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
196 {
197 for( nz=0; nz<numSec/3; ++nz )
198 {
199 if( ++counter < numMul )
200 {
201 nt = np+nm+nz;
202 if( (nt>0) && (nt<=numSec) )
203 {
204 neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
205 neutnorm[nt-1] += neutmul[counter];
206 }
207 }
208 }
209 }
210 }
211 for( i=0; i<numSec; ++i )
212 {
213 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
214 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
215 }
216 //
217 // do the same for annihilation channels
218 //
219 for( i=0; i<numMulA; ++i )protmulA[i] = 0.0;
220 for( i=0; i<numSec; ++i )protnormA[i] = 0.0;
221 counter = -1;
222 for( np=1; np<(numSec/3); ++np )
223 {
224 nm = np-1;
225 for( nz=0; nz<numSec/3; ++nz )
226 {
227 if( ++counter < numMulA )
228 {
229 nt = np+nm+nz;
230 if( nt>1 && nt<=numSec )
231 {
232 protmulA[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
233 protnormA[nt-1] += protmulA[counter];
234 }
235 }
236 }
237 }
238 for( i=0; i<numMulA; ++i )neutmulA[i] = 0.0;
239 for( i=0; i<numSec; ++i )neutnormA[i] = 0.0;
240 counter = -1;
241 for( np=0; np<numSec/3; ++np )
242 {
243 nm = np;
244 for( nz=0; nz<numSec/3; ++nz )
245 {
246 if( ++counter < numMulA )
247 {
248 nt = np+nm+nz;
249 if( nt>1 && nt<=numSec )
250 {
251 neutmulA[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
252 neutnormA[nt-1] += neutmulA[counter];
253 }
254 }
255 }
256 }
257 for( i=0; i<numSec; ++i )
258 {
259 if( protnormA[i] > 0.0 )protnormA[i] = 1.0/protnormA[i];
260 if( neutnormA[i] > 0.0 )neutnormA[i] = 1.0/neutnormA[i];
261 }
262 } // end of initialization
263 const G4double expxu = 82.; // upper bound for arg. of exp
264 const G4double expxl = -expxu; // lower bound for arg. of exp
265 G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
266 G4ParticleDefinition *aProton = G4Proton::Proton();
267 G4ParticleDefinition *anAntiProton = G4AntiProton::AntiProton();
268 G4ParticleDefinition *aPiPlus = G4PionPlus::PionPlus();
269
270 // energetically possible to produce pion(s) --> inelastic scattering
271 // otherwise quasi-elastic scattering
272
273 const G4double anhl[] = {1.00,1.00,1.00,1.00,1.00,1.00,1.00,1.00,0.97,0.88,
274 0.85,0.81,0.75,0.64,0.64,0.55,0.55,0.45,0.47,0.40,
275 0.39,0.36,0.33,0.10,0.01};
276 G4int iplab = G4int( pOriginal/GeV*10.0 );
277 if( iplab > 9 )iplab = G4int( (pOriginal/GeV- 1.0)*5.0 ) + 10;
278 if( iplab > 14 )iplab = G4int( pOriginal/GeV- 2.0 ) + 15;
279 if( iplab > 22 )iplab = G4int( (pOriginal/GeV-10.0)/10.0 ) + 23;
280 if( iplab > 24 )iplab = 24;
281 if( G4UniformRand() > anhl[iplab] )
282 {
283 if( availableEnergy <= aPiPlus->GetPDGMass()/MeV )
284 {
285 quasiElastic = true;
286 return;
287 }
288 G4int ieab = static_cast<G4int>(availableEnergy*5.0/GeV);
289 const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
290 G4double w0, wp, wt, wm;
291 if( (availableEnergy < 2.0*GeV) && (G4UniformRand() >= supp[ieab]) )
292 {
293 // suppress high multiplicity events at low momentum
294 // only one pion will be produced
295 //
296 np = nm = nz = 0;
297 if( targetParticle.GetDefinition() == aProton )
298 {
299 test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
300 w0 = test;
301 wp = test;
302 if( G4UniformRand() < w0/(w0+wp) )
303 nz = 1;
304 else
305 np = 1;
306 }
307 else // target is a neutron
308 {
309 test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[1])*(1.0+b[1])/(2.0*c*c) ) ) );
310 w0 = test;
311 wp = test;
312 test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[1])*(-1.0+b[1])/(2.0*c*c) ) ) );
313 wm = test;
314 wt = w0+wp+wm;
315 wp += w0;
316 G4double ran = G4UniformRand();
317 if( ran < w0/wt )
318 nz = 1;
319 else if( ran < wp/wt )
320 np = 1;
321 else
322 nm = 1;
323 }
324 }
325 else // (availableEnergy >= 2.0*GeV) || (random number < supp[ieab])
326 {
327 G4double n, anpn;
328 GetNormalizationConstant( availableEnergy, n, anpn );
329 G4double ran = G4UniformRand();
330 G4double dum, excs = 0.0;
331 if( targetParticle.GetDefinition() == aProton )
332 {
333 counter = -1;
334 for( np=0; np<numSec/3 && ran>=excs; ++np )
335 {
336 for( nm=std::max(0,np-2); nm<=np && ran>=excs; ++nm )
337 {
338 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
339 {
340 if( ++counter < numMul )
341 {
342 nt = np+nm+nz;
343 if( nt > 0 )
344 {
345 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
346 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
347 if( std::fabs(dum) < 1.0 )
348 {
349 if( test >= 1.0e-10 )excs += dum*test;
350 }
351 else
352 excs += dum*test;
353 }
354 }
355 }
356 }
357 }
358 if( ran >= excs ) // 3 previous loops continued to the end
359 {
360 quasiElastic = true;
361 return;
362 }
363 np--; nm--; nz--;
364 }
365 else // target must be a neutron
366 {
367 counter = -1;
368 for( np=0; np<numSec/3 && ran>=excs; ++np )
369 {
370 for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm )
371 {
372 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
373 {
374 if( ++counter < numMul )
375 {
376 nt = np+nm+nz;
377 if( (nt>=1) && (nt<=numSec) )
378 {
379 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
380 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
381 if( std::fabs(dum) < 1.0 )
382 {
383 if( test >= 1.0e-10 )excs += dum*test;
384 }
385 else
386 excs += dum*test;
387 }
388 }
389 }
390 }
391 }
392 if( ran >= excs ) // 3 previous loops continued to the end
393 {
394 quasiElastic = true;
395 return;
396 }
397 np--; nm--; nz--;
398 }
399 }
400 if( targetParticle.GetDefinition() == aProton )
401 {
402 switch( np-nm )
403 {
404 case 1:
405 if( G4UniformRand() < 0.5 )
406 {
407 currentParticle.SetDefinitionAndUpdateE( anAntiProton );
408 incidentHasChanged = true;
409 }
410 else
411 {
412 targetParticle.SetDefinitionAndUpdateE( aNeutron );
413 targetHasChanged = true;
414 }
415 break;
416 case 2:
417 currentParticle.SetDefinitionAndUpdateE( anAntiProton );
418 targetParticle.SetDefinitionAndUpdateE( aNeutron );
419 incidentHasChanged = true;
420 targetHasChanged = true;
421 break;
422 default:
423 break;
424 }
425 }
426 else // target must be a neutron
427 {
428 switch( np-nm )
429 {
430 case 0:
431 if( G4UniformRand() < 0.33 )
432 {
433 currentParticle.SetDefinitionAndUpdateE( anAntiProton );
434 targetParticle.SetDefinitionAndUpdateE( aProton );
435 incidentHasChanged = true;
436 targetHasChanged = true;
437 }
438 break;
439 case 1:
440 currentParticle.SetDefinitionAndUpdateE( anAntiProton );
441 incidentHasChanged = true;
442 break;
443 default:
444 targetParticle.SetDefinitionAndUpdateE( aProton );
445 targetHasChanged = true;
446 break;
447 }
448 }
449 }
450 else // random number <= anhl[iplab]
451 {
452 if( centerofmassEnergy <= 2*aPiPlus->GetPDGMass()/MeV )
453 {
454 quasiElastic = true;
455 return;
456 }
457 //
458 // annihilation channels
459 //
460 G4double n, anpn;
461 GetNormalizationConstant( -centerofmassEnergy, n, anpn );
462 G4double ran = G4UniformRand();
463 G4double dum, excs = 0.0;
464 if( targetParticle.GetDefinition() == aProton )
465 {
466 counter = -1;
467 for( np=1; (np<numSec/3) && (ran>=excs); ++np )
468 {
469 nm = np-1;
470 for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
471 {
472 if( ++counter < numMulA )
473 {
474 nt = np+nm+nz;
475 if( nt>1 && nt<=numSec )
476 {
477 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
478 dum = (pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
479 if( std::fabs(dum) < 1.0 )
480 {
481 if( test >= 1.0e-10 )excs += dum*test;
482 }
483 else
484 excs += dum*test;
485 }
486 }
487 }
488 }
489 }
490 else // target must be a neutron
491 {
492 counter = -1;
493 for( np=0; (np<numSec/3) && (ran>=excs); ++np )
494 {
495 nm = np;
496 for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
497 {
498 if( ++counter < numMulA )
499 {
500 nt = np+nm+nz;
501 if( (nt>1) && (nt<=numSec) )
502 {
503 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
504 dum = (pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
505 if( std::fabs(dum) < 1.0 )
506 {
507 if( test >= 1.0e-10 )excs += dum*test;
508 }
509 else
510 excs += dum*test;
511 }
512 }
513 }
514 }
515 }
516 if( ran >= excs ) // 3 previous loops continued to the end
517 {
518 quasiElastic = true;
519 return;
520 }
521 np--; nz--;
522 currentParticle.SetMass( 0.0 );
523 targetParticle.SetMass( 0.0 );
524 }
525 while(np+nm+nz<3) nz++;
526
527 SetUpPions( np, nm, nz, vec, vecLen );
528 return;
529}
530
531 /* end of file */
532
Note: See TracBrowser for help on using the repository browser.