Geant4  10.02.p03
G4ChipsKaonMinusElasticXS Class Reference

#include <G4ChipsKaonMinusElasticXS.hh>

Inheritance diagram for G4ChipsKaonMinusElasticXS:
Collaboration 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=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 G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual G4int GetVerboseLevel () const
 
virtual 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 ()
 

Private Member Functions

G4double CalculateCrossSection (G4bool CS, G4int F, G4int I, G4int pPDG, G4int Z, G4int N, G4double pP)
 
G4double GetSlope (G4int tZ, G4int tN, G4int pPDG)
 
G4double GetHMaxT ()
 
G4double GetPTables (G4double lpP, G4double lPm, G4int PDG, G4int tZ, G4int tN)
 
G4double GetTabValues (G4double lp, G4int pPDG, G4int tgZ, G4int tgN)
 
G4double GetQ2max (G4int pPDG, G4int tgZ, G4int tgN, G4double pP)
 

Private Attributes

const G4int nPoints
 
const G4int nLast
 
G4double lPMin
 
G4double lPMax
 
G4double dlnP
 
G4bool onlyCS
 
G4double lastSIG
 
G4double lastLP
 
G4double lastTM
 
G4int lastN
 
G4int lastZ
 
G4double lastP
 
G4double lastTH
 
G4double lastCS
 
G4int lastI
 
G4double theSS
 
G4double theS1
 
G4double theB1
 
G4double theS2
 
G4double theB2
 
G4double theS3
 
G4double theB3
 
G4double theS4
 
G4double theB4
 
G4int lastTZ
 
G4int lastTN
 
G4double lastPIN
 
G4doublelastCST
 
G4doublelastPAR
 
G4doublelastSST
 
G4doublelastS1T
 
G4doublelastB1T
 
G4doublelastS2T
 
G4doublelastB2T
 
G4doublelastS3T
 
G4doublelastB3T
 
G4doublelastS4T
 
G4doublelastB4T
 
std::vector< G4double * > PAR
 
std::vector< G4double * > CST
 
std::vector< G4double * > SST
 
std::vector< G4double * > S1T
 
std::vector< G4double * > B1T
 
std::vector< G4double * > S2T
 
std::vector< G4double * > B2T
 
std::vector< G4double * > S3T
 
std::vector< G4double * > B3T
 
std::vector< G4double * > S4T
 
std::vector< G4double * > B4T
 
std::vector< G4intcolN
 
std::vector< G4intcolZ
 
std::vector< G4doublecolP
 
std::vector< G4doublecolTH
 
std::vector< G4doublecolCS
 
std::vector< G4doublePIN
 

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

Constructor & Destructor Documentation

◆ G4ChipsKaonMinusElasticXS()

G4ChipsKaonMinusElasticXS::G4ChipsKaonMinusElasticXS ( )

Definition at line 76 of file G4ChipsKaonMinusElasticXS.cc.

77 {
78  G4AutoLock l(&initM);
79  mK = G4KaonMinus::KaonMinus()->GetPDGMass()*.001;
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 }
G4VCrossSectionDataSet(const G4String &nam="")
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
Here is the call graph for this function:

◆ ~G4ChipsKaonMinusElasticXS()

G4ChipsKaonMinusElasticXS::~G4ChipsKaonMinusElasticXS ( )

Definition at line 120 of file G4ChipsKaonMinusElasticXS.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 }
static const G4double pos

Member Function Documentation

◆ CalculateCrossSection()

G4double G4ChipsKaonMinusElasticXS::CalculateCrossSection ( G4bool  CS,
G4int  F,
G4int  I,
G4int  pPDG,
G4int  Z,
G4int  N,
G4double  pP 
)
private

Definition at line 257 of file G4ChipsKaonMinusElasticXS.cc.

