Geant4 9.6.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 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)
 
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 GetElementCrossSection (const G4DynamicParticle *, G4int Z, 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 void CrossSectionDescription (std::ostream &) const
 
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 47 of file G4ChipsAntiBaryonElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsAntiBaryonElasticXS()

G4ChipsAntiBaryonElasticXS::G4ChipsAntiBaryonElasticXS ( )

Definition at line 55 of file G4ChipsAntiBaryonElasticXS.cc.

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

◆ ~G4ChipsAntiBaryonElasticXS()

G4ChipsAntiBaryonElasticXS::~G4ChipsAntiBaryonElasticXS ( )

Definition at line 95 of file G4ChipsAntiBaryonElasticXS.cc.

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

Member Function Documentation

◆ Default_Name()

static const char * G4ChipsAntiBaryonElasticXS::Default_Name ( )
inlinestatic

Definition at line 56 of file G4ChipsAntiBaryonElasticXS.hh.

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

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

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

◆ GetExchangeT()

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

Definition at line 640 of file G4ChipsAntiBaryonElasticXS.cc.

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

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

185{
186 G4double pMom=Pt->GetTotalMomentum();
187 G4int tgN = A - tgZ;
188 G4int pdg = Pt->GetDefinition()->GetPDGEncoding();
189
190 return GetChipsCrossSection(pMom, tgZ, tgN, pdg);
191}
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 134 of file G4ChipsAntiBaryonElasticXS.cc.

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

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