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

#include <G4ChipsHyperonElasticXS.hh>

+ Inheritance diagram for G4ChipsHyperonElasticXS:

Public Member Functions

 G4ChipsHyperonElasticXS ()
 
 ~G4ChipsHyperonElasticXS ()
 
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 G4ChipsHyperonElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsHyperonElasticXS()

G4ChipsHyperonElasticXS::G4ChipsHyperonElasticXS ( )

Definition at line 60 of file G4ChipsHyperonElasticXS.cc.

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

◆ ~G4ChipsHyperonElasticXS()

G4ChipsHyperonElasticXS::~G4ChipsHyperonElasticXS ( )

Definition at line 100 of file G4ChipsHyperonElasticXS.cc.

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

Member Function Documentation

◆ Default_Name()

static const char * G4ChipsHyperonElasticXS::Default_Name ( )
inlinestatic

Definition at line 55 of file G4ChipsHyperonElasticXS.hh.

55{return "ChipsHyperonElasticXS";}

Referenced by G4ChipsComponentXS::G4ChipsComponentXS().

◆ GetChipsCrossSection()

G4double G4ChipsHyperonElasticXS::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 188 of file G4ChipsHyperonElasticXS.cc.

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

Referenced by GetIsoCrossSection(), and G4ChipsComponentXS::GetTotalElementCrossSection().

◆ GetExchangeT()

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

Definition at line 634 of file G4ChipsHyperonElasticXS.cc.

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

◆ GetIsoCrossSection()

G4double G4ChipsHyperonElasticXS::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 176 of file G4ChipsHyperonElasticXS.cc.

180{
181 G4double pMom=Pt->GetTotalMomentum();
182 G4int tgN = A - tgZ;
183 G4int pdg = Pt->GetDefinition()->GetPDGEncoding();
184
185 return GetChipsCrossSection(pMom, tgZ, tgN, pdg);
186}
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4ParticleDefinition * GetDefinition() const
G4double GetTotalMomentum() const

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 138 of file G4ChipsHyperonElasticXS.cc.

141{
142 G4ParticleDefinition* particle = Pt->GetDefinition();
143 if (particle == G4Lambda::Lambda())
144 {
145 return true;
146 }
147 else if(particle == G4SigmaPlus::SigmaPlus())
148 {
149 return true;
150 }
151 else if(particle == G4SigmaMinus::SigmaMinus())
152 {
153 return true;
154 }
155 else if(particle == G4SigmaZero::SigmaZero())
156 {
157 return true;
158 }
159 else if(particle == G4XiMinus::XiMinus())
160 {
161 return true;
162 }
163 else if(particle == G4XiZero::XiZero())
164 {
165 return true;
166 }
167 else if(particle == G4OmegaMinus::OmegaMinus())
168 {
169 return true;
170 }
171 return false;
172}
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
static G4OmegaMinus * OmegaMinus()
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:99
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:106
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106

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