Geant4 9.6.0
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 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=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 46 of file G4ChipsProtonElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsProtonElasticXS()

G4ChipsProtonElasticXS::G4ChipsProtonElasticXS ( )

Definition at line 55 of file G4ChipsProtonElasticXS.cc.

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

◆ ~G4ChipsProtonElasticXS()

G4ChipsProtonElasticXS::~G4ChipsProtonElasticXS ( )

Definition at line 97 of file G4ChipsProtonElasticXS.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
134}

Member Function Documentation

◆ 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 160 of file G4ChipsProtonElasticXS.cc.

161{
162 static std::vector <G4int> colN; // Vector of N for calculated nuclei (isotops)
163 static std::vector <G4int> colZ; // Vector of Z for calculated nuclei (isotops)
164 static std::vector <G4double> colP; // Vector of last momenta for the reaction
165 static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
166 static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
167 // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
168
169 G4double pEn=pMom;
170 onlyCS=false;
171
172 G4bool in=false; // By default the isotope must be found in the AMDB
173 lastP = 0.; // New momentum history (nothing to compare with)
174 lastN = tgN; // The last N of the calculated nucleus
175 lastZ = tgZ; // The last Z of the calculated nucleus
176 lastI = colN.size(); // Size of the Associative Memory DB in the heap
177 if(lastI) for(G4int i=0; i<lastI; i++) // Loop over proj/tgZ/tgN lines of DB
178 { // The nucleus with projPDG is found in AMDB
179 if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
180 {
181 lastI=i;
182 lastTH =colTH[i]; // Last THreshold (A-dependent)
183 if(pEn<=lastTH)
184 {
185 return 0.; // Energy is below the Threshold value
186 }
187 lastP =colP [i]; // Last Momentum (A-dependent)
188 lastCS =colCS[i]; // Last CrossSect (A-dependent)
189 if(lastP == pMom) // Do not recalculate
190 {
191 CalculateCrossSection(onlyCS,-1,i,2212,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(onlyCS,-1,i,2212,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(onlyCS,0,lastI,2212,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:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67

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 614 of file G4ChipsProtonElasticXS.cc.

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

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

◆ GetHMaxT()

G4double G4ChipsProtonElasticXS::GetHMaxT ( )

Definition at line 744 of file G4ChipsProtonElasticXS.cc.

745{
746 static const G4double HGeVSQ=gigaelectronvolt*gigaelectronvolt/2.;
747 return lastTM*HGeVSQ;
748}

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 146 of file G4ChipsProtonElasticXS.cc.

150{
151 G4double pMom=Pt->GetTotalMomentum();
152 G4int tgN = A - tgZ;
153
154 return GetChipsCrossSection(pMom, tgZ, tgN, 2212);
155}
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 136 of file G4ChipsProtonElasticXS.cc.

139{
140 G4ParticleDefinition* particle = Pt->GetDefinition();
141 if (particle == G4Proton::Proton() ) return true;
142 return false;
143}
G4ParticleDefinition * GetDefinition() const
static G4Proton * Proton()
Definition: G4Proton.cc:93

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