Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4CrossSectionPairGG Class Reference

#include <G4CrossSectionPairGG.hh>

Inheritance diagram for G4CrossSectionPairGG:
G4VCrossSectionDataSet

Public Member Functions

 G4CrossSectionPairGG (G4VCrossSectionDataSet *low, G4double Etransit)
 
virtual ~G4CrossSectionPairGG ()
 
virtual void CrossSectionDescription (std::ostream &) const
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
- 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 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
 

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 51 of file G4CrossSectionPairGG.hh.

Constructor & Destructor Documentation

G4CrossSectionPairGG::G4CrossSectionPairGG ( G4VCrossSectionDataSet low,
G4double  Etransit 
)

Definition at line 48 of file G4CrossSectionPairGG.cc.

References G4NistManager::Instance(), and G4VCrossSectionDataSet::verboseLevel.

49  :
50  G4VCrossSectionDataSet("G4CrossSectionPairGG"), theLowX(low), ETransition(
51  Etransit) {
52  NistMan = G4NistManager::Instance();
53  theHighX = new G4GlauberGribovCrossSection();
54  verboseLevel = 0;
55 }
G4VCrossSectionDataSet(const G4String &nam="")
static G4NistManager * Instance()
G4CrossSectionPairGG::~G4CrossSectionPairGG ( )
virtual

Definition at line 57 of file G4CrossSectionPairGG.cc.

57  {
58  // The cross section registry wil delete these
59  // delete theLowX;
60  // delete theHighX;
61 }

Member Function Documentation

void G4CrossSectionPairGG::BuildPhysicsTable ( const G4ParticleDefinition pDef)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 127 of file G4CrossSectionPairGG.cc.

References G4VCrossSectionDataSet::BuildPhysicsTable(), G4cout, G4endl, G4lrint(), G4NistManager::GetAtomicMassAmu(), G4VCrossSectionDataSet::GetElementCrossSection(), G4GlauberGribovCrossSection::GetInelasticGlauberGribov(), G4VCrossSectionDataSet::GetName(), G4ParticleDefinition::GetParticleName(), G4VCrossSectionDataSet::IsElementApplicable(), and G4VCrossSectionDataSet::verboseLevel.

Referenced by GetElementCrossSection().

127  {
128  theLowX->BuildPhysicsTable(pDef);
129  theHighX->BuildPhysicsTable(pDef);
130 
131  if (verboseLevel > 0) {
132  G4cout << "G4CrossSectionPairGG::BuildPhysicsTable "
133  << theLowX->GetName() << " " << theHighX->GetName() << G4endl;
134  }
135 
136  G4ParticleDefinition * myDef = const_cast<G4ParticleDefinition*>(&pDef);
137  std::vector<ParticleXScale>::iterator iter;
138  iter = scale_factors.begin();
139  while (iter != scale_factors.end() && (*iter).first != myDef) {
140  ++iter;
141  }
142 
143  // new particle, initialise
144 
145  G4Material* mat = 0;
146 
147  if (iter == scale_factors.end()) {
148  XS_factors factors(93);
149  G4ThreeVector mom(0.0, 0.0, 1.0);
150  G4DynamicParticle DynPart(myDef, mom, ETransition); // last is kinetic Energy
151 
152  if (verboseLevel > 0) {
153  G4cout << "G4CrossSectionPairGG::BuildPhysicsTable for particle "
154  << pDef.GetParticleName() << G4endl;
155  }
156  for (G4int aZ = 1; aZ < 93; ++aZ) {
157  factors[aZ] = 1.; // default, to give reasonable value if only high is applicable
158  G4int AA = G4lrint(NistMan->GetAtomicMassAmu(aZ));
159  G4bool isApplicable = theLowX->IsElementApplicable(&DynPart, aZ,
160  mat) && (aZ > 1);
161 
162  if (isApplicable) {
163  factors[aZ] = theLowX->GetElementCrossSection(&DynPart, aZ, mat)
164  / theHighX->GetInelasticGlauberGribov(&DynPart, aZ, AA);
165 
166  }
167  if (verboseLevel > 0) {
168  G4cout << "Z=" << aZ << ", A=" << AA << ", scale="
169  << factors[aZ];
170  if (verboseLevel == 1) {
171  G4cout << G4endl;
172  } else {
173  if (isApplicable) {
174  G4cout << ", low / high "
175  << theLowX->GetElementCrossSection(&DynPart, aZ,
176  mat) << " "
177  << theHighX->GetInelasticGlauberGribov(&DynPart,
178  aZ, AA) << G4endl;
179  } else {
180  G4cout << ", N/A" << G4endl;
181  }
182  }
183  }
184  }
185  ParticleXScale forPart(myDef, factors);
186  scale_factors.push_back(forPart);
187  }
188 }
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
const G4String & GetName() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetAtomicMassAmu(const G4String &symb) const
#define G4endl
Definition: G4ios.hh:61
void G4CrossSectionPairGG::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 63 of file G4CrossSectionPairGG.cc.

