Geant4  10.02.p03
G4GeneralSpaceNNCrossSection Class Reference

#include <G4GeneralSpaceNNCrossSection.hh>

Inheritance diagram for G4GeneralSpaceNNCrossSection:
Collaboration diagram for G4GeneralSpaceNNCrossSection:

Public Member Functions

 G4GeneralSpaceNNCrossSection ()
 
 ~G4GeneralSpaceNNCrossSection ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *theProjectile, G4int Z, const G4Material *)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *theProjectile, G4int Z, const G4Material *)
 
virtual void CrossSectionDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Private Attributes

G4ProtonInelasticCrossSectionprotonInelastic
 
G4IonProtonCrossSectionionProton
 
G4TripathiLightCrossSectionTripathiLight
 
G4TripathiCrossSectionTripathiGeneral
 
G4IonsShenCrossSectionShen
 
G4ProtontheProton
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 87 of file G4GeneralSpaceNNCrossSection.hh.

Constructor & Destructor Documentation

◆ G4GeneralSpaceNNCrossSection()

G4GeneralSpaceNNCrossSection::G4GeneralSpaceNNCrossSection ( )

Definition at line 78 of file G4GeneralSpaceNNCrossSection.cc.

79  : G4VCrossSectionDataSet("General Space NN")
80 {
87 }
G4VCrossSectionDataSet(const G4String &nam="")
G4ProtonInelasticCrossSection * protonInelastic
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4TripathiLightCrossSection * TripathiLight
Here is the call graph for this function:

◆ ~G4GeneralSpaceNNCrossSection()

G4GeneralSpaceNNCrossSection::~G4GeneralSpaceNNCrossSection ( )

Definition at line 89 of file G4GeneralSpaceNNCrossSection.cc.

90 {
91  delete protonInelastic;
92  delete ionProton;
93  delete TripathiGeneral;
94  delete TripathiLight;
95  delete Shen;
96 }
G4ProtonInelasticCrossSection * protonInelastic
G4TripathiLightCrossSection * TripathiLight
Here is the call graph for this function:

Member Function Documentation

◆ CrossSectionDescription()

void G4GeneralSpaceNNCrossSection::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 197 of file G4GeneralSpaceNNCrossSection.cc.

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 }

◆ GetElementCrossSection()

G4double G4GeneralSpaceNNCrossSection::GetElementCrossSection ( const G4DynamicParticle theProjectile,
G4int  Z,
const G4Material mat 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 111 of file G4GeneralSpaceNNCrossSection.cc.

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 }
int G4int
Definition: G4Types.hh:78
G4ProtonInelasticCrossSection * protonInelastic
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
virtual G4bool IsElementApplicable(const G4DynamicParticle *theProjectile, G4int Z, const G4Material *)
virtual G4double GetElementCrossSection(const G4DynamicParticle *theProjectile, G4int Z, const G4Material *)
G4TripathiLightCrossSection * TripathiLight
int G4lrint(double ad)
Definition: templates.hh:163
virtual G4bool IsElementApplicable(const G4DynamicParticle *aDP, G4int Z, const G4Material *)
static const double millibarn
Definition: G4SIunits.hh:105
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static const double eplus
Definition: G4SIunits.hh:196
G4double GetPDGCharge() const
virtual G4bool IsElementApplicable(const G4DynamicParticle *aPart, G4int Z, const G4Material *)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsElementApplicable()

G4bool G4GeneralSpaceNNCrossSection::IsElementApplicable ( const G4DynamicParticle theProjectile,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 102 of file G4GeneralSpaceNNCrossSection.cc.

103 {
104  return (1 <= theProjectile->GetDefinition()->GetBaryonNumber());
105 }
G4ParticleDefinition * GetDefinition() const
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ ionProton

G4IonProtonCrossSection* G4GeneralSpaceNNCrossSection::ionProton
private

Definition at line 106 of file G4GeneralSpaceNNCrossSection.hh.

◆ protonInelastic

G4ProtonInelasticCrossSection* G4GeneralSpaceNNCrossSection::protonInelastic
private

Definition at line 105 of file G4GeneralSpaceNNCrossSection.hh.

◆ Shen

G4IonsShenCrossSection* G4GeneralSpaceNNCrossSection::Shen
private

Definition at line 109 of file G4GeneralSpaceNNCrossSection.hh.

◆ theProton

G4Proton* G4GeneralSpaceNNCrossSection::theProton
private

Definition at line 110 of file G4GeneralSpaceNNCrossSection.hh.

◆ TripathiGeneral

G4TripathiCrossSection* G4GeneralSpaceNNCrossSection::TripathiGeneral
private

Definition at line 108 of file G4GeneralSpaceNNCrossSection.hh.

◆ TripathiLight

G4TripathiLightCrossSection* G4GeneralSpaceNNCrossSection::TripathiLight
private

Definition at line 107 of file G4GeneralSpaceNNCrossSection.hh.


The documentation for this class was generated from the following files: