Geant4 11.2.2
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 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)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, 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 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
 
void SetName (const G4String &nam)
 
G4VCrossSectionDataSetoperator= (const G4VCrossSectionDataSet &right)=delete
 
 G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete
 

Static Public Member Functions

static const char * Default_Name ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 
G4String name
 

Detailed Description

Definition at line 46 of file G4ChipsHyperonElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsHyperonElasticXS()

G4ChipsHyperonElasticXS::G4ChipsHyperonElasticXS ( )

Definition at line 63 of file G4ChipsHyperonElasticXS.cc.

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

◆ ~G4ChipsHyperonElasticXS()

G4ChipsHyperonElasticXS::~G4ChipsHyperonElasticXS ( )

Definition at line 103 of file G4ChipsHyperonElasticXS.cc.

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

Member Function Documentation

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 142 of file G4ChipsHyperonElasticXS.cc.

143{
144 outFile << "G4ChipsHyperonElasticXS provides the elastic cross\n"
145 << "section for hyperon nucleus scattering as a function of incident\n"
146 << "momentum. The cross section is calculated using M. Kossov's\n"
147 << "CHIPS parameterization of cross section data.\n";
148}

◆ Default_Name()

static const char * G4ChipsHyperonElasticXS::Default_Name ( )
inlinestatic

Definition at line 54 of file G4ChipsHyperonElasticXS.hh.

