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

#include <G4ChipsKaonPlusElasticXS.hh>

+ Inheritance diagram for G4ChipsKaonPlusElasticXS:

Public Member Functions

 G4ChipsKaonPlusElasticXS ()
 
 ~G4ChipsKaonPlusElasticXS ()
 
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 47 of file G4ChipsKaonPlusElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsKaonPlusElasticXS()

G4ChipsKaonPlusElasticXS::G4ChipsKaonPlusElasticXS ( )

Definition at line 54 of file G4ChipsKaonPlusElasticXS.cc.

54 :G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
55{
56 lPMin=-8.; //Min tabulatedLogarithmMomentum/D
57 lPMax= 8.; //Max tabulatedLogarithmMomentum/D
58 dlnP=(lPMax-lPMin)/nLast;// LogStep inTable /D
59 onlyCS=true;//Flag toCalculOnlyCS(not Si/Bi)/L
60 lastSIG=0.; //Last calculated cross section /L
61 lastLP=-10.;//LastLog(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 difrMax/L
65 theB1=0.; //TheLastSlope of 1st difructMax/L
66 theS2=0.; //TheLastMantissa of 2nd difrMax/L
67 theB2=0.; //TheLastSlope of 2nd difructMax/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 difrMax/L
71 theB4=0.; //TheLastSlope of 4th difructMax/L
72 lastTZ=0; // Last atomic number of theTarget
73 lastTN=0; // Last # of neutrons in theTarget
74 lastPIN=0.;// Last initialized max momentum
75 lastCST=0; // Elastic cross-section table
76 lastPAR=0; // ParametersForFunctionCalculation
77 lastSST=0; // E-dep ofSqardSlope 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 inCrossSection 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}
static const char * Default_Name()

◆ ~G4ChipsKaonPlusElasticXS()

G4ChipsKaonPlusElasticXS::~G4ChipsKaonPlusElasticXS ( )

Definition at line 94 of file G4ChipsKaonPlusElasticXS.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 * G4ChipsKaonPlusElasticXS::Default_Name ( )
inlinestatic

◆ GetChipsCrossSection()

G4double G4ChipsKaonPlusElasticXS::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 154 of file G4ChipsKaonPlusElasticXS.cc.

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

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

◆ GetExchangeT()

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

Definition at line 601 of file G4ChipsKaonPlusElasticXS.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!=321) G4cout<<"*Warning*G4ChipsKaonPlusElasticXS::GetExT:PDG="<<PDG<<G4endl;
608 if(onlyCS) G4cout<<"*Warning*G4ChipsKaonPlusElasticXS::GetExT: 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*G4QKaonPlusElasticCS::GetExchT: -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 G4ChipsKaonPlusElasticXS::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 143 of file G4ChipsKaonPlusElasticXS.cc.

147{
148 G4double pMom=Pt->GetTotalMomentum();
149 G4int tgN = A - tgZ;
150
151 return GetChipsCrossSection(pMom, tgZ, tgN, 321);
152}
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4double GetTotalMomentum() const

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 132 of file G4ChipsKaonPlusElasticXS.cc.

135{
136 G4ParticleDefinition* particle = Pt->GetDefinition();
137 if (particle == G4KaonPlus::KaonPlus() ) return true;
138 return false;
139}
G4ParticleDefinition * GetDefinition() const
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113

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