259 {
260  G4double pMom=pIU/GeV; // All calculations are in GeV
261  onlyCS=CS; // Flag to calculate only CS (not Si/Bi)
262  lastLP=std::log(pMom); // Make a logarithm of the momentum for calculation
263  if(F) // This isotope was found in AMDB =>RETRIEVE/UPDATE
264  {
265  if(F<0) // the AMDB must be loded
266  {
267  lastPIN = PIN[I]; // Max log(P) initialised for this table set
268  lastPAR = PAR[I]; // Pointer to the parameter set
269  lastCST = CST[I]; // Pointer to the total sross-section table
270  lastSST = SST[I]; // Pointer to the first squared slope
271  lastS1T = S1T[I]; // Pointer to the first mantissa
272  lastB1T = B1T[I]; // Pointer to the first slope
273  lastS2T = S2T[I]; // Pointer to the second mantissa
274  lastB2T = B2T[I]; // Pointer to the second slope
275  lastS3T = S3T[I]; // Pointer to the third mantissa
276  lastB3T = B3T[I]; // Pointer to the rhird slope
277  lastS4T = S4T[I]; // Pointer to the 4-th mantissa
278  lastB4T = B4T[I]; // Pointer to the 4-th slope
279  }
280  if(lastLP>lastPIN && lastLP<lPMax)
281  {
282  lastPIN=GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);// Can update upper logP-Limit in tabs
283  PIN[I]=lastPIN; // Remember the new P-Limit of the tables
284  }
285  }
286  else // This isotope wasn't initialized => CREATE
287  {
288  lastPAR = new G4double[nPoints]; // Allocate memory for parameters of CS function
289  lastPAR[nLast]=0; // Initialization for VALGRIND
290  lastCST = new G4double[nPoints]; // Allocate memory for Tabulated CS function
291  lastSST = new G4double[nPoints]; // Allocate memory for Tabulated first sqaredSlope
292  lastS1T = new G4double[nPoints]; // Allocate memory for Tabulated first mantissa
293  lastB1T = new G4double[nPoints]; // Allocate memory for Tabulated first slope
294  lastS2T = new G4double[nPoints]; // Allocate memory for Tabulated second mantissa
295  lastB2T = new G4double[nPoints]; // Allocate memory for Tabulated second slope
296  lastS3T = new G4double[nPoints]; // Allocate memory for Tabulated third mantissa
297  lastB3T = new G4double[nPoints]; // Allocate memory for Tabulated third slope
298  lastS4T = new G4double[nPoints]; // Allocate memory for Tabulated 4-th mantissa
299  lastB4T = new G4double[nPoints]; // Allocate memory for Tabulated 4-th slope
300  lastPIN = GetPTables(lastLP,lPMin,PDG,tgZ,tgN); // Returns the new P-limit for tables
301  PIN.push_back(lastPIN); // Fill parameters of CS function to AMDB
302  PAR.push_back(lastPAR); // Fill parameters of CS function to AMDB
303  CST.push_back(lastCST); // Fill Tabulated CS function to AMDB
304  SST.push_back(lastSST); // Fill Tabulated first sq.slope to AMDB
305  S1T.push_back(lastS1T); // Fill Tabulated first mantissa to AMDB
306  B1T.push_back(lastB1T); // Fill Tabulated first slope to AMDB
307  S2T.push_back(lastS2T); // Fill Tabulated second mantissa to AMDB
308  B2T.push_back(lastB2T); // Fill Tabulated second slope to AMDB
309  S3T.push_back(lastS3T); // Fill Tabulated third mantissa to AMDB
310  B3T.push_back(lastB3T); // Fill Tabulated third slope to AMDB
311  S4T.push_back(lastS4T); // Fill Tabulated 4-th mantissa to AMDB
312  B4T.push_back(lastB4T); // Fill Tabulated 4-th slope to AMDB
313  } // End of creation/update of the new set of parameters and tables
314  // =----------= NOW Update (if necessary) and Calculate the Cross Section =-----------=
315  if(lastLP>lastPIN && lastLP<lPMax)
316  {
317  lastPIN = GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);
318  }
319  if(!onlyCS) lastTM=GetQ2max(PDG, tgZ, tgN, pMom); // Calculate (-t)_max=Q2_max (GeV2)
320  if(lastLP>lPMin && lastLP<=lastPIN) // Linear fit is made using precalculated tables
321  {
322  if(lastLP==lastPIN)
323  {
324  G4double shift=(lastLP-lPMin)/dlnP+.000001; // Log distance from lPMin
325  G4int blast=static_cast<int>(shift); // this is a bin number of the lower edge (0)
326  if(blast<0 || blast>=nLast) G4cout<<"G4QKMElCS::CCS:b="<<blast<<",n="<<nLast<<G4endl;
327  lastSIG = lastCST[blast];
328  if(!onlyCS) // Skip the differential cross-section parameters
329  {
330  theSS = lastSST[blast];
331  theS1 = lastS1T[blast];
332  theB1 = lastB1T[blast];
333  theS2 = lastS2T[blast];
334  theB2 = lastB2T[blast];
335  theS3 = lastS3T[blast];
336  theB3 = lastB3T[blast];
337  theS4 = lastS4T[blast];
338  theB4 = lastB4T[blast];
339  }
340  }
341  else
342  {
343  G4double shift=(lastLP-lPMin)/dlnP; // a shift from the beginning of the table
344  G4int blast=static_cast<int>(shift); // the lower bin number
345  if(blast<0) blast=0;
346  if(blast>=nLast) blast=nLast-1; // low edge of the last bin
347  shift-=blast; // step inside the unit bin
348  G4int lastL=blast+1; // the upper bin number
349  G4double SIGL=lastCST[blast]; // the basic value of the cross-section
350  lastSIG= SIGL+shift*(lastCST[lastL]-SIGL); // calculated total elastic cross-section
351  if(!onlyCS) // Skip the differential cross-section parameters
352  {
353  G4double SSTL=lastSST[blast]; // the low bin of the first squared slope
354  theSS=SSTL+shift*(lastSST[lastL]-SSTL); // the basic value of the first sq.slope
355  G4double S1TL=lastS1T[blast]; // the low bin of the first mantissa
356  theS1=S1TL+shift*(lastS1T[lastL]-S1TL); // the basic value of the first mantissa
357  G4double B1TL=lastB1T[blast]; // the low bin of the first slope
358  theB1=B1TL+shift*(lastB1T[lastL]-B1TL); // the basic value of the first slope
359  G4double S2TL=lastS2T[blast]; // the low bin of the second mantissa
360  theS2=S2TL+shift*(lastS2T[lastL]-S2TL); // the basic value of the second mantissa
361  G4double B2TL=lastB2T[blast]; // the low bin of the second slope
362  theB2=B2TL+shift*(lastB2T[lastL]-B2TL); // the basic value of the second slope
363  G4double S3TL=lastS3T[blast]; // the low bin of the third mantissa
364  theS3=S3TL+shift*(lastS3T[lastL]-S3TL); // the basic value of the third mantissa
365  G4double B3TL=lastB3T[blast]; // the low bin of the third slope
366  theB3=B3TL+shift*(lastB3T[lastL]-B3TL); // the basic value of the third slope
367  G4double S4TL=lastS4T[blast]; // the low bin of the 4-th mantissa
368  theS4=S4TL+shift*(lastS4T[lastL]-S4TL); // the basic value of the 4-th mantissa
369  G4double B4TL=lastB4T[blast]; // the low bin of the 4-th slope
370  theB4=B4TL+shift*(lastB4T[lastL]-B4TL); // the basic value of the 4-th slope
371  }
372  }
373  }
374  else lastSIG=GetTabValues(lastLP, PDG, tgZ, tgN); // Direct calculation beyond the table
375  if(lastSIG<0.) lastSIG = 0.; // @@ a Warning print can be added
376  return lastSIG;
377 }
G4double GetPTables(G4double lpP, G4double lPm, G4int PDG, G4int tZ, G4int tN)
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static const double GeV
Definition: G4SIunits.hh:214
G4double GetTabValues(G4double lp, G4int pPDG, G4int tgZ, G4int tgN)
#define G4endl
Definition: G4ios.hh:61
G4double GetQ2max(G4int pPDG, G4int tgZ, G4int tgN, G4double pP)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 159 of file G4ChipsKaonMinusElasticXS.cc.

