source: trunk/source/processes/hadronic/models/high_energy/src/G4HEAntiSigmaPlusInelastic.cc @ 1347

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

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