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

#include <G4ChipsAntiBaryonElasticXS.hh>

Inheritance diagram for G4ChipsAntiBaryonElasticXS:
G4VCrossSectionDataSet

Public Member Functions

 G4ChipsAntiBaryonElasticXS ()
 
 ~G4ChipsAntiBaryonElasticXS ()
 
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)
 
G4double GetExchangeT (G4int tZ, G4int tN, G4int pPDG)
 
- 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 46 of file G4ChipsAntiBaryonElasticXS.hh.

Constructor & Destructor Documentation

G4ChipsAntiBaryonElasticXS::G4ChipsAntiBaryonElasticXS ( )

Definition at line 56 of file G4ChipsAntiBaryonElasticXS.cc.

56  :G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
57 {
58  lPMin=-8.; //Min tabulatedLogarithmMomentum(D)
59  lPMax= 8.; //Max tabulatedLogarithmMomentum(D)
60  dlnP=(lPMax-lPMin)/nLast;// LogStep inTable (D)
61  onlyCS=true;//Flag toCalculOnlyCS(not Si/Bi)(L)
62  lastSIG=0.; //Last calculated cross section (L)
63  lastLP=-10.;//LastLog(mom_of IncidentHadron)(L)
64  lastTM=0.; //Last t_maximum (L)
65  theSS=0.; //TheLastSqSlope of 1st difr.Max(L)
66  theS1=0.; //TheLastMantissa of 1st difrMax(L)
67  theB1=0.; //TheLastSlope of 1st difructMax(L)
68  theS2=0.; //TheLastMantissa of 2nd difrMax(L)
69  theB2=0.; //TheLastSlope of 2nd difructMax(L)
70  theS3=0.; //TheLastMantissa of 3d difr.Max(L)
71  theB3=0.; //TheLastSlope of 3d difruct.Max(L)
72  theS4=0.; //TheLastMantissa of 4th difrMax(L)
73  theB4=0.; //TheLastSlope of 4th difructMax(L)
74  lastTZ=0; // Last atomic number of the target
75  lastTN=0; // Last # of neutrons in the target
76  lastPIN=0.; // Last initialized max momentum
77  lastCST=0; // Elastic cross-section table
78  lastPAR=0; // ParametersForFunctionCalculation
79  lastSST=0; // E-dep ofSqardSlope of 1st difMax
80  lastS1T=0; // E-dep of mantissa of 1st dif.Max
81  lastB1T=0; // E-dep of the slope of 1st difMax
82  lastS2T=0; // E-dep of mantissa of 2nd difrMax
83  lastB2T=0; // E-dep of the slope of 2nd difMax
84  lastS3T=0; // E-dep of mantissa of 3d difr.Max
85  lastB3T=0; // E-dep of the slope of 3d difrMax
86  lastS4T=0; // E-dep of mantissa of 4th difrMax
87  lastB4T=0; // E-dep of the slope of 4th difMax
88  lastN=0; // The last N of calculated nucleus
89  lastZ=0; // The last Z of calculated nucleus
90  lastP=0.; // LastUsed inCrossSection Momentum
91  lastTH=0.; // Last threshold momentum
92  lastCS=0.; // Last value of the Cross Section
93  lastI=0; // The last position in the DAMDB
94 }
G4VCrossSectionDataSet(const G4String &nam="")
G4ChipsAntiBaryonElasticXS::~G4ChipsAntiBaryonElasticXS ( )

Definition at line 96 of file G4ChipsAntiBaryonElasticXS.cc.

97 {
98  std::vector<G4double*>::iterator pos;
99  for (pos=CST.begin(); pos<CST.end(); pos++)
100  { delete [] *pos; }
101  CST.clear();
102  for (pos=PAR.begin(); pos<PAR.end(); pos++)
103  { delete [] *pos; }
104  PAR.clear();
105  for (pos=SST.begin(); pos<SST.end(); pos++)
106  { delete [] *pos; }
107  SST.clear();
108  for (pos=S1T.begin(); pos<S1T.end(); pos++)
109  { delete [] *pos; }
110  S1T.clear();
111  for (pos=B1T.begin(); pos<B1T.end(); pos++)
112  { delete [] *pos; }
113  B1T.clear();
114  for (pos=S2T.begin(); pos<S2T.end(); pos++)
115  { delete [] *pos; }
116  S2T.clear();
117  for (pos=B2T.begin(); pos<B2T.end(); pos++)
118  { delete [] *pos; }
119  B2T.clear();
120  for (pos=S3T.begin(); pos<S3T.end(); pos++)
121  { delete [] *pos; }
122  S3T.clear();
123  for (pos=B3T.begin(); pos<B3T.end(); pos++)
124  { delete [] *pos; }
125  B3T.clear();
126  for (pos=S4T.begin(); pos<S4T.end(); pos++)
127  { delete [] *pos; }
128  S4T.clear();
129  for (pos=B4T.begin(); pos<B4T.end(); pos++)
130  { delete [] *pos; }
131  B4T.clear();
132 }

