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

#include <G4ChipsProtonElasticXS.hh>

+ Inheritance diagram for G4ChipsProtonElasticXS:

Public Member Functions

 G4ChipsProtonElasticXS ()
 
 ~G4ChipsProtonElasticXS ()
 
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)
 
G4double GetHMaxT ()
 
- 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 45 of file G4ChipsProtonElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsProtonElasticXS()

G4ChipsProtonElasticXS::G4ChipsProtonElasticXS ( )

Definition at line 59 of file G4ChipsProtonElasticXS.cc.

59 :G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
60{
61 // Initialization of the parameters
62 lPMin=-8.; // Min tabulated logarithmicMomentum(D)
63 lPMax= 8.; // Max tabulated logarithmicMomentum(D)
64 dlnP=(lPMax-lPMin)/nLast;// LogStep in the table(D)
65 onlyCS=false;// Flag toCalculateOnlyCS(not Si/Bi)(L)
66 lastSIG=0.; // Last calculated cross section (L)
67 lastLP=-10.;// Last log(mom_ofTheIncidentHadron)(L)
68 lastTM=0.; // Last t_maximum (L)
69 theSS=0.; // The Last sq.slope of 1st difr.Max(L)
70 theS1=0.; // The Last mantissa of 1st difr.Max(L)
71 theB1=0.; // The Last slope of 1st difruct.Max(L)
72 theS2=0.; // The Last mantissa of 2nd difr.Max(L)
73 theB2=0.; // The Last slope of 2nd difruct.Max(L)
74 theS3=0.; // The Last mantissa of 3d difr. Max(L)
75 theB3=0.; // The Last slope of 3d difruct. Max(L)
76 theS4=0.; // The Last mantissa of 4th difr.Max(L)
77 theB4=0.; // The Last slope of 4th difruct.Max(L)
78 lastTZ=0; // Last atomic number of the target
79 lastTN=0; // Last # of neutrons in the target
80 lastPIN=0.; // Last initialized max momentum
81 lastCST=0; // Elastic cross-section table
82 lastPAR=0; // Parameters for FunctionalCalculation
83 lastSST=0; // E-dep of sq.slope of the 1st dif.Max
84 lastS1T=0; // E-dep of mantissa of the 1st dif.Max
85 lastB1T=0; // E-dep of the slope of the 1st difMax
86 lastS2T=0; // E-dep of mantissa of the 2nd difrMax
87 lastB2T=0; // E-dep of the slope of the 2nd difMax
88 lastS3T=0; // E-dep of mantissa of the 3d difr.Max
89 lastB3T=0; // E-dep of the slope of the 3d difrMax
90 lastS4T=0; // E-dep of mantissa of the 4th difrMax
91 lastB4T=0; // E-dep of the slope of the 4th difMax
92 lastN=0; // The last N of calculated nucleus
93 lastZ=0; // The last Z of calculated nucleus
94 lastP=0.; // Last used in cross section Momentum
95 lastTH=0.; // Last threshold momentum
96 lastCS=0.; // Last value of the Cross Section
97 lastI=0; // The last position in the DAMDB
98
99 mProt= G4Proton::Proton()->GetPDGMass()*.001; // MeV to GeV
100 mProt2= mProt*mProt;
101
102
103}
static const char * Default_Name()
static G4Proton * Proton()
Definition G4Proton.cc:90
G4VCrossSectionDataSet(const G4String &nam="")

◆ ~G4ChipsProtonElasticXS()

G4ChipsProtonElasticXS::~G4ChipsProtonElasticXS ( )

Definition at line 106 of file G4ChipsProtonElasticXS.cc.

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

Member Function Documentation

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 146 of file G4ChipsProtonElasticXS.cc.

147{
148 outFile << "G4ChipsProtonElasticXS provides the elastic cross\n"
149 << "section for proton nucleus scattering as a function of incident\n"
150 << "momentum. The cross section is calculated using M. Kossov's\n"
151 << "CHIPS parameterization of cross section data.\n";
152}

◆ Default_Name()

static const char * G4ChipsProtonElasticXS::Default_Name ( )
inlinestatic

◆ GetChipsCrossSection()

G4double G4ChipsProtonElasticXS::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 176 of file G4ChipsProtonElasticXS.cc.

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

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

◆ GetExchangeT()

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

Definition at line 621 of file G4ChipsProtonElasticXS.cc.

