Geant4 9.6.0
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 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 48 of file G4ChipsPionMinusElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsPionMinusElasticXS()

G4ChipsPionMinusElasticXS::G4ChipsPionMinusElasticXS ( )

Definition at line 54 of file G4ChipsPionMinusElasticXS.cc.

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

◆ ~G4ChipsPionMinusElasticXS()

G4ChipsPionMinusElasticXS::~G4ChipsPionMinusElasticXS ( )

Definition at line 94 of file G4ChipsPionMinusElasticXS.cc.

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

Member Function Documentation

◆ Default_Name()

static const char * G4ChipsPionMinusElasticXS::Default_Name ( )
inlinestatic

Definition at line 56 of file G4ChipsPionMinusElasticXS.hh.

56{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 155 of file G4ChipsPionMinusElasticXS.cc.

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

Referenced by G4ChipsComponentXS::GetElasticElementCrossSection(), GetIsoCrossSection(), G4ChipsComponentXS::GetTotalElementCrossSection(), and G4ChipsElasticModel::SampleInvariantT().

◆ GetExchangeT()

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

Definition at line 601 of file G4ChipsPionMinusElasticXS.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*G4ChipsPionMinusElasticXS::GetExT:PDG="<<PDG<<G4endl;
608 if(onlyCS)G4cout<<"Warning*G4ChipsPionMinusElasticXS::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.-std::exp(-E1));
615 G4double E2=lastTM*theB2;
616 G4double R2=(1.-std::exp(-E2*E2*E2));
617 G4double E3=lastTM*theB3;
618 G4double R3=(1.-std::exp(-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=-std::log(1.-ran)/theB1;
629 }
630 else if(rand<I12)
631 {
632 G4double ran=R2*G4UniformRand();
633 if(ran>1.) ran=1.;
634 q2=-std::log(1.-ran);
635 if(q2<0.) q2=0.;
636 q2=std::pow(q2,third)/theB2;
637 }
638 else
639 {
640 G4double ran=R3*G4UniformRand();
641 if(ran>1.) ran=1.;
642 q2=-std::log(1.-ran)/theB3;
643 }
644 }
645 else
646 {
647 G4double a=tgZ+tgN;
648 G4double E1=lastTM*(theB1+lastTM*theSS);
649 G4double R1=(1.-std::exp(-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.-std::exp(-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.-std::exp(-E3));
658 G4double E4=lastTM*theB4;
659 G4double R4=(1.-std::exp(-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=-std::log(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=-std::log(1.-ran)/theB2;
679 if(q2<0.) q2=0.;
680 if(a<6.5) q2=std::pow(q2,third);
681 else q2=std::pow(q2,fifth);
682 }
683 else if(rand<I13)
684 {
685 G4double ran=R3*G4UniformRand();
686 if(ran>1.) ran=1.;
687 q2=-std::log(1.-ran)/theB3;
688 if(q2<0.) q2=0.;
689 if(a>6.5) q2=std::pow(q2,sevth);
690 }
691 else
692 {
693 G4double ran=R4*G4UniformRand();
694 if(ran>1.) ran=1.;
695 q2=-std::log(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}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53

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 144 of file G4ChipsPionMinusElasticXS.cc.

148{
149 G4double pMom=Pt->GetTotalMomentum();
150 G4int tgN = A - tgZ;
151
152 return GetChipsCrossSection(pMom, tgZ, tgN, -212);
153}
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 133 of file G4ChipsPionMinusElasticXS.cc.

136{
137 G4ParticleDefinition* particle = Pt->GetDefinition();
138 if (particle == G4PionMinus::PionMinus() ) return true;
139 return false;
140}
G4ParticleDefinition * GetDefinition() const
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98

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