Member Function Documentation

static const char* G4ChipsAntiBaryonElasticXS::Default_Name ( )
inlinestatic

Definition at line 55 of file G4ChipsAntiBaryonElasticXS.hh.

Referenced by G4ChipsComponentXS::G4ChipsComponentXS(), and G4ChipsElasticModel::G4ChipsElasticModel().

55 {return "ChipsAntiBaryonElasticXS";}
G4double G4ChipsAntiBaryonElasticXS::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 194 of file G4ChipsAntiBaryonElasticXS.cc.

References G4ThreadLocal, and python.hepunit::millibarn.

Referenced by GetIsoCrossSection(), and G4ChipsElasticModel::SampleInvariantT().

195 {
196  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)
197  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)
198  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
199  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
200  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
201  // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
202 
203  G4bool fCS = false;
204 
205  G4double pEn=pMom;
206  onlyCS=fCS;
207 
208  G4bool in=false; // By default the isotope must be found in the AMDB
209  lastP = 0.; // New momentum history (nothing to compare with)
210  lastN = tgN; // The last N of the calculated nucleus
211  lastZ = tgZ; // The last Z of the calculated nucleus
212  lastI = colN.size(); // Size of the Associative Memory DB in the heap
213  if(lastI) for(G4int i=0; i<lastI; i++) // Loop over proj/tgZ/tgN lines of DB
214  { // The nucleus with projPDG is found in AMDB
215  if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
216  {
217  lastI=i;
218  lastTH =colTH[i]; // Last THreshold (A-dependent)
219  if(pEn<=lastTH)
220  {
221  return 0.; // Energy is below the Threshold value
222  }
223  lastP =colP [i]; // Last Momentum (A-dependent)
224  lastCS =colCS[i]; // Last CrossSect (A-dependent)
225  // if(std::fabs(lastP/pMom-1.)<tolerance) //VI (do not use tolerance)
226  if(lastP == pMom) // Do not recalculate
227  {
228  CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // Update param's only
229  return lastCS*millibarn; // Use theLastCS
230  }
231  in = true; // This is the case when the isotop is found in DB
232  // Momentum pMom is in IU ! @@ Units
233  lastCS=CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // read & update
234  if(lastCS<=0. && pEn>lastTH) // Correct the threshold
235  {
236  lastTH=pEn;
237  }
238  break; // Go out of the LOOP with found lastI
239  }
240  } // End of attampt to find the nucleus in DB
241  if(!in) // This nucleus has not been calculated previously
242  {
243  //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
244  lastCS=CalculateCrossSection(fCS,0,lastI,pPDG,lastZ,lastN,pMom);//calculate&create
245  if(lastCS<=0.)
246  {
247  lastTH = 0; // ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
248  if(pEn>lastTH)
249  {
250  lastTH=pEn;
251  }
252  }
253  colN.push_back(tgN);
254  colZ.push_back(tgZ);
255  colP.push_back(pMom);
256  colTH.push_back(lastTH);
257  colCS.push_back(lastCS);
258  return lastCS*millibarn;
259  } // End of creation of the new set of parameters
260  else
261  {
262  colP[lastI]=pMom;
263  colCS[lastI]=lastCS;
264  }
265  return lastCS*millibarn;
266 }
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
bool G4bool
Definition: G4Types.hh:79
double G4double
Definition: G4Types.hh:76
G4double G4ChipsAntiBaryonElasticXS::GetExchangeT ( G4int  tZ,
G4int  tN,
G4int  pPDG 
)

Definition at line 641 of file G4ChipsAntiBaryonElasticXS.cc.

References test::a, G4cout, G4endl, G4UniformRand, and python.hepunit::gigaelectronvolt.

Referenced by G4ChipsElasticModel::SampleInvariantT().

