Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/cross_sections/src/G4BGGPionElasticXS.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 Class file
 30 //
 31 //
 32 // File name:     G4BGGPionElasticXS
 33 //
 34 // Author:        Vladimir Ivanchenko
 35 //
 36 // Creation date: 01.10.2003
 37 // Modifications:
 38 //
 39 // -------------------------------------------------------------------
 40 //
 41 
 42 #include "G4BGGPionElasticXS.hh"
 43 #include "G4SystemOfUnits.hh"
 44 #include "G4ComponentGGHadronNucleusXsc.hh"
 45 #include "G4UPiNuclearCrossSection.hh"
 46 #include "G4HadronNucleonXsc.hh"
 47 #include "G4NuclearRadii.hh"
 48 
 49 #include "G4Proton.hh"
 50 #include "G4PionPlus.hh"
 51 #include "G4PionMinus.hh"
 52 #include "G4NistManager.hh"
 53 #include "G4HadronicParameters.hh"
 54 #include "G4Pow.hh"
 55 
 56 G4double G4BGGPionElasticXS::theGlauberFacPiPlus[93] = {0.0};
 57 G4double G4BGGPionElasticXS::theCoulombFacPiPlus[93] = {0.0};
 58 G4double G4BGGPionElasticXS::theGlauberFacPiMinus[93] = {0.0};
 59 G4double G4BGGPionElasticXS::theCoulombFacPiMinus[93] = {0.0};
 60 G4int G4BGGPionElasticXS::theA[93] = {0};
 61 
 62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 63 
 64 G4BGGPionElasticXS::G4BGGPionElasticXS(const G4ParticleDefinition* p) 
 65  : G4VCrossSectionDataSet("BarashenkovGlauberGribov") 
 66 {
 67   verboseLevel = 0;
 68   fGlauberEnergy = 91.*CLHEP::GeV;
 69   fLowEnergy = 20.*CLHEP::MeV;
 70   fLowestEnergy = 1.*CLHEP::MeV;
 71   SetMinKinEnergy(0.0);
 72   SetMaxKinEnergy( G4HadronicParameters::Instance()->GetMaxEnergy() );
 73 
 74   fPion    = new G4UPiNuclearCrossSection();
 75   fGlauber = new G4ComponentGGHadronNucleusXsc();
 76   fHadron  = new G4HadronNucleonXsc();
 77 
 78   fG4pow = G4Pow::GetInstance();
 79 
 80   theProton = G4Proton::Proton();
 81   thePiPlus = G4PionPlus::PionPlus();
 82   isPiplus = (p == thePiPlus);
 83   SetForAllAtomsAndEnergies(true);
 84 
 85   if (0 == theA[0]) { Initialise(); } 
 86 }
 87 
 88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 89 
 90 G4BGGPionElasticXS::~G4BGGPionElasticXS()
 91 {
 92   delete fHadron;
 93 }
 94 
 95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 96 
 97 G4bool 
 98 G4BGGPionElasticXS::IsElementApplicable(const G4DynamicParticle*, G4int,
 99                                         const G4Material*)
