source: trunk/source/processes/hadronic/models/high_energy/src/G4HEAntiNeutronInelastic.cc @ 1348

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

geant4 tag 9.4

File size: 24.3 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: G4HEAntiNeutronInelastic.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 is 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 "G4HEAntiNeutronInelastic.hh"
43
44G4HadFinalState*
45G4HEAntiNeutronInelastic::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 << "GHEAntiNeutronInelastic: incident energy < 1 GeV" << G4endl;;
62
63  if (verboseLevel > 1) {
64    G4cout << "G4HEAntiNeutronInelastic::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  FirstIntInCasAntiNeutron(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  HighEnergyCascading(successful, pv, vecLength,
129                      excitationEnergyGNP, excitationEnergyDTA,
130                      incidentParticle, targetParticle,
131                      atomicWeight, atomicNumber);
132  if (!successful)
133    HighEnergyClusterProduction(successful, pv, vecLength,
134                                excitationEnergyGNP, excitationEnergyDTA,
135                                incidentParticle, targetParticle,
136                                atomicWeight, atomicNumber);
137  if (!successful) 
138    MediumEnergyCascading(successful, pv, vecLength, 
139                          excitationEnergyGNP, excitationEnergyDTA, 
140                          incidentParticle, targetParticle,
141                          atomicWeight, atomicNumber);
142
143  if (!successful)
144    MediumEnergyClusterProduction(successful, pv, vecLength,
145                                  excitationEnergyGNP, excitationEnergyDTA,       
146                                  incidentParticle, targetParticle,
147                                  atomicWeight, atomicNumber);
148  if (!successful)
149    QuasiElasticScattering(successful, pv, vecLength,
150                           excitationEnergyGNP, excitationEnergyDTA,
151                           incidentParticle, targetParticle, 
152                           atomicWeight, atomicNumber);
153  if (!successful)
154    ElasticScattering(successful, pv, vecLength,
155                      incidentParticle,   
156                      atomicWeight, atomicNumber);
157
158  if (!successful)
159    G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles" 
160           << G4endl;
161
162  FillParticleChange(pv,  vecLength);
163  delete [] pv;
164  theParticleChange.SetStatusChange(stopAndKill);
165  return &theParticleChange;
166}
167
168
169void
170G4HEAntiNeutronInelastic::FirstIntInCasAntiNeutron(G4bool& inElastic,
171                                                   const G4double availableEnergy,
172                                                   G4HEVector pv[],
173                                                   G4int& vecLen,
174                                                   const G4HEVector& incidentParticle,
175                                                   const G4HEVector& targetParticle,
176                                                   const G4double atomicWeight)
177
178// AntiNeutron undergoes interaction with nucleon within a nucleus.  Check if
179// it is energetically possible to produce pions/kaons.  If not, assume
180// nuclear excitation occurs and input particle is degraded in energy.  No
181// 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 neutronCode = Neutron.getCode();
199  G4int protonCode = Proton.getCode();
200
201  G4int targetCode = targetParticle.getCode();
202  G4double incidentTotalMomentum = incidentParticle.getTotalMomentum();
203
204  static G4bool first = true;
205  static G4double protmul[numMul],  protnorm[numSec];   // proton constants
206  static G4double protmulAn[numMulAn],protnormAn[numSec]; 
207  static G4double neutmul[numMul],  neutnorm[numSec];   // neutron constants
208  static G4double neutmulAn[numMulAn],neutnormAn[numSec];
209
210  //  misc. local variables
211  //  np = number of pi+,  nm = number of pi-,  nz = number of pi0
212
213  G4int i, counter, nt, np, nm, nz;
214
215   if( first ) 
216     {                 // compute normalization constants, this will only be done once
217       first = false;
218       for( i=0; i<numMul  ; i++ ) protmul[i]  = 0.0;
219       for( i=0; i<numSec  ; i++ ) protnorm[i] = 0.0;
220       counter = -1;
221       for( np=0; np<(numSec/3); np++ ) 
222          {
223            for( nm=std::max(0,np-2); nm<=np; nm++ ) 
224               {
225                 for( nz=0; nz<numSec/3; nz++ ) 
226                    {
227                      if( ++counter < numMul ) 
228                        {
229                          nt = np+nm+nz;
230                          if( (nt>0) && (nt<=numSec) ) 
231                            {
232                              protmul[counter] = pmltpc(np,nm,nz,nt,protb,c);
233                              protnorm[nt-1] += protmul[counter];
234                            }
235                        }
236                    }
237               }
238          }
239       for( i=0; i<numMul; i++ )neutmul[i]  = 0.0;
240       for( i=0; i<numSec; i++ )neutnorm[i] = 0.0;
241       counter = -1;
242       for( np=0; np<numSec/3; np++ ) 
243          {
244            for( nm=std::max(0,np-1); nm<=(np+1); nm++ ) 
245               {
246                 for( nz=0; nz<numSec/3; nz++ ) 
247                    {
248                      if( ++counter < numMul ) 
249                        {
250                          nt = np+nm+nz;
251                          if( (nt>0) && (nt<=numSec) ) 
252                            {
253                               neutmul[counter] = pmltpc(np,nm,nz,nt,neutb,c);
254                               neutnorm[nt-1] += neutmul[counter];
255                            }
256                        }
257                    }
258               }
259          }
260       for( i=0; i<numSec; i++ ) 
261          {
262            if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
263            if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
264          }
265                                                                   // annihilation
266       for( i=0; i<numMulAn  ; i++ ) protmulAn[i]  = 0.0;
267       for( i=0; i<numSec    ; i++ ) protnormAn[i] = 0.0;
268       counter = -1;
269       for( np=1; np<(numSec/3); np++ ) 
270          {
271            nm = std::max(0,np-1); 
272            for( nz=0; nz<numSec/3; nz++ ) 
273               {
274                 if( ++counter < numMulAn ) 
275                   {
276                     nt = np+nm+nz;
277                     if( (nt>1) && (nt<=numSec) ) 
278                       {
279                         protmulAn[counter] = pmltpc(np,nm,nz,nt,protb,c);
280                         protnormAn[nt-1] += protmulAn[counter];
281                       }
282                   }
283               }
284          }
285       for( i=0; i<numMulAn; i++ ) neutmulAn[i]  = 0.0;
286       for( i=0; i<numSec;   i++ ) neutnormAn[i] = 0.0;
287       counter = -1;
288       for( np=0; np<numSec/3; np++ ) 
289          {
290            nm = np; 
291            for( nz=0; nz<numSec/3; nz++ ) 
292               {
293                 if( ++counter < numMulAn ) 
294                   {
295                     nt = np+nm+nz;
296                     if( (nt>1) && (nt<=numSec) ) 
297                       {
298                          neutmulAn[counter] = pmltpc(np,nm,nz,nt,neutb,c);
299                          neutnormAn[nt-1] += neutmulAn[counter];
300                       }
301                   }
302               }
303          }
304       for( i=0; i<numSec; i++ ) 
305          {
306            if( protnormAn[i] > 0.0 )protnormAn[i] = 1.0/protnormAn[i];
307            if( neutnormAn[i] > 0.0 )neutnormAn[i] = 1.0/neutnormAn[i];
308          }
309     }                                          // end of initialization
310
311         
312                                              // initialize the first two places
313                                              // the same as beam and target                                   
314   pv[0] = incidentParticle;
315   pv[1] = targetParticle;
316   vecLen = 2;
317
318   if( !inElastic ) 
319     {                                        // nb n --> pb p
320       if( targetCode == neutronCode ) 
321         {
322           G4double cech[] = {0.50, 0.45, 0.40, 0.35, 0.30, 0.25, 0.06, 0.04, 0.005, 0.};
323
324           G4int iplab = std::min(9, G4int( incidentTotalMomentum*2.5));
325           if( G4UniformRand() < cech[iplab]/std::pow(atomicWeight,0.42) ) 
326             {
327               pv[0] = AntiProton;
328               pv[1] = Proton;
329             }
330         }
331       return;
332     }
333   else if (availableEnergy <= PionPlus.getMass())
334       return;
335
336                                                  //   inelastic scattering
337
338   np = 0, nm = 0, nz = 0;
339   G4double anhl[] = {1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 0.97, 0.88, 
340                      0.85, 0.81, 0.75, 0.64, 0.64, 0.55, 0.55, 0.45, 0.47, 0.40, 
341                      0.39, 0.36, 0.33, 0.10, 0.01};
342   G4int            iplab =      G4int( incidentTotalMomentum*10.);
343   if ( iplab >  9) iplab = 10 + G4int( (incidentTotalMomentum  -1.)*5. );         
344   if ( iplab > 14) iplab = 15 + G4int(  incidentTotalMomentum  -2.     );
345   if ( iplab > 22) iplab = 23 + G4int( (incidentTotalMomentum -10.)/10.); 
346                    iplab = std::min(24, iplab);
347
348   if ( G4UniformRand() > anhl[iplab] )
349     {
350
351       G4double eab = availableEnergy;
352       G4int ieab = G4int( eab*5.0 );
353   
354       G4double supp[] = {0., 0.4, 0.55, 0.65, 0.75, 0.82, 0.86, 0.90, 0.94, 0.98};
355       if( (ieab <= 9) && (G4UniformRand() >= supp[ieab]) ) 
356         {
357                                              //   suppress high multiplicity events at low momentum
358                                              //   only one additional pion will be produced
359           G4double w0, wp, wm, wt, ran;
360           if( targetCode == protonCode )                    // target is a proton
361             {
362               w0 = - sqr(1.+protb)/(2.*c*c);
363               w0 = wp = std::exp(w0);
364               if( G4UniformRand() < w0/(w0+wp) ) 
365                 { np = 0; nm = 0; nz = 1; }
366               else 
367                 { np = 1; nm = 0; nz = 0; }       
368             } 
369           else 
370             {                                               // target is a neutron
371               w0 = -sqr(1.+neutb)/(2.*c*c);
372               w0 = wp = std::exp(w0);
373               wm = -sqr(-1.+neutb)/(2.*c*c);
374               wm = std::exp(wm);
375               wt = w0+wp+wm;
376               wp = w0+wp;
377               ran = G4UniformRand();
378               if( ran < w0/wt)
379                 { np = 0; nm = 0; nz = 1; }       
380               else if( ran < wp/wt)
381                 { np = 1; nm = 0; nz = 0; }       
382               else
383                 { np = 0; nm = 1; nz = 0; }       
384             }
385         }
386       else
387         {
388                         //  number of total particles vs. centre of mass Energy - 2*proton mass
389   
390           G4double aleab = std::log(availableEnergy);
391           G4double n     = 3.62567+aleab*(0.665843+aleab*(0.336514
392                            + aleab*(0.117712+0.0136912*aleab))) - 2.0;
393   
394                         // normalization constant for kno-distribution.
395                         // calculate first the sum of all constants, check for numerical problems.   
396           G4double test, dum, anpn = 0.0;
397
398           for (nt=1; nt<=numSec; nt++) {
399               test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
400               dum = pi*nt/(2.0*n*n);
401
402               if (std::fabs(dum) < 1.0) { 
403                 if( test >= 1.0e-10 )anpn += dum*test;
404               } else { 
405                 anpn += dum*test;
406               }
407           }
408   
409           G4double ran = G4UniformRand();
410           G4double excs = 0.0;
411           if (targetCode == protonCode) {
412             counter = -1;
413             for (np=0; np<numSec/3; np++) {
414               for (nm=std::max(0,np-2); nm<=np; nm++) {
415                 for (nz=0; nz<numSec/3; nz++) {
416                   if (++counter < numMul) {
417                     nt = np+nm+nz;
418                     if ( (nt>0) && (nt<=numSec) ) {
419                       test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
420                       dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
421
422                       if (std::fabs(dum) < 1.0) { 
423                         if( test >= 1.0e-10 )excs += dum*test;
424                       } else { 
425                         excs += dum*test;
426                       }
427
428                       if (ran < excs) goto outOfLoop;      //----------------------->
429                     } 
430                   }   
431                 }     
432               }                                                                                 
433             }
434       
435                                           // 3 previous loops continued to the end
436             inElastic = false;            // quasi-elastic scattering   
437             return;
438
439           } else {                        // target must be a neutron
440             counter = -1;
441             for (np=0; np<numSec/3; np++) {
442               for (nm=std::max(0,np-1); nm<=(np+1); nm++) {
443                 for (nz=0; nz<numSec/3; nz++) {
444                   if (++counter < numMul) {
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       if( targetCode == protonCode)
469         {
470           if( np == nm)
471             {
472             }
473           else if (np == (nm+1))
474             {
475               if( G4UniformRand() < 0.5)
476                 {
477                   pv[1] = Neutron;
478                 }
479               else
480                 {
481                   pv[0] = AntiProton;
482                 }
483             }
484           else     
485             {
486               pv[0] = AntiProton;
487               pv[1] = Neutron;
488             } 
489         } 
490       else
491         {
492           if( np == nm)
493             {
494               if( G4UniformRand() < 0.25)
495                 {
496                   pv[0] = AntiProton;
497                   pv[1] = Proton;
498                 }
499               else
500                 {
501                 }
502             } 
503           else if ( np == (nm-1))
504             {
505               pv[1] = Proton;
506             }
507           else
508             {
509               pv[0] = AntiProton;
510             }
511         }     
512     
513     }
514   else                                                               // annihilation
515     { 
516       if ( availableEnergy > 2. * PionPlus.getMass() )
517         {
518
519           G4double aleab = std::log(availableEnergy);
520           G4double n     = 3.62567+aleab*(0.665843+aleab*(0.336514
521                            + aleab*(0.117712+0.0136912*aleab))) - 2.0;
522   
523                      //   normalization constant for kno-distribution.
524                      //   calculate first the sum of all constants, check for numerical problems.   
525           G4double test, dum, anpn = 0.0;
526
527           for (nt=2; nt<=numSec; nt++) {
528             test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
529             dum = pi*nt/(2.0*n*n);
530             if (std::fabs(dum) < 1.0) { 
531               if( test >= 1.0e-10 )anpn += dum*test;
532             } else { 
533               anpn += dum*test;
534             }
535           }
536   
537           G4double ran = G4UniformRand();
538           G4double excs = 0.0;
539           if (targetCode == protonCode) {
540             counter = -1;
541             for (np=1; np<numSec/3; np++) {
542               nm = np-1; 
543               for (nz=0; nz<numSec/3; nz++) {
544                 if (++counter < numMulAn) {
545                   nt = np+nm+nz;
546                   if ( (nt>1) && (nt<=numSec) ) {
547                     test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
548                     dum = (pi/anpn)*nt*protmulAn[counter]*protnormAn[nt-1]/(2.0*n*n);
549
550                     if (std::fabs(dum) < 1.0) { 
551                       if( test >= 1.0e-10 )excs += dum*test;
552                     } else { 
553                       excs += dum*test;
554                     }
555
556                     if (ran < excs) goto outOfLoopAn;      //----------------------->
557                   }
558                 }
559               }
560             }
561                                             // 3 previous loops continued to the end
562             inElastic = false;            // quasi-elastic scattering   
563             return;
564
565           } else {                          // target must be a neutron
566             counter = -1;
567             for (np=0; np<numSec/3; np++) { 
568               nm = np; 
569               for (nz=0; nz<numSec/3; nz++) {
570                 if (++counter < numMulAn) {
571                   nt = np+nm+nz;
572                   if ( (nt>1) && (nt<=numSec) ) {
573                     test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
574                     dum = (pi/anpn)*nt*neutmulAn[counter]*neutnormAn[nt-1]/(2.0*n*n);
575
576                     if (std::fabs(dum) < 1.0) { 
577                       if( test >= 1.0e-10 )excs += dum*test;
578                     } else { 
579                       excs += dum*test;
580                     }
581
582                     if (ran < excs) goto outOfLoopAn;       // -------------------------->
583                   }
584                 }
585               }
586             }
587             inElastic = false;                     // quasi-elastic scattering.
588             return;
589           }
590       outOfLoopAn:           //  <------------------------------------------------------------------   
591       vecLen = 0;
592         }
593     }
594
595   nt = np + nm + nz;
596   while ( nt > 0)
597       {
598         G4double ran = G4UniformRand();
599         if ( ran < (G4double)np/nt)
600            { 
601              if( np > 0 ) 
602                { pv[vecLen++] = PionPlus;
603                  np--;
604                }
605            }
606         else if ( ran < (G4double)(np+nm)/nt)
607            {   
608              if( nm > 0 )
609                { 
610                  pv[vecLen++] = PionMinus;
611                  nm--;
612                }
613            }
614         else
615            {
616              if( nz > 0 )
617                { 
618                  pv[vecLen++] = PionZero;
619                  nz--;
620                }
621            }
622         nt = np + nm + nz;
623       } 
624   if (verboseLevel > 1)
625      {
626        G4cout << "Particles produced: " ;
627        G4cout << pv[0].getName() << " " ;
628        G4cout << pv[1].getName() << " " ;
629        for (i=2; i < vecLen; i++)   
630            { 
631              G4cout << pv[i].getName() << " " ;
632            }
633         G4cout << G4endl;
634      }
635   return;
636 }
637
638
639
640
641
642
643
644
645
Note: See TracBrowser for help on using the repository browser.