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