642 {
643  static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
644  static const G4double third=1./3.;
645  static const G4double fifth=1./5.;
646  static const G4double sevth=1./7.;
647 
648  if(PDG<-3334 || PDG>-1111)G4cout<<"*Warning*G4QAntiBaryonElCS::GetExT:PDG="<<PDG<<G4endl;
649  if(onlyCS)G4cout<<"WarningG4ChipsAntiBaryonElasticXS::GetExchanT:onlyCS=1"<<G4endl;
650  if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
651  G4double q2=0.;
652  if(tgZ==1 && tgN==0) // ===> p+p=p+p
653  {
654  G4double E1=lastTM*theB1;
655  G4double R1=(1.-std::exp(-E1));
656  G4double E2=lastTM*theB2;
657  G4double R2=(1.-std::exp(-E2*E2*E2));
658  G4double E3=lastTM*theB3;
659  G4double R3=(1.-std::exp(-E3));
660  G4double I1=R1*theS1/theB1;
661  G4double I2=R2*theS2;
662  G4double I3=R3*theS3;
663  G4double I12=I1+I2;
664  G4double rand=(I12+I3)*G4UniformRand();
665  if (rand<I1 )
666  {
667  G4double ran=R1*G4UniformRand();
668  if(ran>1.) ran=1.;
669  q2=-std::log(1.-ran)/theB1;
670  }
671  else if(rand<I12)
672  {
673  G4double ran=R2*G4UniformRand();
674  if(ran>1.) ran=1.;
675  q2=-std::log(1.-ran);
676  if(q2<0.) q2=0.;
677  q2=std::pow(q2,third)/theB2;
678  }
679  else
680  {
681  G4double ran=R3*G4UniformRand();
682  if(ran>1.) ran=1.;
683  q2=-std::log(1.-ran)/theB3;
684  }
685  }
686  else
687  {
688  G4double a=tgZ+tgN;
689  G4double E1=lastTM*(theB1+lastTM*theSS);
690  G4double R1=(1.-std::exp(-E1));
691  G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
692  G4double tm2=lastTM*lastTM;
693  G4double E2=lastTM*tm2*theB2; // power 3 for lowA, 5 for HighA (1st)
694  if(a>6.5)E2*=tm2; // for heavy nuclei
695  G4double R2=(1.-std::exp(-E2));
696  G4double E3=lastTM*theB3;
697  if(a>6.5)E3*=tm2*tm2*tm2; // power 1 for lowA, 7 (2nd) for HighA
698  G4double R3=(1.-std::exp(-E3));
699  G4double E4=lastTM*theB4;
700  G4double R4=(1.-std::exp(-E4));
701  G4double I1=R1*theS1;
702  G4double I2=R2*theS2;
703  G4double I3=R3*theS3;
704  G4double I4=R4*theS4;
705  G4double I12=I1+I2;
706  G4double I13=I12+I3;
707  G4double rand=(I13+I4)*G4UniformRand();
708  if(rand<I1)
709  {
710  G4double ran=R1*G4UniformRand();
711  if(ran>1.) ran=1.;
712  q2=-std::log(1.-ran)/theB1;
713  if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
714  }
715  else if(rand<I12)
716  {
717  G4double ran=R2*G4UniformRand();
718  if(ran>1.) ran=1.;
719  q2=-std::log(1.-ran)/theB2;
720  if(q2<0.) q2=0.;
721  if(a<6.5) q2=std::pow(q2,third);
722  else q2=std::pow(q2,fifth);
723  }
724  else if(rand<I13)
725  {
726  G4double ran=R3*G4UniformRand();
727  if(ran>1.) ran=1.;
728  q2=-std::log(1.-ran)/theB3;
729  if(q2<0.) q2=0.;
730  if(a>6.5) q2=std::pow(q2,sevth);
731  }
732  else
733  {
734  G4double ran=R4*G4UniformRand();
735  if(ran>1.) ran=1.;
736  q2=-std::log(1.-ran)/theB4;
737  if(a<6.5) q2=lastTM-q2; // u reduced for lightA (starts from 0)
738  }
739  }
740  if(q2<0.) q2=0.;
741  if(!(q2>=-1.||q2<=1.))G4cout<<"*NAN*G4QaBElasticCrossSect::GetExchangeT:-t="<<q2<<G4endl;
742  if(q2>lastTM)
743  {
744  q2=lastTM;
745  }
746  return q2*GeVSQ;
747 }
int gigaelectronvolt
Definition: hepunit.py:110
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4ChipsAntiBaryonElasticXS::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 182 of file G4ChipsAntiBaryonElasticXS.cc.

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

186 {
187  G4double pMom=Pt->GetTotalMomentum();
188  G4int tgN = A - tgZ;
189  G4int pdg = Pt->GetDefinition()->GetPDGEncoding();
190 
191  return GetChipsCrossSection(pMom, tgZ, tgN, pdg);
192 }
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 G4ChipsAntiBaryonElasticXS::IsIsoApplicable ( const G4DynamicParticle Pt,
G4int  Z,
G4int  A,
const G4Element elm,
const G4Material mat 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 135 of file G4ChipsAntiBaryonElasticXS.cc.

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

138 {
139  G4ParticleDefinition* particle = Pt->GetDefinition();
140 
141  if(particle == G4AntiNeutron::AntiNeutron())
142  {
143  return true;
144  }
145  else if(particle == G4AntiProton::AntiProton())
146  {
147  return true;
148  }
149  else if(particle == G4AntiLambda::AntiLambda())
150  {
151  return true;
152  }
153  else if(particle == G4AntiSigmaPlus::AntiSigmaPlus())
154  {
155  return true;
156  }
157  else if(particle == G4AntiSigmaMinus::AntiSigmaMinus())
158  {
159  return true;
160  }
161  else if(particle == G4AntiSigmaZero::AntiSigmaZero())
162  {
163  return true;
164  }
165  else if(particle == G4AntiXiMinus::AntiXiMinus())
166  {
167  return true;
168  }
169  else if(particle == G4AntiXiZero::AntiXiZero())
170  {
171  return true;
172  }
173  else if(particle == G4AntiOmegaMinus::AntiOmegaMinus())
174  {
175  return true;
176  }
177  return false;
178 }
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: