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

#include <G4ChipsPionPlusElasticXS.hh>

+ Inheritance diagram for G4ChipsPionPlusElasticXS:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4ChipsPionPlusElasticXS()

G4ChipsPionPlusElasticXS::G4ChipsPionPlusElasticXS ( )

Definition at line 57 of file G4ChipsPionPlusElasticXS.cc.

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

◆ ~G4ChipsPionPlusElasticXS()

G4ChipsPionPlusElasticXS::~G4ChipsPionPlusElasticXS ( )

Definition at line 97 of file G4ChipsPionPlusElasticXS.cc.

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

Member Function Documentation

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 136 of file G4ChipsPionPlusElasticXS.cc.

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

◆ Default_Name()

static const char * G4ChipsPionPlusElasticXS::Default_Name ( )
inlinestatic

Definition at line 54 of file G4ChipsPionPlusElasticXS.hh.

54{return "ChipsPionPlusElasticXS";}

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

◆ GetChipsCrossSection()

G4double G4ChipsPionPlusElasticXS::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 164 of file G4ChipsPionPlusElasticXS.cc.

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

Definition at line 600 of file G4ChipsPionPlusElasticXS.cc.

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

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

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 144 of file G4ChipsPionPlusElasticXS.cc.

147{
148 return true;
149}

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