source: trunk/source/processes/hadronic/models/low_energy/src/G4LEAntiProtonInelastic.cc@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

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