source: trunk/source/processes/hadronic/models/rpg/src/G4RPGKZeroInelastic.cc@ 1038

Last change on this file since 1038 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 14.2 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: G4RPGKZeroInelastic.cc,v 1.1 2007/07/18 21:04:20 dennis Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29
30#include "G4RPGKZeroInelastic.hh"
31#include "Randomize.hh"
32
33G4HadFinalState*
34G4RPGKZeroInelastic::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 << "G4RPGKZeroInelastic::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 G4ReactionProduct currentParticle = modifiedOriginal;
87 G4ReactionProduct targetParticle;
88 targetParticle = *originalTarget;
89 currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
90 targetParticle.SetSide( -1 ); // target always goes in backward hemisphere
91 G4bool incidentHasChanged = false;
92 G4bool targetHasChanged = false;
93 G4bool quasiElastic = false;
94 G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec; // vec will contain the secondary particles
95 G4int vecLen = 0;
96 vec.Initialize( 0 );
97
98 const G4double cutOff = 0.1;
99 if( currentParticle.GetKineticEnergy()/MeV > cutOff )
100 Cascade( vec, vecLen,
101 originalIncident, currentParticle, targetParticle,
102 incidentHasChanged, targetHasChanged, quasiElastic );
103
104 CalculateMomenta( vec, vecLen,
105 originalIncident, originalTarget, modifiedOriginal,
106 targetNucleus, currentParticle, targetParticle,
107 incidentHasChanged, targetHasChanged, quasiElastic );
108
109 SetUpChange( vec, vecLen,
110 currentParticle, targetParticle,
111 incidentHasChanged );
112
113 delete originalTarget;
114 return &theParticleChange;
115}
116
117void G4RPGKZeroInelastic::Cascade(
118 G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
119 G4int& vecLen,
120 const G4HadProjectile *originalIncident,
121 G4ReactionProduct &currentParticle,
122 G4ReactionProduct &targetParticle,
123 G4bool &incidentHasChanged,
124 G4bool &targetHasChanged,
125 G4bool &quasiElastic )
126{
127 // Derived from H. Fesefeldt's original FORTRAN code CASK0
128 //
129 // K0Short undergoes interaction with nucleon within a nucleus. Check if it is
130 // energetically possible to produce pions/kaons. In not, assume nuclear excitation
131 // occurs and input particle is degraded in energy. No other particles are produced.
132 // If reaction is possible, find the correct number of pions/protons/neutrons
133 // produced using an interpolation to multiplicity data. Replace some pions or
134 // protons/neutrons by kaons or strange baryons according to the average
135 // multiplicity per Inelastic reaction.
136
137 const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
138 const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
139 const G4double targetMass = targetParticle.GetMass()/MeV;
140 G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
141 targetMass*targetMass +
142 2.0*targetMass*etOriginal );
143 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
144 if( availableEnergy <= G4PionPlus::PionPlus()->GetPDGMass()/MeV )
145 {
146 quasiElastic = true;
147 return;
148 }
149 static G4bool first = true;
150 const G4int numMul = 1200;
151 const G4int numSec = 60;
152 static G4double protmul[numMul], protnorm[numSec]; // proton constants
153 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
154 // np = number of pi+, nm = number of pi-, nz = number of pi0
155 G4int counter, nt=0, np=0, nm=0, nz=0;
156 const G4double c = 1.25;
157 const G4double b[] = { 0.7, 0.7 };
158 if( first ) // compute normalization constants, this will only be Done once
159 {
160 first = false;
161 G4int i;
162 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
163 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
164 counter = -1;
165 for( np=0; np<(numSec/3); ++np )
166 {
167 for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
168 {
169 for( nz=0; nz<numSec/3; ++nz )
170 {
171 if( ++counter < numMul )
172 {
173 nt = np+nm+nz;
174 if( nt>0 && nt<=numSec )
175 {
176 protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
177 protnorm[nt-1] += protmul[counter];
178 }
179 }
180 }
181 }
182 }
183 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
184 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
185 counter = -1;
186 for( np=0; np<numSec/3; ++np )
187 {
188 for( nm=std::max(0,np-2); nm<=np; ++nm )
189 {
190 for( nz=0; nz<numSec/3; ++nz )
191 {
192 if( ++counter < numMul )
193 {
194 nt = np+nm+nz;
195 if( nt>0 && nt<=numSec )
196 {
197 neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
198 neutnorm[nt-1] += neutmul[counter];
199 }
200 }
201 }
202 }
203 }
204 for( i=0; i<numSec; ++i )
205 {
206 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
207 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
208 }
209 } // end of initialization
210
211 const G4double expxu = 82.; // upper bound for arg. of exp
212 const G4double expxl = -expxu; // lower bound for arg. of exp
213 G4ParticleDefinition *aKaonPlus = G4KaonPlus::KaonPlus();
214 G4ParticleDefinition *aKaonZL = G4KaonZeroLong::KaonZeroLong();
215 G4ParticleDefinition *aKaonZS = G4KaonZeroShort::KaonZeroShort();
216 G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
217 G4ParticleDefinition *aProton = G4Proton::Proton();
218 G4int ieab = static_cast<G4int>(5.0*availableEnergy*MeV/GeV);
219 const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
220 G4double test, w0, wp, wt, wm;
221 if( (availableEnergy*MeV/GeV < 2.0) && (G4UniformRand() >= supp[ieab]) )
222 {
223 //
224 // suppress high multiplicity events at low momentum
225 // only one pion will be produced
226 //
227 nm = np = nz = 0;
228 if( targetParticle.GetDefinition() == aNeutron )
229 {
230 test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
231 w0 = test/2.0;
232 test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
233 wm = test*1.5;
234 if( G4UniformRand() < w0/(w0+wm) )
235 nz = 1;
236 else
237 nm = 1;
238 }
239 else // target is a proton
240 {
241 test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[1])*(1.0+b[1])/(2.0*c*c) ) ) );
242 w0 = test;
243 wp = test;
244 test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[1])*(-1.0+b[1])/(2.0*c*c) ) ) );
245 wm = test;
246 wt = w0+wp+wm;
247 wp += w0;
248 G4double ran = G4UniformRand();
249 if( ran < w0/wt )
250 nz = 1;
251 else if( ran < wp/wt )
252 np = 1;
253 else
254 nm = 1;
255 }
256 }
257 else // (availableEnergy*MeV/GeV >= 2.0) || (G4UniformRand() < supp[ieab])
258 {
259 G4double n, anpn;
260 GetNormalizationConstant( availableEnergy, n, anpn );
261 G4double ran = G4UniformRand();
262 G4double dum, excs = 0.0;
263 if( targetParticle.GetDefinition() == aProton )
264 {
265 counter = -1;
266 for( np=0; np<numSec/3 && ran>=excs; ++np )
267 {
268 for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm )
269 {
270 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
271 {
272 if( ++counter < numMul )
273 {
274 nt = np+nm+nz;
275 if( nt>0 && nt<=numSec )
276 {
277 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
278 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
279 if( std::fabs(dum) < 1.0 )
280 {
281 if( test >= 1.0e-10 )excs += dum*test;
282 }
283 else
284 excs += dum*test;
285 }
286 }
287 }
288 }
289 }
290 if( ran >= excs ) // 3 previous loops continued to the end
291 {
292 quasiElastic = true;
293 return;
294 }
295 np--; nm--; nz--;
296 }
297 else // target must be a neutron
298 {
299 counter = -1;
300 for( np=0; np<numSec/3 && ran>=excs; ++np )
301 {
302 for( nm=std::max(0,np-2); nm<=np && ran>=excs; ++nm )
303 {
304 for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
305 {
306 if( ++counter < numMul )
307 {
308 nt = np+nm+nz;
309 if( nt>0 && nt<=numSec )
310 {
311 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
312 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
313 if( std::fabs(dum) < 1.0 )
314 {
315 if( test >= 1.0e-10 )excs += dum*test;
316 }
317 else
318 excs += dum*test;
319 }
320 }
321 }
322 }
323 }
324 if( ran >= excs ) // 3 previous loops continued to the end
325 {
326 quasiElastic = true;
327 return;
328 }
329 np--; nm--; nz--;
330 }
331 }
332 if( targetParticle.GetDefinition() == aProton )
333 {
334 switch( np-nm )
335 {
336 case 0:
337 if( G4UniformRand() < 0.25 )
338 {
339 currentParticle.SetDefinitionAndUpdateE( aKaonPlus );
340 targetParticle.SetDefinitionAndUpdateE( aNeutron );
341 incidentHasChanged = true;
342 targetHasChanged = true;
343 }
344 break;
345 case 1:
346 targetParticle.SetDefinitionAndUpdateE( aNeutron );
347 targetHasChanged = true;
348 break;
349 default:
350 targetParticle.SetDefinitionAndUpdateE( aNeutron );
351 targetHasChanged = true;
352 break;
353 }
354 }
355 else // targetParticle is a neutron
356 {
357 switch( np-nm ) // seems wrong, charge not conserved
358 {
359 case 1:
360 if( G4UniformRand() < 0.5 )
361 {
362 currentParticle.SetDefinitionAndUpdateE( aKaonPlus );
363 incidentHasChanged = true;
364 }
365 else
366 {
367 targetParticle.SetDefinitionAndUpdateE( aProton );
368 targetHasChanged = true;
369 }
370 break;
371 case 2:
372 currentParticle.SetDefinitionAndUpdateE( aKaonPlus );
373 incidentHasChanged = true;
374 targetParticle.SetDefinitionAndUpdateE( aProton );
375 targetHasChanged = true;
376 break;
377 default:
378 break;
379 }
380 }
381 if( currentParticle.GetDefinition() == aKaonZS )
382 {
383 if( G4UniformRand() >= 0.5 )
384 {
385 currentParticle.SetDefinitionAndUpdateE( aKaonZL);
386 incidentHasChanged = true;
387 }
388 }
389 if( targetParticle.GetDefinition() == aKaonZS )
390 {
391 if( G4UniformRand() >= 0.5 )
392 {
393 targetParticle.SetDefinitionAndUpdateE( aKaonZL );
394 targetHasChanged = true;
395 }
396 }
397 SetUpPions( np, nm, nz, vec, vecLen );
398 return;
399}
400
401 /* end of file */
402
Note: See TracBrowser for help on using the repository browser.