source: trunk/source/processes/hadronic/models/rpg/src/G4RPGLambdaInelastic.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: 13.5 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: G4RPGLambdaInelastic.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 "G4RPGLambdaInelastic.hh"
31#include "Randomize.hh"
32
33G4HadFinalState*
34G4RPGLambdaInelastic::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 << "G4RPGLambdaInelastic::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 if( currentParticle.GetKineticEnergy()/MeV > cutOff )
101 Cascade( vec, vecLen,
102 originalIncident, currentParticle, targetParticle,
103 incidentHasChanged, targetHasChanged, quasiElastic );
104
105 CalculateMomenta( vec, vecLen,
106 originalIncident, originalTarget, modifiedOriginal,
107 targetNucleus, currentParticle, targetParticle,
108 incidentHasChanged, targetHasChanged, quasiElastic );
109
110 SetUpChange( vec, vecLen,
111 currentParticle, targetParticle,
112 incidentHasChanged );
113
114 delete originalTarget;
115 return &theParticleChange;
116}
117
118
119void G4RPGLambdaInelastic::Cascade(
120 G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
121 G4int& vecLen,
122 const G4HadProjectile *originalIncident,
123 G4ReactionProduct &currentParticle,
124 G4ReactionProduct &targetParticle,
125 G4bool &incidentHasChanged,
126 G4bool &targetHasChanged,
127 G4bool &quasiElastic )
128{
129 // Derived from H. Fesefeldt's original FORTRAN code CASL0
130 //
131 // Lambda 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 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
143 targetMass*targetMass +
144 2.0*targetMass*etOriginal );
145 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
146 if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
147 {
148 quasiElastic = true;
149 return;
150 }
151 static G4bool first = true;
152 const G4int numMul = 1200;
153 const G4int numSec = 60;
154 static G4double protmul[numMul], protnorm[numSec]; // proton constants
155 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
156
157 // np = number of pi+, nm = number of pi-, nz = number of pi0
158
159 G4int counter, nt=0, np=0, nm=0, nz=0;
160 G4double test;
161 const G4double c = 1.25;
162 const G4double b[] = { 0.70, 0.35 };
163 if( first ) { // compute normalization constants, this will only be Done once
164 first = false;
165 G4int i;
166 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
167 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
168 counter = -1;
169 for( np=0; np<(numSec/3); ++np ) {
170 for( nm=std::max(0,np-2); nm<=(np+1); ++nm ) {
171 for( nz=0; nz<numSec/3; ++nz ) {
172 if( ++counter < numMul ) {
173 nt = np+nm+nz;
174 if( nt>0 && nt<=numSec ) {
175 protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
176 protnorm[nt-1] += protmul[counter];
177 }
178 }
179 }
180 }
181 }
182 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
183 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
184 counter = -1;
185 for( np=0; np<numSec/3; ++np ) {
186 for( nm=std::max(0,np-1); nm<=(np+2); ++nm ) {
187 for( nz=0; nz<numSec/3; ++nz ) {
188 if( ++counter < numMul ) {
189 nt = np+nm+nz;
190 if( nt>0 && nt<=numSec ) {
191 neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
192 neutnorm[nt-1] += neutmul[counter];
193 }
194 }
195 }
196 }
197 }
198 for( i=0; i<numSec; ++i ) {
199 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
200 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
201 }
202 } // end of initialization
203
204 const G4double expxu = 82.; // upper bound for arg. of exp
205 const G4double expxl = -expxu; // lower bound for arg. of exp
206 G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
207 G4ParticleDefinition *aProton = G4Proton::Proton();
208 G4ParticleDefinition *aSigmaPlus = G4SigmaPlus::SigmaPlus();
209 G4ParticleDefinition *aSigmaMinus = G4SigmaMinus::SigmaMinus();
210 G4ParticleDefinition *aSigmaZero = G4SigmaZero::SigmaZero();
211
212 // energetically possible to produce pion(s) --> inelastic scattering
213 // otherwise quasi-elastic scattering
214
215 G4double n, anpn;
216 GetNormalizationConstant( availableEnergy, n, anpn );
217 G4double ran = G4UniformRand();
218 G4double dum, excs = 0.0;
219 if( targetParticle.GetDefinition() == aProton ) {
220 counter = -1;
221 for( np=0; np<numSec/3 && ran>=excs; ++np ) {
222 for( nm=std::max(0,np-2); nm<=(np+1) && ran>=excs; ++nm ) {
223 for( nz=0; nz<numSec/3 && ran>=excs; ++nz ) {
224 if( ++counter < numMul ) {
225 nt = np+nm+nz;
226 if( nt>0 && nt<=numSec ) {
227 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
228 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
229 if( std::fabs(dum) < 1.0 ) {
230 if( test >= 1.0e-10 )excs += dum*test;
231 } else {
232 excs += dum*test;
233 }
234 }
235 }
236 }
237 }
238 }
239 if( ran >= excs ) // 3 previous loops continued to the end
240 {
241 quasiElastic = true;
242 return;
243 }
244 np--; nm--; nz--;
245 G4int ncht = std::max( 1, np-nm );
246 switch( ncht ) {
247 case 1:
248 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
249 incidentHasChanged = true;
250 break;
251 case 2:
252 if( G4UniformRand() < 0.5 ) {
253 if( G4UniformRand() < 0.5 ) {
254 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
255 incidentHasChanged = true;
256 }
257 } else {
258 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
259 incidentHasChanged = true;
260 targetParticle.SetDefinitionAndUpdateE( aNeutron );
261 targetHasChanged = true;
262 }
263 break;
264 case 3:
265 if( G4UniformRand() < 0.5 ) {
266 if( G4UniformRand() < 0.5 ) {
267 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
268 incidentHasChanged = true;
269 }
270 targetParticle.SetDefinitionAndUpdateE( aNeutron );
271 targetHasChanged = true;
272 } else {
273 currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
274 incidentHasChanged = true;
275 }
276 break;
277 default:
278 currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
279 incidentHasChanged = true;
280 targetParticle.SetDefinitionAndUpdateE( aNeutron );
281 targetHasChanged = true;
282 break;
283 }
284 }
285 else // target must be a neutron
286 {
287 counter = -1;
288 for( np=0; np<numSec/3 && ran>=excs; ++np ) {
289 for( nm=std::max(0,np-1); nm<=(np+2) && ran>=excs; ++nm ) {
290 for( nz=0; nz<numSec/3 && ran>=excs; ++nz ) {
291 if( ++counter < numMul ) {
292 nt = np+nm+nz;
293 if( nt>0 && nt<=numSec ) {
294 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
295 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
296 if( std::fabs(dum) < 1.0 ) {
297 if( test >= 1.0e-10 )excs += dum*test;
298 } else {
299 excs += dum*test;
300 }
301 }
302 }
303 }
304 }
305 }
306 if( ran >= excs ) // 3 previous loops continued to the end
307 {
308 quasiElastic = true;
309 return;
310 }
311 np--; nm--; nz--;
312 G4int ncht = std::max( 1, np-nm+3 );
313 switch( ncht ) {
314 case 1:
315 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
316 incidentHasChanged = true;
317 targetParticle.SetDefinitionAndUpdateE( aProton );
318 targetHasChanged = true;
319 break;
320 case 2:
321 if( G4UniformRand() < 0.5 ) {
322 if( G4UniformRand() < 0.5 ) {
323 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
324 incidentHasChanged = true;
325 }
326 targetParticle.SetDefinitionAndUpdateE( aProton );
327 targetHasChanged = true;
328 } else {
329 currentParticle.SetDefinitionAndUpdateE( aSigmaPlus );
330 incidentHasChanged = true;
331 }
332 break;
333 case 3:
334 if( G4UniformRand() < 0.5 ) {
335 if( G4UniformRand() < 0.5 ) {
336 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
337 incidentHasChanged = true;
338 }
339 } else {
340 currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
341 incidentHasChanged = true;
342 targetParticle.SetDefinitionAndUpdateE( aProton );
343 targetHasChanged = true;
344 }
345 break;
346 default:
347 currentParticle.SetDefinitionAndUpdateE( aSigmaMinus );
348 incidentHasChanged = true;
349 break;
350 }
351 }
352 SetUpPions( np, nm, nz, vec, vecLen );
353 return;
354}
355
356 /* end of file */
357
Note: See TracBrowser for help on using the repository browser.