source: trunk/source/processes/hadronic/models/high_energy/src/G4HEAntiLambdaInelastic.cc@ 1201

Last change on this file since 1201 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 28.0 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: G4HEAntiLambdaInelastic.cc,v 1.15 2008/03/17 20:49:17 dennis Exp $
28// GEANT4 tag $Name: geant4-09-03-cand-01 $
29//
30//
31
32#include "globals.hh"
33#include "G4ios.hh"
34
35//
36// G4 Process: Gheisha High Energy Collision model.
37// This includes the high energy cascading model, the two-body-resonance model
38// and the low energy two-body model. Not included are the low energy stuff like
39// nuclear reactions, nuclear fission without any cascading and all processes for
40// particles at rest.
41// First work done by J.L.Chuma and F.W.Jones, TRIUMF, June 96.
42// H. Fesefeldt, RWTH-Aachen, 23-October-1996
43// Last modified: 29-July-1998
44
45#include "G4HEAntiLambdaInelastic.hh"
46
47G4HadFinalState * G4HEAntiLambdaInelastic::
48ApplyYourself( const G4HadProjectile &aTrack, G4Nucleus &targetNucleus )
49 {
50 G4HEVector * pv = new G4HEVector[MAXPART];
51 const G4HadProjectile *aParticle = &aTrack;
52// G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
53 const G4double atomicWeight = targetNucleus.GetN();
54 const G4double atomicNumber = targetNucleus.GetZ();
55 G4HEVector incidentParticle(aParticle);
56
57 G4int incidentCode = incidentParticle.getCode();
58 G4double incidentMass = incidentParticle.getMass();
59 G4double incidentTotalEnergy = incidentParticle.getEnergy();
60 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
61 G4double incidentKineticEnergy = incidentTotalEnergy - incidentMass;
62
63 if(incidentKineticEnergy < 1.)
64 {
65 G4cout << "GHEAntiLambdaInelastic: incident energy < 1 GeV" << G4endl;
66 }
67 if(verboseLevel > 1)
68 {
69 G4cout << "G4HEAntiLambdaInelastic::ApplyYourself" << G4endl;
70 G4cout << "incident particle " << incidentParticle.getName()
71 << "mass " << incidentMass
72 << "kinetic energy " << incidentKineticEnergy
73 << G4endl;
74 G4cout << "target material with (A,Z) = ("
75 << atomicWeight << "," << atomicNumber << ")" << G4endl;
76 }
77
78 G4double inelasticity = NuclearInelasticity(incidentKineticEnergy,
79 atomicWeight, atomicNumber);
80 if(verboseLevel > 1)
81 G4cout << "nuclear inelasticity = " << inelasticity << G4endl;
82
83 incidentKineticEnergy -= inelasticity;
84
85 G4double excitationEnergyGNP = 0.;
86 G4double excitationEnergyDTA = 0.;
87
88 G4double excitation = NuclearExcitation(incidentKineticEnergy,
89 atomicWeight, atomicNumber,
90 excitationEnergyGNP,
91 excitationEnergyDTA);
92 if(verboseLevel > 1)
93 G4cout << "nuclear excitation = " << excitation << excitationEnergyGNP
94 << excitationEnergyDTA << G4endl;
95
96
97 incidentKineticEnergy -= excitation;
98 incidentTotalEnergy = incidentKineticEnergy + incidentMass;
99 incidentTotalMomentum = std::sqrt( (incidentTotalEnergy-incidentMass)
100 *(incidentTotalEnergy+incidentMass));
101
102
103 G4HEVector targetParticle;
104 if(G4UniformRand() < atomicNumber/atomicWeight)
105 {
106 targetParticle.setDefinition("Proton");
107 }
108 else
109 {
110 targetParticle.setDefinition("Neutron");
111 }
112
113 G4double targetMass = targetParticle.getMass();
114 G4double centerOfMassEnergy = std::sqrt( incidentMass*incidentMass + targetMass*targetMass
115 + 2.0*targetMass*incidentTotalEnergy);
116 G4double availableEnergy = centerOfMassEnergy - targetMass - incidentMass;
117
118 // this was the meaning of inElastic in the
119 // original Gheisha stand-alone version.
120// G4bool inElastic = InElasticCrossSectionInFirstInt
121// (availableEnergy, incidentCode, incidentTotalMomentum);
122 // by unknown reasons, it has been replaced
123 // to the following code in Geant???
124 G4bool inElastic = true;
125// if (G4UniformRand() < elasticCrossSection/totalCrossSection) inElastic = false;
126
127 vecLength = 0;
128
129 if(verboseLevel > 1)
130 G4cout << "ApplyYourself: CallFirstIntInCascade for particle "
131 << incidentCode << G4endl;
132
133 G4bool successful = false;
134
135 if(inElastic || (!inElastic && atomicWeight < 1.5))
136 {
137 FirstIntInCasAntiLambda(inElastic, availableEnergy, pv, vecLength,
138 incidentParticle, targetParticle, atomicWeight);
139
140 if(verboseLevel > 1)
141 G4cout << "ApplyYourself::StrangeParticlePairProduction" << G4endl;
142
143
144 if ((vecLength > 0) && (availableEnergy > 1.))
145 StrangeParticlePairProduction( availableEnergy, centerOfMassEnergy,
146 pv, vecLength,
147 incidentParticle, targetParticle);
148 HighEnergyCascading( successful, pv, vecLength,
149 excitationEnergyGNP, excitationEnergyDTA,
150 incidentParticle, targetParticle,
151 atomicWeight, atomicNumber);
152 if (!successful)
153 HighEnergyClusterProduction( successful, pv, vecLength,
154 excitationEnergyGNP, excitationEnergyDTA,
155 incidentParticle, targetParticle,
156 atomicWeight, atomicNumber);
157 if (!successful)
158 MediumEnergyCascading( successful, pv, vecLength,
159 excitationEnergyGNP, excitationEnergyDTA,
160 incidentParticle, targetParticle,
161 atomicWeight, atomicNumber);
162
163 if (!successful)
164 MediumEnergyClusterProduction( successful, pv, vecLength,
165 excitationEnergyGNP, excitationEnergyDTA,
166 incidentParticle, targetParticle,
167 atomicWeight, atomicNumber);
168 if (!successful)
169 QuasiElasticScattering( successful, pv, vecLength,
170 excitationEnergyGNP, excitationEnergyDTA,
171 incidentParticle, targetParticle,
172 atomicWeight, atomicNumber);
173 }
174 if (!successful)
175 {
176 ElasticScattering( successful, pv, vecLength,
177 incidentParticle,
178 atomicWeight, atomicNumber);
179 }
180
181 if (!successful)
182 {
183 G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles" << G4endl;
184 }
185 FillParticleChange(pv, vecLength);
186 delete [] pv;
187 theParticleChange.SetStatusChange(stopAndKill);
188 return & theParticleChange;
189 }
190
191void
192G4HEAntiLambdaInelastic::FirstIntInCasAntiLambda( G4bool &inElastic,
193 const G4double availableEnergy,
194 G4HEVector pv[],
195 G4int &vecLen,
196 G4HEVector incidentParticle,
197 G4HEVector targetParticle,
198 const G4double atomicWeight)
199
200// AntiLambda undergoes interaction with nucleon within a nucleus. Check if it is
201// energetically possible to produce pions/kaons. In not, assume nuclear excitation
202// occurs and input particle is degraded in energy. No other particles are produced.
203// If reaction is possible, find the correct number of pions/protons/neutrons
204// produced using an interpolation to multiplicity data. Replace some pions or
205// protons/neutrons by kaons or strange baryons according to the average
206// multiplicity per inelastic reaction.
207
208 {
209 static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
210 static const G4double expxl = -expxu; // lower bound for arg. of exp
211
212 static const G4double protb = 0.7;
213 static const G4double neutb = 0.7;
214 static const G4double c = 1.25;
215
216 static const G4int numMul = 1200;
217 static const G4int numMulAn = 400;
218 static const G4int numSec = 60;
219
220// G4int neutronCode = Neutron.getCode();
221 G4int protonCode = Proton.getCode();
222
223 G4int targetCode = targetParticle.getCode();
224// G4double incidentMass = incidentParticle.getMass();
225// G4double incidentEnergy = incidentParticle.getEnergy();
226 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
227
228 static G4bool first = true;
229 static G4double protmul[numMul], protnorm[numSec]; // proton constants
230 static G4double protmulAn[numMulAn],protnormAn[numSec];
231 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
232 static G4double neutmulAn[numMulAn],neutnormAn[numSec];
233
234 // misc. local variables
235 // np = number of pi+, nm = number of pi-, nz = number of pi0
236
237 G4int i, counter, nt, np, nm, nz;
238
239 if( first )
240 { // compute normalization constants, this will only be done once
241 first = false;
242 for( i=0; i<numMul ; i++ ) protmul[i] = 0.0;
243 for( i=0; i<numSec ; i++ ) protnorm[i] = 0.0;
244 counter = -1;
245 for( np=0; np<(numSec/3); np++ )
246 {
247 for( nm=std::max(0,np-2); nm<=(np+1); nm++ )
248 {
249 for( nz=0; nz<numSec/3; nz++ )
250 {
251 if( ++counter < numMul )
252 {
253 nt = np+nm+nz;
254 if( (nt>0) && (nt<=numSec) )
255 {
256 protmul[counter] = pmltpc(np,nm,nz,nt,protb,c);
257 protnorm[nt-1] += protmul[counter];
258 }
259 }
260 }
261 }
262 }
263 for( i=0; i<numMul; i++ )neutmul[i] = 0.0;
264 for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
265 counter = -1;
266 for( np=0; np<numSec/3; np++ )
267 {
268 for( nm=std::max(0,np-1); nm<=(np+2); nm++ )
269 {
270 for( nz=0; nz<numSec/3; nz++ )
271 {
272 if( ++counter < numMul )
273 {
274 nt = np+nm+nz;
275 if( (nt>0) && (nt<=numSec) )
276 {
277 neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c);
278 neutnorm[nt-1] += neutmul[counter];
279 }
280 }
281 }
282 }
283 }
284 for( i=0; i<numSec; i++ )
285 {
286 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
287 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
288 }
289 // annihilation
290 for( i=0; i<numMulAn ; i++ ) protmulAn[i] = 0.0;
291 for( i=0; i<numSec ; i++ ) protnormAn[i] = 0.0;
292 counter = -1;
293 for( np=1; np<(numSec/3); np++ )
294 {
295 nm = std::max(0,np-1);
296 for( nz=0; nz<numSec/3; nz++ )
297 {
298 if( ++counter < numMulAn )
299 {
300 nt = np+nm+nz;
301 if( (nt>1) && (nt<=numSec) )
302 {
303 protmulAn[counter] = pmltpc(np,nm,nz,nt,protb,c);
304 protnormAn[nt-1] += protmulAn[counter];
305 }
306 }
307 }
308 }
309 for( i=0; i<numMulAn; i++ ) neutmulAn[i] = 0.0;
310 for( i=0; i<numSec; i++ ) neutnormAn[i] = 0.0;
311 counter = -1;
312 for( np=0; np<numSec/3; np++ )
313 {
314 nm = np;
315 for( nz=0; nz<numSec/3; nz++ )
316 {
317 if( ++counter < numMulAn )
318 {
319 nt = np+nm+nz;
320 if( (nt>1) && (nt<=numSec) )
321 {
322 neutmulAn[counter] = pmltpc(np,nm,nz,nt,neutb,c);
323 neutnormAn[nt-1] += neutmulAn[counter];
324 }
325 }
326 }
327 }
328 for( i=0; i<numSec; i++ )
329 {
330 if( protnormAn[i] > 0.0 )protnormAn[i] = 1.0/protnormAn[i];
331 if( neutnormAn[i] > 0.0 )neutnormAn[i] = 1.0/neutnormAn[i];
332 }
333 } // end of initialization
334
335
336 // initialize the first two places
337 // the same as beam and target
338 pv[0] = incidentParticle;
339 pv[1] = targetParticle;
340 vecLen = 2;
341
342 if( !inElastic )
343 { // some two-body reactions
344 G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
345
346 G4int iplab = std::min(9, G4int( incidentTotalMomentum*2.5));
347 if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) )
348 {
349 G4double ran = G4UniformRand();
350
351 if ( targetCode == protonCode)
352 {
353 if(ran < 0.2)
354 {
355 pv[0] = AntiSigmaZero;
356 }
357 else if (ran < 0.4)
358 {
359 pv[0] = AntiSigmaMinus;
360 pv[1] = Neutron;
361 }
362 else if (ran < 0.6)
363 {
364 pv[0] = Proton;
365 pv[1] = AntiLambda;
366 }
367 else if (ran < 0.8)
368 {
369 pv[0] = Proton;
370 pv[1] = AntiSigmaZero;
371 }
372 else
373 {
374 pv[0] = Neutron;
375 pv[1] = AntiSigmaMinus;
376 }
377 }
378 else
379 {
380 if (ran < 0.2)
381 {
382 pv[0] = AntiSigmaZero;
383 }
384 else if (ran < 0.4)
385 {
386 pv[0] = AntiSigmaPlus;
387 pv[1] = Proton;
388 }
389 else if (ran < 0.6)
390 {
391 pv[0] = Neutron;
392 pv[1] = AntiLambda;
393 }
394 else if (ran < 0.8)
395 {
396 pv[0] = Neutron;
397 pv[1] = AntiSigmaZero;
398 }
399 else
400 {
401 pv[0] = Proton;
402 pv[1] = AntiSigmaPlus;
403 }
404 }
405 }
406 return;
407 }
408 else if (availableEnergy <= PionPlus.getMass())
409 return;
410
411 // inelastic scattering
412
413 np = 0; nm = 0; nz = 0;
414 G4double anhl[] = {1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 0.97, 0.88,
415 0.85, 0.81, 0.75, 0.64, 0.64, 0.55, 0.55, 0.45, 0.47, 0.40,
416 0.39, 0.36, 0.33, 0.10, 0.01};
417 G4int iplab = G4int( incidentTotalMomentum*10.);
418 if ( iplab > 9) iplab = 10 + G4int( (incidentTotalMomentum -1.)*5. );
419 if ( iplab > 14) iplab = 15 + G4int( incidentTotalMomentum -2. );
420 if ( iplab > 22) iplab = 23 + G4int( (incidentTotalMomentum -10.)/10.);
421 iplab = std::min(24, iplab);
422
423 if ( G4UniformRand() > anhl[iplab] )
424 { // non- annihilation channels
425
426 // number of total particles vs. centre of mass Energy - 2*proton mass
427
428 G4double aleab = std::log(availableEnergy);
429 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
430 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
431
432 // normalization constant for kno-distribution.
433 // calculate first the sum of all constants, check for numerical problems.
434 G4double test, dum, anpn = 0.0;
435
436 for (nt=1; nt<=numSec; nt++) {
437 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
438 dum = pi*nt/(2.0*n*n);
439 if (std::fabs(dum) < 1.0) {
440 if( test >= 1.0e-10 )anpn += dum*test;
441 } else {
442 anpn += dum*test;
443 }
444 }
445
446 G4double ran = G4UniformRand();
447 G4double excs = 0.0;
448 if( targetCode == protonCode )
449 {
450 counter = -1;
451 for( np=0; np<numSec/3; np++ )
452 {
453 for( nm=std::max(0,np-2); nm<=(np+1); nm++ )
454 {
455 for( nz=0; nz<numSec/3; nz++ )
456 {
457 if( ++counter < numMul )
458 {
459 nt = np+nm+nz;
460 if( (nt>0) && (nt<=numSec) )
461 {
462 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
463 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
464
465 if (std::fabs(dum) < 1.0) {
466 if( test >= 1.0e-10 )excs += dum*test;
467 } else {
468 excs += dum*test;
469 }
470
471 if (ran < excs) goto outOfLoop; //----------------------->
472 }
473 }
474 }
475 }
476 }
477
478 // 3 previous loops continued to the end
479 inElastic = false; // quasi-elastic scattering
480 return;
481 }
482 else
483 { // target must be a neutron
484 counter = -1;
485 for( np=0; np<numSec/3; np++ )
486 {
487 for( nm=std::max(0,np-1); nm<=(np+2); nm++ )
488 {
489 for( nz=0; nz<numSec/3; nz++ )
490 {
491 if( ++counter < numMul )
492 {
493 nt = np+nm+nz;
494 if( (nt>0) && (nt<=numSec) )
495 {
496 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
497 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
498 if (std::fabs(dum) < 1.0) {
499 if( test >= 1.0e-10 )excs += dum*test;
500 } else {
501 excs += dum*test;
502 }
503
504 if (ran < excs) goto outOfLoop; // -------------------------->
505 }
506 }
507 }
508 }
509 }
510 // 3 previous loops continued to the end
511 inElastic = false; // quasi-elastic scattering.
512 return;
513 }
514
515 outOfLoop: // <------------------------------------------------------------------------
516
517 ran = G4UniformRand();
518
519 if( targetCode == protonCode)
520 {
521 if( np == nm)
522 {
523 if (ran < 0.40)
524 {
525 }
526 else if (ran < 0.8)
527 {
528 pv[0] = AntiSigmaZero;
529 }
530 else
531 {
532 pv[0] = AntiSigmaMinus;
533 pv[1] = Neutron;
534 }
535 }
536 else if (np == (nm+1))
537 {
538 if( ran < 0.25)
539 {
540 pv[1] = Neutron;
541 }
542 else if (ran < 0.5)
543 {
544 pv[0] = AntiSigmaZero;
545 pv[1] = Neutron;
546 }
547 else
548 {
549 pv[0] = AntiSigmaPlus;
550 }
551 }
552 else if (np == (nm-1))
553 {
554 pv[0] = AntiSigmaMinus;
555 }
556 else
557 {
558 pv[0] = AntiSigmaPlus;
559 pv[1] = Neutron;
560 }
561 }
562 else
563 {
564 if( np == nm)
565 {
566 if (ran < 0.4)
567 {
568 }
569 else if(ran < 0.8)
570 {
571 pv[0] = AntiSigmaZero;
572 }
573 else
574 {
575 pv[0] = AntiSigmaPlus;
576 pv[1] = Proton;
577 }
578 }
579 else if ( np == (nm-1))
580 {
581 if (ran < 0.5)
582 {
583 pv[0] = AntiSigmaMinus;
584 }
585 else if (ran < 0.75)
586 {
587 pv[1] = Proton;
588 }
589 else
590 {
591 pv[0] = AntiSigmaZero;
592 pv[1] = Proton;
593 }
594 }
595 else if (np == (nm+1))
596 {
597 pv[0] = AntiSigmaPlus;
598 }
599 else
600 {
601 pv[0] = AntiSigmaMinus;
602 pv[1] = Proton;
603 }
604 }
605
606 }
607 else // annihilation
608 {
609 if ( availableEnergy > 2. * PionPlus.getMass() )
610 {
611
612 G4double aleab = std::log(availableEnergy);
613 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
614 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
615
616 // normalization constant for kno-distribution.
617 // calculate first the sum of all constants, check for numerical problems.
618 G4double test, dum, anpn = 0.0;
619
620 for (nt=2; nt<=numSec; nt++) {
621 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
622 dum = pi*nt/(2.0*n*n);
623
624 if (std::fabs(dum) < 1.0) {
625 if( test >= 1.0e-10 )anpn += dum*test;
626 } else {
627 anpn += dum*test;
628 }
629 }
630
631 G4double ran = G4UniformRand();
632 G4double excs = 0.0;
633 if (targetCode == protonCode) {
634 counter = -1;
635 for (np=1; np<numSec/3; np++) {
636 nm = np-1;
637 for( nz=0; nz<numSec/3; nz++ )
638 {
639 if( ++counter < numMulAn )
640 {
641 nt = np+nm+nz;
642 if( (nt>1) && (nt<=numSec) ) {
643 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
644 dum = (pi/anpn)*nt*protmulAn[counter]*protnormAn[nt-1]/(2.0*n*n);
645
646 if (std::fabs(dum) < 1.0) {
647 if( test >= 1.0e-10 )excs += dum*test;
648 } else {
649 excs += dum*test;
650 }
651
652 if (ran < excs) goto outOfLoopAn; //----------------------->
653 }
654 }
655 }
656 }
657 // 3 previous loops continued to the end
658 inElastic = false; // quasi-elastic scattering
659 return;
660
661 } else { // target must be a neutron
662 counter = -1;
663 for (np=0; np<numSec/3; np++) {
664 nm = np;
665 for( nz=0; nz<numSec/3; nz++ ) {
666 if (++counter < numMulAn) {
667 nt = np+nm+nz;
668 if( (nt>1) && (nt<=numSec) ) {
669 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
670 dum = (pi/anpn)*nt*neutmulAn[counter]*neutnormAn[nt-1]/(2.0*n*n);
671
672 if (std::fabs(dum) < 1.0) {
673 if( test >= 1.0e-10 )excs += dum*test;
674 } else {
675 excs += dum*test;
676 }
677
678 if (ran < excs) goto outOfLoopAn; // -------------------------->
679 }
680 }
681 }
682 }
683
684 inElastic = false; // quasi-elastic scattering.
685 return;
686 }
687 outOfLoopAn: // <---------------------------------------------------------
688 vecLen = 0;
689 }
690 }
691
692 nt = np + nm + nz;
693 while ( nt > 0)
694 {
695 G4double ran = G4UniformRand();
696 if ( ran < (G4double)np/nt)
697 {
698 if( np > 0 )
699 { pv[vecLen++] = PionPlus;
700 np--;
701 }
702 }
703 else if ( ran < (G4double)(np+nm)/nt)
704 {
705 if( nm > 0 )
706 {
707 pv[vecLen++] = PionMinus;
708 nm--;
709 }
710 }
711 else
712 {
713 if( nz > 0 )
714 {
715 pv[vecLen++] = PionZero;
716 nz--;
717 }
718 }
719 nt = np + nm + nz;
720 }
721 if (verboseLevel > 1)
722 {
723 G4cout << "Particles produced: " ;
724 G4cout << pv[0].getName() << " " ;
725 G4cout << pv[1].getName() << " " ;
726 for (i=2; i < vecLen; i++)
727 {
728 G4cout << pv[i].getName() << " " ;
729 }
730 G4cout << G4endl;
731 }
732 return;
733 }
734
735
736
737
738
739
740
741
742
743
Note: See TracBrowser for help on using the repository browser.