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