source: trunk/source/processes/hadronic/models/high_energy/src/G4HEAntiXiZeroInelastic.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.1 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: G4HEAntiXiZeroInelastic.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 "G4HEAntiXiZeroInelastic.hh"
46
47G4HadFinalState * G4HEAntiXiZeroInelastic::
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 A = targetNucleus.GetN();
54 const G4double Z = targetNucleus.GetZ();
55 G4HEVector incidentParticle(aParticle);
56
57 G4double atomicNumber = Z;
58 G4double atomicWeight = A;
59
60 G4int incidentCode = incidentParticle.getCode();
61 G4double incidentMass = incidentParticle.getMass();
62 G4double incidentTotalEnergy = incidentParticle.getEnergy();
63 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
64 G4double incidentKineticEnergy = incidentTotalEnergy - incidentMass;
65
66 if(incidentKineticEnergy < 1.)
67 {
68 G4cout << "GHEAntiXiZeroInelastic: incident energy < 1 GeV" << G4endl;
69 }
70 if(verboseLevel > 1)
71 {
72 G4cout << "G4HEAntiXiZeroInelastic::ApplyYourself" << G4endl;
73 G4cout << "incident particle " << incidentParticle.getName()
74 << "mass " << incidentMass
75 << "kinetic energy " << incidentKineticEnergy
76 << G4endl;
77 G4cout << "target material with (A,Z) = ("
78 << atomicWeight << "," << atomicNumber << ")" << G4endl;
79 }
80
81 G4double inelasticity = NuclearInelasticity(incidentKineticEnergy,
82 atomicWeight, atomicNumber);
83 if(verboseLevel > 1)
84 G4cout << "nuclear inelasticity = " << inelasticity << G4endl;
85
86 incidentKineticEnergy -= inelasticity;
87
88 G4double excitationEnergyGNP = 0.;
89 G4double excitationEnergyDTA = 0.;
90
91 G4double excitation = NuclearExcitation(incidentKineticEnergy,
92 atomicWeight, atomicNumber,
93 excitationEnergyGNP,
94 excitationEnergyDTA);
95 if(verboseLevel > 1)
96 G4cout << "nuclear excitation = " << excitation << excitationEnergyGNP
97 << excitationEnergyDTA << G4endl;
98
99
100 incidentKineticEnergy -= excitation;
101 incidentTotalEnergy = incidentKineticEnergy + incidentMass;
102 incidentTotalMomentum = std::sqrt( (incidentTotalEnergy-incidentMass)
103 *(incidentTotalEnergy+incidentMass));
104
105
106 G4HEVector targetParticle;
107 if(G4UniformRand() < atomicNumber/atomicWeight)
108 {
109 targetParticle.setDefinition("Proton");
110 }
111 else
112 {
113 targetParticle.setDefinition("Neutron");
114 }
115
116 G4double targetMass = targetParticle.getMass();
117 G4double centerOfMassEnergy = std::sqrt( incidentMass*incidentMass + targetMass*targetMass
118 + 2.0*targetMass*incidentTotalEnergy);
119 G4double availableEnergy = centerOfMassEnergy - targetMass - incidentMass;
120
121 // this was the meaning of inElastic in the
122 // original Gheisha stand-alone version.
123// G4bool inElastic = InElasticCrossSectionInFirstInt
124// (availableEnergy, incidentCode, incidentTotalMomentum);
125 // by unknown reasons, it has been replaced
126 // to the following code in Geant???
127 G4bool inElastic = true;
128// if (G4UniformRand() < elasticCrossSection/totalCrossSection) inElastic = false;
129
130 vecLength = 0;
131
132 if(verboseLevel > 1)
133 G4cout << "ApplyYourself: CallFirstIntInCascade for particle "
134 << incidentCode << G4endl;
135
136 G4bool successful = false;
137
138 if(inElastic || (!inElastic && atomicWeight < 1.5))
139 {
140 FirstIntInCasAntiXiZero(inElastic, availableEnergy, pv, vecLength,
141 incidentParticle, targetParticle, atomicWeight);
142
143 if(verboseLevel > 1)
144 G4cout << "ApplyYourself::StrangeParticlePairProduction" << G4endl;
145
146
147 if ((vecLength > 0) && (availableEnergy > 1.))
148 StrangeParticlePairProduction( availableEnergy, centerOfMassEnergy,
149 pv, vecLength,
150 incidentParticle, targetParticle);
151
152 HighEnergyCascading( successful, pv, vecLength,
153 excitationEnergyGNP, excitationEnergyDTA,
154 incidentParticle, targetParticle,
155 atomicWeight, atomicNumber);
156 if (!successful)
157 HighEnergyClusterProduction( successful, pv, vecLength,
158 excitationEnergyGNP, excitationEnergyDTA,
159 incidentParticle, targetParticle,
160 atomicWeight, atomicNumber);
161 if (!successful)
162 MediumEnergyCascading( successful, pv, vecLength,
163 excitationEnergyGNP, excitationEnergyDTA,
164 incidentParticle, targetParticle,
165 atomicWeight, atomicNumber);
166
167 if (!successful)
168 MediumEnergyClusterProduction( successful, pv, vecLength,
169 excitationEnergyGNP, excitationEnergyDTA,
170 incidentParticle, targetParticle,
171 atomicWeight, atomicNumber);
172 if (!successful)
173 QuasiElasticScattering( successful, pv, vecLength,
174 excitationEnergyGNP, excitationEnergyDTA,
175 incidentParticle, targetParticle,
176 atomicWeight, atomicNumber);
177 }
178 if (!successful)
179 {
180 ElasticScattering( successful, pv, vecLength,
181 incidentParticle,
182 atomicWeight, atomicNumber);
183 }
184
185 if (!successful)
186 {
187 G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles" << G4endl;
188 }
189 FillParticleChange(pv, vecLength);
190 delete [] pv;
191 theParticleChange.SetStatusChange(stopAndKill);
192 return & theParticleChange;
193 }
194
195void
196G4HEAntiXiZeroInelastic::FirstIntInCasAntiXiZero( G4bool &inElastic,
197 const G4double availableEnergy,
198 G4HEVector pv[],
199 G4int &vecLen,
200 G4HEVector incidentParticle,
201 G4HEVector targetParticle,
202 const G4double atomicWeight)
203
204// AntiXi0 undergoes interaction with nucleon within a nucleus.
205// As in Geant3, we think that this routine has absolutely no influence
206// on the whole performance of the program. Take AntiLambda instaed.
207// ( decay Xi0 -> L Pi > 99 % )
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 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
462 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
463 if (std::fabs(dum) < 1.0) {
464 if( test >= 1.0e-10 )excs += dum*test;
465 } else {
466 excs += dum*test;
467 }
468
469 if (ran < excs) goto outOfLoop; //----------------------->
470 }
471 }
472 }
473 }
474 }
475
476 // 3 previous loops continued to the end
477 inElastic = false; // quasi-elastic scattering
478 return;
479 }
480 else
481 { // target must be a neutron
482 counter = -1;
483 for( np=0; np<numSec/3; np++ )
484 {
485 for( nm=std::max(0,np-1); nm<=(np+2); nm++ )
486 {
487 for( nz=0; nz<numSec/3; nz++ )
488 {
489 if( ++counter < numMul )
490 {
491 nt = np+nm+nz;
492 if ( (nt>0) && (nt<=numSec) ) {
493 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
494 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
495 if (std::fabs(dum) < 1.0) {
496 if( test >= 1.0e-10 )excs += dum*test;
497 } else {
498 excs += dum*test;
499 }
500
501 if (ran < excs) goto outOfLoop; // -------------------------->
502 }
503 }
504 }
505 }
506 }
507 // 3 previous loops continued to the end
508 inElastic = false; // quasi-elastic scattering.
509 return;
510 }
511
512 outOfLoop: // <------------------------------------------------------------------------
513
514 ran = G4UniformRand();
515
516 if( targetCode == protonCode)
517 {
518 if( np == nm)
519 {
520 if (ran < 0.40)
521 {
522 }
523 else if (ran < 0.8)
524 {
525 pv[0] = AntiSigmaZero;
526 }
527 else
528 {
529 pv[0] = AntiSigmaMinus;
530 pv[1] = Neutron;
531 }
532 }
533 else if (np == (nm+1))
534 {
535 if( ran < 0.25)
536 {
537 pv[1] = Neutron;
538 }
539 else if (ran < 0.5)
540 {
541 pv[0] = AntiSigmaZero;
542 pv[1] = Neutron;
543 }
544 else
545 {
546 pv[0] = AntiSigmaPlus;
547 }
548 }
549 else if (np == (nm-1))
550 {
551 pv[0] = AntiSigmaMinus;
552 }
553 else
554 {
555 pv[0] = AntiSigmaPlus;
556 pv[1] = Neutron;
557 }
558 }
559 else
560 {
561 if( np == nm)
562 {
563 if (ran < 0.4)
564 {
565 }
566 else if(ran < 0.8)
567 {
568 pv[0] = AntiSigmaZero;
569 }
570 else
571 {
572 pv[0] = AntiSigmaPlus;
573 pv[1] = Proton;
574 }
575 }
576 else if ( np == (nm-1))
577 {
578 if (ran < 0.5)
579 {
580 pv[0] = AntiSigmaMinus;
581 }
582 else if (ran < 0.75)
583 {
584 pv[1] = Proton;
585 }
586 else
587 {
588 pv[0] = AntiSigmaZero;
589 pv[1] = Proton;
590 }
591 }
592 else if (np == (nm+1))
593 {
594 pv[0] = AntiSigmaPlus;
595 }
596 else
597 {
598 pv[0] = AntiSigmaMinus;
599 pv[1] = Proton;
600 }
601 }
602
603 }
604 else // annihilation
605 {
606 if ( availableEnergy > 2. * PionPlus.getMass() )
607 {
608
609 G4double aleab = std::log(availableEnergy);
610 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
611 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
612
613 // normalization constant for kno-distribution.
614 // calculate first the sum of all constants, check for numerical problems.
615 G4double test, dum, anpn = 0.0;
616
617 for (nt=2; nt<=numSec; nt++) {
618 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
619 dum = pi*nt/(2.0*n*n);
620 if (std::fabs(dum) < 1.0) {
621 if( test >= 1.0e-10 )anpn += dum*test;
622 } else {
623 anpn += dum*test;
624 }
625 }
626
627 G4double ran = G4UniformRand();
628 G4double excs = 0.0;
629 if( targetCode == protonCode )
630 {
631 counter = -1;
632 for( np=1; np<numSec/3; np++ )
633 {
634 nm = np-1;
635 for( nz=0; nz<numSec/3; nz++ )
636 {
637 if( ++counter < numMulAn )
638 {
639 nt = np+nm+nz;
640 if ( (nt>1) && (nt<=numSec) ) {
641 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
642 dum = (pi/anpn)*nt*protmulAn[counter]*protnormAn[nt-1]/(2.0*n*n);
643 if (std::fabs(dum) < 1.0) {
644 if( test >= 1.0e-10 )excs += dum*test;
645 } else {
646 excs += dum*test;
647 }
648
649 if (ran < excs) goto outOfLoopAn; //----------------------->
650 }
651 }
652 }
653 }
654 // 3 previous loops continued to the end
655 inElastic = false; // quasi-elastic scattering
656 return;
657 }
658 else
659 { // target must be a neutron
660 counter = -1;
661 for( np=0; np<numSec/3; np++ )
662 {
663 nm = np;
664 for( nz=0; nz<numSec/3; nz++ )
665 {
666 if( ++counter < numMulAn )
667 {
668 nt = np+nm+nz;
669 if ( (nt>1) && (nt<=numSec) ) {
670 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
671 dum = (pi/anpn)*nt*neutmulAn[counter]*neutnormAn[nt-1]/(2.0*n*n);
672 if (std::fabs(dum) < 1.0) {
673 if( test >= 1.0e-10 )excs += dum*test;
674 } else {
675 excs += dum*test;
676 }
677 if (ran < excs) goto outOfLoopAn; // -------------------------->
678 }
679 }
680 }
681 }
682 inElastic = false; // quasi-elastic scattering.
683 return;
684 }
685 outOfLoopAn: // <------------------------------------------------------------------
686 vecLen = 0;
687 }
688 }
689
690 nt = np + nm + nz;
691 while ( nt > 0)
692 {
693 G4double ran = G4UniformRand();
694 if ( ran < (G4double)np/nt)
695 {
696 if( np > 0 )
697 { pv[vecLen++] = PionPlus;
698 np--;
699 }
700 }
701 else if ( ran < (G4double)(np+nm)/nt)
702 {
703 if( nm > 0 )
704 {
705 pv[vecLen++] = PionMinus;
706 nm--;
707 }
708 }
709 else
710 {
711 if( nz > 0 )
712 {
713 pv[vecLen++] = PionZero;
714 nz--;
715 }
716 }
717 nt = np + nm + nz;
718 }
719 if (verboseLevel > 1)
720 {
721 G4cout << "Particles produced: " ;
722 G4cout << pv[0].getCode() << " " ;
723 G4cout << pv[1].getCode() << " " ;
724 for (i=2; i < vecLen; i++)
725 {
726 G4cout << pv[i].getCode() << " " ;
727 }
728 G4cout << G4endl;
729 }
730 return;
731 }
732
733
734
735
736
737
738
739
740
741
Note: See TracBrowser for help on using the repository browser.