100 {
101   return true;
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
106 G4bool G4BGGPionElasticXS::IsIsoApplicable(const G4DynamicParticle*, 
107                                            G4int Z, G4int, 
108                                            const G4Element*, const G4Material*)
109 {
110   return (1 == Z);
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
114 
115 G4double
116 G4BGGPionElasticXS::GetElementCrossSection(const G4DynamicParticle* dp,
117                                            G4int ZZ, const G4Material*)
118 {
119   // this method should be called only for Z > 1
120   G4double cross = 0.0;
121   G4double ekin = std::max(dp->GetKineticEnergy(), fLowestEnergy);
122   G4int Z = std::min(ZZ, 92);
123   if(1 == Z) {
124     cross = 1.0115*GetIsoCrossSection(dp,1,1);
125   } else {
126     if(ekin <= fLowEnergy) {
127       cross = (isPiplus) ? theCoulombFacPiPlus[Z]*CoulombFactorPiPlus(ekin, Z)
128         : theCoulombFacPiMinus[Z]*FactorPiMinus(ekin);
129     } else if(ekin > fGlauberEnergy) {
130       cross = (isPiplus) ? theGlauberFacPiPlus[Z] : theGlauberFacPiMinus[Z];
131       cross *= fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
132     } else {
133       cross = fPion->GetElasticCrossSection(dp, Z, theA[Z]);
134     }
135   }
136 #ifdef G4VERBOSE
137   if(verboseLevel > 1) {
138     G4cout << "G4BGGPionElasticXS::GetElementCrossSection  for "
139            << dp->GetDefinition()->GetParticleName()
140            << "  Ekin(GeV)= " << dp->GetKineticEnergy()
141            << " in nucleus Z= " << Z << "  A= " << theA[Z]
142            << " XS(b)= " << cross/barn 
143            << G4endl;
144   }
145 #endif
146   return cross;
147 }
148 
149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150 
151 G4double
152 G4BGGPionElasticXS::GetIsoCrossSection(const G4DynamicParticle* dp, 
153                                        G4int, G4int A, 
154                                        const G4Isotope*,
155                                        const G4Element*,
156                                        const G4Material*)
157 {
158   // this method should be called only for Z = 1
159   fHadron->HadronNucleonXscNS(dp->GetDefinition(), theProton, 
160                               dp->GetKineticEnergy());
161   G4double cross = A*fHadron->GetElasticHadronNucleonXsc();
162 
163 #ifdef G4VERBOSE
164   if(verboseLevel > 1) {
165     G4cout << "G4BGGPionElasticXS::GetIsoCrossSection  for "
166            << dp->GetDefinition()->GetParticleName()
167            << "  Ekin(GeV)= " << dp->GetKineticEnergy()
168            << " in nucleus Z=1  A=" << A
169            << " XS(b)= " << cross/barn 
170            << G4endl;
171   }
172 #endif
173   return cross;
174 }
175 
176 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
177 
178 void G4BGGPionElasticXS::BuildPhysicsTable(const G4ParticleDefinition& p)
179 {
180   if(verboseLevel > 1) {
181     G4cout << "G4BGGPionElasticXS::BuildPhysicsTable for " 
182            << p.GetParticleName() << G4endl;
183   } 
184   if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
185     isPiplus = (&p == G4PionPlus::PionPlus());
186   } else {
187     G4ExceptionDescription ed;
188     ed << "This BGG cross section is applicable only to pions and not to " 
189        << p.GetParticleName() << G4endl; 
190     G4Exception("G4BGGPionElasticXS::BuildPhysicsTable", "had001", 
191                 FatalException, ed);
192   }
193 }
194 
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196 
197 void G4BGGPionElasticXS::Initialise()
198 {
199   theA[0] = theA[1] = 1;
200   G4ThreeVector mom(0.0,0.0,1.0);
201   G4DynamicParticle dp(thePiPlus, mom, fGlauberEnergy);
202 
203   G4NistManager* nist = G4NistManager::Instance();
204 
205   G4double csup, csdn;
206   for (G4int iz=2; iz<93; ++iz) {
207 
208     G4int A = G4lrint(nist->GetAtomicMassAmu(iz));
209     theA[iz] = A;
210 
211     csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
212     csdn = fPion->GetElasticCrossSection(&dp, iz, A);
213     theGlauberFacPiPlus[iz] = csdn/csup;
214   }
215 
216   dp.SetDefinition(G4PionMinus::PionMinus());
217   for (G4int iz=2; iz<93; ++iz) {
218     csup = fGlauber->GetElasticGlauberGribov(&dp, iz, theA[iz]);
219     csdn = fPion->GetElasticCrossSection(&dp, iz, theA[iz]);
220     theGlauberFacPiMinus[iz] = csdn/csup;
221 
222     if (verboseLevel > 1) {
223       G4cout << "Z= " << iz <<  "  A= " << theA[iz] 
224        << " factorPiPlus= " << theGlauberFacPiPlus[iz] 
225        << " factorPiMinus= " << theGlauberFacPiMinus[iz] 
226        << G4endl; 
227     }
228   }
229   theCoulombFacPiPlus[1] = 1.0;
230   theCoulombFacPiMinus[1]= 1.0;
231   dp.SetKineticEnergy(fLowEnergy);
232   dp.SetDefinition(thePiPlus);
233   for (G4int iz=2; iz<93; ++iz) {
234     theCoulombFacPiPlus[iz] = fPion->GetElasticCrossSection(&dp, iz, theA[iz])
235       /CoulombFactorPiPlus(fLowEnergy, iz);
236   }
237   dp.SetDefinition(G4PionMinus::PionMinus());
238   for(G4int iz=2; iz<93; ++iz) {
239     theCoulombFacPiMinus[iz] = fPion->GetElasticCrossSection(&dp, iz, theA[iz])
240       /FactorPiMinus(fLowEnergy);
241 
242     if(verboseLevel > 1) {
243       G4cout << "Z= " << iz <<  "  A= " << theA[iz] 
244        << " CoulombFactorPiPlus= " << theCoulombFacPiPlus[iz] 
245        << " CoulombFactorPiMinus= " << theCoulombFacPiMinus[iz] 
246        << G4endl; 
247     }
248   }
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
252 
253 G4double G4BGGPionElasticXS::CoulombFactorPiPlus(G4double kinEnergy, G4int Z)
254 {
255   return (kinEnergy > 0.0) ? 
256     G4NuclearRadii::CoulombFactor(Z, theA[Z], thePiPlus, kinEnergy) : 0.0;
257 }
258 
259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
260 
261 G4double G4BGGPionElasticXS::FactorPiMinus(G4double kinEnergy)
262 {
263   return 1.0/std::sqrt(kinEnergy);
264 }
265 
266 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
267 
268 void
269 G4BGGPionElasticXS::CrossSectionDescription(std::ostream& outFile) const 
270 {
271   outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
272           << "scattering of pions from nuclei at all energies. The\n"
273           << "Barashenkov parameterization is used below 91 GeV and the\n"
274           << "Glauber-Gribov parameterization is used above 91 GeV.\n";
275 }
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
278