160 {
161  outFile << "G4ChipsKaonMinusElasticXS 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 }
Here is the caller graph for this function:

◆ Default_Name()

static const char* G4ChipsKaonMinusElasticXS::Default_Name ( )
inlinestatic

Definition at line 55 of file G4ChipsKaonMinusElasticXS.hh.

55 {return "ChipsKaonMinusElasticXS";}
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 187 of file G4ChipsKaonMinusElasticXS.cc.

188 {
189 
190  G4bool fCS = false;
191 
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  {
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 }
ifstream in
Definition: comparison.C:7
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
G4double CalculateCrossSection(G4bool CS, G4int F, G4int I, G4int pPDG, G4int Z, G4int N, G4double pP)
static const double millibarn
Definition: G4SIunits.hh:105
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetExchangeT()

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

Definition at line 617 of file G4ChipsKaonMinusElasticXS.cc.

618 {
619  if(PDG==310 || PDG==130) PDG=-321;
620  if(PDG!=-321)G4cout<<"*Warning*G4ChipsKaonMinusElasticXS::GetET:PDG="<<PDG<<G4endl;
621  if(onlyCS) G4cout<<"*Warning*G4ChipsKaonMinusElasticXS::GetExT: 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*G4QKaonMinusElasticCS::GetExchT: -t="<<q2<<G4endl;
714  if(q2>lastTM)
715  {
716  q2=lastTM;
717  }
718  return q2*GeVSQ;
719 }
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static const G4double third
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ GetHMaxT()

G4double G4ChipsKaonMinusElasticXS::GetHMaxT ( )
private

Definition at line 741 of file G4ChipsKaonMinusElasticXS.cc.

742 {
743  return lastTM*HGeVSQ;
744 }
Here is the caller graph for this function:

◆ 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 176 of file G4ChipsKaonMinusElasticXS.cc.

180 {
181  G4double pMom=Pt->GetTotalMomentum();
182  G4int tgN = A - tgZ;
183 
184  return GetChipsCrossSection(pMom, tgZ, tgN, -321);
185 }
G4double GetTotalMomentum() const
int G4int
Definition: G4Types.hh:78
double A(double temperature)
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetPTables()

G4double G4ChipsKaonMinusElasticXS::GetPTables ( G4double  lpP,
G4double  lPm,
G4int  PDG,
G4int  tZ,
G4int  tN 
)
private

Definition at line 380 of file G4ChipsKaonMinusElasticXS.cc.

382 {
383  // @@ At present all nA==pA ---------> Each neucleus can have not more than 51 parameters
384  if(PDG == -321)
385  {
386  // -- Total pp elastic cross section cs & s1/b1 (main), s2/b2 (tail1), s3/b3 (tail2) --
387  //p2=p*p;p3=p2*p;sp=sqrt(p);p2s=p2*sp;lp=log(p);dl1=lp-(3.=par(3));p4=p2*p2; p=|3-mom|
388  //CS=2.865/p2s/(1+.0022/p2s)+(18.9+.6461*dl1*dl1+9./p)/(1.+.425*lp)/(1.+.4276/p4);
389  // par(0) par(7) par(1) par(2) par(4) par(5) par(6)
390  //dl2=lp-5., s1=(74.+3.*dl2*dl2)/(1+3.4/p4/p)+(.2/p2+17.*p)/(p4+.001*sp),
391  // par(8) par(9) par(10) par(11) par(12)par(13) par(14)
392  // b1=8.*p**.055/(1.+3.64/p3); s2=5.e-5+4000./(p4+1500.*p); b2=.46+1.2e6/(p4+3.5e6/sp);
393  // par(15) par(16) par(17) par(18) par(19) par(20) par(21) par(22) par(23)
394  // s3=5.e-5+1.e10/(p4*p4+8.5e8*p2+1.e10); b3=1.1+3.4e6/(p4+6.8e6); ss=0.
395  // par(24) par(25) par(26) par(27) par(28) par(29) par(30) par(31)
396  //
397  if(lastPAR[nLast]!=pwd) // A unique flag to avoid the repeatable definition
398  {
399  if ( tgZ == 1 && tgN == 0 )
400  {
401  for (G4int ip=0; ip<n_kmpel; ip++) lastPAR[ip]=kmp_el[ip]; // PiMinus+P
402  }
403  else
404  {
405  G4double a=tgZ+tgN;
406  G4double sa=std::sqrt(a);
407  G4double ssa=std::sqrt(sa);
408  G4double asa=a*sa;
409  G4double a2=a*a;
410  G4double a3=a2*a;
411  G4double a4=a3*a;
412  G4double a5=a4*a;
413  G4double a6=a4*a2;
414  G4double a7=a6*a;
415  G4double a8=a7*a;
416  G4double a9=a8*a;
417  G4double a10=a5*a5;
418  G4double a12=a6*a6;
419  G4double a14=a7*a7;
420  G4double a16=a8*a8;
421  G4double a17=a16*a;
422  //G4double a20=a16*a4;
423  G4double a32=a16*a16;
424  // Reaction cross-section parameters (kmael_fit.f)
425  lastPAR[0]=.06*asa/(1.+a*(.01+.1/ssa)); // p1
426  lastPAR[1]=.75*asa/(1.+.009*a); // p2
427  lastPAR[2]=.1*a2*ssa/(1.+.0015*a2/ssa); // p3
428  lastPAR[3]=1./(1.+500./a2); // p4
429  lastPAR[4]=4.2; // p5
430  lastPAR[5]=0.; // p6 not used
431  lastPAR[6]=0.; // p7 not used
432  lastPAR[7]=0.; // p8 not used
433  lastPAR[8]=0.; // p9 not used
434  // @@ the differential cross-section is parameterized separately for A>6 & A<7
435  if(a<6.5)
436  {
437  G4double a28=a16*a12;
438  // The main pre-exponent (pel_sg)
439  lastPAR[ 9]=4000*a; // p1
440  lastPAR[10]=1.2e7*a8+380*a17; // p2
441  lastPAR[11]=.7/(1.+4.e-12*a16); // p3
442  lastPAR[12]=2.5/a8/(a4+1.e-16*a32); // p4
443  lastPAR[13]=.28*a; // p5
444  lastPAR[14]=1.2*a2+2.3; // p6
445  lastPAR[15]=3.8/a; // p7
446  // The main slope (pel_sl)
447  lastPAR[16]=.01/(1.+.0024*a5); // p1
448  lastPAR[17]=.2*a; // p2
449  lastPAR[18]=9.e-7/(1.+.035*a5); // p3
450  lastPAR[19]=(42.+2.7e-11*a16)/(1.+.14*a); // p4
451  // The main quadratic (pel_sh)
452  lastPAR[20]=2.25*a3; // p1
453  lastPAR[21]=18.; // p2
454  lastPAR[22]=2.4e-3*a8/(1.+2.6e-4*a7); // p3
455  lastPAR[23]=3.5e-36*a32*a8/(1.+5.e-15*a32/a); // p4
456  // The 1st max pre-exponent (pel_qq)
457  lastPAR[24]=1.e5/(a8+2.5e12/a16); // p1
458  lastPAR[25]=8.e7/(a12+1.e-27*a28*a28); // p2
459  lastPAR[26]=.0006*a3; // p3
460  // The 1st max slope (pel_qs)
461  lastPAR[27]=10.+4.e-8*a12*a; // p1
462  lastPAR[28]=.114; // p2
463  lastPAR[29]=.003; // p3
464  lastPAR[30]=2.e-23; // p4
465  // The effective pre-exponent (pel_ss)
466  lastPAR[31]=1./(1.+.0001*a8); // p1
467  lastPAR[32]=1.5e-4/(1.+5.e-6*a12); // p2
468  lastPAR[33]=.03; // p3
469  // The effective slope (pel_sb)
470  lastPAR[34]=a/2; // p1
471  lastPAR[35]=2.e-7*a4; // p2
472  lastPAR[36]=4.; // p3
473  lastPAR[37]=64./a3; // p4
474  // The gloria pre-exponent (pel_us)
475  lastPAR[38]=1.e8*std::exp(.32*asa); // p1
476  lastPAR[39]=20.*std::exp(.45*asa); // p2
477  lastPAR[40]=7.e3+2.4e6/a5; // p3
478  lastPAR[41]=2.5e5*std::exp(.085*a3); // p4
479  lastPAR[42]=2.5*a; // p5
480  // The gloria slope (pel_ub)
481  lastPAR[43]=920.+.03*a8*a3; // p1
482  lastPAR[44]=93.+.0023*a12; // p2
483  }
484  else
485  {
486  G4double p1a10=2.2e-28*a10;
487  G4double r4a16=6.e14/a16;
488  G4double s4a16=r4a16*r4a16;
489  // a24
490  // a36
491  // The main pre-exponent (peh_sg)
492  lastPAR[ 9]=4.5*std::pow(a,1.15); // p1
493  lastPAR[10]=.06*std::pow(a,.6); // p2
494  lastPAR[11]=.6*a/(1.+2.e15/a16); // p3
495  lastPAR[12]=.17/(a+9.e5/a3+1.5e33/a32); // p4
496  lastPAR[13]=(.001+7.e-11*a5)/(1.+4.4e-11*a5); // p5
497  lastPAR[14]=(p1a10*p1a10+2.e-29)/(1.+2.e-22*a12); // p6
498  // The main slope (peh_sl)
499  lastPAR[15]=400./a12+2.e-22*a9; // p1
500  lastPAR[16]=1.e-32*a12/(1.+5.e22/a14); // p2
501  lastPAR[17]=1000./a2+9.5*sa*ssa; // p3
502  lastPAR[18]=4.e-6*a*asa+1.e11/a16; // p4
503  lastPAR[19]=(120./a+.002*a2)/(1.+2.e14/a16); // p5
504  lastPAR[20]=9.+100./a; // p6
505  // The main quadratic (peh_sh)
506  lastPAR[21]=.002*a3+3.e7/a6; // p1
507  lastPAR[22]=7.e-15*a4*asa; // p2
508  lastPAR[23]=9000./a4; // p3
509  // The 1st max pre-exponent (peh_qq)
510  lastPAR[24]=.0011*asa/(1.+3.e34/a32/a4); // p1
511  lastPAR[25]=1.e-5*a2+2.e14/a16; // p2
512  lastPAR[26]=1.2e-11*a2/(1.+1.5e19/a12); // p3
513  lastPAR[27]=.016*asa/(1.+5.e16/a16); // p4
514  // The 1st max slope (peh_qs)
515  lastPAR[28]=.002*a4/(1.+7.e7/std::pow(a-6.83,14)); // p1
516  lastPAR[29]=2.e6/a6+7.2/std::pow(a,.11); // p2
517  lastPAR[30]=11.*a3/(1.+7.e23/a16/a8); // p3
518  lastPAR[31]=100./asa; // p4
519  // The 2nd max pre-exponent (peh_ss)
520  lastPAR[32]=(.1+4.4e-5*a2)/(1.+5.e5/a4); // p1
521  lastPAR[33]=3.5e-4*a2/(1.+1.e8/a8); // p2
522  lastPAR[34]=1.3+3.e5/a4; // p3
523  lastPAR[35]=500./(a2+50.)+3; // p4
524  lastPAR[36]=1.e-9/a+s4a16*s4a16; // p5
525  // The 2nd max slope (peh_sb)
526  lastPAR[37]=.4*asa+3.e-9*a6; // p1
527  lastPAR[38]=.0005*a5; // p2
528  lastPAR[39]=.002*a5; // p3
529  lastPAR[40]=10.; // p4
530  // The effective pre-exponent (peh_us)
531  lastPAR[41]=.05+.005*a; // p1
532  lastPAR[42]=7.e-8/sa; // p2
533  lastPAR[43]=.8*sa; // p3
534  lastPAR[44]=.02*sa; // p4
535  lastPAR[45]=1.e8/a3; // p5
536  lastPAR[46]=3.e32/(a32+1.e32); // p6
537  // The effective slope (peh_ub)
538  lastPAR[47]=24.; // p1
539  lastPAR[48]=20./sa; // p2
540  lastPAR[49]=7.e3*a/(sa+1.); // p3
541  lastPAR[50]=900.*sa/(1.+500./a3); // p4
542  }
543  // Parameter for lowEnergyNeutrons
544  lastPAR[51]=1.e15+2.e27/a4/(1.+2.e-18*a16);
545  }
546  lastPAR[nLast]=pwd;
547  // and initialize the zero element of the table
548  G4double lp=lPMin; // ln(momentum)
549  G4bool memCS=onlyCS; // ??
550  onlyCS=false;
551  lastCST[0]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables
552  onlyCS=memCS;
553  lastSST[0]=theSS;
554  lastS1T[0]=theS1;
555  lastB1T[0]=theB1;
556  lastS2T[0]=theS2;
557  lastB2T[0]=theB2;
558  lastS3T[0]=theS3;
559  lastB3T[0]=theB3;
560  lastS4T[0]=theS4;
561  lastB4T[0]=theB4;
562  }
563  if(LP>ILP)
564  {
565  G4int ini = static_cast<int>((ILP-lPMin+.000001)/dlnP)+1; // already inited till this
566  if(ini<0) ini=0;
567  if(ini<nPoints)
568  {
569  G4int fin = static_cast<int>((LP-lPMin)/dlnP)+1; // final bin of initialization
570  if(fin>=nPoints) fin=nLast; // Limit of the tabular initialization
571  if(fin>=ini)
572  {
573  G4double lp=0.;
574  for(G4int ip=ini; ip<=fin; ip++) // Calculate tabular CS,S1,B1,S2,B2,S3,B3
575  {
576  lp=lPMin+ip*dlnP; // ln(momentum)
577  G4bool memCS=onlyCS;
578  onlyCS=false;
579  lastCST[ip]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables (ret CS)
580  onlyCS=memCS;
581  lastSST[ip]=theSS;
582  lastS1T[ip]=theS1;
583  lastB1T[ip]=theB1;
584  lastS2T[ip]=theS2;
585  lastB2T[ip]=theB2;
586  lastS3T[ip]=theS3;
587  lastB3T[ip]=theB3;
588  lastS4T[ip]=theS4;
589  lastB4T[ip]=theB4;
590  }
591  return lp;
592  }
593  else G4cout<<"*Warning*G4ChipsKaonMinusElasticXS::GetPTables: PDG="<<PDG
594  <<", Z="<<tgZ<<", N="<<tgN<<", i="<<ini<<" > fin="<<fin<<", LP="<<LP
595  <<" > ILP="<<ILP<<" nothing is done!"<<G4endl;
596  }
597  else G4cout<<"*Warning*G4ChipsKaonMinusElasticXS::GetPTables: PDG="<<PDG
598  <<", Z="<<tgZ<<", N="<<tgN<<", i="<<ini<<">= max="<<nPoints<<", LP="<<LP
599  <<" > ILP="<<ILP<<", lPMax="<<lPMax<<" nothing is done!"<<G4endl;
600  }
601  }
602  else
603  {
604  // G4cout<<"*Error*G4ChipsKaonMinusElasticXS::GetPTables: PDG="<<PDG<<", Z="<<tgZ
605  // <<", N="<<tgN<<", while it is defined only for PDG=-321"<<G4endl;
606  // throw G4QException("G4ChipsKaonMinusElasticXS::GetPTables:onlyK-'s implemented");
608  ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
609  << ", while it is defined only for PDG=-321 (K-) " << G4endl;
610  G4Exception("G4ChipsKaonMinusElasticXS::GetPTables()", "HAD_CHPS_0000",
611  FatalException, ed);
612  }
613  return ILP;
614 }
TString fin
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
static const G4double a4
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4double GetTabValues(G4double lp, G4int pPDG, G4int tgZ, G4int tgN)
static const G4double a3
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
static const G4double a5
double G4double
Definition: G4Types.hh:76
static const G4double a2
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetQ2max()

G4double G4ChipsKaonMinusElasticXS::GetQ2max ( G4int  pPDG,
G4int  tgZ,
G4int  tgN,
G4double  pP 
)
private

Definition at line 840 of file G4ChipsKaonMinusElasticXS.cc.

842 {
843  G4double pP2=pP*pP; // squared momentum of the projectile
844  if(tgZ || tgN>-1) // ---> pipA
845  {
846  G4double mt=G4ParticleTable::GetParticleTable()->GetIonTable()->GetIon(tgZ,tgZ+tgN,0)->GetPDGMass()*.001; // Target mass in GeV
847 
848  G4double dmt=mt+mt;
849  G4double mds=dmt*std::sqrt(pP2+mK2)+mK2+mt*mt; // Mondelstam mds
850  return dmt*dmt*pP2/mds;
851  }
852  else
853  {
855  ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
856  << ", while it is defined only for p projectiles & Z_target>0" << G4endl;
857  G4Exception("G4ChipsKaonMinusElasticXS::GetQ2max()", "HAD_CHPS_0000",
858  FatalException, ed);
859  return 0;
860  }
861 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
G4IonTable * GetIonTable() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4ParticleTable * GetParticleTable()
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetSlope()

G4double G4ChipsKaonMinusElasticXS::GetSlope ( G4int  tZ,
G4int  tN,
G4int  pPDG 
)
private

Definition at line 722 of file G4ChipsKaonMinusElasticXS.cc.

723 {
724  if(onlyCS)G4cout<<"*Warning*G4ChipsKaonMinusElasticXS::GetSl:onlCS=true"<<G4endl;
725  if(lastLP<-4.3) return 0.; // S-wave for p<14 MeV/c (kinE<.1MeV)
726  if(PDG != -321)
727  {
728  // G4cout<<"*Error*G4ChipsKaonMinusElasticXS::GetSlope: PDG="<<PDG<<", Z="<<tgZ
729  // <<", N="<<tgN<<", while it is defined only for PDG=-321"<<G4endl;
730  // throw G4QException("G4ChipsKaonMinusElasticXS::GetSlope:Only K- is implemented");
732  ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
733  << ", while it is defined only for PDG=-321 (K-)" << G4endl;
734  }
735  if(theB1<0.) theB1=0.;
736  if(!(theB1>=-1.||theB1<=1.))G4cout<<"*NAN*G4QKaonMinusElCS::GetSlope:B1="<<theB1<<G4endl;
737  return theB1/GeVSQ;
738 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the caller graph for this function:

◆ GetTabValues()

G4double G4ChipsKaonMinusElasticXS::GetTabValues ( G4double  lp,
G4int  pPDG,
G4int  tgZ,
G4int  tgN 
)
private

Definition at line 747 of file G4ChipsKaonMinusElasticXS.cc.

749 {
750  if(PDG!=-321)G4cout<<"*Warning*G4ChipsKaonMinusElasticXS::GetTV:PDG="<<PDG<<G4endl;
751  if(tgZ<0 || tgZ>92)
752  {
753  G4cout<<"*Warning*G4QKaonMinusElasticCS::GetTabV:(1-92)NoIsotopes for Z="<<tgZ<<G4endl;
754  return 0.;
755  }
756  G4int iZ=tgZ-1; // Z index
757  if(iZ<0)
758  {
759  iZ=0; // conversion of the neutron target to the proton target
760  tgZ=1;
761  tgN=0;
762  }
763  G4double p=std::exp(lp); // momentum
764  G4double sp=std::sqrt(p); // sqrt(p)
765  G4double psp=p*sp; // p*sqrt(p)
766  G4double p2=p*p;
767  G4double p3=p2*p;
768  G4double p4=p3*p;
769  if ( tgZ == 1 && tgN == 0 ) // KaonMinus+P
770  {
771  G4double dl2=lp-lastPAR[12];
772  theSS=lastPAR[35];
773  theS1=(lastPAR[13]+lastPAR[14]*dl2*dl2)/(1.+lastPAR[15]/p4/p)+
774  (lastPAR[16]/p2+lastPAR[17]*p)/(p4+lastPAR[18]*sp);
775  theB1=lastPAR[19]*std::pow(p,lastPAR[20])/(1.+lastPAR[21]/p3);
776  theS2=lastPAR[22]+lastPAR[23]/(p4+lastPAR[24]*p);
777  theB2=lastPAR[25]+lastPAR[26]/(p4+lastPAR[27]/sp);
778  theS3=lastPAR[28]+lastPAR[29]/(p4*p4+lastPAR[30]*p2+lastPAR[31]);
779  theB3=lastPAR[32]+lastPAR[33]/(p4+lastPAR[34]);
780  theS4=0.;
781  theB4=0.;
782  // Returns the total elastic pim-p cross-section (to avoid spoiling lastSIG)
783  G4double dp=lp-lastPAR[2];
784  return lastPAR[0]/psp+(lastPAR[1]*dp*dp+lastPAR[3])/(1.-lastPAR[4]/sp+lastPAR[5]/p4)+
785  lastPAR[6]/(sqr(p-lastPAR[7])+lastPAR[8])+lastPAR[9]/(sqr(p-lastPAR[10])+lastPAR[11]);
786  }
787  else
788  {
789  G4double p5=p4*p;
790  G4double p6=p5*p;
791  G4double p8=p6*p2;
792  G4double p10=p8*p2;
793  G4double p12=p10*p2;
794  G4double p16=p8*p8;
795  //G4double p24=p16*p8;
796  G4double dl=lp-5.;
797  G4double a=tgZ+tgN;
798  G4double pah=std::pow(p,a/2);
799  G4double pa=pah*pah;
800  G4double pa2=pa*pa;
801  if(a<6.5)
802  {
803  theS1=lastPAR[9]/(1.+lastPAR[10]*p4*pa)+lastPAR[11]/(p4+lastPAR[12]*p4/pa2)+
804  (lastPAR[13]*dl*dl+lastPAR[14])/(1.+lastPAR[15]/p2);
805  theB1=(lastPAR[16]+lastPAR[17]*p2)/(p4+lastPAR[18]/pah)+lastPAR[19];
806  theSS=lastPAR[20]/(1.+lastPAR[21]/p2)+lastPAR[22]/(p6/pa+lastPAR[23]/p16);
807  theS2=lastPAR[24]/(pa/p2+lastPAR[25]/p4)+lastPAR[26];
808  theB2=lastPAR[27]*std::pow(p,lastPAR[28])+lastPAR[29]/(p8+lastPAR[30]/p16);
809  theS3=lastPAR[31]/(pa*p+lastPAR[32]/pa)+lastPAR[33];
810  theB3=lastPAR[34]/(p3+lastPAR[35]/p6)+lastPAR[36]/(1.+lastPAR[37]/p2);
811  theS4=p2*(pah*lastPAR[38]*std::exp(-pah*lastPAR[39])+
812  lastPAR[40]/(1.+lastPAR[41]*std::pow(p,lastPAR[42])));
813  theB4=lastPAR[43]*pa/p2/(1.+pa*lastPAR[44]);
814  }
815  else
816  {
817  theS1=lastPAR[9]/(1.+lastPAR[10]/p4)+lastPAR[11]/(p4+lastPAR[12]/p2)+
818  lastPAR[13]/(p5+lastPAR[14]/p16);
819  theB1=(lastPAR[15]/p8+lastPAR[19])/(p+lastPAR[16]/std::pow(p,lastPAR[20]))+
820  lastPAR[17]/(1.+lastPAR[18]/p4);
821  theSS=lastPAR[21]/(p4/std::pow(p,lastPAR[23])+lastPAR[22]/p4);
822  theS2=lastPAR[24]/p4/(std::pow(p,lastPAR[25])+lastPAR[26]/p12)+lastPAR[27];
823  theB2=lastPAR[28]/std::pow(p,lastPAR[29])+lastPAR[30]/std::pow(p,lastPAR[31]);
824  theS3=lastPAR[32]/std::pow(p,lastPAR[35])/(1.+lastPAR[36]/p12)+
825  lastPAR[33]/(1.+lastPAR[34]/p6);
826  theB3=lastPAR[37]/p8+lastPAR[38]/p2+lastPAR[39]/(1.+lastPAR[40]/p8);
827  theS4=(lastPAR[41]/p4+lastPAR[46]/p)/(1.+lastPAR[42]/p10)+
828  (lastPAR[43]+lastPAR[44]*dl*dl)/(1.+lastPAR[45]/p12);
829  theB4=lastPAR[47]/(1.+lastPAR[48]/p)+lastPAR[49]*p4/(1.+lastPAR[50]*p5);
830  }
831  // Returns the total elastic (n/p)A cross-section (to avoid spoiling lastSIG)
832  G4double dlp=lp-lastPAR[4]; // ax
833  // p1 p2 p3 p4
834  return (lastPAR[0]*dlp*dlp+lastPAR[1]+lastPAR[2]/p3)/(1.+lastPAR[3]/p2/sp);
835  }
836  return 0.;
837 } // End of GetTableValues
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 167 of file G4ChipsKaonMinusElasticXS.cc.

170 {
171  return true;
172 }
Here is the caller graph for this function:

Member Data Documentation

◆ B1T

std::vector<G4double*> G4ChipsKaonMinusElasticXS::B1T
private

Definition at line 132 of file G4ChipsKaonMinusElasticXS.hh.

◆ B2T

std::vector<G4double*> G4ChipsKaonMinusElasticXS::B2T
private

Definition at line 134 of file G4ChipsKaonMinusElasticXS.hh.

◆ B3T

std::vector<G4double*> G4ChipsKaonMinusElasticXS::B3T
private

Definition at line 136 of file G4ChipsKaonMinusElasticXS.hh.

◆ B4T

std::vector<G4double*> G4ChipsKaonMinusElasticXS::B4T
private

Definition at line 138 of file G4ChipsKaonMinusElasticXS.hh.

◆ colCS

std::vector<G4double> G4ChipsKaonMinusElasticXS::colCS
private

Definition at line 144 of file G4ChipsKaonMinusElasticXS.hh.

◆ colN

std::vector<G4int> G4ChipsKaonMinusElasticXS::colN
private

Definition at line 140 of file G4ChipsKaonMinusElasticXS.hh.

◆ colP

std::vector<G4double> G4ChipsKaonMinusElasticXS::colP
private

Definition at line 142 of file G4ChipsKaonMinusElasticXS.hh.

◆ colTH

std::vector<G4double> G4ChipsKaonMinusElasticXS::colTH
private

Definition at line 143 of file G4ChipsKaonMinusElasticXS.hh.

◆ colZ

std::vector<G4int> G4ChipsKaonMinusElasticXS::colZ
private

Definition at line 141 of file G4ChipsKaonMinusElasticXS.hh.

◆ CST

std::vector<G4double*> G4ChipsKaonMinusElasticXS::CST
private

Definition at line 129 of file G4ChipsKaonMinusElasticXS.hh.

◆ dlnP

G4double G4ChipsKaonMinusElasticXS::dlnP
private

Definition at line 91 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastB1T

G4double* G4ChipsKaonMinusElasticXS::lastB1T
private

Definition at line 120 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastB2T

G4double* G4ChipsKaonMinusElasticXS::lastB2T
private

Definition at line 122 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastB3T

G4double* G4ChipsKaonMinusElasticXS::lastB3T
private

Definition at line 124 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastB4T

G4double* G4ChipsKaonMinusElasticXS::lastB4T
private

Definition at line 126 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastCS

G4double G4ChipsKaonMinusElasticXS::lastCS
private

Definition at line 101 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastCST

G4double* G4ChipsKaonMinusElasticXS::lastCST
private

Definition at line 116 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastI

G4int G4ChipsKaonMinusElasticXS::lastI
private

Definition at line 102 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastLP

G4double G4ChipsKaonMinusElasticXS::lastLP
private

Definition at line 95 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastN

G4int G4ChipsKaonMinusElasticXS::lastN
private

Definition at line 97 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastP

G4double G4ChipsKaonMinusElasticXS::lastP
private

Definition at line 99 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastPAR

G4double* G4ChipsKaonMinusElasticXS::lastPAR
private

Definition at line 117 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastPIN

G4double G4ChipsKaonMinusElasticXS::lastPIN
private

Definition at line 115 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastS1T

G4double* G4ChipsKaonMinusElasticXS::lastS1T
private

Definition at line 119 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastS2T

G4double* G4ChipsKaonMinusElasticXS::lastS2T
private

Definition at line 121 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastS3T

G4double* G4ChipsKaonMinusElasticXS::lastS3T
private

Definition at line 123 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastS4T

G4double* G4ChipsKaonMinusElasticXS::lastS4T
private

Definition at line 125 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastSIG

G4double G4ChipsKaonMinusElasticXS::lastSIG
private

Definition at line 94 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastSST

G4double* G4ChipsKaonMinusElasticXS::lastSST
private

Definition at line 118 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastTH

G4double G4ChipsKaonMinusElasticXS::lastTH
private

Definition at line 100 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastTM

G4double G4ChipsKaonMinusElasticXS::lastTM
private

Definition at line 96 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastTN

G4int G4ChipsKaonMinusElasticXS::lastTN
private

Definition at line 114 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastTZ

G4int G4ChipsKaonMinusElasticXS::lastTZ
private

Definition at line 113 of file G4ChipsKaonMinusElasticXS.hh.

◆ lastZ

G4int G4ChipsKaonMinusElasticXS::lastZ
private

Definition at line 98 of file G4ChipsKaonMinusElasticXS.hh.

◆ lPMax

G4double G4ChipsKaonMinusElasticXS::lPMax
private

Definition at line 90 of file G4ChipsKaonMinusElasticXS.hh.

◆ lPMin

G4double G4ChipsKaonMinusElasticXS::lPMin
private

Definition at line 89 of file G4ChipsKaonMinusElasticXS.hh.

◆ nLast

const G4int G4ChipsKaonMinusElasticXS::nLast
private

Definition at line 88 of file G4ChipsKaonMinusElasticXS.hh.

◆ nPoints

const G4int G4ChipsKaonMinusElasticXS::nPoints
private

Definition at line 87 of file G4ChipsKaonMinusElasticXS.hh.

◆ onlyCS

G4bool G4ChipsKaonMinusElasticXS::onlyCS
private

Definition at line 93 of file G4ChipsKaonMinusElasticXS.hh.

◆ PAR

std::vector<G4double*> G4ChipsKaonMinusElasticXS::PAR
private

Definition at line 128 of file G4ChipsKaonMinusElasticXS.hh.

◆ PIN

std::vector<G4double> G4ChipsKaonMinusElasticXS::PIN
private

Definition at line 146 of file G4ChipsKaonMinusElasticXS.hh.

◆ S1T

std::vector<G4double*> G4ChipsKaonMinusElasticXS::S1T
private

Definition at line 131 of file G4ChipsKaonMinusElasticXS.hh.

◆ S2T

std::vector<G4double*> G4ChipsKaonMinusElasticXS::S2T
private

Definition at line 133 of file G4ChipsKaonMinusElasticXS.hh.

◆ S3T

std::vector<G4double*> G4ChipsKaonMinusElasticXS::S3T
private

Definition at line 135 of file G4ChipsKaonMinusElasticXS.hh.

◆ S4T

std::vector<G4double*> G4ChipsKaonMinusElasticXS::S4T
private

Definition at line 137 of file G4ChipsKaonMinusElasticXS.hh.

◆ SST

std::vector<G4double*> G4ChipsKaonMinusElasticXS::SST
private

Definition at line 130 of file G4ChipsKaonMinusElasticXS.hh.

◆ theB1

G4double G4ChipsKaonMinusElasticXS::theB1
private

Definition at line 105 of file G4ChipsKaonMinusElasticXS.hh.

◆ theB2

G4double G4ChipsKaonMinusElasticXS::theB2
private

Definition at line 107 of file G4ChipsKaonMinusElasticXS.hh.

◆ theB3

G4double G4ChipsKaonMinusElasticXS::theB3
private

Definition at line 109 of file G4ChipsKaonMinusElasticXS.hh.

◆ theB4

G4double G4ChipsKaonMinusElasticXS::theB4
private

Definition at line 111 of file G4ChipsKaonMinusElasticXS.hh.

◆ theS1

G4double G4ChipsKaonMinusElasticXS::theS1
private

Definition at line 104 of file G4ChipsKaonMinusElasticXS.hh.

◆ theS2

G4double G4ChipsKaonMinusElasticXS::theS2
private

Definition at line 106 of file G4ChipsKaonMinusElasticXS.hh.

◆ theS3

G4double G4ChipsKaonMinusElasticXS::theS3
private

Definition at line 108 of file G4ChipsKaonMinusElasticXS.hh.

◆ theS4

G4double G4ChipsKaonMinusElasticXS::theS4
private

Definition at line 110 of file G4ChipsKaonMinusElasticXS.hh.

◆ theSS

G4double G4ChipsKaonMinusElasticXS::theSS
private

Definition at line 103 of file G4ChipsKaonMinusElasticXS.hh.


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