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