source: trunk/source/processes/hadronic/models/rpg/src/G4RPGAntiLambdaInelastic.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: 21.8 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: G4RPGAntiLambdaInelastic.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 "G4RPGAntiLambdaInelastic.hh"
31#include "Randomize.hh"
32
33G4HadFinalState*
34G4RPGAntiLambdaInelastic::ApplyYourself( const G4HadProjectile &aTrack,
35 G4Nucleus &targetNucleus )
36{
37 const G4HadProjectile *originalIncident = &aTrack;
38
39 // Choose the target particle
40
41 G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
42
43 if( verboseLevel > 1 )
44 {
45 const G4Material *targetMaterial = aTrack.GetMaterial();
46 G4cout << "G4RPGAntiLambdaInelastic::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;
100 const G4double anni = std::min( 1.3*currentParticle.GetTotalMomentum()/GeV, 0.4 );
101 if( (originalIncident->GetKineticEnergy()/MeV > cutOff) || (G4UniformRand() > anni) )
102 Cascade( vec, vecLen,
103 originalIncident, currentParticle, targetParticle,
104 incidentHasChanged, targetHasChanged, quasiElastic );
105
106 CalculateMomenta( vec, vecLen,
107 originalIncident, originalTarget, modifiedOriginal,
108 targetNucleus, currentParticle, targetParticle,
109 incidentHasChanged, targetHasChanged, quasiElastic );
110
111 SetUpChange( vec, vecLen,
112 currentParticle, targetParticle,
113 incidentHasChanged );
114
115 delete originalTarget;
116 return &theParticleChange;
117}
118
119
120void G4RPGAntiLambdaInelastic::Cascade(
121 G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
122 G4int &vecLen,
123 const G4HadProjectile *originalIncident,
124 G4ReactionProduct &currentParticle,
125 G4ReactionProduct &targetParticle,
126 G4bool &incidentHasChanged,
127 G4bool &targetHasChanged,
128 G4bool &quasiElastic )
129{
130 // Derived from H. Fesefeldt's original FORTRAN code CASAL0
131 // AntiLambda undergoes interaction with nucleon within a nucleus. Check if it is
132 // energetically possible to produce pions/kaons. In not, assume nuclear excitation
133 // occurs and input particle is degraded in energy. No other particles are produced.
134 // If reaction is possible, find the correct number of pions/protons/neutrons
135 // produced using an interpolation to multiplicity data. Replace some pions or
136 // protons/neutrons by kaons or strange baryons according to the average
137 // multiplicity per Inelastic reaction.
138
139 const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
140 const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
141 const G4double targetMass = targetParticle.GetMass()/MeV;
142 const G4double pOriginal = originalIncident->GetTotalMomentum()/GeV;
143 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
144 targetMass*targetMass +
145 2.0*targetMass*etOriginal );
146 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
147
148 static G4bool first = true;
149 const G4int numMul = 1200;
150 const G4int numMulA = 400;
151 const G4int numSec = 60;
152 static G4double protmul[numMul], protnorm[numSec]; // proton constants
153 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
154 static G4double protmulA[numMulA], protnormA[numSec]; // proton constants
155 static G4double neutmulA[numMulA], neutnormA[numSec]; // neutron constants
156 // np = number of pi+, nm = number of pi-, nz = number of pi0
157 G4int nt=0, np=0, nm=0, nz=0;
158 G4double test;
159 const G4double c = 1.25;
160 const G4double b[] = { 0.7, 0.7 };
161 if( first ) // compute normalization constants, this will only be Done once
162 {
163 first = false;
164 G4int i;
165 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
166 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
167 G4int counter = -1;
168 for( np=0; np<(numSec/3); ++np )
169 {
170 for( nm=std::max(0,np-2); nm<=(np+1); ++nm )
171 {
172 for( nz=0; nz<numSec/3; ++nz )
173 {
174 if( ++counter < numMul )
175 {
176 nt = np+nm+nz;
177 if( nt>0 && nt<=numSec )
178 {
179 protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
180 protnorm[nt-1] += protmul[counter];
181 }
182 }
183 }
184 }
185 }
186 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
187 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
188 counter = -1;
189 for( np=0; np<numSec/3; ++np )
190 {
191 for( nm=std::max(0,np-1); nm<=(np+2); ++nm )
192 {
193 for( nz=0; nz<numSec/3; ++nz )
194 {
195 if( ++counter < numMul )
196 {
197 nt = np+nm+nz;
198 if( nt>0 && nt<=numSec )
199 {
200 neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
201 neutnorm[nt-1] += neutmul[counter];
202 }
203 }
204 }
205 }
206 }
207 for( i=0; i<numSec; ++i )
208 {
209 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
210 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
211 }
212 //
213 // do the same for annihilation channels
214 //
215 for( i=0; i<numMulA; ++i )protmulA[i] = 0.0;
216 for( i=0; i<numSec; ++i )protnormA[i] = 0.0;
217 counter = -1;
218 for( np=1; np<(numSec/3); ++np )
219 {
220 nm = np-1;
221 for( nz=0; nz<numSec/3; ++nz )
222 {
223 if( ++counter < numMulA )
224 {
225 nt = np+nm+nz;
226 if( nt>1 && nt<=numSec )
227 {
228 protmulA[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
229 protnormA[nt-1] += protmulA[counter];
230 }
231 }
232 }
233 }
234 for( i=0; i<numMulA; ++i )neutmulA[i] = 0.0;
235 for( i=0; i<numSec; ++i )neutnormA[i] = 0.0;
236 counter = -1;
237 for( np=0; np<numSec/3; ++np )
238 {
239 nm = np;
240 for( nz=0; nz<numSec/3; ++nz )
241 {
242 if( ++counter < numMulA )
243 {
244 nt = np+nm+nz;
245 if( nt>1 && nt<=numSec )
246 {
247 neutmulA[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
248 neutnormA[nt-1] += neutmulA[counter];
249 }
250 }
251 }
252 }
253 for( i=0; i<numSec; ++i )
254 {
255 if( protnormA[i] > 0.0 )protnormA[i] = 1.0/protnormA[i];
256 if( neutnormA[i] > 0.0 )neutnormA[i] = 1.0/neutnormA[i];
257 }
258 } // end of initialization
259 const G4double expxu = 82.; // upper bound for arg. of exp
260 const G4double expxl = -expxu; // lower bound for arg. of exp
261
262 G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
263 G4ParticleDefinition *aProton = G4Proton::Proton();
264 G4ParticleDefinition *aPiPlus = G4PionPlus::PionPlus();
265 G4ParticleDefinition *aPiMinus = G4PionMinus::PionMinus();
266 G4ParticleDefinition *aPiZero = G4PionZero::PionZero();
267 G4ParticleDefinition *aKaonPlus = G4KaonPlus::KaonPlus();
268 G4ParticleDefinition *aKaonMinus = G4KaonMinus::KaonMinus();
269 G4ParticleDefinition *aKaonZL = G4KaonZeroLong::KaonZeroLong();
270 G4ParticleDefinition *anAntiSigmaZero = G4AntiSigmaZero::AntiSigmaZero();
271 G4ParticleDefinition *anAntiSigmaPlus = G4AntiSigmaPlus::AntiSigmaPlus();
272 G4ParticleDefinition *anAntiSigmaMinus = G4AntiSigmaMinus::AntiSigmaMinus();
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*10.0 );
277 if( iplab > 9 )iplab = G4int( (pOriginal- 1.0)*5.0 ) + 10;
278 if( iplab > 14 )iplab = G4int( pOriginal- 2.0 ) + 15;
279 if( iplab > 22 )iplab = G4int( (pOriginal-10.0)/10.0 ) + 23;
280 if( iplab > 24 )iplab = 24;
281 if( G4UniformRand() > anhl[iplab] )
282 {
283 if( availableEnergy <= aPiPlus->GetPDGMass()/MeV )
284 { // not energetically possible to produce pion(s)
285 quasiElastic = true;
286 return;
287 }
288 G4double n, anpn;
289 GetNormalizationConstant( availableEnergy, n, anpn );
290 G4double ran = G4UniformRand();
291 G4double dum, excs = 0.0;
292 if( targetParticle.GetDefinition() == aProton )
293 {
294 G4int counter = -1;
295 for( np=0; np<numSec/3 && ran>=excs; ++np )
296 {
297 for( nm=std::max(0,np-2); nm<=(np+1) && ran>=excs; ++nm )
298 {
299 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
300 {
301 if( ++counter < numMul )
302 {
303 nt = np+nm+nz;
304 if( nt>0 && nt<=numSec )
305 {
306 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
307 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
308 if( std::fabs(dum) < 1.0 )
309 {
310 if( test >= 1.0e-10 )excs += dum*test;
311 }
312 else
313 excs += dum*test;
314 }
315 }
316 }
317 }
318 }
319 if( ran >= excs ) // 3 previous loops continued to the end
320 {
321 quasiElastic = true;
322 return;
323 }
324 np--; nm--; nz--;
325 G4int ncht = std::min( 4, std::max( 1, np-nm+2 ) );
326 switch( ncht )
327 {
328 case 1:
329 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
330 incidentHasChanged = true;
331 break;
332 case 2:
333 if( G4UniformRand() < 0.5 )
334 {
335 if( G4UniformRand() < 0.5 )
336 {
337 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
338 incidentHasChanged = true;
339 }
340 else
341 {
342 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
343 incidentHasChanged = true;
344 targetParticle.SetDefinitionAndUpdateE( aNeutron );
345 targetHasChanged = true;
346 }
347 }
348 else
349 {
350 if( G4UniformRand() >= 0.5 )
351 {
352 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
353 incidentHasChanged = true;
354 targetParticle.SetDefinitionAndUpdateE( aNeutron );
355 targetHasChanged = true;
356 }
357 }
358 break;
359 case 3:
360 if( G4UniformRand() < 0.5 )
361 {
362 if( G4UniformRand() < 0.5 )
363 {
364 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
365 incidentHasChanged = true;
366 targetParticle.SetDefinitionAndUpdateE( aNeutron );
367 targetHasChanged = true;
368 }
369 else
370 {
371 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
372 incidentHasChanged = true;
373 }
374 }
375 else
376 {
377 if( G4UniformRand() < 0.5 )
378 {
379 targetParticle.SetDefinitionAndUpdateE( aNeutron );
380 targetHasChanged = true;
381 }
382 else
383 {
384 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
385 incidentHasChanged = true;
386 }
387 }
388 break;
389 case 4:
390 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
391 incidentHasChanged = true;
392 targetParticle.SetDefinitionAndUpdateE( aNeutron );
393 targetHasChanged = true;
394 break;
395 }
396 }
397 else // target must be a neutron
398 {
399 G4int counter = -1;
400 for( np=0; np<numSec/3 && ran>=excs; ++np )
401 {
402 for( nm=std::max(0,np-1); nm<=(np+2) && ran>=excs; ++nm )
403 {
404 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
405 {
406 if( ++counter < numMul )
407 {
408 nt = np+nm+nz;
409 if( nt>0 && nt<=numSec )
410 {
411 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
412 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
413 if( std::fabs(dum) < 1.0 )
414 {
415 if( test >= 1.0e-10 )excs += dum*test;
416 }
417 else
418 excs += dum*test;
419 }
420 }
421 }
422 }
423 }
424 if( ran >= excs ) // 3 previous loops continued to the end
425 {
426 quasiElastic = true;
427 return;
428 }
429 np--; nm--; nz--;
430 G4int ncht = std::min( 4, std::max( 1, np-nm+3 ) );
431 switch( ncht )
432 {
433 case 1:
434 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
435 incidentHasChanged = true;
436 targetParticle.SetDefinitionAndUpdateE( aProton );
437 targetHasChanged = true;
438 break;
439 case 2:
440 if( G4UniformRand() < 0.5 )
441 {
442 if( G4UniformRand() < 0.5 )
443 {
444 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
445 incidentHasChanged = true;
446 targetParticle.SetDefinitionAndUpdateE( aProton );
447 targetHasChanged = true;
448 }
449 else
450 {
451 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
452 incidentHasChanged = true;
453 }
454 }
455 else
456 {
457 if( G4UniformRand() < 0.5 )
458 {
459 targetParticle.SetDefinitionAndUpdateE( aProton );
460 targetHasChanged = true;
461 }
462 else
463 {
464 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaMinus );
465 incidentHasChanged = true;
466 }
467 }
468 break;
469 case 3:
470 if( G4UniformRand() < 0.5 )
471 {
472 if( G4UniformRand() < 0.5 )
473 {
474 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
475 incidentHasChanged = true;
476 }
477 else
478 {
479 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
480 incidentHasChanged = true;
481 targetParticle.SetDefinitionAndUpdateE( aProton );
482 targetHasChanged = true;
483 }
484 }
485 else
486 {
487 if( G4UniformRand() >= 0.5 )
488 {
489 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
490 incidentHasChanged = true;
491 targetParticle.SetDefinitionAndUpdateE( aProton );
492 targetHasChanged = true;
493 }
494 }
495 break;
496 default:
497 currentParticle.SetDefinitionAndUpdateE( anAntiSigmaPlus );
498 incidentHasChanged = true;
499 break;
500 }
501 }
502 }
503 else // random number <= anhl[iplab]
504 {
505 if( centerofmassEnergy <= aPiPlus->GetPDGMass()/MeV+aKaonPlus->GetPDGMass()/MeV )
506 {
507 quasiElastic = true;
508 return;
509 }
510 //
511 // annihilation channels
512 //
513 G4double n, anpn;
514 GetNormalizationConstant( -centerofmassEnergy, n, anpn );
515 G4double ran = G4UniformRand();
516 G4double dum, excs = 0.0;
517 if( targetParticle.GetDefinition() == aProton )
518 {
519 G4int counter = -1;
520 for( np=1; np<numSec/3 && ran>=excs; ++np )
521 {
522 nm = np-1;
523 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
524 {
525 if( ++counter < numMulA )
526 {
527 nt = np+nm+nz;
528 if( nt>1 && nt<=numSec )
529 {
530 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
531 dum = (pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
532 if( std::fabs(dum) < 1.0 )
533 {
534 if( test >= 1.0e-10 )excs += dum*test;
535 }
536 else
537 excs += dum*test;
538 }
539 }
540 }
541 }
542 }
543 else // target must be a neutron
544 {
545 G4int counter = -1;
546 for( np=0; np<numSec/3 && ran>=excs; ++np )
547 {
548 nm = np;
549 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
550 {
551 if( ++counter < numMulA )
552 {
553 nt = np+nm+nz;
554 if( nt>1 && nt<=numSec )
555 {
556 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
557 dum = (pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
558 if( std::fabs(dum) < 1.0 )
559 {
560 if( test >= 1.0e-10 )excs += dum*test;
561 }
562 else
563 excs += dum*test;
564 }
565 }
566 }
567 }
568 }
569 if( ran >= excs ) // 3 previous loops continued to the end
570 {
571 quasiElastic = true;
572 return;
573 }
574 np--; nz--;
575 currentParticle.SetMass( 0.0 );
576 targetParticle.SetMass( 0.0 );
577 }
578 SetUpPions( np, nm, nz, vec, vecLen );
579 if( currentParticle.GetMass() == 0.0 )
580 {
581 if( nz == 0 )
582 {
583 if( nm > 0 )
584 {
585 for( G4int i=0; i<vecLen; ++i )
586 {
587 if( vec[i]->GetDefinition() == aPiMinus )
588 {
589 vec[i]->SetDefinitionAndUpdateE( aKaonMinus );
590 break;
591 }
592 }
593 }
594 }
595 else // nz > 0
596 {
597 if( nm == 0 )
598 {
599 for( G4int i=0; i<vecLen; ++i )
600 {
601 if( vec[i]->GetDefinition() == aPiZero )
602 {
603 vec[i]->SetDefinitionAndUpdateE( aKaonZL );
604 break;
605 }
606 }
607 }
608 else // nm > 0
609 {
610 if( G4UniformRand() < 0.5 )
611 {
612 if( nm > 0 )
613 {
614 for( G4int i=0; i<vecLen; ++i )
615 {
616 if( vec[i]->GetDefinition() == aPiMinus )
617 {
618 vec[i]->SetDefinitionAndUpdateE( aKaonMinus );
619 break;
620 }
621 }
622 }
623 }
624 else // random number >= 0.5
625 {
626 for( G4int i=0; i<vecLen; ++i )
627 {
628 if( vec[i]->GetDefinition() == aPiZero )
629 {
630 vec[i]->SetDefinitionAndUpdateE( aKaonZL );
631 break;
632 }
633 }
634 }
635 }
636 }
637 }
638 return;
639}
640
641 /* end of file */
642
Note: See TracBrowser for help on using the repository browser.