source: trunk/source/processes/hadronic/models/low_energy/src/G4LESigmaMinusInelastic.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: 13.3 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: G4LESigmaMinusInelastic.cc,v 1.11 2006/06/29 20:45:17 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02 $
29//
30 // Hadronic Process: SigmaMinus 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 "G4LESigmaMinusInelastic.hh"
36#include "Randomize.hh"
37
38 G4HadFinalState *
39 G4LESigmaMinusInelastic::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 << "G4LESigmaMinusInelastic::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 if( originalIncident->GetKineticEnergy()/MeV > cutOff )
112 Cascade( vec, vecLen,
113 originalIncident, currentParticle, targetParticle,
114 incidentHasChanged, targetHasChanged, quasiElastic );
115
116 CalculateMomenta( vec, vecLen,
117 originalIncident, originalTarget, modifiedOriginal,
118 targetNucleus, currentParticle, targetParticle,
119 incidentHasChanged, targetHasChanged, quasiElastic );
120
121 SetUpChange( vec, vecLen,
122 currentParticle, targetParticle,
123 incidentHasChanged );
124
125 delete originalTarget;
126 return &theParticleChange;
127 }
128
129 void
130 G4LESigmaMinusInelastic::Cascade(
131 G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
132 G4int& vecLen,
133 const G4HadProjectile *originalIncident,
134 G4ReactionProduct &currentParticle,
135 G4ReactionProduct &targetParticle,
136 G4bool &incidentHasChanged,
137 G4bool &targetHasChanged,
138 G4bool &quasiElastic )
139 {
140 // derived from original FORTRAN code CASSM by H. Fesefeldt (13-Sep-1987)
141 //
142 // SigmaMinus 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 targetMass = targetParticle.GetMass()/MeV;
153 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
154 targetMass*targetMass +
155 2.0*targetMass*etOriginal );
156 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
157 if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
158 {
159 quasiElastic = true;
160 return;
161 }
162 static G4bool first = true;
163 const G4int numMul = 1200;
164 const G4int numSec = 60;
165 static G4double protmul[numMul], protnorm[numSec]; // proton constants
166 static G4double neutmul[numMul], neutnorm[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.70, 0.70 };
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 } // end of initialization
224
225 const G4double expxu = 82.; // upper bound for arg. of exp
226 const G4double expxl = -expxu; // lower bound for arg. of exp
227 G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
228 G4ParticleDefinition *aProton = G4Proton::Proton();
229 G4ParticleDefinition *aLambda = G4Lambda::Lambda();
230 G4ParticleDefinition *aSigmaZero = G4SigmaZero::SigmaZero();
231
232 // energetically possible to produce pion(s) --> inelastic scattering
233
234 G4double n, anpn;
235 GetNormalizationConstant( availableEnergy, n, anpn );
236 G4double ran = G4UniformRand();
237 G4double dum, excs = 0.0;
238 if( targetParticle.GetDefinition() == aProton )
239 {
240 counter = -1;
241 for( np=0; np<numSec/3 && ran>=excs; ++np )
242 {
243 for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm )
244 {
245 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
246 {
247 if( ++counter < numMul )
248 {
249 nt = np+nm+nz;
250 if( nt>0 && nt<=numSec )
251 {
252 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
253 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
254 if( std::fabs(dum) < 1.0 )
255 {
256 if( test >= 1.0e-10 )excs += dum*test;
257 }
258 else
259 excs += dum*test;
260 }
261 }
262 }
263 }
264 }
265 if( ran >= excs ) // 3 previous loops continued to the end
266 {
267 quasiElastic = true;
268 return;
269 }
270 np--; nm--; nz--;
271 G4int ncht = std::max( 1, np-nm+2 );
272 switch( ncht )
273 {
274 case 1:
275 if( G4UniformRand() < 0.5 )
276 currentParticle.SetDefinitionAndUpdateE( aLambda );
277 else
278 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
279 incidentHasChanged = true;
280 break;
281 case 2:
282 if( G4UniformRand() >= 0.5 )
283 {
284 if( G4UniformRand() < 0.5 )
285 currentParticle.SetDefinitionAndUpdateE( aLambda );
286 else
287 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
288 incidentHasChanged = true;
289 targetParticle.SetDefinitionAndUpdateE( aNeutron );
290 targetHasChanged = true;
291 }
292 break;
293 default:
294 targetParticle.SetDefinitionAndUpdateE( aNeutron );
295 targetHasChanged = true;
296 break;
297 }
298 }
299 else // target must be a neutron
300 {
301 counter = -1;
302 for( np=0; np<numSec/3 && ran>=excs; ++np )
303 {
304 for( nm=np; nm<=(np+2) && ran>=excs; ++nm )
305 {
306 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
307 {
308 if( ++counter < numMul )
309 {
310 nt = np+nm+nz;
311 if( nt>0 && nt<=numSec )
312 {
313 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
314 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
315 if( std::fabs(dum) < 1.0 )
316 {
317 if( test >= 1.0e-10 )excs += dum*test;
318 }
319 else
320 excs += dum*test;
321 }
322 }
323 }
324 }
325 }
326 if( ran >= excs ) // 3 previous loops continued to the end
327 {
328 quasiElastic = true;
329 return;
330 }
331 np--; nm--; nz--;
332 G4int ncht = std::max( 1, np-nm+3 );
333 switch( ncht )
334 {
335 case 1:
336 if( G4UniformRand() < 0.5 )
337 currentParticle.SetDefinitionAndUpdateE( aLambda );
338 else
339 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
340 incidentHasChanged = true;
341 targetParticle.SetDefinitionAndUpdateE( aProton );
342 targetHasChanged = true;
343 break;
344 case 2:
345 if( G4UniformRand() < 0.5 )
346 {
347 if( G4UniformRand() < 0.5 )
348 currentParticle.SetDefinitionAndUpdateE( aLambda );
349 else
350 currentParticle.SetDefinitionAndUpdateE( aSigmaZero );
351 incidentHasChanged = true;
352 }
353 else
354 {
355 targetParticle.SetDefinitionAndUpdateE( aProton );
356 targetHasChanged = true;
357 }
358 break;
359 default:
360 break;
361 }
362 }
363 SetUpPions( np, nm, nz, vec, vecLen );
364 return;
365 }
366
367 /* end of file */
368
Note: See TracBrowser for help on using the repository browser.