622{
623 static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
624 static const G4double third=1./3.;
625 static const G4double fifth=1./5.;
626 static const G4double sevth=1./7.;
627 if(PDG!=2212) G4cout<<"**Warning*G4ChipsProtonElasticXS::GetExT:PDG="<<PDG<<G4endl;
628 if(onlyCS) G4cout<<"**Warning*G4ChipsProtonElasticXS::GetExchanT:onlyCS=1"<<G4endl;
629 if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
630 G4double q2=0.;
631 if(tgZ==1 && tgN==0) // ===> p+p=p+p
632 {
633 G4double E1=lastTM*theB1;
634 G4double R1=(1.-std::exp(-E1));
635 G4double E2=lastTM*theB2;
636 G4double R2=(1.-std::exp(-E2*E2*E2));
637 G4double E3=lastTM*theB3;
638 G4double R3=(1.-std::exp(-E3));
639 G4double I1=R1*theS1/theB1;
640 G4double I2=R2*theS2;
641 G4double I3=R3*theS3;
642 G4double I12=I1+I2;
643 G4double rand=(I12+I3)*G4UniformRand();
644 if (rand<I1 )
645 {
646 G4double ran=R1*G4UniformRand();
647 if(ran>1.) ran=1.;
648 q2=-std::log(1.-ran)/theB1;
649 }
650 else if(rand<I12)
651 {
652 G4double ran=R2*G4UniformRand();
653 if(ran>1.) ran=1.;
654 q2=-std::log(1.-ran);
655 if(q2<0.) q2=0.;
656 q2=std::pow(q2,third)/theB2;
657 }
658 else
659 {
660 G4double ran=R3*G4UniformRand();
661 if(ran>1.) ran=1.;
662 q2=-std::log(1.-ran)/theB3;
663 }
664 }
665 else
666 {
667 G4double a=tgZ+tgN;
668 G4double E1=lastTM*(theB1+lastTM*theSS);
669 G4double R1=(1.-std::exp(-E1));
670 G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
671 G4double tm2=lastTM*lastTM;
672 G4double E2=lastTM*tm2*theB2; // power 3 for lowA, 5 for HighA (1st)
673 if(a>6.5)E2*=tm2; // for heavy nuclei
674 G4double R2=(1.-std::exp(-E2));
675 G4double E3=lastTM*theB3;
676 if(a>6.5)E3*=tm2*tm2*tm2; // power 1 for lowA, 7 (2nd) for HighA
677 G4double R3=(1.-std::exp(-E3));
678 G4double E4=lastTM*theB4;
679 G4double R4=(1.-std::exp(-E4));
680 G4double I1=R1*theS1;
681 G4double I2=R2*theS2;
682 G4double I3=R3*theS3;
683 G4double I4=R4*theS4;
684 G4double I12=I1+I2;
685 G4double I13=I12+I3;
686 G4double rand=(I13+I4)*G4UniformRand();
687 if(rand<I1)
688 {
689 G4double ran=R1*G4UniformRand();
690 if(ran>1.) ran=1.;
691 q2=-std::log(1.-ran)/theB1;
692 if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
693 }
694 else if(rand<I12)
695 {
696 G4double ran=R2*G4UniformRand();
697 if(ran>1.) ran=1.;
698 q2=-std::log(1.-ran)/theB2;
699 if(q2<0.) q2=0.;
700 if(a<6.5) q2=std::pow(q2,third);
701 else q2=std::pow(q2,fifth);
702 }
703 else if(rand<I13)
704 {
705 G4double ran=R3*G4UniformRand();
706 if(ran>1.) ran=1.;
707 q2=-std::log(1.-ran)/theB3;
708 if(q2<0.) q2=0.;
709 if(a>6.5) q2=std::pow(q2,sevth);
710 }
711 else
712 {
713 G4double ran=R4*G4UniformRand();
714 if(ran>1.) ran=1.;
715 q2=-std::log(1.-ran)/theB4;
716 if(a<6.5) q2=lastTM-q2; // u reduced for lightA (starts from 0)
717 }
718 }
719 if(q2<0.) q2=0.;
720 if(!(q2>=-1.||q2<=1.)) G4cout<<"*NAN*G4QElasticCrossSect::GetExchangeT: -t="<<q2<<G4endl;
721 if(q2>lastTM)
722 {
723 q2=lastTM;
724 }
725 return q2*GeVSQ;
726}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52

Referenced by G4QuasiElRatios::ChExer(), G4ChipsElasticModel::SampleInvariantT(), and G4QuasiElRatios::Scatter().

◆ GetHMaxT()

G4double G4ChipsProtonElasticXS::GetHMaxT ( )

Definition at line 751 of file G4ChipsProtonElasticXS.cc.

752{
753 static const G4double HGeVSQ=gigaelectronvolt*gigaelectronvolt/2.;
754 return lastTM*HGeVSQ;
755}

Referenced by G4QuasiElRatios::ChExer(), and G4QuasiElRatios::Scatter().

◆ GetIsoCrossSection()

G4double G4ChipsProtonElasticXS::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 162 of file G4ChipsProtonElasticXS.cc.

166{
167 G4double pMom=Pt->GetTotalMomentum();
168 G4int tgN = A - tgZ;
169
170 return GetChipsCrossSection(pMom, tgZ, tgN, 2212);
171}
const G4double A[17]
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4double GetTotalMomentum() const

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 154 of file G4ChipsProtonElasticXS.cc.

157{
158 return true;
159}

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