54{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 171 of file G4ChipsHyperonElasticXS.cc.

172{
173
174 G4bool fCS = false;
175 G4double pEn=pMom;
176
177 onlyCS=fCS;
178
179 G4bool in=false; // By default the isotope must be found in the AMDB
180 lastP = 0.; // New momentum history (nothing to compare with)
181 lastN = tgN; // The last N of the calculated nucleus
182 lastZ = tgZ; // The last Z of the calculated nucleus
183 lastI = (G4int)colN.size(); // Size of the Associative Memory DB in the heap
184 if(lastI) for(G4int i=0; i<lastI; ++i) // Loop over proj/tgZ/tgN lines of DB
185 { // The nucleus with projPDG is found in AMDB
186 if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
187 {
188 lastI=i;
189 lastTH =colTH[i]; // Last THreshold (A-dependent)
190 if(pEn<=lastTH)
191 {
192 return 0.; // Energy is below the Threshold value
193 }
194 lastP =colP [i]; // Last Momentum (A-dependent)
195 lastCS =colCS[i]; // Last CrossSect (A-dependent)
196 // if(std::fabs(lastP/pMom-1.)<tolerance) //VI (do not use tolerance)
197 if(lastP == pMom) // Do not recalculate
198 {
199 CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // Update param's only
200 return lastCS*millibarn; // Use theLastCS
201 }
202 in = true; // This is the case when the isotop is found in DB
203 // Momentum pMom is in IU ! @@ Units
204 lastCS=CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // read & update
205 if(lastCS<=0. && pEn>lastTH) // Correct the threshold
206 {
207 lastTH=pEn;
208 }
209 break; // Go out of the LOOP with found lastI
210 }
211 } // End of attampt to find the nucleus in DB
212 if(!in) // This nucleus has not been calculated previously
213 {
214 //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
215 lastCS=CalculateCrossSection(fCS,0,lastI,pPDG,lastZ,lastN,pMom);//calculate&create
216 if(lastCS<=0.)
217 {
218 lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
219 if(pEn>lastTH)
220 {
221 lastTH=pEn;
222 }
223 }
224 colN.push_back(tgN);
225 colZ.push_back(tgZ);
226 colP.push_back(pMom);
227 colTH.push_back(lastTH);
228 colCS.push_back(lastCS);
229 return lastCS*millibarn;
230 } // End of creation of the new set of parameters
231 else
232 {
233 colP[lastI]=pMom;
234 colCS[lastI]=lastCS;
235 }
236 return lastCS*millibarn;
237}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85

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

◆ GetExchangeT()

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

Definition at line 609 of file G4ChipsHyperonElasticXS.cc.

610{
611 static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
612 static const G4double third=1./3.;
613 static const G4double fifth=1./5.;
614 static const G4double sevth=1./7.;
615 //AR-04Jun2014 if(PDG==3222 || PDG<3000 || PDG>3334)G4cout<<"*Warning*G4QHyElCS::GET:PDG="<<PDG<<G4endl;
616 if(PDG<3000 || PDG>3334)G4cout<<"*Warning*G4QHyElCS::GET:PDG="<<PDG<<G4endl;
617 if(onlyCS)G4cout<<"*Warning*G4ChipsHyperonElasticXS::GetExchanT: onlyCS=1"<<G4endl;
618 if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
619 G4double q2=0.;
620 if(tgZ==1 && tgN==0) // ===> p+p=p+p
621 {
622 G4double E1=lastTM*theB1;
623 G4double R1=(1.-G4Exp(-E1));
624 G4double E2=lastTM*theB2;
625 G4double R2=(1.-G4Exp(-E2*E2*E2));
626 G4double E3=lastTM*theB3;
627 G4double R3=(1.-G4Exp(-E3));
628 G4double I1=R1*theS1/theB1;
629 G4double I2=R2*theS2;
630 G4double I3=R3*theS3;
631 G4double I12=I1+I2;
632 G4double rand=(I12+I3)*G4UniformRand();
633 if (rand<I1 )
634 {
635 G4double ran=R1*G4UniformRand();
636 if(ran>1.) ran=1.;
637 q2=-G4Log(1.-ran)/theB1;
638 }
639 else if(rand<I12)
640 {
641 G4double ran=R2*G4UniformRand();
642 if(ran>1.) ran=1.;
643 q2=-G4Log(1.-ran);
644 if(q2<0.) q2=0.;
645 q2=G4Pow::GetInstance()->powA(q2,third)/theB2;
646 }
647 else
648 {
649 G4double ran=R3*G4UniformRand();
650 if(ran>1.) ran=1.;
651 q2=-G4Log(1.-ran)/theB3;
652 }
653 }
654 else
655 {
656 G4double a=tgZ+tgN;
657 G4double E1=lastTM*(theB1+lastTM*theSS);
658 G4double R1=(1.-G4Exp(-E1));
659 G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
660 G4double tm2=lastTM*lastTM;
661 G4double E2=lastTM*tm2*theB2; // power 3 for lowA, 5 for HighA (1st)
662 if(a>6.5)E2*=tm2; // for heavy nuclei
663 G4double R2=(1.-G4Exp(-E2));
664 G4double E3=lastTM*theB3;
665 if(a>6.5)E3*=tm2*tm2*tm2; // power 1 for lowA, 7 (2nd) for HighA
666 G4double R3=(1.-G4Exp(-E3));
667 G4double E4=lastTM*theB4;
668 G4double R4=(1.-G4Exp(-E4));
669 G4double I1=R1*theS1;
670 G4double I2=R2*theS2;
671 G4double I3=R3*theS3;
672 G4double I4=R4*theS4;
673 G4double I12=I1+I2;
674 G4double I13=I12+I3;
675 G4double rand=(I13+I4)*G4UniformRand();
676 if(rand<I1)
677 {
678 G4double ran=R1*G4UniformRand();
679 if(ran>1.) ran=1.;
680 q2=-G4Log(1.-ran)/theB1;
681 if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
682 }
683 else if(rand<I12)
684 {
685 G4double ran=R2*G4UniformRand();
686 if(ran>1.) ran=1.;
687 q2=-G4Log(1.-ran)/theB2;
688 if(q2<0.) q2=0.;
689 if(a<6.5) q2=G4Pow::GetInstance()->powA(q2,third);
690 else q2=G4Pow::GetInstance()->powA(q2,fifth);
691 }
692 else if(rand<I13)
693 {
694 G4double ran=R3*G4UniformRand();
695 if(ran>1.) ran=1.;
696 q2=-G4Log(1.-ran)/theB3;
697 if(q2<0.) q2=0.;
698 if(a>6.5) q2=G4Pow::GetInstance()->powA(q2,sevth);
699 }
700 else
701 {
702 G4double ran=R4*G4UniformRand();
703 if(ran>1.) ran=1.;
704 q2=-G4Log(1.-ran)/theB4;
705 if(a<6.5) q2=lastTM-q2; // u reduced for lightA (starts from 0)
706 }
707 }
708 if(q2<0.) q2=0.;
709 if(!(q2>=-1.||q2<=1.))G4cout<<"*NAN*G4QHyElasticCrossSect::GetExchangeT:-t="<<q2<<G4endl;
710 if(q2>lastTM)
711 {
712 q2=lastTM;
713 }
714 return q2*GeVSQ;
715}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double G4Log(G4double x)
Definition G4Log.hh:227
#define G4endl
Definition G4ios.hh:67
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

◆ 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 159 of file G4ChipsHyperonElasticXS.cc.

163{
164 G4double pMom=Pt->GetTotalMomentum();
165 G4int tgN = A - tgZ;
166 G4int pdg = Pt->GetDefinition()->GetPDGEncoding();
167
168 return GetChipsCrossSection(pMom, tgZ, tgN, pdg);
169}
const G4double A[17]
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 150 of file G4ChipsHyperonElasticXS.cc.

153{
154 return true;
155}

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