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

#include <G4ChipsAntiBaryonInelasticXS.hh>

Inheritance diagram for G4ChipsAntiBaryonInelasticXS:
G4VCrossSectionDataSet

Public Member Functions

 G4ChipsAntiBaryonInelasticXS ()
 
 ~G4ChipsAntiBaryonInelasticXS ()
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4double GetChipsCrossSection (G4double momentum, G4int Z, G4int N, G4int pdg)
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, 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 GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
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
 

Static Public Member Functions

static const char * Default_Name ()
 

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 44 of file G4ChipsAntiBaryonInelasticXS.hh.

Constructor & Destructor Documentation

G4ChipsAntiBaryonInelasticXS::G4ChipsAntiBaryonInelasticXS ( )

Definition at line 59 of file G4ChipsAntiBaryonInelasticXS.cc.

60 {
61  lastLEN=0; // Pointer to lastArray of LowEn CS
62  lastHEN=0; // Pointer to lastArray of HighEn CS
63  lastN=0; // The last N of calculated nucleus
64  lastZ=0; // The last Z of calculated nucleus
65  lastP=0.; // Last used Cross Section Momentum
66  lastTH=0.; // Last threshold momentum
67  lastCS=0.; // Last value of the Cross Section
68  lastI=0; // The last position in the DAMDB
69  LEN = new std::vector<G4double*>;
70  HEN = new std::vector<G4double*>;
71 }
G4VCrossSectionDataSet(const G4String &nam="")
Definition: inflate.h:41
G4ChipsAntiBaryonInelasticXS::~G4ChipsAntiBaryonInelasticXS ( )

Definition at line 73 of file G4ChipsAntiBaryonInelasticXS.cc.

74 {
75  G4int lens=LEN->size();
76  for(G4int i=0; i<lens; ++i) delete[] (*LEN)[i];
77  delete LEN;
78  G4int hens=HEN->size();
79  for(G4int i=0; i<hens; ++i) delete[] (*HEN)[i];
80  delete HEN;
81 }
int G4int
Definition: G4Types.hh:78
Definition: inflate.h:41

Member Function Documentation

static const char* G4ChipsAntiBaryonInelasticXS::Default_Name ( )
inlinestatic

Definition at line 52 of file G4ChipsAntiBaryonInelasticXS.hh.

Referenced by G4ChipsComponentXS::G4ChipsComponentXS().

52 {return "ChipsAntiBaryonInelasticXS";}
G4double G4ChipsAntiBaryonInelasticXS::GetChipsCrossSection ( G4double  momentum,
G4int  Z,
G4int  N,
G4int  pdg 
)
virtual

!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)

Definition at line 142 of file G4ChipsAntiBaryonInelasticXS.cc.

References G4ThreadLocal, and python.hepunit::millibarn.

Referenced by GetIsoCrossSection().

143 {
144  static G4ThreadLocal G4int j; // A#0f Z/N-records already tested in AMDB
145  static G4ThreadLocal std::vector <G4int> *colN_G4MT_TLS_ = 0 ; if (!colN_G4MT_TLS_) colN_G4MT_TLS_ = new std::vector <G4int> ; std::vector <G4int> &colN = *colN_G4MT_TLS_; // Vector of N for calculated nuclei (isotops)
146  static G4ThreadLocal std::vector <G4int> *colZ_G4MT_TLS_ = 0 ; if (!colZ_G4MT_TLS_) colZ_G4MT_TLS_ = new std::vector <G4int> ; std::vector <G4int> &colZ = *colZ_G4MT_TLS_; // Vector of Z for calculated nuclei (isotops)
147  static G4ThreadLocal std::vector <G4double> *colP_G4MT_TLS_ = 0 ; if (!colP_G4MT_TLS_) colP_G4MT_TLS_ = new std::vector <G4double> ; std::vector <G4double> &colP = *colP_G4MT_TLS_; // Vector of last momenta for the reaction
148  static G4ThreadLocal std::vector <G4double> *colTH_G4MT_TLS_ = 0 ; if (!colTH_G4MT_TLS_) colTH_G4MT_TLS_ = new std::vector <G4double> ; std::vector <G4double> &colTH = *colTH_G4MT_TLS_; // Vector of energy thresholds for the reaction
149  static G4ThreadLocal std::vector <G4double> *colCS_G4MT_TLS_ = 0 ; if (!colCS_G4MT_TLS_) colCS_G4MT_TLS_ = new std::vector <G4double> ; std::vector <G4double> &colCS = *colCS_G4MT_TLS_; // Vector of last cross sections for the reaction
150  // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
151 
152  G4bool in=false; // By default the isotope must be found in the AMDB
153  if(tgN!=lastN || tgZ!=lastZ) // The nucleus was not the last used isotope
154  {
155  in = false; // By default the isotope haven't be found in AMDB
156  lastP = 0.; // New momentum history (nothing to compare with)
157  lastN = tgN; // The last N of the calculated nucleus
158  lastZ = tgZ; // The last Z of the calculated nucleus
159  lastI = colN.size(); // Size of the Associative Memory DB in the heap
160  j = 0; // A#0f records found in DB for this projectile
161  if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
162  {
163  if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
164  {
165  lastI=i; // Remember the index for future fast/last use
166  lastTH =colTH[i]; // The last THreshold (A-dependent)
167  if(pMom<=lastTH)
168  {
169  return 0.; // Energy is below the Threshold value
170  }
171  lastP =colP [i]; // Last Momentum (A-dependent)
172  lastCS =colCS[i]; // Last CrossSect (A-dependent)
173  in = true; // This is the case when the isotop is found in DB
174  // Momentum pMom is in IU ! @@ Units
175  lastCS=CalculateCrossSection(-1,j,cPDG,lastZ,lastN,pMom); // read & update
176  if(lastCS<=0. && pMom>lastTH) // Correct the threshold (@@ No intermediate Zeros)
177  {
178  lastCS=0.;
179  lastTH=pMom;
180  }
181  break; // Go out of the LOOP
182  }
183  j++; // Increment a#0f records found in DB
184  }
185  if(!in) // This isotope has not been calculated previously
186  {
187  //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
188  lastCS=CalculateCrossSection(0,j,cPDG,lastZ,lastN,pMom); //calculate & create
189  //if(lastCS>0.) // It means that the AMBD was initialized
190  //{
191 
192  lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
193  colN.push_back(tgN);
194  colZ.push_back(tgZ);
195  colP.push_back(pMom);
196  colTH.push_back(lastTH);
197  colCS.push_back(lastCS);
198  //} // M.K. Presence of H1 with high threshold breaks the syncronization
199  return lastCS*millibarn;
200  } // End of creation of the new set of parameters
201  else
202  {
203  colP[lastI]=pMom;
204  colCS[lastI]=lastCS;
205  }
206  } // End of parameters udate
207  else if(pMom<=lastTH)
208  {
209  return 0.; // Momentum is below the Threshold Value -> CS=0
210  }
211  else // It is the last used -> use the current tables
212  {
213  lastCS=CalculateCrossSection(1,j,cPDG,lastZ,lastN,pMom); // Only read and UpdateDB
214  lastP=pMom;
215  }
216  return lastCS*millibarn;
217 }
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
bool G4bool
Definition: G4Types.hh:79
G4double G4ChipsAntiBaryonInelasticXS::GetIsoCrossSection ( const G4DynamicParticle Pt,
G4int  tgZ,
G4int  A,
const G4Isotope iso = 0,
const G4Element elm = 0,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 130 of file G4ChipsAntiBaryonInelasticXS.cc.

References GetChipsCrossSection(), G4DynamicParticle::GetDefinition(), G4ParticleDefinition::GetPDGEncoding(), and G4DynamicParticle::GetTotalMomentum().

134 {
135  G4double pMom=Pt->GetTotalMomentum();
136  G4int tgN = A - tgZ;
137  G4int pdg = Pt->GetDefinition()->GetPDGEncoding();
138 
139  return GetChipsCrossSection(pMom, tgZ, tgN, pdg);
140 }
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double GetTotalMomentum() const
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
double G4double
Definition: G4Types.hh:76
G4bool G4ChipsAntiBaryonInelasticXS::IsIsoApplicable ( const G4DynamicParticle Pt,
G4int  Z,
G4int  A,
const G4Element elm,
const G4Material mat 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 83 of file G4ChipsAntiBaryonInelasticXS.cc.

References G4AntiLambda::AntiLambda(), G4AntiNeutron::AntiNeutron(), G4AntiOmegaMinus::AntiOmegaMinus(), G4AntiProton::AntiProton(), G4AntiSigmaMinus::AntiSigmaMinus(), G4AntiSigmaPlus::AntiSigmaPlus(), G4AntiSigmaZero::AntiSigmaZero(), G4AntiXiMinus::AntiXiMinus(), G4AntiXiZero::AntiXiZero(), and G4DynamicParticle::GetDefinition().

86 {
87  G4ParticleDefinition* particle = Pt->GetDefinition();
88 
89  if(particle == G4AntiNeutron::AntiNeutron())
90  {
91  return true;
92  }
93  else if(particle == G4AntiProton::AntiProton())
94  {
95  return true;
96  }
97  else if(particle == G4AntiLambda::AntiLambda())
98  {
99  return true;
100  }
101  else if(particle == G4AntiSigmaPlus::AntiSigmaPlus())
102  {
103  return true;
104  }
105  else if(particle == G4AntiSigmaMinus::AntiSigmaMinus())
106  {
107  return true;
108  }
109  else if(particle == G4AntiSigmaZero::AntiSigmaZero())
110  {
111  return true;
112  }
113  else if(particle == G4AntiXiMinus::AntiXiMinus())
114  {
115  return true;
116  }
117  else if(particle == G4AntiXiZero::AntiXiZero())
118  {
119  return true;
120  }
121  else if(particle == G4AntiOmegaMinus::AntiOmegaMinus())
122  {
123  return true;
124  }
125  return false;
126 }
static G4AntiOmegaMinus * AntiOmegaMinus()
G4ParticleDefinition * GetDefinition() const
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4AntiXiMinus * AntiXiMinus()
static G4AntiLambda * AntiLambda()
static G4AntiSigmaZero * AntiSigmaZero()
static G4AntiXiZero * AntiXiZero()
static G4AntiNeutron * AntiNeutron()

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