Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChipsAntiBaryonElasticXS Class Reference

#include <G4ChipsAntiBaryonElasticXS.hh>

+ Inheritance diagram for G4ChipsAntiBaryonElasticXS:

Public Member Functions

 G4ChipsAntiBaryonElasticXS ()
 
 ~G4ChipsAntiBaryonElasticXS ()
 
virtual void CrossSectionDescription (std::ostream &) const
 
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=nullptr)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
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)
 
bool ForAllAtomsAndEnergies () const
 
void SetForAllAtomsAndEnergies (G4bool val)
 
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::G4ChipsAntiBaryonElasticXS ( )

Definition at line 58 of file G4ChipsAntiBaryonElasticXS.cc.

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

◆ ~G4ChipsAntiBaryonElasticXS()

G4ChipsAntiBaryonElasticXS::~G4ChipsAntiBaryonElasticXS ( )

Definition at line 98 of file G4ChipsAntiBaryonElasticXS.cc.

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

Member Function Documentation

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 137 of file G4ChipsAntiBaryonElasticXS.cc.

138{
139 outFile << "G4ChipsAntiBaryonElasticXS provides the elastic cross\n"
140 << "section for anti-baryon nucleus scattering as a function of incident\n"
141 << "momentum. The cross section is calculated using M. Kossov's\n"
142 << "CHIPS parameterization of cross section data.\n";
143}

◆ Default_Name()

static const char * G4ChipsAntiBaryonElasticXS::Default_Name ( )
inlinestatic

Definition at line 55 of file G4ChipsAntiBaryonElasticXS.hh.

55{return "ChipsAntiBaryonElasticXS";}

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

◆ GetChipsCrossSection()

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 205 of file G4ChipsAntiBaryonElasticXS.cc.

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

Referenced by G4ChipsComponentXS::GetElasticElementCrossSection(), GetIsoCrossSection(), G4ChipsComponentXS::GetTotalElementCrossSection(), and G4ChipsElasticModel::SampleInvariantT().

◆ GetExchangeT()

G4double G4ChipsAntiBaryonElasticXS::GetExchangeT ( G4int  tZ,
G4int  tN,
G4int  pPDG 
)

Definition at line 643 of file G4ChipsAntiBaryonElasticXS.cc.

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

Referenced by G4ChipsElasticModel::SampleInvariantT().

◆ GetIsoCrossSection()

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 193 of file G4ChipsAntiBaryonElasticXS.cc.

197{
198 G4double pMom=Pt->GetTotalMomentum();
199 G4int tgN = A - tgZ;
200 G4int pdg = Pt->GetDefinition()->GetPDGEncoding();
201
202 return GetChipsCrossSection(pMom, tgZ, tgN, pdg);
203}
double A(double temperature)
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4ParticleDefinition * GetDefinition() const
G4double GetTotalMomentum() const

◆ IsIsoApplicable()

G4bool G4ChipsAntiBaryonElasticXS::IsIsoApplicable ( const G4DynamicParticle Pt,
G4int  Z,
G4int  A,
const G4Element elm,
const G4Material mat 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 145 of file G4ChipsAntiBaryonElasticXS.cc.

148{
149
150 /*
151 if(particle == G4AntiNeutron::AntiNeutron())
152 {
153 return true;
154 }
155 else if(particle == G4AntiProton::AntiProton())
156 {
157 return true;
158 }
159 else if(particle == G4AntiLambda::AntiLambda())
160 {
161 return true;
162 }
163 else if(particle == G4AntiSigmaPlus::AntiSigmaPlus())
164 {
165 return true;
166 }
167 else if(particle == G4AntiSigmaMinus::AntiSigmaMinus())
168 {
169 return true;
170 }
171 else if(particle == G4AntiSigmaZero::AntiSigmaZero())
172 {
173 return true;
174 }
175 else if(particle == G4AntiXiMinus::AntiXiMinus())
176 {
177 return true;
178 }
179 else if(particle == G4AntiXiZero::AntiXiZero())
180 {
181 return true;
182 }
183 else if(particle == G4AntiOmegaMinus::AntiOmegaMinus())
184 {
185 return true;
186 }
187 */
188 return true;
189}

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