Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4GeneralSpaceNNCrossSection.cc
Go to the documentation of this file.
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 // * *
21 // * Parts of this code which have been developed by QinetiQ Ltd *
22 // * under contract to the European Space Agency (ESA) are the *
23 // * intellectual property of ESA. Rights to use, copy, modify and *
24 // * redistribute this software for general public use are granted *
25 // * in compliance with any licensing, distribution and development *
26 // * policy adopted by the Geant4 Collaboration. This code has been *
27 // * written by QinetiQ Ltd for the European Space Agency, under ESA *
28 // * contract 17191/03/NL/LvH (Aurora Programme). *
29 // * *
30 // * By using, copying, modifying or distributing the software (or *
31 // * any work based on the software) you agree to acknowledge its *
32 // * use in resulting scientific publications, and indicate your *
33 // * acceptance of all terms of the Geant4 Software license. *
34 // ********************************************************************
35 //
36 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
37 //
38 // MODULE: G4GeneralSpaceNNCrossSection.cc
39 //
40 // Version: B.1
41 // Date: 15/04/04
42 // Author: P R Truscott
43 // Organisation: QinetiQ Ltd, UK
44 // Customer: ESA/ESTEC, NOORDWIJK
45 // Contract: 17191/03/NL/LvH
46 //
47 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
48 //
49 // CHANGE HISTORY
50 // --------------
51 //
52 // 6 October 2003, P R Truscott, QinetiQ Ltd, UK
53 // Created.
54 //
55 // 15 March 2004, P R Truscott, QinetiQ Ltd, UK
56 // Beta release
57 //
58 // 19 Aug 2011, V.Ivanchenko move to new design and make x-section per element
59 //
60 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
62 //
63 #include <iomanip>
64 
66 #include "G4SystemOfUnits.hh"
72 #include "G4DynamicParticle.hh"
73 #include "G4Element.hh"
74 #include "G4ParticleDefinition.hh"
75 #include "G4HadTmpUtil.hh"
76 #include "G4Proton.hh"
77 
79  : G4VCrossSectionDataSet("General Space NN")
80 {
81  protonInelastic = new G4ProtonInelasticCrossSection();
82  ionProton = new G4IonProtonCrossSection();
83  TripathiGeneral = new G4TripathiCrossSection();
84  TripathiLight = new G4TripathiLightCrossSection();
85  Shen = new G4IonsShenCrossSection();
86  theProton = G4Proton::Proton();
87 }
88 
90 {
91  delete protonInelastic;
92  delete ionProton;
93  delete TripathiGeneral;
94  delete TripathiLight;
95  delete Shen;
96 }
97 
99 //
100 
102  (const G4DynamicParticle* theProjectile, G4int, const G4Material*)
103 {
104  return (1 <= theProjectile->GetDefinition()->GetBaryonNumber());
105 }
106 
108 //
109 
111  (const G4DynamicParticle* theProjectile, G4int ZT, const G4Material* mat)
112 {
113  G4double result = 0.0;
114  G4int ZP = G4lrint(theProjectile->GetDefinition()->GetPDGCharge()/eplus);
115 
116  if (verboseLevel >= 2)
117  {
118  G4int AP = theProjectile->GetDefinition()->GetBaryonNumber();
119  G4cout <<"In G4GeneralSpaceNNCrossSection::GetCrossSection" <<G4endl;
120  G4cout <<"Projectile A = " <<std::setw(8) <<AP
121  <<" Z = " <<std::setw(8) <<ZP
122  <<" Energy = " <<theProjectile->GetKineticEnergy()/AP
123  <<" MeV/nuc" <<G4endl;
124  G4cout <<"Target Z = " <<std::setw(8) <<ZT
125  <<G4endl;
126  }
127  if (theProjectile->GetDefinition()==theProton)
128  {
129  if (ZT>5)
130  {
131  result = protonInelastic->
132  GetElementCrossSection(theProjectile, ZT, mat);
133  if (verboseLevel >= 2)
134  G4cout <<"Selecting G4ProtonInelasticCrossSection" <<G4endl;
135  }
136  else
137  {
138  result = TripathiLight->
139  GetElementCrossSection(theProjectile, ZT, mat);
140  if (verboseLevel >= 2)
141  G4cout <<"Selecting G4TripathiLightCrossSection" <<G4endl;
142  }
143  }
144  else if (ZT==1)
145  {
146  if (ZP>5)
147  {
148  result = ionProton->
149  GetElementCrossSection(theProjectile, ZT, mat);
150  if (verboseLevel >= 2)
151  G4cout <<"Selecting G4IonProtonCrossSection" <<G4endl;
152  }
153  else
154  {
155  result = TripathiLight->
156  GetElementCrossSection(theProjectile, ZT, mat);
157  if (verboseLevel >= 2)
158  G4cout <<"Selecting G4TripathiLightCrossSection" <<G4endl;
159  }
160  }
161  else
162  {
163  if (TripathiLight->IsElementApplicable(theProjectile, ZT, mat))
164  {
165  result = TripathiLight->
166  GetElementCrossSection(theProjectile, ZT, mat);
167  if (verboseLevel >= 2)
168  G4cout <<"Selecting G4TripathiLightCrossSection" <<G4endl;
169  }
170  else if (TripathiGeneral->IsElementApplicable(theProjectile, ZT, mat))
171  {
172  result = TripathiGeneral->
173  GetElementCrossSection(theProjectile, ZT, mat);
174  if (verboseLevel >= 2)
175  G4cout <<"Selecting G4TripathiCrossSection" <<G4endl;
176  }
177  else if (Shen->IsElementApplicable(theProjectile, ZT, mat))
178  {
179  result = Shen->
180  GetElementCrossSection(theProjectile, ZT, mat);
181  if (verboseLevel >= 2)
182  G4cout <<"Selecting G4IonsShenCrossSection" <<G4endl;
183  }
184  }
185  if (verboseLevel >= 2)
186  {
187  G4cout <<"Cross-section = " <<result/millibarn <<" mbarn" <<G4endl;
188  G4cout <<G4endl;
189  }
190 
191  return result;
192 }
193 
195 //
196 
198 {
199  outFile << "G4GeneralSpaceNNCrossSection calculates hadronic inelastic\n"
200  << "cross sections of interest in space science, by using the\n"
201  << "following cross sections:\n"
202  << "- G4ProtonInelasticCrossSection : for proton projectile\n"
203  << " on targets with Z > 5;\n"
204  << "- G4TripathiLightCrossSection : for proton projectile\n"
205  << " on targets with Z <= 5;\n"
206  << " for targets with Z = 1 and projectile Z <= 5;\n"
207  << " for neutron, or deuteron, or 3He, or alpha projectile\n"
208  << " with kinetic energy less than 10 GeV per nucleon,\n"
209  << " in any target;\n"
210  << " for 3He and 4He targets, for any projectile with\n"
211  << " kinetic energy less than 10 GeV per nucleon;\n"
212  << "- G4IonProtonCrossSection : for projectile with Z > 5\n"
213  << " on hydrogen target;\n"
214  << "- G4TripathiCrossSection : for any projectile with A >=3\n"
215  << " and kinetic energy less than 1 GeV per nucleon,\n"
216  << " for any target, if the previous cross section is\n"
217  << " not applicable;\n"
218  << "- G4IonsShenCrossSection : in all remaining cases, up to\n"
219  << " projectile kinetic energy of 1 TeV per nucleon.\n";
220 }
221 
G4double G4ParticleHPJENDLHEData::G4double result
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
virtual void CrossSectionDescription(std::ostream &outFile) const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
static constexpr double eplus
Definition: G4SIunits.hh:199
virtual G4double GetElementCrossSection(const G4DynamicParticle *theProjectile, G4int Z, const G4Material *)
static G4Proton * Proton()
Definition: G4Proton.cc:93
int G4lrint(double ad)
Definition: templates.hh:163
virtual G4bool IsElementApplicable(const G4DynamicParticle *theProjectile, G4int Z, const G4Material *)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
static constexpr double millibarn
Definition: G4SIunits.hh:106