Geant4 10.7.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 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 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 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 G4int GetVerboseLevel () 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
 

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 G4ChipsKaonPlusElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsKaonPlusElasticXS()

G4ChipsKaonPlusElasticXS::G4ChipsKaonPlusElasticXS ( )

Definition at line 76 of file G4ChipsKaonPlusElasticXS.cc.

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

◆ ~G4ChipsKaonPlusElasticXS()

G4ChipsKaonPlusElasticXS::~G4ChipsKaonPlusElasticXS ( )

Definition at line 120 of file G4ChipsKaonPlusElasticXS.cc.

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

Member Function Documentation

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 159 of file G4ChipsKaonPlusElasticXS.cc.

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

◆ 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 188 of file G4ChipsKaonPlusElasticXS.cc.

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

Definition at line 616 of file G4ChipsKaonPlusElasticXS.cc.

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

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 177 of file G4ChipsKaonPlusElasticXS.cc.

181{
182 G4double pMom=Pt->GetTotalMomentum();
183 G4int tgN = A - tgZ;
184
185 return GetChipsCrossSection(pMom, tgZ, tgN, 321);
186}
double A(double temperature)
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 168 of file G4ChipsKaonPlusElasticXS.cc.

171{
172 return true;
173}

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