Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4ChipsKaonPlusInelasticXS.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// The lust update: M.V. Kossov, CERN/ITEP(Moscow) 17-June-02
28// GEANT4 tag $Name: not supported by cvs2svn $
29//
30//
31// G4 Physics class: G4QKaonPlusNuclearCrossSection for gamma+A cross sections
32// Created: M.V. Kossov, CERN/ITEP(Moscow), 20-Dec-03
33// The last update: M.V. Kossov, CERN/ITEP (Moscow) 15-Feb-04
34//
35// --------------------------------------------------------------------------------
36// Short description: Cross-sections extracted from the CHIPS package for
37// kaon(minus)-nuclear interactions. Author: M. Kossov
38// -------------------------------------------------------------------------------------
39//
40
42#include "G4SystemOfUnits.hh"
43#include "G4DynamicParticle.hh"
45#include "G4KaonPlus.hh"
46#include "G4Proton.hh"
47#include "G4PionPlus.hh"
48
49// factory
51//
53
54
56{
57 // Initialization of the
58 lastLEN=0; // Pointer to the lastArray of LowEn CS
59 lastHEN=0; // Pointer to the lastArray of HighEn CS
60 lastN=0; // The last N of calculated nucleus
61 lastZ=0; // The last Z of calculated nucleus
62 lastP=0.; // Last used in cross section Momentum
63 lastTH=0.; // Last threshold momentum
64 lastCS=0.; // Last value of the Cross Section
65 lastI=0; // The last position in the DAMDB
66 LEN = new std::vector<G4double*>;
67 HEN = new std::vector<G4double*>;
68}
69
71{
72 G4int lens=LEN->size();
73 for(G4int i=0; i<lens; ++i) delete[] (*LEN)[i];
74 delete LEN;
75
76 G4int hens=HEN->size();
77 for(G4int i=0; i<hens; ++i) delete[] (*HEN)[i];
78 delete HEN;
79}
80
81
83 const G4Element*,
84 const G4Material*)
85{
86 G4ParticleDefinition* particle = Pt->GetDefinition();
87 if (particle == G4KaonPlus::KaonPlus() ) return true;
88 return false;
89}
90
91
92// The main member function giving the collision cross section (P is in IU, CS is in mb)
93// Make pMom in independent units ! (Now it is MeV)
95 const G4Isotope*,
96 const G4Element*,
97 const G4Material*)
98{
99 G4double pMom=Pt->GetTotalMomentum();
100 G4int tgN = A - tgZ;
101
102 return GetChipsCrossSection(pMom, tgZ, tgN, 321);
103}
104
106{
107 static G4int j; // A#0f Z/N-records already tested in AMDB
108 static std::vector <G4int> colN; // Vector of N for calculated nuclei (isotops)
109 static std::vector <G4int> colZ; // Vector of Z for calculated nuclei (isotops)
110 static std::vector <G4double> colP; // Vector of last momenta for the reaction
111 static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
112 static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
113 // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
114
115 G4bool in=false; // By default the isotope must be found in the AMDB
116 if(tgN!=lastN || tgZ!=lastZ) // The nucleus was not the last used isotope
117 {
118 in = false; // By default the isotope haven't be found in AMDB
119 lastP = 0.; // New momentum history (nothing to compare with)
120 lastN = tgN; // The last N of the calculated nucleus
121 lastZ = tgZ; // The last Z of the calculated nucleus
122 lastI = colN.size(); // Size of the Associative Memory DB in the heap
123 j = 0; // A#0f records found in DB for this projectile
124
125 if(lastI) for(G4int i=0; i<lastI; i++) // AMDB exists, try to find the (Z,N) isotope
126 {
127 if(colN[i]==tgN && colZ[i]==tgZ) // Try the record "i" in the AMDB
128 {
129 lastI=i; // Remember the index for future fast/last use
130 lastTH =colTH[i]; // The last THreshold (A-dependent)
131
132 if(pMom<=lastTH)
133 {
134 return 0.; // Energy is below the Threshold value
135 }
136 lastP =colP [i]; // Last Momentum (A-dependent)
137 lastCS =colCS[i]; // Last CrossSect (A-dependent)
138 in = true; // This is the case when the isotop is found in DB
139 // Momentum pMom is in IU ! @@ Units
140 lastCS=CalculateCrossSection(-1,j,321,lastZ,lastN,pMom); // read & update
141
142 if(lastCS<=0. && pMom>lastTH) // Correct the threshold (@@ No intermediate Zeros)
143 {
144 lastCS=0.;
145 lastTH=pMom;
146 }
147 break; // Go out of the LOOP
148 }
149 j++; // Increment a#0f records found in DB
150 }
151 if(!in) // This isotope has not been calculated previously
152 {
153 //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
154 lastCS=CalculateCrossSection(0,j,321,lastZ,lastN,pMom); //calculate & create
155
156 //if(lastCS>0.) // It means that the AMBD was initialized
157 //{
158
159 lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
160 colN.push_back(tgN);
161 colZ.push_back(tgZ);
162 colP.push_back(pMom);
163 colTH.push_back(lastTH);
164 colCS.push_back(lastCS);
165 //} // M.K. Presence of H1 with high threshold breaks the syncronization
166 return lastCS*millibarn;
167 } // End of creation of the new set of parameters
168 else
169 {
170 colP[lastI]=pMom;
171 colCS[lastI]=lastCS;
172 }
173 } // End of parameters udate
174 else if(pMom<=lastTH)
175 {
176 return 0.; // Momentum is below the Threshold Value -> CS=0
177 }
178 else // It is the last used -> use the current tables
179 {
180 lastCS=CalculateCrossSection(1,j,321,lastZ,lastN,pMom); // Only read and UpdateDB
181 lastP=pMom;
182 }
183 return lastCS*millibarn;
184}
185
186// The main member function giving the gamma-A cross section (E in GeV, CS in mb)
187G4double G4ChipsKaonPlusInelasticXS::CalculateCrossSection(G4int F, G4int I,
188 G4int, G4int targZ, G4int targN, G4double Momentum)
189{
190 static const G4double THmin=27.; // default minimum Momentum (MeV/c) Threshold
191 static const G4double THmiG=THmin*.001; // minimum Momentum (GeV/c) Threshold
192 static const G4double dP=10.; // step for the LEN (Low ENergy) table MeV/c
193 static const G4double dPG=dP*.001; // step for the LEN (Low ENergy) table GeV/c
194 static const G4int nL=105; // A#of LEN points in E (step 10 MeV/c)
195 static const G4double Pmin=THmin+(nL-1)*dP; // minP for the HighE part with safety
196 static const G4double Pmax=227000.; // maxP for the HEN (High ENergy) part 227 GeV
197 static const G4int nH=224; // A#of HEN points in lnE
198 static const G4double milP=std::log(Pmin);// Low logarithm energy for the HEN part
199 static const G4double malP=std::log(Pmax);// High logarithm energy (each 2.75 percent)
200 static const G4double dlP=(malP-milP)/(nH-1); // Step in log energy in the HEN part
201 static const G4double milPG=std::log(.001*Pmin);// Low logarithmEnergy for HEN part GeV/c
202
203 G4double sigma=0.;
204 if(F&&I) sigma=0.; // @@ *!* Fake line *!* to use F & I !!!Temporary!!!
205 G4double A=targN+targZ; // A of the target
206
207 if(F<=0) // This isotope was not the last used isotop
208 {
209 if(F<0) // This isotope was found in DAMDB =-----=> RETRIEVE
210 {
211 G4int sync=LEN->size();
212 if(sync<=I) G4cerr<<"*!*G4ChipsKPlusNuclCS::CalcCrosSect:Sync="<<sync<<"<="<<I<<G4endl;
213 lastLEN=(*LEN)[I]; // Pointer to prepared LowEnergy cross sections
214 lastHEN=(*HEN)[I]; // Pointer to prepared High Energy cross sections
215 }
216 else // This isotope wasn't calculated before => CREATE
217 {
218 lastLEN = new G4double[nL]; // Allocate memory for the new LEN cross sections
219 lastHEN = new G4double[nH]; // Allocate memory for the new HEN cross sections
220 // --- Instead of making a separate function ---
221 G4double P=THmiG; // Table threshold in GeV/c
222 for(G4int k=0; k<nL; k++)
223 {
224 lastLEN[k] = CrossSectionLin(targZ, targN, P);
225 P+=dPG;
226 }
227 G4double lP=milPG;
228 for(G4int n=0; n<nH; n++)
229 {
230 lastHEN[n] = CrossSectionLog(targZ, targN, lP);
231 lP+=dlP;
232 }
233 // --- End of possible separate function
234 // *** The synchronization check ***
235 G4int sync=LEN->size();
236 if(sync!=I)
237 {
238 G4cerr<<"***G4ChipsKPlusNuclCS::CalcCrossSect: Sinc="<<sync<<"#"<<I<<", Z=" <<targZ
239 <<", N="<<targN<<", F="<<F<<G4endl;
240 //G4Exception("G4PiMinusNuclearCS::CalculateCS:","39",FatalException,"DBoverflow");
241 }
242 LEN->push_back(lastLEN); // remember the Low Energy Table
243 HEN->push_back(lastHEN); // remember the High Energy Table
244 } // End of creation of the new set of parameters
245 } // End of parameters udate
246 // =--------------------------= NOW the Magic Formula =---------------------------------=
247
248 if (Momentum<lastTH) return 0.; // It must be already checked in the interface class
249 else if (Momentum<Pmin) // Low Energy region
250 {
251 if(A<=1. && Momentum < 600.) sigma=0.; // Approximation tot/el uncertainty
252 else sigma=EquLinearFit(Momentum,nL,THmin,dP,lastLEN);
253 }
254 else if (Momentum<Pmax) // High Energy region
255 {
256 G4double lP=std::log(Momentum);
257 sigma=EquLinearFit(lP,nH,milP,dlP,lastHEN);
258 }
259 else // UHE region (calculation, not frequent)
260 {
261 G4double P=0.001*Momentum; // Approximation formula is for P in GeV/c
262 sigma=CrossSectionFormula(targZ, targN, P, std::log(P));
263 }
264 if(sigma<0.) return 0.;
265 return sigma;
266}
267
268// Electromagnetic momentum-threshold (in MeV/c)
269G4double G4ChipsKaonPlusInelasticXS::ThresholdMomentum(G4int tZ, G4int tN)
270{
271 static const G4double third=1./3.;
272 static const G4double prM = G4Proton::Proton()->GetPDGMass(); // Proton mass in MeV
273 static const G4double piM = G4PionPlus::PionPlus()->GetPDGMass()+.1; // Pion mass in MeV+Safety (WP)??
274 static const G4double pM = G4KaonPlus::KaonPlus()->GetPDGMass(); // Projectile mass in MeV
275 static const G4double tpM= pM+pM; // Doubled projectile mass (MeV)
276 G4double tA=tZ+tN;
277 if(tZ<.99 || tN<0.) return 0.;
278 G4double tM=931.5*tA;
279 G4double dE=piM; // At least one Pi0 must be created
280 if(tZ==1 && tN==0) tM=prM; // A threshold on the free proton
281 else dE=tZ/(1.+std::pow(tA,third)); // Safety for diffused edge of the nucleus (QE)
282 //G4double dE=1.263*tZ/(1.+std::pow(tA,third));
283 G4double T=dE+dE*(dE/2+pM)/tM;
284 return std::sqrt(T*(tpM+T));
285}
286
287// Calculation formula for piMinus-nuclear inelastic cross-section (mb) (P in GeV/c)
288G4double G4ChipsKaonPlusInelasticXS::CrossSectionLin(G4int tZ, G4int tN, G4double P)
289{
290 G4double lP=std::log(P);
291 return CrossSectionFormula(tZ, tN, P, lP);
292}
293
294// Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
295G4double G4ChipsKaonPlusInelasticXS::CrossSectionLog(G4int tZ, G4int tN, G4double lP)
296{
297 G4double P=std::exp(lP);
298 return CrossSectionFormula(tZ, tN, P, lP);
299}
300// Calculation formula for piMinus-nuclear inelastic cross-section (mb) log(P in GeV/c)
301G4double G4ChipsKaonPlusInelasticXS::CrossSectionFormula(G4int tZ, G4int tN,
302 G4double P, G4double lP)
303{
304 G4double sigma=0.;
305 if(tZ==1 && !tN) // KPlus-Proton interaction from G4QuasiElRatios
306 {
307 G4double ld=lP-3.5;
308 G4double ld2=ld*ld;
309 G4double sp=std::sqrt(P);
310 G4double p2=P*P;
311 G4double p4=p2*p2;
312 G4double lm=P-1.;
313 G4double md=lm*lm+.372;
314 G4double El=(.0557*ld2+2.23)/(1.-.7/sp+.1/p4);
315 G4double To=(.3*ld2+19.5)/(1.+.46/sp+1.6/p4);
316 sigma=(To-El)+.6/md;
317 }
318 else if(tZ<97 && tN<152) // General solution
319 {
320 G4double p2=P*P;
321 G4double p4=p2*p2;
322 G4double a=tN+tZ; // A of the target
323 G4double al=std::log(a);
324 G4double sa=std::sqrt(a);
325 G4double asa=a*sa;
326 G4double a2=a*a;
327 G4double a3=a2*a;
328 G4double a4=a2*a2;
329 G4double a8=a4*a4;
330 G4double a12=a8*a4;
331 G4double f=.6; // Default values for deutrons
332 G4double r=.5;
333 G4double gg=3.7;
334 G4double c=36.;
335 G4double ss=3.5;
336 G4double t=3.;
337 G4double u=.44;
338 G4double v=5.E-9;
339 if(tZ>1 && tN>1) // More than deuteron
340 {
341 f=1.;
342 r=1./(1.+.007*a2);
343 gg=4.2;
344 c=52.*std::exp(al*.6)*(1.+95./a2)/(1.+9./a)/(1.+46./a2);
345 ss=(40.+.14*a)/(1.+12./a);
346 G4double y=std::exp(al*1.7);
347 t=.185*y/(1.+.00012*y);
348 u=(1.+80./asa)/(1.+200./asa);
349 v=(1.+3.E-6*a4*(1.+6.E-7*a3+4.E10/a12))/a3/20000.;
350 }
351 G4double d=lP-gg;
352 G4double w=P-1.;
353 G4double rD=ss/(w*w+.36);
354 G4double h=P-.44;
355 G4double rR=t/(h*h+u*u);
356 sigma=(f*d*d+c)/(1.+r/std::sqrt(P)+1./p4)+(rD+rR)/(1+v/p4/p4);
357 }
358 else
359 {
360 G4cerr<<"-Warning-G4ChipsKaonPlusNuclearCroSect::CSForm:Bad A, Z="<<tZ<<", N="<<tN<<G4endl;
361 sigma=0.;
362 }
363 if(sigma<0.) return 0.;
364 return sigma;
365}
366
367G4double G4ChipsKaonPlusInelasticXS::EquLinearFit(G4double X, G4int N, G4double X0, G4double DX, G4double* Y)
368{
369 if(DX<=0. || N<2)
370 {
371 G4cerr<<"***G4ChipsKaonPlusInelasticXS::EquLinearFit: DX="<<DX<<", N="<<N<<G4endl;
372 return Y[0];
373 }
374
375 G4int N2=N-2;
376 G4double d=(X-X0)/DX;
377 G4int j=static_cast<int>(d);
378 if (j<0) j=0;
379 else if(j>N2) j=N2;
380 d-=j; // excess
381 G4double yi=Y[j];
382 G4double sigma=yi+(Y[j+1]-yi)*d;
383
384 return sigma;
385}
#define G4_DECLARE_XS_FACTORY(cross_section)
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
G4ParticleDefinition * GetDefinition() const
G4double GetTotalMomentum() const
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Proton * Proton()
Definition: G4Proton.cc:93