Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/coherent_elastic/src/G4ChargeExchangeProcess.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 28 //
 29 // Geant4 Hadron Charge Exchange Process -- source file
 30 //
 31 // Created 21 April 2006 V.Ivanchenko
 32 //
 33 // Modified:
 34 // 24-Apr-06 V.Ivanchenko add neutron scattering on hydrogen from CHIPS
 35 // 07-Jun-06 V.Ivanchenko fix problem of rotation of final state
 36 // 25-Jul-06 V.Ivanchenko add 19 MeV low energy for CHIPS
 37 // 23-Jan-07 V.Ivanchenko add cross section interfaces with Z and A
 38 //                        and do not use CHIPS for cross sections
 39 // 14-Sep-12 M.Kelsey -- Pass subType code to base ctor
 40 // 06-Aug-15 A.Ribon migrating to G4Pow
 41 
 42 #include "G4ChargeExchangeProcess.hh"
 43 #include "globals.hh"
 44 #include "G4SystemOfUnits.hh"
 45 #include "G4CrossSectionDataStore.hh"
 46 #include "G4CrossSectionElastic.hh"
 47 #include "G4ComponentGGHadronNucleusXsc.hh"
 48 #include "G4Element.hh"
 49 #include "G4ElementVector.hh"
 50 #include "G4IsotopeVector.hh"
 51 #include "G4Neutron.hh"
 52 #include "G4Proton.hh"
 53 #include "G4PhysicsLinearVector.hh"
 54 
 55 #include "G4Pow.hh"
 56 
 57 
 58 G4ChargeExchangeProcess::G4ChargeExchangeProcess(const G4String& procName)
 59   : G4HadronicProcess(procName,fChargeExchange), first(true)
 60 {
 61   G4cout << "###=== The class G4ChargeExchangeProcess is obsolete!!!" << G4endl;
 62   G4cout << "###=== It will be removed at the next public release" << G4endl;
 63   thEnergy = 20.*MeV;
 64   pPDG = 0;
 65   verboseLevel= 1;
 66   AddDataSet( new G4CrossSectionElastic( new G4ComponentGGHadronNucleusXsc ) );
 67   theProton   = G4Proton::Proton();
 68   theNeutron  = G4Neutron::Neutron();
 69   theAProton  = G4AntiProton::AntiProton();
 70   theANeutron = G4AntiNeutron::AntiNeutron();
 71   thePiPlus   = G4PionPlus::PionPlus();
 72   thePiMinus  = G4PionMinus::PionMinus();
 73   thePiZero   = G4PionZero::PionZero();
 74   theKPlus    = G4KaonPlus::KaonPlus();
 75   theKMinus   = G4KaonMinus::KaonMinus();
 76   theK0S      = G4KaonZeroShort::KaonZeroShort();
 77   theK0L      = G4KaonZeroLong::KaonZeroLong();
 78   theL        = G4Lambda::Lambda();
 79   theAntiL    = G4AntiLambda::AntiLambda();
 80   theSPlus    = G4SigmaPlus::SigmaPlus();
 81   theASPlus   = G4AntiSigmaPlus::AntiSigmaPlus();
 82   theSMinus   = G4SigmaMinus::SigmaMinus();
 83   theASMinus  = G4AntiSigmaMinus::AntiSigmaMinus();
 84   theS0       = G4SigmaZero::SigmaZero();
 85   theAS0      = G4AntiSigmaZero::AntiSigmaZero();
 86   theXiMinus  = G4XiMinus::XiMinus();
 87   theXi0      = G4XiZero::XiZero();
 88   theAXiMinus = G4AntiXiMinus::AntiXiMinus();
 89   theAXi0     = G4AntiXiZero::AntiXiZero();
 90   theOmega    = G4OmegaMinus::OmegaMinus();
 91   theAOmega   = G4AntiOmegaMinus::AntiOmegaMinus();
 92   theD        = G4Deuteron::Deuteron();
 93   theT        = G4Triton::Triton();
 94   theA        = G4Alpha::Alpha();
 95   theHe3      = G4He3::He3();
 96 }
 97 
 98 G4ChargeExchangeProcess::~G4ChargeExchangeProcess()
 99 {
100   if (factors) delete factors;
101 }
102 
103 void G4ChargeExchangeProcess::
104 BuildPhysicsTable(const G4ParticleDefinition& aParticleType)
105 {
106   if(first) {
107     first = false;
108     theParticle = &aParticleType;
109     pPDG = theParticle->GetPDGEncoding();
110 
111     store = G4HadronicProcess::GetCrossSectionDataStore();
112 
113     const size_t n = 10;
114     if(theParticle == thePiPlus || theParticle == thePiMinus ||
115        theParticle == theKPlus  || theParticle == theKMinus ||
116        theParticle == theK0S    || theParticle == theK0L) {
117 
118       G4double F[n] = {0.33,0.27,0.29,0.31,0.27,0.18,0.13,0.1,0.09,0.07};
119       factors = new G4PhysicsLinearVector(0.0,2.0*GeV,n);
120       for(size_t i=0; i<n; i++) {factors->PutValue(i,F[i]);}
121 
122     } else {
123 
124       G4double F[n] = {0.50,0.45,0.40,0.35,0.30,0.25,0.06,0.04,0.005,0.0};
125       factors = new G4PhysicsLinearVector(0.0,4.0*GeV,n);
126       for(size_t i=0; i<n; i++) {factors->PutValue(i,F[i]);}
127     }
128 
129     if(verboseLevel>1)
130       G4cout << "G4ChargeExchangeProcess for "
131        << theParticle->GetParticleName()
132        << G4endl;
133   }
134   G4HadronicProcess::BuildPhysicsTable(aParticleType);
135 }
136 
137 G4double G4ChargeExchangeProcess::GetElementCrossSection(
138                                   const G4DynamicParticle* dp,
139           const G4Element* elm,
140           const G4Material* mat)
141 {
142   // gives the microscopic cross section in GEANT4 internal units
143   G4double Z = elm->GetZ();
144   G4int iz = G4int(Z);
145   G4double x = 0.0;
146 
147   // The process is effective only above the threshold
148   if(iz == 1 || dp->GetKineticEnergy() < thEnergy) return x;
149 
150   if(verboseLevel>1)
151     G4cout << "G4ChargeExchangeProcess compute GHAD CS for element "
152      << elm->GetName()
153      << G4endl;
154   x = store->GetCrossSection(dp, elm, mat);
155 
156   if(verboseLevel>1)
157     G4cout << "G4ChargeExchangeProcess cross(mb)= " << x/millibarn
158            << "  E(MeV)= " << dp->GetKineticEnergy()
159      << "  " << theParticle->GetParticleName()
160            << "  in Z= " << iz
161      << G4endl;
162   G4bool b;
163   G4double A = elm->GetN();
164   G4double ptot = dp->GetTotalMomentum();
165   x *= factors->GetValue(ptot, b)/G4Pow::GetInstance()->powA(A, 0.42);
166   if(theParticle == thePiPlus || theParticle == theProton ||
167      theParticle == theKPlus  || theParticle == theANeutron)
168     { x *= (1.0 - Z/A); }
169 
170   else if(theParticle == thePiMinus || theParticle == theNeutron ||
171           theParticle == theKMinus  || theParticle == theAProton)
172     { x *= Z/A; }
173 
174   if(theParticle->GetPDGMass() < GeV) {
175     if(ptot > 2.*GeV) x *= 4.0*GeV*GeV/(ptot*ptot);
176   }
177 
178   if(verboseLevel>1) 
179     G4cout << "Corrected cross(mb)= " << x/millibarn << G4endl;
180 
181   return x;
182 }
183 
184 G4bool G4ChargeExchangeProcess::
185 IsApplicable(const G4ParticleDefinition& aParticleType)
186 {
187   const G4ParticleDefinition* p = &aParticleType;
188   return (p == thePiPlus || p == thePiMinus ||
189           p == theProton || p == theNeutron ||
190           p == theAProton|| p == theANeutron||
191     p == theKPlus  || p == theKMinus  ||
192     p == theK0S    || p == theK0L     ||
193     p == theL);
194 }
195 
196 void G4ChargeExchangeProcess::
197 DumpPhysicsTable(const G4ParticleDefinition& aParticleType)
198 {
199   store->DumpPhysicsTable(aParticleType);
200 }
201