source: trunk/source/processes/hadronic/models/low_energy/src/G4LEAntiSigmaMinusInelastic.cc

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

update ti head

File size: 20.0 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: G4LEAntiSigmaMinusInelastic.cc,v 1.11 2006/06/29 20:44:49 gunter Exp $
28// GEANT4 tag $Name: geant4-09-03-ref-09 $
29//
30 // Hadronic Process: AntiSigmaMinus Inelastic Process
31 // J.L. Chuma, TRIUMF, 19-Feb-1997
32 // Last modified: 27-Mar-1997
33 // J.P. Wellisch: 25.Apr-97: counter errors removed lines 426, 447
34 // Modified by J.L.Chuma 30-Apr-97: added originalTarget for CalculateMomenta
35 
36#include "G4LEAntiSigmaMinusInelastic.hh"
37#include "Randomize.hh"
38
39 G4HadFinalState *
40  G4LEAntiSigmaMinusInelastic::ApplyYourself( const G4HadProjectile &aTrack,
41                                                   G4Nucleus &targetNucleus )
42  { 
43    const G4HadProjectile *originalIncident = &aTrack;
44    if (originalIncident->GetKineticEnergy()<= 0.1*MeV) 
45    {
46      theParticleChange.SetStatusChange(isAlive);
47      theParticleChange.SetEnergyChange(aTrack.GetKineticEnergy());
48      theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit()); 
49      return &theParticleChange;     
50    }
51    //
52    // create the target particle
53    //
54    G4DynamicParticle *originalTarget = targetNucleus.ReturnTargetParticle();
55   
56    if( verboseLevel > 1 )
57    {
58      const G4Material *targetMaterial = aTrack.GetMaterial();
59      G4cout << "G4LEAntiSigmaMinusInelastic::ApplyYourself called" << G4endl;
60      G4cout << "kinetic energy = " << originalIncident->GetKineticEnergy()/MeV << "MeV, ";
61      G4cout << "target material = " << targetMaterial->GetName() << ", ";
62      G4cout << "target particle = " << originalTarget->GetDefinition()->GetParticleName()
63           << G4endl;
64    }
65    //
66    // Fermi motion and evaporation
67    // As of Geant3, the Fermi energy calculation had not been Done
68    //
69    G4double ek = originalIncident->GetKineticEnergy()/MeV;
70    G4double amas = originalIncident->GetDefinition()->GetPDGMass()/MeV;
71    G4ReactionProduct modifiedOriginal;
72    modifiedOriginal = *originalIncident;
73   
74    G4double tkin = targetNucleus.Cinema( ek );
75    ek += tkin;
76    modifiedOriginal.SetKineticEnergy( ek*MeV );
77    G4double et = ek + amas;
78    G4double p = std::sqrt( std::abs((et-amas)*(et+amas)) );
79    G4double pp = modifiedOriginal.GetMomentum().mag()/MeV;
80    if( pp > 0.0 )
81    {
82      G4ThreeVector momentum = modifiedOriginal.GetMomentum();
83      modifiedOriginal.SetMomentum( momentum * (p/pp) );
84    }
85    //
86    // calculate black track energies
87    //
88    tkin = targetNucleus.EvaporationEffects( ek );
89    ek -= tkin;
90    modifiedOriginal.SetKineticEnergy( ek*MeV );
91    et = ek + amas;
92    p = std::sqrt( std::abs((et-amas)*(et+amas)) );
93    pp = modifiedOriginal.GetMomentum().mag()/MeV;
94    if( pp > 0.0 )
95    {
96      G4ThreeVector momentum = modifiedOriginal.GetMomentum();
97      modifiedOriginal.SetMomentum( momentum * (p/pp) );
98    }
99    G4ReactionProduct currentParticle = modifiedOriginal;
100    G4ReactionProduct targetParticle;
101    targetParticle = *originalTarget;
102    currentParticle.SetSide( 1 ); // incident always goes in forward hemisphere
103    targetParticle.SetSide( -1 );  // target always goes in backward hemisphere
104    G4bool incidentHasChanged = false;
105    G4bool targetHasChanged = false;
106    G4bool quasiElastic = false;
107    G4FastVector<G4ReactionProduct,GHADLISTSIZE> vec;  // vec will contain the secondary particles
108    G4int vecLen = 0;
109    vec.Initialize( 0 );
110   
111    const G4double cutOff = 0.1;
112    const G4double anni = std::min( 1.3*currentParticle.GetTotalMomentum()/GeV, 0.4 );
113    if( (currentParticle.GetKineticEnergy()/MeV > cutOff) || (G4UniformRand() > anni) )
114      Cascade( vec, vecLen,
115               originalIncident, currentParticle, targetParticle,
116               incidentHasChanged, targetHasChanged, quasiElastic );
117   
118    CalculateMomenta( vec, vecLen,
119                      originalIncident, originalTarget, modifiedOriginal,
120                      targetNucleus, currentParticle, targetParticle,
121                      incidentHasChanged, targetHasChanged, quasiElastic );
122   
123    SetUpChange( vec, vecLen,
124                 currentParticle, targetParticle,
125                 incidentHasChanged );
126   
127    delete originalTarget;
128    return &theParticleChange;
129  }
130 
131 void
132  G4LEAntiSigmaMinusInelastic::Cascade(
133   G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec,
134   G4int& vecLen,
135   const G4HadProjectile *originalIncident,
136   G4ReactionProduct &currentParticle,
137   G4ReactionProduct &targetParticle,
138   G4bool &incidentHasChanged,
139   G4bool &targetHasChanged,
140   G4bool &quasiElastic )
141  {
142    // derived from original FORTRAN code CASASM by H. Fesefeldt (13-Sep-1987)
143    //
144    // AntiSigmaMinus undergoes interaction with nucleon within a nucleus.  Check if it is
145    // energetically possible to produce pions/kaons.  In not, assume nuclear excitation
146    // occurs and input particle is degraded in energy. No other particles are produced.
147    // If reaction is possible, find the correct number of pions/protons/neutrons
148    // produced using an interpolation to multiplicity data.  Replace some pions or
149    // protons/neutrons by kaons or strange baryons according to the average
150    // multiplicity per Inelastic reaction.
151    //
152    const G4double mOriginal = originalIncident->GetDefinition()->GetPDGMass()/MeV;
153    const G4double etOriginal = originalIncident->GetTotalEnergy()/MeV;
154    const G4double pOriginal = originalIncident->GetTotalMomentum()/MeV;
155    const G4double targetMass = targetParticle.GetMass()/MeV;
156    G4double centerofmassEnergy = std::sqrt( mOriginal*mOriginal +
157                                        targetMass*targetMass +
158                                        2.0*targetMass*etOriginal );
159    G4double availableEnergy = centerofmassEnergy-(targetMass+mOriginal);
160   
161    static G4bool first = true;
162    const G4int numMul = 1200;
163    const G4int numMulA = 400;
164    const G4int numSec = 60;
165    static G4double protmul[numMul], protnorm[numSec]; // proton constants
166    static G4double neutmul[numMul], neutnorm[numSec]; // neutron constants
167    static G4double protmulA[numMulA], protnormA[numSec]; // proton constants
168    static G4double neutmulA[numMulA], neutnormA[numSec]; // neutron constants
169    // np = number of pi+, nm = number of pi-, nz = number of pi0
170    G4int counter, nt=0, np=0, nm=0, nz=0;
171    G4double test;
172    const G4double c = 1.25;   
173    const G4double b[2] = { 0.7, 0.7 };
174    if( first )       // compute normalization constants, this will only be Done once
175    {
176      first = false;
177      G4int i;
178      for( i=0; i<numMul; ++i )protmul[i] = 0.0;
179      for( i=0; i<numSec; ++i )protnorm[i] = 0.0;
180      counter = -1;
181      for( np=0; np<(numSec/3); ++np )
182      {
183        for( nm=std::max(0,np-2); nm<=np; ++nm )
184        {
185          for( nz=0; nz<numSec/3; ++nz )
186          {
187            if( ++counter < numMul )
188            {
189              nt = np+nm+nz;
190              if( nt>0 && nt<=numSec )
191              {
192                protmul[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
193                protnorm[nt-1] += protmul[counter];
194              }
195            }
196          }
197        }
198      }
199      for( i=0; i<numMul; ++i )neutmul[i] = 0.0;
200      for( i=0; i<numSec; ++i )neutnorm[i] = 0.0;
201      counter = -1;
202      for( np=0; np<numSec/3; ++np )
203      {
204        for( nm=std::max(0,np-1); nm<=(np+1); ++nm )
205        {
206          for( nz=0; nz<numSec/3; ++nz )
207          {
208            if( ++counter < numMul )
209            {
210              nt = np+nm+nz;
211              if( nt>0 && nt<=numSec )
212              {
213                neutmul[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
214                neutnorm[nt-1] += neutmul[counter];
215              }
216            }
217          }
218        }
219      }
220      for( i=0; i<numSec; ++i )
221      {
222        if( protnorm[i] > 0.0 )protnorm[i] = 1.0/protnorm[i];
223        if( neutnorm[i] > 0.0 )neutnorm[i] = 1.0/neutnorm[i];
224      }
225      //
226      // do the same for annihilation channels
227      //
228      for( i=0; i<numMulA; ++i )protmulA[i] = 0.0;
229      for( i=0; i<numSec; ++i )protnormA[i] = 0.0;
230      counter = -1;
231      for( np=2; np<(numSec/3); ++np )
232      {
233        nm = np-2;
234        for( nz=0; nz<numSec/3; ++nz )
235        {
236          if( ++counter < numMulA )
237          {
238            nt = np+nm+nz;
239            if( nt>1 && nt<=numSec )
240            {
241              protmulA[counter] = Pmltpc(np,nm,nz,nt,b[0],c);
242              protnormA[nt-1] += protmulA[counter];
243            }
244          }
245        }
246      }
247      for( i=0; i<numMulA; ++i )neutmulA[i] = 0.0;
248      for( i=0; i<numSec; ++i )neutnormA[i] = 0.0;
249      counter = -1;
250       for( np=1; np<numSec/3; ++np )
251       {
252        nm = np-1;
253        for( nz=0; nz<numSec/3; ++nz )
254        {
255          if( ++counter < numMulA )
256          {
257            nt = np+nm+nz;
258            if( nt>1 && nt<=numSec )
259            {
260              neutmulA[counter] = Pmltpc(np,nm,nz,nt,b[1],c);
261              neutnormA[nt-1] += neutmulA[counter];
262            }
263          }
264        }
265      }
266      for( i=0; i<numSec; ++i )
267      {
268        if( protnormA[i] > 0.0 )protnormA[i] = 1.0/protnormA[i];
269        if( neutnormA[i] > 0.0 )neutnormA[i] = 1.0/neutnormA[i];
270      }
271    }   // end of initialization
272    const G4double expxu = 82.;           // upper bound for arg. of exp
273    const G4double expxl = -expxu;        // lower bound for arg. of exp
274    G4ParticleDefinition *aNeutron = G4Neutron::Neutron();
275    G4ParticleDefinition *aProton = G4Proton::Proton();
276    G4ParticleDefinition *aPiPlus = G4PionPlus::PionPlus();
277    G4ParticleDefinition *aKaonMinus = G4KaonMinus::KaonMinus();
278    G4ParticleDefinition *aKaonZL = G4KaonZeroLong::KaonZeroLong();
279    G4ParticleDefinition *aKaonPlus = G4KaonPlus::KaonPlus();
280    G4ParticleDefinition *anAntiLambda = G4AntiLambda::AntiLambda();
281    G4ParticleDefinition *anAntiSigmaZero = G4AntiSigmaZero::AntiSigmaZero();
282    const G4double anhl[] = {1.00,1.00,1.00,1.00,1.00,1.00,1.00,1.00,0.97,0.88,
283                             0.85,0.81,0.75,0.64,0.64,0.55,0.55,0.45,0.47,0.40,
284                             0.39,0.36,0.33,0.10,0.01};
285    G4int iplab = G4int( pOriginal/GeV*10.0 );
286    if( iplab >  9 )iplab = G4int( (pOriginal/GeV- 1.0)*5.0  ) + 10;
287    if( iplab > 14 )iplab = G4int(  pOriginal/GeV- 2.0       ) + 15;
288    if( iplab > 23 )iplab = G4int( (pOriginal/GeV-10.0)/10.0 ) + 23;
289    if( iplab > 24 )iplab = 24;
290    if( G4UniformRand() > anhl[iplab] )
291    {
292      if( availableEnergy <= aPiPlus->GetPDGMass()/MeV )
293      {
294        quasiElastic = true;
295        return;
296      }
297      G4double n, anpn;
298      GetNormalizationConstant( availableEnergy, n, anpn );
299      G4double ran = G4UniformRand();
300      G4double dum, excs = 0.0;
301      if( targetParticle.GetDefinition() == aProton )
302      {
303        counter = -1;
304        for( np=0; np<numSec/3 && ran>=excs; ++np )
305        {
306          for( nm=std::max(0,np-2); nm<=np && ran>=excs; ++nm )
307          {
308            for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
309            {
310              if( ++counter < numMul )
311              {
312                nt = np+nm+nz;
313                if( nt>0 && nt<=numSec )
314                {
315                  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
316                  dum = (pi/anpn)*nt*protmul[counter]*protnorm[nt-1]/(2.0*n*n);
317                  if( std::fabs(dum) < 1.0 )
318                  {
319                    if( test >= 1.0e-10 )excs += dum*test;
320                  }
321                  else
322                    excs += dum*test;
323                }
324              }
325            }
326          }
327        }
328        if( ran >= excs )  // 3 previous loops continued to the end
329        {
330          quasiElastic = true;
331          return;
332        }
333        np--; nm--; nz--;
334        G4int ncht = std::min( 3, std::max( 1, np-nm+1 ) );
335        switch( ncht )
336        {
337         case 1:
338           break;
339         case 2:
340           if( G4UniformRand() < 0.5 )
341           {
342             targetParticle.SetDefinitionAndUpdateE( aNeutron );
343             targetHasChanged = true;
344           }
345           else
346           {
347             if( G4UniformRand() < 0.5 )
348               currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
349             else
350               currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
351             incidentHasChanged = true;
352           }             
353           break;
354         case 3:
355           if( G4UniformRand() < 0.5 )
356             currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
357           else
358             currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
359           incidentHasChanged = true;
360           targetParticle.SetDefinitionAndUpdateE( aNeutron );
361           targetHasChanged = true;
362           break;
363        }
364      }
365      else  // target must be a neutron
366      {
367        counter = -1;
368        for( np=0; np<numSec/3 && ran>=excs; ++np )
369        {
370          for( nm=std::max(0,np-1); nm<=(np+1) && ran>=excs; ++nm )
371          {
372            for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
373            {
374              if( ++counter < numMul )
375              {
376                nt = np+nm+nz;
377                if( nt>0 && nt<=numSec )
378                {
379                  test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
380                  dum = (pi/anpn)*nt*neutmul[counter]*neutnorm[nt-1]/(2.0*n*n);
381                  if( std::fabs(dum) < 1.0 )
382                  {
383                    if( test >= 1.0e-10 )excs += dum*test;
384                  }
385                  else
386                    excs += dum*test;
387                }
388              }
389            }
390          }
391        }
392        if( ran >= excs )  // 3 previous loops continued to the end
393        {
394          quasiElastic = true;
395          return;
396        }
397        np--; nm--; nz--;
398        G4int ncht = std::min( 3, std::max( 1, np-nm+2 ) );
399        switch( ncht )
400        {
401         case 1:
402           {
403             targetParticle.SetDefinitionAndUpdateE( aProton );
404             targetHasChanged = true;
405           }
406           break;
407         case 2:
408           if( G4UniformRand() < 0.5 )
409           {
410             if( G4UniformRand() < 0.5 )
411             {
412               currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
413               incidentHasChanged = true;
414               targetParticle.SetDefinitionAndUpdateE( aProton );
415               targetHasChanged = true;
416             }
417           }
418           else
419           {
420             if( G4UniformRand() < 0.5 )
421             {
422               currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
423               incidentHasChanged = true;
424               targetParticle.SetDefinitionAndUpdateE( aProton );
425               targetHasChanged = true;
426             }
427           }
428           break;
429         case 3:
430           if( G4UniformRand() < 0.5 )
431             currentParticle.SetDefinitionAndUpdateE( anAntiLambda );
432           else
433             currentParticle.SetDefinitionAndUpdateE( anAntiSigmaZero );
434           incidentHasChanged = true;
435           break;
436        }
437      }
438    }
439    else  // random number <= anhl[iplab]
440    {
441      if( centerofmassEnergy <= aPiPlus->GetPDGMass()/MeV+aKaonPlus->GetPDGMass()/MeV )
442      {
443        quasiElastic = true;
444        return;
445      }
446      G4double n, anpn;
447      GetNormalizationConstant( -centerofmassEnergy, n, anpn );
448      G4double ran = G4UniformRand();
449      G4double dum, excs = 0.0;
450      if( targetParticle.GetDefinition() == aProton )
451      {
452        counter = -1;
453        for( np=2; np<numSec/3 && ran>=excs; ++np )
454        {
455          nm=np-2;
456          for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
457          {
458            if( ++counter < numMulA )
459            {
460              nt = np+nm+nz;
461              if( nt>1 && nt<=numSec )
462              {
463                test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
464                dum = (pi/anpn)*nt*protmulA[counter]*protnormA[nt-1]/(2.0*n*n);
465                if( std::fabs(dum) < 1.0 )
466                {
467                  if( test >= 1.0e-10 )excs += dum*test;
468                }
469                else
470                  excs += dum*test;
471              }
472            }
473          }
474        }
475        if( ran >= excs )  // 3 previous loops continued to the end
476        {
477          quasiElastic = true;
478          return;
479        }
480        np--; nz--;
481      }
482      else  // target must be a neutron
483      {
484        counter = -1;
485        for( np=1; np<numSec/3 && ran>=excs; ++np )
486        {
487          nm = np-1;
488          for( nz=0; nz<numSec/3 && ran>=excs; ++nz )
489          {
490            if( ++counter < numMulA )
491            {
492              nt = np+nm+nz;
493              if( nt>1 && nt<=numSec )
494              {
495                test = std::exp( std::min( expxu, std::max( expxl, -(pi/4.0)*(nt*nt)/(n*n) ) ) );
496                dum = (pi/anpn)*nt*neutmulA[counter]*neutnormA[nt-1]/(2.0*n*n);
497                if( std::fabs(dum) < 1.0 )
498                {
499                  if( test >= 1.0e-10 )excs += dum*test;
500                }
501                else
502                  excs += dum*test;
503              }
504            }
505          }
506        }
507        if( ran >= excs )  // 3 previous loops continued to the end
508        {
509          quasiElastic = true;
510          return;
511        }
512        np--; nz--;
513      }
514      if( nz > 0 )
515      {
516        if( nm > 0 )
517        {
518          if( G4UniformRand() < 0.5 )
519          {
520            vec.Initialize( 1 );
521            G4ReactionProduct *p = new G4ReactionProduct;
522            p->SetDefinition( aKaonMinus );
523            (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
524            vec.SetElement( vecLen++, p );
525            --nm;
526          }
527          else   // random number >= 0.5
528          {
529            vec.Initialize( 1 );
530            G4ReactionProduct *p = new G4ReactionProduct;
531            p->SetDefinition( aKaonZL );
532            (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
533            vec.SetElement( vecLen++, p );
534            --nz;
535          }
536        }
537        else   // nm == 0
538        {
539          vec.Initialize( 1 );
540          G4ReactionProduct *p = new G4ReactionProduct;
541          p->SetDefinition( aKaonZL );
542          (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
543          vec.SetElement( vecLen++, p );
544          --nz;
545        }
546      }
547      else    //  nz == 0
548      {
549        if( nm > 0 )
550        {
551          vec.Initialize( 1 );
552          G4ReactionProduct *p = new G4ReactionProduct;
553          p->SetDefinition( aKaonMinus );
554          (G4UniformRand() < 0.5) ? p->SetSide( -1 ) : p->SetSide( 1 );
555          vec.SetElement( vecLen++, p );
556          --nm;
557        }
558      }
559      currentParticle.SetMass( 0.0 );
560      targetParticle.SetMass( 0.0 );
561    }
562    SetUpPions( np, nm, nz, vec, vecLen );
563    return;
564  }
565
566 /* end of file */
567 
Note: See TracBrowser for help on using the repository browser.