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

#include <G4ChipsKaonMinusElasticXS.hh>

+ Inheritance diagram for G4ChipsKaonMinusElasticXS:

Public Member Functions

 G4ChipsKaonMinusElasticXS ()
 
 ~G4ChipsKaonMinusElasticXS ()
 
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 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 GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, 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 CrossSectionDescription (std::ostream &) 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
 
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 46 of file G4ChipsKaonMinusElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsKaonMinusElasticXS()

G4ChipsKaonMinusElasticXS::G4ChipsKaonMinusElasticXS ( )

Definition at line 75 of file G4ChipsKaonMinusElasticXS.cc.

75 :G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
76{
77 G4AutoLock l(&initM);
79 mK2 = mK*mK;
80 l.unlock();
81 lPMin=-8.; //Min tabulatedLogarithmMomentum/D
82 lPMax= 8.; //Max tabulatedLogarithmMomentum/D
83 dlnP=(lPMax-lPMin)/nLast;// LogStep inTable /D
84 onlyCS=true;//Flag toCalculOnlyCS(not Si/Bi)/L
85 lastSIG=0.; //Last calculated cross section /L
86 lastLP=-10.;//LastLog(mom_of IncidentHadron)/L
87 lastTM=0.; //Last t_maximum /L
88 theSS=0.; //TheLastSqSlope of 1st difr.Max/L
89 theS1=0.; //TheLastMantissa of 1st difrMax/L
90 theB1=0.; //TheLastSlope of 1st difructMax/L
91 theS2=0.; //TheLastMantissa of 2nd difrMax/L
92 theB2=0.; //TheLastSlope of 2nd difructMax/L
93 theS3=0.; //TheLastMantissa of 3d difr.Max/L
94 theB3=0.; //TheLastSlope of 3d difruct.Max/L
95 theS4=0.; //TheLastMantissa of 4th difrMax/L
96 theB4=0.; //TheLastSlope of 4th difructMax/L
97 lastTZ=0; // Last atomic number of theTarget
98 lastTN=0; // Last # of neutrons in theTarget
99 lastPIN=0.;// Last initialized max momentum
100 lastCST=0; // Elastic cross-section table
101 lastPAR=0; // ParametersForFunctionCalculation
102 lastSST=0; // E-dep ofSqardSlope of 1st difMax
103 lastS1T=0; // E-dep of mantissa of 1st dif.Max
104 lastB1T=0; // E-dep of the slope of 1st difMax
105 lastS2T=0; // E-dep of mantissa of 2nd difrMax
106 lastB2T=0; // E-dep of the slope of 2nd difMax
107 lastS3T=0; // E-dep of mantissa of 3d difr.Max
108 lastB3T=0; // E-dep of the slope of 3d difrMax
109 lastS4T=0; // E-dep of mantissa of 4th difrMax
110 lastB4T=0; // E-dep of the slope of 4th difMax
111 lastN=0; // The last N of calculated nucleus
112 lastZ=0; // The last Z of calculated nucleus
113 lastP=0.; // LastUsed inCrossSection Momentum
114 lastTH=0.; // Last threshold momentum
115 lastCS=0.; // Last value of the Cross Section
116 lastI=0; // The last position in the DAMDB
117}
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:112

◆ ~G4ChipsKaonMinusElasticXS()

G4ChipsKaonMinusElasticXS::~G4ChipsKaonMinusElasticXS ( )

Definition at line 119 of file G4ChipsKaonMinusElasticXS.cc.

120{
121 std::vector<G4double*>::iterator pos;
122 for (pos=CST.begin(); pos<CST.end(); pos++)
123 { delete [] *pos; }
124 CST.clear();
125 for (pos=PAR.begin(); pos<PAR.end(); pos++)
126 { delete [] *pos; }
127 PAR.clear();
128 for (pos=SST.begin(); pos<SST.end(); pos++)
129 { delete [] *pos; }
130 SST.clear();
131 for (pos=S1T.begin(); pos<S1T.end(); pos++)
132 { delete [] *pos; }
133 S1T.clear();
134 for (pos=B1T.begin(); pos<B1T.end(); pos++)
135 { delete [] *pos; }
136 B1T.clear();
137 for (pos=S2T.begin(); pos<S2T.end(); pos++)
138 { delete [] *pos; }
139 S2T.clear();
140 for (pos=B2T.begin(); pos<B2T.end(); pos++)
141 { delete [] *pos; }
142 B2T.clear();
143 for (pos=S3T.begin(); pos<S3T.end(); pos++)
144 { delete [] *pos; }
145 S3T.clear();
146 for (pos=B3T.begin(); pos<B3T.end(); pos++)
147 { delete [] *pos; }
148 B3T.clear();
149 for (pos=S4T.begin(); pos<S4T.end(); pos++)
150 { delete [] *pos; }
151 S4T.clear();
152 for (pos=B4T.begin(); pos<B4T.end(); pos++)
153 { delete [] *pos; }
154 B4T.clear();
155}

Member Function Documentation

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 158 of file G4ChipsKaonMinusElasticXS.cc.

159{
160 outFile << "G4ChipsKaonMinusElasticXS provides the elastic cross\n"
161 << "section for K- nucleus scattering as a function of incident\n"
162 << "momentum. The cross section is calculated using M. Kossov's\n"
163 << "CHIPS parameterization of cross section data.\n";
164}

◆ Default_Name()

static const char * G4ChipsKaonMinusElasticXS::Default_Name ( )
inlinestatic

◆ GetChipsCrossSection()

G4double G4ChipsKaonMinusElasticXS::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 186 of file G4ChipsKaonMinusElasticXS.cc.

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

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

◆ GetExchangeT()

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

Definition at line 616 of file G4ChipsKaonMinusElasticXS.cc.

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

Referenced by G4ChipsElasticModel::SampleInvariantT().

◆ GetIsoCrossSection()

G4double G4ChipsKaonMinusElasticXS::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 175 of file G4ChipsKaonMinusElasticXS.cc.

179{
180 G4double pMom=Pt->GetTotalMomentum();
181 G4int tgN = A - tgZ;
182
183 return GetChipsCrossSection(pMom, tgZ, tgN, -321);
184}
const G4double A[17]
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4double GetTotalMomentum() const

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 166 of file G4ChipsKaonMinusElasticXS.cc.

169{
170 return true;
171}

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