source: trunk/source/processes/hadronic/models/high_energy/src/G4HEAntiSigmaMinusInelastic.cc @ 1340

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

update ti head

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: G4HEAntiSigmaMinusInelastic.cc,v 1.15 2008/03/17 20:49:17 dennis Exp $
28// GEANT4 tag $Name: geant4-09-03-ref-09 $
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 "G4HEAntiSigmaMinusInelastic.hh"
46
47G4HadFinalState *  G4HEAntiSigmaMinusInelastic::
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 << "GHEAntiSigmaMinusInelastic: incident energy < 1 GeV" << G4endl;
66      }
67    if(verboseLevel > 1)
68      {
69        G4cout << "G4HEAntiSigmaMinusInelastic::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    G4HEVector targetParticle;
103    if(G4UniformRand() < atomicNumber/atomicWeight)
104      { 
105        targetParticle.setDefinition("Proton");
106      }
107    else
108      { 
109        targetParticle.setDefinition("Neutron");
110      }
111
112    G4double targetMass         = targetParticle.getMass();
113    G4double centerOfMassEnergy = std::sqrt( incidentMass*incidentMass + targetMass*targetMass
114                                       + 2.0*targetMass*incidentTotalEnergy);
115    G4double availableEnergy    = centerOfMassEnergy - targetMass - incidentMass;
116
117                                                                // this was the meaning of inElastic in the
118                                                                // original Gheisha stand-alone version.
119//    G4bool   inElastic          = InElasticCrossSectionInFirstInt
120//                                    (availableEnergy, incidentCode, incidentTotalMomentum); 
121                                                                // by unknown reasons, it has been replaced
122                                                                // to the following code in Geant???
123    G4bool inElastic = true;
124//    if (G4UniformRand() < elasticCrossSection/totalCrossSection) inElastic = false;   
125
126    vecLength  = 0;           
127       
128    if(verboseLevel > 1)
129      G4cout << "ApplyYourself: CallFirstIntInCascade for particle "
130           << incidentCode << G4endl;
131
132    G4bool successful = false; 
133   
134    if(inElastic || (!inElastic && atomicWeight < 1.5))
135      { 
136        FirstIntInCasAntiSigmaMinus(inElastic, availableEnergy, pv, vecLength,
137                                    incidentParticle, targetParticle, atomicWeight);
138
139        if(verboseLevel > 1)
140           G4cout << "ApplyYourself::StrangeParticlePairProduction" << G4endl; 
141
142
143        if ((vecLength > 0) && (availableEnergy > 1.)) 
144                   StrangeParticlePairProduction( availableEnergy, centerOfMassEnergy,
145                                                  pv, vecLength,
146                                                  incidentParticle, targetParticle);
147            HighEnergyCascading( successful, pv, vecLength,
148                                 excitationEnergyGNP, excitationEnergyDTA,
149                                 incidentParticle, targetParticle,
150                                 atomicWeight, atomicNumber);
151        if (!successful)
152            HighEnergyClusterProduction( successful, pv, vecLength,
153                                         excitationEnergyGNP, excitationEnergyDTA,
154                                         incidentParticle, targetParticle,
155                                         atomicWeight, atomicNumber);
156        if (!successful) 
157            MediumEnergyCascading( successful, pv, vecLength, 
158                                   excitationEnergyGNP, excitationEnergyDTA, 
159                                   incidentParticle, targetParticle,
160                                   atomicWeight, atomicNumber);
161
162        if (!successful)
163            MediumEnergyClusterProduction( successful, pv, vecLength,
164                                           excitationEnergyGNP, excitationEnergyDTA,       
165                                           incidentParticle, targetParticle,
166                                           atomicWeight, atomicNumber);
167        if (!successful)
168            QuasiElasticScattering( successful, pv, vecLength,
169                                    excitationEnergyGNP, excitationEnergyDTA,
170                                    incidentParticle, targetParticle, 
171                                    atomicWeight, atomicNumber);
172      }
173    if (!successful)
174      { 
175            ElasticScattering( successful, pv, vecLength,
176                               incidentParticle,   
177                               atomicWeight, atomicNumber);
178      }
179
180    if (!successful)
181      { 
182         G4cout << "GHEInelasticInteraction::ApplyYourself fails to produce final state particles" << G4endl;
183      }
184      FillParticleChange(pv,  vecLength);
185      delete [] pv;
186      theParticleChange.SetStatusChange(stopAndKill);
187      return & theParticleChange;
188  }
189
190void
191G4HEAntiSigmaMinusInelastic::FirstIntInCasAntiSigmaMinus( G4bool &inElastic,
192                                                          const G4double availableEnergy,
193                                                          G4HEVector pv[],
194                                                          G4int &vecLen,
195                                                          G4HEVector incidentParticle,
196                                                          G4HEVector targetParticle,
197                                                          const G4double atomicWeight)
198
199// AntiSigma- undergoes interaction with nucleon within a nucleus.  Check if it is
200// energetically possible to produce pions/kaons.  In not, assume nuclear excitation
201// occurs and input particle is degraded in energy. No other particles are produced.
202// If reaction is possible, find the correct number of pions/protons/neutrons
203// produced using an interpolation to multiplicity data.  Replace some pions or
204// protons/neutrons by kaons or strange baryons according to the average
205// multiplicity per inelastic reaction.
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; 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+1); 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-2); 
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-1; 
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 == neutronCode)
351             {
352               if(ran < 0.2)
353                 {
354                   pv[0] = AntiSigmaZero;
355                   pv[1] = Proton;
356                 }
357               else if (ran < 0.4)
358                 {
359                   pv[0] = AntiLambda;
360                   pv[1] = Proton;
361                 }
362               else if (ran < 0.6)
363                 {
364                   pv[0] = Proton;
365                   pv[1] = AntiLambda;
366                 }
367               else if (ran < 0.8)
368                 {
369                   pv[0] = Proton;
370                   pv[1] = AntiSigmaZero;
371                 }
372               else
373                 {
374                   pv[0] = Neutron;
375                   pv[1] = AntiSigmaMinus;
376                 }     
377             }
378           else
379             {
380               pv[0] = Proton;
381               pv[1] = AntiSigmaMinus;
382             } 
383         }   
384       return;
385     }
386   else if (availableEnergy <= PionPlus.getMass())
387       return;
388
389                                                  //   inelastic scattering
390
391   np = 0; nm = 0; nz = 0;
392   G4double anhl[] = {1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 1.00, 0.97, 0.88, 
393                      0.85, 0.81, 0.75, 0.64, 0.64, 0.55, 0.55, 0.45, 0.47, 0.40, 
394                      0.39, 0.36, 0.33, 0.10, 0.01};
395   G4int            iplab =      G4int( incidentTotalMomentum*10.);
396   if ( iplab >  9) iplab = 10 + G4int( (incidentTotalMomentum  -1.)*5. );         
397   if ( iplab > 14) iplab = 15 + G4int(  incidentTotalMomentum  -2.     );
398   if ( iplab > 22) iplab = 23 + G4int( (incidentTotalMomentum -10.)/10.); 
399                    iplab = std::min(24, iplab);
400
401   if ( G4UniformRand() > anhl[iplab] )
402     {                                           // non- annihilation channels
403
404                         //  number of total particles vs. centre of mass Energy - 2*proton mass
405   
406           G4double aleab = std::log(availableEnergy);
407           G4double n     = 3.62567+aleab*(0.665843+aleab*(0.336514
408                            + aleab*(0.117712+0.0136912*aleab))) - 2.0;
409   
410                         // normalization constant for kno-distribution.
411                         // calculate first the sum of all constants, check for numerical problems.   
412           G4double test, dum, anpn = 0.0;
413
414           for (nt=1; nt<=numSec; nt++) {
415             test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
416             dum = pi*nt/(2.0*n*n);
417             if (std::fabs(dum) < 1.0) { 
418               if( test >= 1.0e-10 )anpn += dum*test;
419             } else { 
420               anpn += dum*test;
421             }
422           }
423   
424           G4double ran = G4UniformRand();
425           G4double excs = 0.0;
426           if( targetCode == protonCode ) 
427             {
428               counter = -1;
429               for( np=0; np<numSec/3; np++ ) 
430                  {
431                    for( nm=std::max(0,np-2); nm<=np; nm++ ) 
432                       {
433                         for( nz=0; nz<numSec/3; nz++ ) 
434                            {
435                              if( ++counter < numMul ) 
436                                {
437                                  nt = np+nm+nz;
438                                  if ( (nt>0) && (nt<=numSec) ) {
439                                    test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
440                                    dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
441                                    if (std::fabs(dum) < 1.0) { 
442                                      if( test >= 1.0e-10 )excs += dum*test;
443                                    } else { 
444                                      excs += dum*test;
445                                    }
446
447                                    if (ran < excs) goto outOfLoop;      //----------------------->
448                                  }   
449                                }   
450                            }     
451                       }                                                                                 
452                  }
453       
454                                              // 3 previous loops continued to the end
455               inElastic = false;                 // quasi-elastic scattering   
456               return;
457             }
458           else   
459             {                                         // target must be a neutron
460               counter = -1;
461               for( np=0; np<numSec/3; np++ ) 
462                  {
463                    for( nm=std::max(0,np-1); nm<=(np+1); nm++ ) 
464                       {
465                         for( nz=0; nz<numSec/3; nz++ ) 
466                            {
467                              if( ++counter < numMul ) 
468                                {
469                                  nt = np+nm+nz;
470                                  if ( (nt>0) && (nt<=numSec) ) {
471                                    test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
472                                    dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
473                                    if (std::fabs(dum) < 1.0) { 
474                                      if( test >= 1.0e-10 )excs += dum*test;
475                                    } else { 
476                                      excs += dum*test;
477                                    }
478
479                                    if (ran < excs) goto outOfLoop;       // -------------------------->
480                                  }
481                                }
482                            }
483                       }
484                  }
485                                                      // 3 previous loops continued to the end
486               inElastic = false;                     // quasi-elastic scattering.
487               return;
488             }
489         
490       outOfLoop:           //  <------------------------------------------------------------------------   
491   
492       ran = G4UniformRand();
493
494       if( targetCode == protonCode)
495         {
496           if( np == nm)
497             {
498             }
499           else if (np == (nm+1))
500             {
501               if( ran < 0.50)
502                 {
503                   pv[1] = Neutron;
504                 }
505               else if (ran < 0.75)
506                 {
507                   pv[0] = AntiSigmaZero;
508                 }
509               else
510                 {
511                   pv[0] = AntiLambda;
512                 }
513             }
514           else     
515             {
516               if (ran < 0.5)
517                 {
518                   pv[0] = AntiSigmaZero;
519                   pv[1] = Neutron;
520                 }
521               else
522                 {
523                   pv[0] = AntiLambda;
524                   pv[1] = Neutron;
525                 } 
526             } 
527         } 
528       else
529         {
530           if( np == nm)
531             {
532               if (ran < 0.5)
533                 {
534                 }
535               else if(ran < 0.75)
536                 {
537                   pv[0] = AntiSigmaZero;
538                   pv[1] = Proton;
539                 }
540               else
541                 {
542                   pv[0] = AntiLambda;
543                   pv[1] = Proton;
544                 }
545             } 
546           else if ( np == (nm+1))
547             {
548               if (ran < 0.5)
549                 {
550                   pv[0] = AntiSigmaZero;
551                 }
552               else 
553                 {
554                   pv[0] = AntiLambda;
555                 }
556             }
557           else 
558             {
559               pv[1] = Proton;
560             }
561         }     
562     
563     }
564   else                                                               // annihilation
565     { 
566       if ( availableEnergy > 2. * PionPlus.getMass() )
567         {
568
569           G4double aleab = std::log(availableEnergy);
570           G4double n     = 3.62567+aleab*(0.665843+aleab*(0.336514
571                            + aleab*(0.117712+0.0136912*aleab))) - 2.0;
572   
573                      //   normalization constant for kno-distribution.
574                      //   calculate first the sum of all constants, check for numerical problems.   
575           G4double test, dum, anpn = 0.0;
576
577           for (nt=2; nt<=numSec; nt++) {
578             test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
579             dum = pi*nt/(2.0*n*n);
580             if (std::fabs(dum) < 1.0) { 
581               if( test >= 1.0e-10 )anpn += dum*test;
582             } else { 
583               anpn += dum*test;
584             }
585           }
586   
587           G4double ran = G4UniformRand();
588           G4double excs = 0.0;
589           if( targetCode == protonCode ) 
590             {
591               counter = -1;
592               for( np=2; np<numSec/3; np++ ) 
593                  {
594                    nm = np-2; 
595                    for( nz=0; nz<numSec/3; nz++ ) 
596                      {
597                        if( ++counter < numMulAn ) 
598                          {
599                            nt = np+nm+nz;
600                            if ( (nt>1) && (nt<=numSec) ) {
601                              test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
602                              dum = (pi/anpn)*nt*protmulAn[counter]*protnormAn[nt-1]/(2.0*n*n);
603                              if (std::fabs(dum) < 1.0) { 
604                                if( test >= 1.0e-10 )excs += dum*test;
605                              } else { 
606                                excs += dum*test;
607                              }
608
609                              if (ran < excs) goto outOfLoopAn;      //----------------------->
610                            }   
611                          }   
612                      }     
613                 }                                                                                 
614                                              // 3 previous loops continued to the end
615               inElastic = false;                 // quasi-elastic scattering   
616               return;
617             }
618           else   
619             {                                         // target must be a neutron
620               counter = -1;
621               for( np=1; np<numSec/3; np++ ) 
622                 { 
623                   nm = np-1; 
624                   for( nz=0; nz<numSec/3; nz++ ) 
625                      {
626                        if (++counter < numMulAn) {
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*neutmulAn[counter]*neutnormAn[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               inElastic = false;                     // quasi-elastic scattering.
643               return;
644             }
645       outOfLoopAn:           //  <------------------------------------------------------------------   
646       vecLen = 0;
647         }
648     }
649
650   nt = np + nm + nz;
651   while ( nt > 0)
652       {
653         G4double ran = G4UniformRand();
654         if ( ran < (G4double)np/nt)
655            { 
656              if( np > 0 ) 
657                { pv[vecLen++] = PionPlus;
658                  np--;
659                }
660            }
661         else if ( ran < (G4double)(np+nm)/nt)
662            {   
663              if( nm > 0 )
664                { 
665                  pv[vecLen++] = PionMinus;
666                  nm--;
667                }
668            }
669         else
670            {
671              if( nz > 0 )
672                { 
673                  pv[vecLen++] = PionZero;
674                  nz--;
675                }
676            }
677         nt = np + nm + nz;
678       } 
679   if (verboseLevel > 1)
680      {
681        G4cout << "Particles produced: " ;
682        G4cout << pv[0].getName() << " " ;
683        G4cout << pv[1].getName() << " " ;
684        for (i=2; i < vecLen; i++)   
685            { 
686              G4cout << pv[i].getName() << " " ;
687            }
688         G4cout << G4endl;
689      }
690   return;
691 }
692
693
694
695
696
697
698
699
700
701
Note: See TracBrowser for help on using the repository browser.