Geant4 10.7.0
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)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, 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 GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, 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 CrossSectionDescription (std::ostream &) const
 
virtual G4int GetVerboseLevel () const
 
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
 

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 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}
static const char * Default_Name()

◆ ~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 165 of file G4ChipsPionPlusElasticXS.cc.

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

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

158{
159 G4double pMom=Pt->GetTotalMomentum();
160 G4int tgN = A - tgZ;
161
162 return GetChipsCrossSection(pMom, tgZ, tgN, 211);
163}
double A(double temperature)
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: