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