source: trunk/source/processes/hadronic/models/rpg/src/G4RPGSigmaMinusInelastic.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: 12.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// $Id: G4RPGSigmaMinusInelastic.cc,v 1.1 2007/07/18 21:04:21 dennis Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29
30#include "G4RPGSigmaMinusInelastic.hh"
31#include "Randomize.hh"
32
33G4HadFinalState*
34G4RPGSigmaMinusInelastic::ApplyYourself( const G4HadProjectile &aTrack,
35 G4Nucleus &targetNucleus )
36{
37 const G4HadProjectile *originalIncident = &aTrack;
38 if (originalIncident->GetKineticEnergy()<= 0.1*MeV)
39 {
40 theParticleChange.SetStatusChange(isAlive);
41 theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
42 theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
43 return &theParticleChange;
44 }
45
46 // create the target particle
47
48 G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
49
50 if( verboseLevel > 1 )
51 {
52 const G4Material *targetMaterial = aTrack.GetMaterial();
53 G4cout << "G4RPGSigmaMinusInelastic::ApplyYourself called" << G4endl;
54 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
55 G4cout << "target material = " << targetMaterial->GetName() << ", ";
56 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
57 << G4endl;
58 }
59
60 // Fermi motion and evaporation
61 // As of Geant3, the Fermi energy calculation had not been Done
62
63 G4double ek = originalIncident->GetKineticEnergy()/MeV;
64 G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
65 G4ReactionProduct modifiedOriginal;
66 modifiedOriginal = *originalIncident;
67
68 G4double tkin = targetNucleus.Cinema( ek );
69 ek += tkin;
70 modifiedOriginal.SetKineticEnergy( ek*MeV );
71 G4double et = ek + amas;
72 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
73 G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
74 if( pp > 0.0 )
75 {
76 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
77 modifiedOriginal.SetMomentum( momentum * (p/pp) );
78 }
79 //
80 // calculate black track energies
81 //
82 tkin = targetNucleus.EvaporationEffects( ek );
83 ek -= tkin;
84 modifiedOriginal.SetKineticEnergy( ek*MeV );
85 et = ek + amas;
86 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
87 pp = modifiedOriginal.GetMomentum().mag()/MeV;
88 if( pp > 0.0 )
89 {
90 G4ThreeVector momentum = modifiedOriginal.GetMomentum();
91 modifiedOriginal.SetMomentum( momentum * (p/pp) );
92 }
93 G4ReactionProduct currentParticle = modifiedOriginal;
94 G4ReactionProduct targetParticle;
95 targetParticle = *originalTarget;
96 currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
97 targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
98 G4bool incidentHasChanged = false;
99 G4bool targetHasChanged = false;
100 G4bool quasiElastic = false;
101 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
102 G4int vecLen = 0;
103 vec.Initialize( 0 );
104
105 const G4double cutOff = 0.1;
106 if( originalIncident->GetKineticEnergy()/MeV > cutOff )
107 Cascade( vec, vecLen,
108 originalIncident, currentParticle, targetParticle,
109 incidentHasChanged, targetHasChanged, quasiElastic );
110
111 CalculateMomenta( vec, vecLen,
112 originalIncident, originalTarget, modifiedOriginal,
113 targetNucleus, currentParticle, targetParticle,
114 incidentHasChanged, targetHasChanged, quasiElastic );
115
116 SetUpChange( vec, vecLen,
117 currentParticle, targetParticle,
118 incidentHasChanged );
119
120 delete originalTarget;
121 return &theParticleChange;
122}
123
124
125void G4RPGSigmaMinusInelastic::Cascade(
126 G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
127 G4int& vecLen,
128 const G4HadProjectile *originalIncident,
129 G4ReactionProduct &currentParticle,
130 G4ReactionProduct &targetParticle,
131 G4bool &incidentHasChanged,
132 G4bool &targetHasChanged,
133 G4bool &quasiElastic )
134{
135 // Derived from H. Fesefeldt's original FORTRAN code CASSM
136 //
137 // SigmaMinus undergoes interaction with nucleon within a nucleus. Check if it is
138 // energetically possible to produce pions/kaons. In not, assume nuclear excitation
139 // occurs and input particle is degraded in energy. No other particles are produced.
140 // If reaction is possible, find the correct number of pions/protons/neutrons
141 // produced using an interpolation to multiplicity data. Replace some pions or
142 // protons/neutrons by kaons or strange baryons according to the average
143 // multiplicity per Inelastic reaction.
144
145 const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
146 const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
147 const G4double targetMass = targetParticle.GetMass()/MeV;
148 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
149 targetMass*targetMass +
150 2.0*targetMass*etOriginal );
151 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
152 if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
153 {
154 quasiElastic = true;
155 return;
156 }
157 static G4bool first = true;
158 const G4int numMul = 1200;
159 const G4int numSec = 60;
160 static G4double protmul[numMul], protnorm[numSec]; // proton constants
161 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
162 // np = number of pi+, nm = number of pi-, nz = number of pi0
163 G4int counter, nt=0, np=0, nm=0, nz=0;
164 G4double test;
165 const G4double c = 1.25;
166 const G4double b[] = { 0.70, 0.70 };
167 if( first ) // compute normalization constants, this will only be Done once
168 {
169 first = false;
170 G4int i;
171 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
172 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
173 counter = -1;
174 for( np=0; np<(numSec/3); ++np )
175 {
176 for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
177 {
178 for( nz=0; nz<numSec/3; ++nz )
179 {
180 if( ++counter < numMul )
181 {
182 nt = np+nm+nz;
183 if( nt>0 && nt<=numSec )
184 {
185 protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
186 protnorm[nt-1] += protmul[counter];
187 }
188 }
189 }
190 }
191 }
192 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
193 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
194 counter = -1;
195 for( np=0; np<numSec/3; ++np )
196 {
197 for( nm=np; nm<=(np+2); ++nm )
198 {
199 for( nz=0; nz<numSec/3; ++nz )
200 {
201 if( ++counter < numMul )
202 {
203 nt = np+nm+nz;
204 if( nt>0 && nt<=numSec )
205 {
206 neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
207 neutnorm[nt-1] += neutmul[counter];
208 }
209 }
210 }
211 }
212 }
213 for( i=0; i<numSec; ++i )
214 {
215 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
216 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
217 }
218 } // end of initialization
219
220 const G4double expxu = 82.; // upper bound for arg. of exp
221 const G4double expxl = -expxu; // lower bound for arg. of exp
222 G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
223 G4ParticleDefinition *aProton = G4Proton::Proton();
224 G4ParticleDefinition *aLambda = G4Lambda::Lambda();
225 G4ParticleDefinition *aSigmaZero = G4SigmaZero::SigmaZero();
226
227 // energetically possible to produce pion(s) --> inelastic scattering
228
229 G4double n, anpn;
230 GetNormalizationConstant( availableEnergy, n, anpn );
231 G4double ran = G4UniformRand();
232 G4double dum, excs = 0.0;
233 if( targetParticle.GetDefinition() == aProton )
234 {
235 counter = -1;
236 for( np=0; np<numSec/3 && ran>=excs; ++np )
237 {
238 for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm )
239 {
240 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
241 {
242 if( ++counter < numMul )
243 {
244 nt = np+nm+nz;
245 if( nt>0 && nt<=numSec )
246 {
247 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
248 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
249 if( std::fabs(dum) < 1.0 )
250 {
251 if( test >= 1.0e-10 )excs += dum*test;
252 }
253 else
254 excs += dum*test;
255 }
256 }
257 }
258 }
259 }
260 if( ran >= excs ) // 3 previous loops continued to the end
261 {
262 quasiElastic = true;
263 return;
264 }
265 np--; nm--; nz--;
266 G4int ncht = std::max( 1, np-nm+2 );
267 switch( ncht )
268 {
269 case 1:
270 if( G4UniformRand() < 0.5 )
271 currentParticle.SetDefinitionAndUpdateE( aLambda );
272 else
273 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
274 incidentHasChanged = true;
275 break;
276 case 2:
277 if( G4UniformRand() >= 0.5 )
278 {
279 if( G4UniformRand() < 0.5 )
280 currentParticle.SetDefinitionAndUpdateE( aLambda );
281 else
282 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
283 incidentHasChanged = true;
284 targetParticle.SetDefinitionAndUpdateE( aNeutron );
285 targetHasChanged = true;
286 }
287 break;
288 default:
289 targetParticle.SetDefinitionAndUpdateE( aNeutron );
290 targetHasChanged = true;
291 break;
292 }
293 }
294 else // target must be a neutron
295 {
296 counter = -1;
297 for( np=0; np<numSec/3 && ran>=excs; ++np )
298 {
299 for( nm=np; nm<=(np+2) && ran>=excs; ++nm )
300 {
301 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
302 {
303 if( ++counter < numMul )
304 {
305 nt = np+nm+nz;
306 if( nt>0 && nt<=numSec )
307 {
308 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
309 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
310 if( std::fabs(dum) < 1.0 )
311 {
312 if( test >= 1.0e-10 )excs += dum*test;
313 }
314 else
315 excs += dum*test;
316 }
317 }
318 }
319 }
320 }
321 if( ran >= excs ) // 3 previous loops continued to the end
322 {
323 quasiElastic = true;
324 return;
325 }
326 np--; nm--; nz--;
327 G4int ncht = std::max( 1, np-nm+3 );
328 switch( ncht )
329 {
330 case 1:
331 if( G4UniformRand() < 0.5 )
332 currentParticle.SetDefinitionAndUpdateE( aLambda );
333 else
334 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
335 incidentHasChanged = true;
336 targetParticle.SetDefinitionAndUpdateE( aProton );
337 targetHasChanged = true;
338 break;
339 case 2:
340 if( G4UniformRand() < 0.5 )
341 {
342 if( G4UniformRand() < 0.5 )
343 currentParticle.SetDefinitionAndUpdateE( aLambda );
344 else
345 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
346 incidentHasChanged = true;
347 }
348 else
349 {
350 targetParticle.SetDefinitionAndUpdateE( aProton );
351 targetHasChanged = true;
352 }
353 break;
354 default:
355 break;
356 }
357 }
358 SetUpPions( np, nm, nz, vec, vecLen );
359 return;
360}
361
362 /* end of file */
363
Note: See TracBrowser for help on using the repository browser.