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

#include <G4ChipsPionMinusElasticXS.hh>

+ Inheritance diagram for G4ChipsPionMinusElasticXS:

Public Member Functions

 G4ChipsPionMinusElasticXS ()
 
 ~G4ChipsPionMinusElasticXS ()
 
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 47 of file G4ChipsPionMinusElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsPionMinusElasticXS()

G4ChipsPionMinusElasticXS::G4ChipsPionMinusElasticXS ( )

Definition at line 58 of file G4ChipsPionMinusElasticXS.cc.

58 :G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
59{
60 lPMin=-8.; //Min tabulated logarithmMomentum(D)
61 lPMax= 8.; //Max tabulated logarithmMomentum(D)
62 dlnP=(lPMax-lPMin)/nLast;//LogStep inTheTable(D)
63 onlyCS=true;//Flag toCalcul OnlyCS(not Si/Bi)(L)
64 lastSIG=0.; //Last calculated cross section (L)
65 lastLP=-10.;//Last log(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 difr.Max(L)
69 theB1=0.; //TheLastSlope of 1st difruct.Max(L)
70 theS2=0.; //TheLastMantissa of 2nd difr.Max(L)
71 theB2=0.; //TheLastSlope of 2nd difruct.Max(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 difr.Max(L)
75 theB4=0.; //TheLastSlope of 4th difruct.Max(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; // Parameters ForFunctionCalculation
81 lastSST=0; // E-dep of SqardSlope 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 in CrossSection 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}
G4VCrossSectionDataSet(const G4String &nam="")

◆ ~G4ChipsPionMinusElasticXS()

G4ChipsPionMinusElasticXS::~G4ChipsPionMinusElasticXS ( )

Definition at line 98 of file G4ChipsPionMinusElasticXS.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 G4ChipsPionMinusElasticXS::CrossSectionDescription ( std::ostream & outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 137 of file G4ChipsPionMinusElasticXS.cc.

138{
139 outFile << "G4ChipsPionMinusElasticXS provides the elastic cross\n"
140 << "section for pion- 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 * G4ChipsPionMinusElasticXS::Default_Name ( )
inlinestatic

Definition at line 55 of file G4ChipsPionMinusElasticXS.hh.

55{return "ChipsPionMinusElasticXS";}

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

◆ GetChipsCrossSection()

G4double G4ChipsPionMinusElasticXS::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 166 of file G4ChipsPionMinusElasticXS.cc.

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

Definition at line 603 of file G4ChipsPionMinusElasticXS.cc.

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

Referenced by G4ChipsElasticModel::SampleInvariantT().

◆ GetIsoCrossSection()

G4double G4ChipsPionMinusElasticXS::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 155 of file G4ChipsPionMinusElasticXS.cc.

159{
160 G4double pMom=Pt->GetTotalMomentum();
161 G4int tgN = A - tgZ;
162
163 return GetChipsCrossSection(pMom, tgZ, tgN, -212);
164}
const G4double A[17]
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4double GetTotalMomentum() const

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 146 of file G4ChipsPionMinusElasticXS.cc.

149{
150 return true;
151}

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