source: trunk/source/processes/hadronic/models/rpg/src/G4RPGPiMinusInelastic.cc@ 819

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

import all except CVS

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: G4RPGPiMinusInelastic.cc,v 1.1 2007/07/18 21:04:20 dennis Exp $
27// GEANT4 tag $Name: geant4-09-01-patch-02 $
28//
29
30#include "G4RPGPiMinusInelastic.hh"
31#include "Randomize.hh"
32
33
34G4HadFinalState*
35G4RPGPiMinusInelastic::ApplyYourself(const G4HadProjectile& aTrack,
36 G4Nucleus& targetNucleus)
37{
38 const G4HadProjectile* originalIncident = &aTrack;
39 if (originalIncident->GetKineticEnergy() <= 0.1*MeV)
40 {
41 theParticleChange.SetStatusChange(isAlive);
42 theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
43 theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
44 return &theParticleChange;
45 }
46
47 // create the target particle
48
49 G4DynamicParticle* originalTarget = targetNucleus.ReturnTargetParticle();
50 G4ReactionProduct targetParticle( originalTarget->GetDefinition() );
51
52 if( verboseLevel > 1 )
53 {
54 const G4Material *targetMaterial = aTrack.GetMaterial();
55 G4cout << "G4PionMinusInelastic::ApplyYourself called" << G4endl;
56 G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy() << "MeV, ";
57 G4cout << "target material = " << targetMaterial->GetName() << ", ";
58 G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
59 << G4endl;
60 }
61
62 G4ReactionProduct currentParticle(
63 const_cast<G4ParticleDefinition *>(originalIncident->GetDefinition() ) );
64 currentParticle.SetMomentum( originalIncident->Get4Momentum().vect() );
65 currentParticle.SetKineticEnergy( originalIncident->GetKineticEnergy() );
66
67 // Fermi motion and evaporation
68 // As of Geant3, the Fermi energy calculation had not been Done
69
70 G4double ek = originalIncident->GetKineticEnergy();
71 G4double amas = originalIncident->GetDefinition()->GetPDGMass();
72
73 G4double tkin = targetNucleus.Cinema( ek );
74 ek += tkin;
75 currentParticle.SetKineticEnergy( ek );
76 G4double et = ek + amas;
77 G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
78 G4double pp = currentParticle.GetMomentum().mag();
79 if( pp > 0.0 )
80 {
81 G4ThreeVector momentum = currentParticle.GetMomentum();
82 currentParticle.SetMomentum( momentum * (p/pp) );
83 }
84
85 // calculate black track energies
86
87 tkin = targetNucleus.EvaporationEffects( ek );
88 ek -= tkin;
89 currentParticle.SetKineticEnergy( ek );
90 et = ek + amas;
91 p = std::sqrt( std::abs((et-amas)*(et+amas)) );
92 pp = currentParticle.GetMomentum().mag();
93 if( pp > 0.0 )
94 {
95 G4ThreeVector momentum = currentParticle.GetMomentum();
96 currentParticle.SetMomentum( momentum * (p/pp) );
97 }
98
99 G4ReactionProduct modifiedOriginal = currentParticle;
100
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,256> vec;
107 G4int vecLen = 0;
108 vec.Initialize( 0 );
109
110 const G4double cutOff = 0.1*MeV;
111 if( currentParticle.GetKineticEnergy() > 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, currentParticle, targetParticle,
122 incidentHasChanged);
123
124 delete originalTarget;
125 return &theParticleChange;
126}
127
128void G4RPGPiMinusInelastic::Cascade(
129 G4FastVector<G4ReactionProduct,256> &vec,
130 G4int& vecLen,
131 const G4HadProjectile *originalIncident,
132 G4ReactionProduct &currentParticle,
133 G4ReactionProduct &targetParticle,
134 G4bool &incidentHasChanged,
135 G4bool &targetHasChanged,
136 G4bool &quasiElastic)
137{
138 // Derived from H. Fesefeldt's original FORTRAN code CASPIM
139 //
140 // pi- undergoes interaction with nucleon within nucleus.
141 // Check if energetically possible to produce pions/kaons.
142 // If not assume nuclear excitation occurs and input particle
143 // is degraded in energy. No other particles produced.
144 // If reaction is possible find correct number of pions/protons/neutrons
145 // produced using an interpolation to multiplicity data.
146 // Replace some pions or protons/neutrons by kaons or strange baryons
147 // according to average multiplicity per inelastic reactions.
148 //
149 const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass();
150 const G4double etOriginal = originalIncident->GetTotalEnergy();
151 const G4double pOriginal = originalIncident->GetTotalMomentum();
152 const G4double targetMass = targetParticle.GetMass();
153 G4double centerofmassEnergy = std::sqrt(mOriginal*mOriginal +
154 targetMass*targetMass +
155 2.0*targetMass*etOriginal);
156 G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
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 const G4double c = 1.25;
165 const G4double b[] = { 0.70, 0.70 };
166 if( first ) // compute normalization constants, this will only be Done once
167 {
168 first = false;
169 G4int i;
170 for( i=0; i<numMul; ++i )protmul[i] = 0.0;
171 for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
172 counter = -1;
173 for( np=0; np<(numSec/3); ++np )
174 {
175 for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
176 {
177 for( nz=0; nz<numSec/3; ++nz )
178 {
179 if( ++counter < numMul )
180 {
181 nt = np+nm+nz;
182 if( nt > 0 )
183 {
184 protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
185 protnorm[nt-1] += protmul[counter];
186 }
187 }
188 }
189 }
190 }
191 for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
192 for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
193 counter = -1;
194 for( np=0; np<numSec/3; ++np )
195 {
196 for( nm=np; nm<=(np+2); ++nm )
197 {
198 for( nz=0; nz<numSec/3; ++nz )
199 {
200 if( ++counter < numMul )
201 {
202 nt = np+nm+nz;
203 if( (nt>0) && (nt<=numSec) )
204 {
205 neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
206 neutnorm[nt-1] += neutmul[counter];
207 }
208 }
209 }
210 }
211 }
212 for( i=0; i<numSec; ++i ) {
213 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
214 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
215 }
216 } // end of initialization
217
218 const G4double expxu = 82.; // upper bound for arg. of exp
219 const G4double expxl = -expxu; // lower bound for arg. of exp
220 G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
221 G4ParticleDefinition *aProton = G4Proton::Proton();
222 G4ParticleDefinition *aPiZero = G4PionZero::PionZero();
223 G4int ieab = G4int(availableEnergy*5.0/GeV);
224 const G4double supp[] = {0.,0.4,0.55,0.65,0.75,0.82,0.86,0.90,0.94,0.98};
225 G4double test, w0, wp, wt, wm;
226 if( (availableEnergy<2.0*GeV) && (G4UniformRand()>=supp[ieab]) )
227 {
228 // suppress high multiplicity events at low momentum
229 // only one pion will be produced
230
231 // charge exchange reaction is included in inelastic cross section
232
233 const G4double cech[] = {1.,0.95,0.79,0.32,0.19,0.16,0.14,0.12,0.10,0.08};
234 G4int iplab = G4int(std::min( 9.0, pOriginal/GeV*5.0 ));
235 if( G4UniformRand() <= cech[iplab] )
236 {
237 if( targetParticle.GetDefinition() == aProton )
238 {
239 currentParticle.SetDefinitionAndUpdateE( aPiZero ); // charge exchange
240 targetParticle.SetDefinitionAndUpdateE( aNeutron );
241 incidentHasChanged = true;
242 targetHasChanged = true;
243 }
244 }
245
246 if( availableEnergy <= G4PionMinus::PionMinus()->GetPDGMass() )
247 {
248 quasiElastic = true;
249 return;
250 }
251
252 nm = np = nz = 0;
253 if( targetParticle.GetDefinition() == aProton )
254 {
255 test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[0])*(1.0+b[0])/(2.0*c*c) ) ) );
256 w0 = test;
257 wp = 10.0*test;
258 test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[0])*(-1.0+b[0])/(2.0*c*c) ) ) );
259 wm = test;
260 wt = w0+wp+wm;
261 wp += w0;
262 G4double ran = G4UniformRand();
263 if( ran < w0/wt )
264 nz = 1;
265 else if( ran < wp/wt )
266 np = 1;
267 else
268 nm = 1;
269 }
270 else // target is a neutron
271 {
272 test = std::exp( std::min( expxu, std::max( expxl, -(1.0+b[1])*(1.0+b[1])/(2.0*c*c) ) ) );
273 w0 = test;
274 test = std::exp( std::min( expxu, std::max( expxl, -(-1.0+b[1])*(-1.0+b[1])/(2.0*c*c) ) ) );
275 wm = test;
276 G4double ran = G4UniformRand();
277 if( ran < w0/(w0+wm) )
278 nz = 1;
279 else
280 nm = 1;
281 }
282 }
283 else
284 {
285 if( availableEnergy <= G4PionMinus::PionMinus()->GetPDGMass() )
286 {
287 quasiElastic = true;
288 return;
289 }
290 G4double n, anpn;
291 GetNormalizationConstant( availableEnergy, n, anpn );
292 G4double ran = G4UniformRand();
293 G4double dum, excs = 0.0;
294 if( targetParticle.GetDefinition() == aProton )
295 {
296 counter = -1;
297 for( np=0; (np<numSec/3) && (ran>=excs); ++np )
298 {
299 for( nm=std::max(0,np-1); (nm<=(np+1)) && (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 )
307 {
308 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
309 dum = (pi/anpn)*nt*protmul[counter]*protnorm[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 }
328 else // target must be a neutron
329 {
330 counter = -1;
331 for( np=0; (np<numSec/3) && (ran>=excs); ++np )
332 {
333 for( nm=np; (nm<=(np+2)) && (ran>=excs); ++nm )
334 {
335 for( nz=0; (nz<numSec/3) && (ran>=excs); ++nz )
336 {
337 if( ++counter < numMul )
338 {
339 nt = np+nm+nz;
340 if( (nt>=1) && (nt<=numSec) )
341 {
342 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
343 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
344 if( std::fabs(dum) < 1.0 )
345 {
346 if( test >= 1.0e-10 )excs += dum*test;
347 }
348 else
349 excs += dum*test;
350 }
351 }
352 }
353 }
354 }
355 if( ran >= excs ) // 3 previous loops continued to the end
356 {
357 quasiElastic = true;
358 return;
359 }
360 np--; nm--; nz--;
361 }
362 }
363 if( targetParticle.GetDefinition() == aProton )
364 {
365 switch( np-nm )
366 {
367 case 0:
368 if( G4UniformRand() >= 0.75 )
369 {
370 currentParticle.SetDefinitionAndUpdateE( aPiZero );
371 targetParticle.SetDefinitionAndUpdateE( aNeutron );
372 incidentHasChanged = true;
373 targetHasChanged = true;
374 }
375 break;
376 case 1:
377 targetParticle.SetDefinitionAndUpdateE( aNeutron );
378 targetHasChanged = true;
379 break;
380 default:
381 currentParticle.SetDefinitionAndUpdateE( aPiZero );
382 incidentHasChanged = true;
383 break;
384 }
385 }
386 else
387 {
388 switch( np-nm )
389 {
390 case -1:
391 if( G4UniformRand() < 0.5 )
392 {
393 targetParticle.SetDefinitionAndUpdateE( aProton );
394 targetHasChanged = true;
395 } else {
396 currentParticle.SetDefinitionAndUpdateE( aPiZero );
397 incidentHasChanged = true;
398 }
399 break;
400 case 0:
401 break;
402 default:
403 currentParticle.SetDefinitionAndUpdateE( aPiZero );
404 incidentHasChanged = true;
405 break;
406 }
407 }
408
409 SetUpPions( np, nm, nz, vec, vecLen );
410 return;
411}
412
413 /* end of file */
414
Note: See TracBrowser for help on using the repository browser.