64  {
65  outFile << "G4CrossSectionPairGG is used to add the relativistic rise to\n"
66  << "hadronic cross section data sets above a given energy. In this\n"
67  << "case, the Glauber-Gribov cross section is used above 91 GeV.\n"
68  << "At this energy the low energy cross section is smoothly joined\n"
69  << "to the high energy cross section. Below 91 GeV the Barashenkov\n"
70  << "cross section is used for pions (G4PiNuclearCrossSection), the\n"
71  << "Axen-Wellisch cross section is used for protons\n"
72  << "(G4ProtonInelasticCrossSection), and the Wellisch-Laidlaw cross\n"
73  << "section is used for neutrons (G4NeutronInelasticCrossSection).\n";
74 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
void G4CrossSectionPairGG::DumpPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 200 of file G4CrossSectionPairGG.cc.

References G4cout, G4endl, G4VCrossSectionDataSet::GetName(), and python.hepunit::GeV.

200  {
201  G4cout << std::setw(24) << " " << " G4CrossSectionPairGG: "
202  << theLowX->GetName() << " cross sections " << G4endl;
203  G4cout << std::setw(27) << " " << "below " << ETransition / GeV
204  << " GeV, Glauber-Gribov above " << G4endl;
205 }
const G4String & GetName() const
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
G4double G4CrossSectionPairGG::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 88 of file G4CrossSectionPairGG.cc.

References BuildPhysicsTable(), G4cout, G4endl, G4lrint(), G4NistManager::GetAtomicMassAmu(), G4DynamicParticle::GetDefinition(), G4VCrossSectionDataSet::GetElementCrossSection(), G4GlauberGribovCrossSection::GetInelasticGlauberGribov(), G4DynamicParticle::GetKineticEnergy(), and G4VCrossSectionDataSet::verboseLevel.

90 {
91  G4double Xsec(0.);
92 
93  if (aParticle->GetKineticEnergy() < ETransition)
94  {
95  Xsec = theLowX->GetElementCrossSection(aParticle, ZZ, mat);
96  } else {
97 
98  std::vector<ParticleXScale>::iterator iter = scale_factors.begin();
99  G4ParticleDefinition * pDef = aParticle->GetDefinition();
100  while (iter != scale_factors.end() && (*iter).first != pDef)
101  {
102  ++iter;
103  }
104  if (iter != scale_factors.end() )
105  {
106  G4int AA = G4lrint(NistMan->GetAtomicMassAmu(ZZ));
107  Xsec = theHighX->GetInelasticGlauberGribov(aParticle, ZZ, AA)
108  * (*iter).second[ZZ];
109  if (verboseLevel > 2)
110  {
111  G4cout << " scaling .." << ZZ << " " << AA << " "
112  << (*iter).second[ZZ] << " "
113  << theHighX->GetInelasticGlauberGribov(aParticle, ZZ, AA)
114  << " " << Xsec << G4endl;
115  }
116  } else {
117  // BuildPhysicsTable not done for pDef=aParticle->GetDefinition
118  // build table, and recurse
119  BuildPhysicsTable(*pDef);
120  Xsec=GetElementCrossSection(aParticle, ZZ, mat);
121  }
122  }
123 
124  return Xsec;
125 }
G4double GetKineticEnergy() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetAtomicMassAmu(const G4String &symb) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
G4bool G4CrossSectionPairGG::IsElementApplicable ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 76 of file G4CrossSectionPairGG.cc.

References G4DynamicParticle::GetKineticEnergy(), and G4VCrossSectionDataSet::IsElementApplicable().

77  {
78  G4bool isApplicable(false);
79  G4double Ekin = aParticle->GetKineticEnergy();
80  if (Ekin <= ETransition) {
81  isApplicable = theLowX->IsElementApplicable(aParticle, Z, mat);
82  } else if (Z > 1) {
83  isApplicable = true;
84  }
85  return isApplicable;
86 }
G4double GetKineticEnergy() const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
bool G4bool
Definition: G4Types.hh:79
double G4double
Definition: G4Types.hh:76

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