source: trunk/source/processes/hadronic/models/high_energy/src/G4HEAntiXiMinusInelastic.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: 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: G4HEAntiXiMinusInelastic.cc,v 1.15 2008/03/17 20:49:17 dennis Exp $
28// GEANT4 tag $Name: geant4-09-02 $
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 "G4HEAntiXiMinusInelastic.hh"
46
47G4HadFinalState * G4HEAntiXiMinusInelastic::
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 << "GHEAntiXiMinusInelastic: incident energy < 1 GeV" << G4endl;
69 }
70 if(verboseLevel > 1)
71 {
72 G4cout << "G4HEAntiXiMinusInelastic::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 FirstIntInCasAntiXiMinus(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 HighEnergyCascading( successful, pv, vecLength,
152 excitationEnergyGNP, excitationEnergyDTA,
153 incidentParticle, targetParticle,
154 atomicWeight, atomicNumber);
155 if (!successful)
156 HighEnergyClusterProduction( successful, pv, vecLength,
157 excitationEnergyGNP, excitationEnergyDTA,
158 incidentParticle, targetParticle,
159 atomicWeight, atomicNumber);
160 if (!successful)
161 MediumEnergyCascading( successful, pv, vecLength,
162 excitationEnergyGNP, excitationEnergyDTA,
163 incidentParticle, targetParticle,
164 atomicWeight, atomicNumber);
165
166 if (!successful)
167 MediumEnergyClusterProduction( successful, pv, vecLength,
168 excitationEnergyGNP, excitationEnergyDTA,
169 incidentParticle, targetParticle,
170 atomicWeight, atomicNumber);
171 if (!successful)
172 QuasiElasticScattering( successful, pv, vecLength,
173 excitationEnergyGNP, excitationEnergyDTA,
174 incidentParticle, targetParticle,
175 atomicWeight, atomicNumber);
176 }
177 if (!successful)
178 {
179 ElasticScattering( successful, pv, vecLength,
180 incidentParticle,
181 atomicWeight, atomicNumber);
182 }
183
184 if (!successful)
185 {
186 G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles";
187 }
188 FillParticleChange(pv, vecLength);
189 delete [] pv;
190 theParticleChange.SetStatusChange(stopAndKill);
191 return & theParticleChange;
192 }
193
194void
195G4HEAntiXiMinusInelastic::FirstIntInCasAntiXiMinus( G4bool &inElastic,
196 const G4double availableEnergy,
197 G4HEVector pv[],
198 G4int &vecLen,
199 G4HEVector incidentParticle,
200 G4HEVector targetParticle,
201 const G4double atomicWeight)
202
203// AntiXi- undergoes interaction with nucleon within a nucleus.
204// As in Geant3, we think that this routine has absolutely no influence
205// on the whole performance of the program. Take AntiLambda instaed.
206
207 {
208 static const G4double expxu = std::log(MAXFLOAT); // upper bound for arg. of exp
209 static const G4double expxl = -expxu; // lower bound for arg. of exp
210
211 static const G4double protb = 0.7;
212 static const G4double neutb = 0.7;
213 static const G4double c = 1.25;
214
215 static const G4int numMul = 1200;
216 static const G4int numMulAn = 400;
217 static const G4int numSec = 60;
218
219// G4int neutronCode = Neutron.getCode();
220 G4int protonCode = Proton.getCode();
221
222 G4int targetCode = targetParticle.getCode();
223// G4double incidentMass = incidentParticle.getMass();
224// G4double incidentEnergy = incidentParticle.getEnergy();
225 G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
226
227 static G4bool first = true;
228 static G4double protmul[numMul], protnorm[numSec]; // proton constants
229 static G4double protmulAn[numMulAn],protnormAn[numSec];
230 static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
231 static G4double neutmulAn[numMulAn],neutnormAn[numSec];
232
233 // misc. local variables
234 // np = number of pi+, nm = number of pi-, nz = number of pi0
235
236 G4int i, counter, nt, np, nm, nz;
237
238 if( first )
239 { // compute normalization constants, this will only be done once
240 first = false;
241 for( i=0; i<numMul ; i++ ) protmul[i] = 0.0;
242 for( i=0; i<numSec ; i++ ) protnorm[i] = 0.0;
243 counter = -1;
244 for( np=0; np<(numSec/3); np++ )
245 {
246 for( nm=std::max(0,np-2); nm<=(np+1); nm++ )
247 {
248 for( nz=0; nz<numSec/3; nz++ )
249 {
250 if( ++counter < numMul )
251 {
252 nt = np+nm+nz;
253 if( (nt>0) && (nt<=numSec) )
254 {
255 protmul[counter] = pmltpc(np,nm,nz,nt,protb,c);
256 protnorm[nt-1] += protmul[counter];
257 }
258 }
259 }
260 }
261 }
262 for( i=0; i<numMul; i++ )neutmul[i] = 0.0;
263 for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
264 counter = -1;
265 for( np=0; np<numSec/3; np++ )
266 {
267 for( nm=std::max(0,np-1); nm<=(np+2); nm++ )
268 {
269 for( nz=0; nz<numSec/3; nz++ )
270 {
271 if( ++counter < numMul )
272 {
273 nt = np+nm+nz;
274 if( (nt>0) && (nt<=numSec) )
275 {
276 neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c);
277 neutnorm[nt-1] += neutmul[counter];
278 }
279 }
280 }
281 }
282 }
283 for( i=0; i<numSec; i++ )
284 {
285 if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
286 if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
287 }
288 // annihilation
289 for( i=0; i<numMulAn ; i++ ) protmulAn[i] = 0.0;
290 for( i=0; i<numSec ; i++ ) protnormAn[i] = 0.0;
291 counter = -1;
292 for( np=1; np<(numSec/3); np++ )
293 {
294 nm = std::max(0,np-1);
295 for( nz=0; nz<numSec/3; nz++ )
296 {
297 if( ++counter < numMulAn )
298 {
299 nt = np+nm+nz;
300 if( (nt>1) && (nt<=numSec) )
301 {
302 protmulAn[counter] = pmltpc(np,nm,nz,nt,protb,c);
303 protnormAn[nt-1] += protmulAn[counter];
304 }
305 }
306 }
307 }
308 for( i=0; i<numMulAn; i++ ) neutmulAn[i] = 0.0;
309 for( i=0; i<numSec; i++ ) neutnormAn[i] = 0.0;
310 counter = -1;
311 for( np=0; np<numSec/3; np++ )
312 {
313 nm = np;
314 for( nz=0; nz<numSec/3; nz++ )
315 {
316 if( ++counter < numMulAn )
317 {
318 nt = np+nm+nz;
319 if( (nt>1) && (nt<=numSec) )
320 {
321 neutmulAn[counter] = pmltpc(np,nm,nz,nt,neutb,c);
322 neutnormAn[nt-1] += neutmulAn[counter];
323 }
324 }
325 }
326 }
327 for( i=0; i<numSec; i++ )
328 {
329 if( protnormAn[i] > 0.0 )protnormAn[i] = 1.0/protnormAn[i];
330 if( neutnormAn[i] > 0.0 )neutnormAn[i] = 1.0/neutnormAn[i];
331 }
332 } // end of initialization
333
334
335 // initialize the first two places
336 // the same as beam and target
337 pv[0] = incidentParticle;
338 pv[1] = targetParticle;
339 vecLen = 2;
340
341 if( !inElastic )
342 { // some two-body reactions
343 G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
344
345 G4int iplab = std::min(9, G4int( incidentTotalMomentum*2.5));
346 if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) )
347 {
348 G4double ran = G4UniformRand();
349
350 if ( targetCode == protonCode)
351 {
352 if(ran < 0.2)
353 {
354 pv[0] = AntiSigmaZero;
355 }
356 else if (ran < 0.4)
357 {
358 pv[0] = AntiSigmaMinus;
359 pv[1] = Neutron;
360 }
361 else if (ran < 0.6)
362 {
363 pv[0] = Proton;
364 pv[1] = AntiLambda;
365 }
366 else if (ran < 0.8)
367 {
368 pv[0] = Proton;
369 pv[1] = AntiSigmaZero;
370 }
371 else
372 {
373 pv[0] = Neutron;
374 pv[1] = AntiSigmaMinus;
375 }
376 }
377 else
378 {
379 if (ran < 0.2)
380 {
381 pv[0] = AntiSigmaZero;
382 }
383 else if (ran < 0.4)
384 {
385 pv[0] = AntiSigmaPlus;
386 pv[1] = Proton;
387 }
388 else if (ran < 0.6)
389 {
390 pv[0] = Neutron;
391 pv[1] = AntiLambda;
392 }
393 else if (ran < 0.8)
394 {
395 pv[0] = Neutron;
396 pv[1] = AntiSigmaZero;
397 }
398 else
399 {
400 pv[0] = Proton;
401 pv[1] = AntiSigmaPlus;
402 }
403 }
404 }
405 return;
406 }
407 else if (availableEnergy <= PionPlus.getMass())
408 return;
409
410 // inelastic scattering
411
412 np = 0; nm = 0; nz = 0;
413 G4double anhl[] = {1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 0.97, 0.88,
414 0.85, 0.81, 0.75, 0.64, 0.64, 0.55, 0.55, 0.45, 0.47, 0.40,
415 0.39, 0.36, 0.33, 0.10, 0.01};
416 G4int iplab = G4int( incidentTotalMomentum*10.);
417 if ( iplab > 9) iplab = 10 + G4int( (incidentTotalMomentum -1.)*5. );
418 if ( iplab > 14) iplab = 15 + G4int( incidentTotalMomentum -2. );
419 if ( iplab > 22) iplab = 23 + G4int( (incidentTotalMomentum -10.)/10.);
420 iplab = std::min(24, iplab);
421
422 if ( G4UniformRand() > anhl[iplab] )
423 { // non- annihilation channels
424
425 // number of total particles vs. centre of mass Energy - 2*proton mass
426
427 G4double aleab = std::log(availableEnergy);
428 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
429 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
430
431 // normalization constant for kno-distribution.
432 // calculate first the sum of all constants, check for numerical problems.
433 G4double test, dum, anpn = 0.0;
434
435 for (nt=1; nt<=numSec; nt++) {
436 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
437 dum = pi*nt/(2.0*n*n);
438 if (std::fabs(dum) < 1.0) {
439 if( test >= 1.0e-10 )anpn += dum*test;
440 } else {
441 anpn += dum*test;
442 }
443 }
444
445 G4double ran = G4UniformRand();
446 G4double excs = 0.0;
447 if( targetCode == protonCode )
448 {
449 counter = -1;
450 for( np=0; np<numSec/3; np++ )
451 {
452 for( nm=std::max(0,np-2); nm<=(np+1); nm++ )
453 {
454 for( nz=0; nz<numSec/3; nz++ )
455 {
456 if( ++counter < numMul )
457 {
458 nt = np+nm+nz;
459 if ( (nt>0) && (nt<=numSec) ) {
460 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
461 dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
462 if (std::fabs(dum) < 1.0) {
463 if( test >= 1.0e-10 )excs += dum*test;
464 } else {
465 excs += dum*test;
466 }
467
468 if (ran < excs) goto outOfLoop; //----------------------->
469 }
470 }
471 }
472 }
473 }
474
475 // 3 previous loops continued to the end
476 inElastic = false; // quasi-elastic scattering
477 return;
478 }
479 else
480 { // target must be a neutron
481 counter = -1;
482 for( np=0; np<numSec/3; np++ )
483 {
484 for( nm=std::max(0,np-1); nm<=(np+2); nm++ )
485 {
486 for( nz=0; nz<numSec/3; nz++ )
487 {
488 if( ++counter < numMul )
489 {
490 nt = np+nm+nz;
491 if ( (nt>0) && (nt<=numSec) ) {
492 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
493 dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
494 if (std::fabs(dum) < 1.0) {
495 if( test >= 1.0e-10 )excs += dum*test;
496 } else {
497 excs += dum*test;
498 }
499
500 if (ran < excs) goto outOfLoop; // -------------------------->
501 }
502 }
503 }
504 }
505 }
506 // 3 previous loops continued to the end
507 inElastic = false; // quasi-elastic scattering.
508 return;
509 }
510
511 outOfLoop: // <------------------------------------------------------------------------
512
513 ran = G4UniformRand();
514
515 if( targetCode == protonCode)
516 {
517 if( np == nm)
518 {
519 if (ran < 0.40)
520 {
521 }
522 else if (ran < 0.8)
523 {
524 pv[0] = AntiSigmaZero;
525 }
526 else
527 {
528 pv[0] = AntiSigmaMinus;
529 pv[1] = Neutron;
530 }
531 }
532 else if (np == (nm+1))
533 {
534 if( ran < 0.25)
535 {
536 pv[1] = Neutron;
537 }
538 else if (ran < 0.5)
539 {
540 pv[0] = AntiSigmaZero;
541 pv[1] = Neutron;
542 }
543 else
544 {
545 pv[0] = AntiSigmaPlus;
546 }
547 }
548 else if (np == (nm-1))
549 {
550 pv[0] = AntiSigmaMinus;
551 }
552 else
553 {
554 pv[0] = AntiSigmaPlus;
555 pv[1] = Neutron;
556 }
557 }
558 else
559 {
560 if( np == nm)
561 {
562 if (ran < 0.4)
563 {
564 }
565 else if(ran < 0.8)
566 {
567 pv[0] = AntiSigmaZero;
568 }
569 else
570 {
571 pv[0] = AntiSigmaPlus;
572 pv[1] = Proton;
573 }
574 }
575 else if ( np == (nm-1))
576 {
577 if (ran < 0.5)
578 {
579 pv[0] = AntiSigmaMinus;
580 }
581 else if (ran < 0.75)
582 {
583 pv[1] = Proton;
584 }
585 else
586 {
587 pv[0] = AntiSigmaZero;
588 pv[1] = Proton;
589 }
590 }
591 else if (np == (nm+1))
592 {
593 pv[0] = AntiSigmaPlus;
594 }
595 else
596 {
597 pv[0] = AntiSigmaMinus;
598 pv[1] = Proton;
599 }
600 }
601
602 }
603 else // annihilation
604 {
605 if ( availableEnergy > 2. * PionPlus.getMass() )
606 {
607
608 G4double aleab = std::log(availableEnergy);
609 G4double n = 3.62567+aleab*(0.665843+aleab*(0.336514
610 + aleab*(0.117712+0.0136912*aleab))) - 2.0;
611
612 // normalization constant for kno-distribution.
613 // calculate first the sum of all constants, check for numerical problems.
614 G4double test, dum, anpn = 0.0;
615
616 for (nt=2; nt<=numSec; nt++) {
617 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
618 dum = pi*nt/(2.0*n*n);
619 if (std::fabs(dum) < 1.0) {
620 if( test >= 1.0e-10 )anpn += dum*test;
621 } else {
622 anpn += dum*test;
623 }
624 }
625
626 G4double ran = G4UniformRand();
627 G4double excs = 0.0;
628 if( targetCode == protonCode )
629 {
630 counter = -1;
631 for( np=1; np<numSec/3; np++ )
632 {
633 nm = np-1;
634 for( nz=0; nz<numSec/3; nz++ )
635 {
636 if( ++counter < numMulAn )
637 {
638 nt = np+nm+nz;
639 if ( (nt>1) && (nt<=numSec) ) {
640 test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
641 dum = (pi/anpn)*nt*protmulAn[counter]*protnormAn[nt-1]/(2.0*n*n);
642 if (std::fabs(dum) < 1.0) {
643 if( test >= 1.0e-10 )excs += dum*test;
644 } else {
645 excs += dum*test;
646 }
647
648 if (ran < excs) goto outOfLoopAn; //----------------------->
649 }
650 }
651 }
652 }
653 // 3 previous loops continued to the end
654 inElastic = false; // quasi-elastic scattering
655 return;
656 }
657 else
658 { // target must be a neutron
659 counter = -1;
660 for( np=0; np<numSec/3; np++ )
661 {
662 nm = np;
663 for( nz=0; nz<numSec/3; nz++ )
664 {
665 if( ++counter < numMulAn )
666 {
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 if (std::fabs(dum) < 1.0) {
672 if( test >= 1.0e-10 )excs += dum*test;
673 } else {
674 excs += dum*test;
675 }
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].getName() << " " ;
723 G4cout << pv[1].getName() << " " ;
724 for (i=2; i < vecLen; i++)
725 {
726 G4cout << pv[i].getName() << " " ;
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.