Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChipsAntiBaryonElasticXS.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// $Id$
28//
29//
30// G4 Physics class: G4ChipsAntiBaryonElasticXS for pA elastic cross sections
31// Created: M.V. Kossov, CERN/ITEP(Moscow), 5-Feb-2010
32// The last update: M.V. Kossov, CERN/ITEP (Moscow) 5-Feb-2010
33//
34//
35// -------------------------------------------------------------------------------
36// Short description: Interaction cross-sections for the elastic process.
37// Class extracted from CHIPS and integrated in Geant4 by W.Pokorski
38// -------------------------------------------------------------------------------
39
40
42#include "G4SystemOfUnits.hh"
43#include "G4DynamicParticle.hh"
45#include "G4AntiProton.hh"
46#include "G4Nucleus.hh"
47#include "G4ParticleTable.hh"
48#include "G4NucleiProperties.hh"
49
50// factory
52//
54
56{
57 lPMin=-8.; //Min tabulatedLogarithmMomentum(D)
58 lPMax= 8.; //Max tabulatedLogarithmMomentum(D)
59 dlnP=(lPMax-lPMin)/nLast;// LogStep inTable (D)
60 onlyCS=true;//Flag toCalculOnlyCS(not Si/Bi)(L)
61 lastSIG=0.; //Last calculated cross section (L)
62 lastLP=-10.;//LastLog(mom_of IncidentHadron)(L)
63 lastTM=0.; //Last t_maximum (L)
64 theSS=0.; //TheLastSqSlope of 1st difr.Max(L)
65 theS1=0.; //TheLastMantissa of 1st difrMax(L)
66 theB1=0.; //TheLastSlope of 1st difructMax(L)
67 theS2=0.; //TheLastMantissa of 2nd difrMax(L)
68 theB2=0.; //TheLastSlope of 2nd difructMax(L)
69 theS3=0.; //TheLastMantissa of 3d difr.Max(L)
70 theB3=0.; //TheLastSlope of 3d difruct.Max(L)
71 theS4=0.; //TheLastMantissa of 4th difrMax(L)
72 theB4=0.; //TheLastSlope of 4th difructMax(L)
73 lastTZ=0; // Last atomic number of the target
74 lastTN=0; // Last # of neutrons in the target
75 lastPIN=0.; // Last initialized max momentum
76 lastCST=0; // Elastic cross-section table
77 lastPAR=0; // ParametersForFunctionCalculation
78 lastSST=0; // E-dep ofSqardSlope of 1st difMax
79 lastS1T=0; // E-dep of mantissa of 1st dif.Max
80 lastB1T=0; // E-dep of the slope of 1st difMax
81 lastS2T=0; // E-dep of mantissa of 2nd difrMax
82 lastB2T=0; // E-dep of the slope of 2nd difMax
83 lastS3T=0; // E-dep of mantissa of 3d difr.Max
84 lastB3T=0; // E-dep of the slope of 3d difrMax
85 lastS4T=0; // E-dep of mantissa of 4th difrMax
86 lastB4T=0; // E-dep of the slope of 4th difMax
87 lastN=0; // The last N of calculated nucleus
88 lastZ=0; // The last Z of calculated nucleus
89 lastP=0.; // LastUsed inCrossSection Momentum
90 lastTH=0.; // Last threshold momentum
91 lastCS=0.; // Last value of the Cross Section
92 lastI=0; // The last position in the DAMDB
93}
94
96{
97 std::vector<G4double*>::iterator pos;
98 for (pos=CST.begin(); pos<CST.end(); pos++)
99 { delete [] *pos; }
100 CST.clear();
101 for (pos=PAR.begin(); pos<PAR.end(); pos++)
102 { delete [] *pos; }
103 PAR.clear();
104 for (pos=SST.begin(); pos<SST.end(); pos++)
105 { delete [] *pos; }
106 SST.clear();
107 for (pos=S1T.begin(); pos<S1T.end(); pos++)
108 { delete [] *pos; }
109 S1T.clear();
110 for (pos=B1T.begin(); pos<B1T.end(); pos++)
111 { delete [] *pos; }
112 B1T.clear();
113 for (pos=S2T.begin(); pos<S2T.end(); pos++)
114 { delete [] *pos; }
115 S2T.clear();
116 for (pos=B2T.begin(); pos<B2T.end(); pos++)
117 { delete [] *pos; }
118 B2T.clear();
119 for (pos=S3T.begin(); pos<S3T.end(); pos++)
120 { delete [] *pos; }
121 S3T.clear();
122 for (pos=B3T.begin(); pos<B3T.end(); pos++)
123 { delete [] *pos; }
124 B3T.clear();
125 for (pos=S4T.begin(); pos<S4T.end(); pos++)
126 { delete [] *pos; }
127 S4T.clear();
128 for (pos=B4T.begin(); pos<B4T.end(); pos++)
129 { delete [] *pos; }
130 B4T.clear();
131}
132
133
135 const G4Element*,
136 const G4Material*)
137{
138 G4ParticleDefinition* particle = Pt->GetDefinition();
139
140 if(particle == G4AntiNeutron::AntiNeutron())
141 {
142 return true;
143 }
144 else if(particle == G4AntiProton::AntiProton())
145 {
146 return true;
147 }
148 else if(particle == G4AntiLambda::AntiLambda())
149 {
150 return true;
151 }
152 else if(particle == G4AntiSigmaPlus::AntiSigmaPlus())
153 {
154 return true;
155 }
156 else if(particle == G4AntiSigmaMinus::AntiSigmaMinus())
157 {
158 return true;
159 }
160 else if(particle == G4AntiSigmaZero::AntiSigmaZero())
161 {
162 return true;
163 }
164 else if(particle == G4AntiXiMinus::AntiXiMinus())
165 {
166 return true;
167 }
168 else if(particle == G4AntiXiZero::AntiXiZero())
169 {
170 return true;
171 }
172 else if(particle == G4AntiOmegaMinus::AntiOmegaMinus())
173 {
174 return true;
175 }
176 return false;
177}
178
179// The main member function giving the collision cross section (P is in IU, CS is in mb)
180// Make pMom in independent units ! (Now it is MeV)
182 const G4Isotope*,
183 const G4Element*,
184 const G4Material*)
185{
186 G4double pMom=Pt->GetTotalMomentum();
187 G4int tgN = A - tgZ;
188 G4int pdg = Pt->GetDefinition()->GetPDGEncoding();
189
190 return GetChipsCrossSection(pMom, tgZ, tgN, pdg);
191}
192
194{
195 static std::vector <G4int> colN; // Vector of N for calculated nuclei (isotops)
196 static std::vector <G4int> colZ; // Vector of Z for calculated nuclei (isotops)
197 static std::vector <G4double> colP; // Vector of last momenta for the reaction
198 static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
199 static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
200 // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
201
202 G4bool fCS = false;
203
204 G4double pEn=pMom;
205 onlyCS=fCS;
206
207 G4bool in=false; // By default the isotope must be found in the AMDB
208 lastP = 0.; // New momentum history (nothing to compare with)
209 lastN = tgN; // The last N of the calculated nucleus
210 lastZ = tgZ; // The last Z of the calculated nucleus
211 lastI = colN.size(); // Size of the Associative Memory DB in the heap
212 if(lastI) for(G4int i=0; i<lastI; i++) // Loop over proj/tgZ/tgN lines of DB
213 { // The nucleus with projPDG is found in AMDB
214 if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
215 {
216 lastI=i;
217 lastTH =colTH[i]; // Last THreshold (A-dependent)
218 if(pEn<=lastTH)
219 {
220 return 0.; // Energy is below the Threshold value
221 }
222 lastP =colP [i]; // Last Momentum (A-dependent)
223 lastCS =colCS[i]; // Last CrossSect (A-dependent)
224 // if(std::fabs(lastP/pMom-1.)<tolerance) //VI (do not use tolerance)
225 if(lastP == pMom) // Do not recalculate
226 {
227 CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // Update param's only
228 return lastCS*millibarn; // Use theLastCS
229 }
230 in = true; // This is the case when the isotop is found in DB
231 // Momentum pMom is in IU ! @@ Units
232 lastCS=CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // read & update
233 if(lastCS<=0. && pEn>lastTH) // Correct the threshold
234 {
235 lastTH=pEn;
236 }
237 break; // Go out of the LOOP with found lastI
238 }
239 } // End of attampt to find the nucleus in DB
240 if(!in) // This nucleus has not been calculated previously
241 {
242 //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
243 lastCS=CalculateCrossSection(fCS,0,lastI,pPDG,lastZ,lastN,pMom);//calculate&create
244 if(lastCS<=0.)
245 {
246 lastTH = 0; // ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
247 if(pEn>lastTH)
248 {
249 lastTH=pEn;
250 }
251 }
252 colN.push_back(tgN);
253 colZ.push_back(tgZ);
254 colP.push_back(pMom);
255 colTH.push_back(lastTH);
256 colCS.push_back(lastCS);
257 return lastCS*millibarn;
258 } // End of creation of the new set of parameters
259 else
260 {
261 colP[lastI]=pMom;
262 colCS[lastI]=lastCS;
263 }
264 return lastCS*millibarn;
265}
266
267// Calculation of total elastic cross section (p in IU, CS in mb) @@ Units (?)
268// F=0 - create AMDB, F=-1 - read&update AMDB, F=1 - update AMDB (sinchro with higher AMDB)
269G4double G4ChipsAntiBaryonElasticXS::CalculateCrossSection(G4bool CS,G4int F,G4int I,
270 G4int PDG, G4int tgZ, G4int tgN, G4double pIU)
271{
272 // *** Begin of Associative Memory DB for acceleration of the cross section calculations
273 static std::vector <G4double> PIN; // Vector of max initialized log(P) in the table
274 // *** End of Static Definitions (Associative Memory Data Base) ***
275 G4double pMom=pIU/GeV; // All calculations are in GeV
276 onlyCS=CS; // Flag to calculate only CS (not Si/Bi)
277 lastLP=std::log(pMom); // Make a logarithm of the momentum for calculation
278 if(F) // This isotope was found in AMDB =>RETRIEVE/UPDATE
279 {
280 if(F<0) // the AMDB must be loded
281 {
282 lastPIN = PIN[I]; // Max log(P) initialised for this table set
283 lastPAR = PAR[I]; // Pointer to the parameter set
284 lastCST = CST[I]; // Pointer to the total sross-section table
285 lastSST = SST[I]; // Pointer to the first squared slope
286 lastS1T = S1T[I]; // Pointer to the first mantissa
287 lastB1T = B1T[I]; // Pointer to the first slope
288 lastS2T = S2T[I]; // Pointer to the second mantissa
289 lastB2T = B2T[I]; // Pointer to the second slope
290 lastS3T = S3T[I]; // Pointer to the third mantissa
291 lastB3T = B3T[I]; // Pointer to the rhird slope
292 lastS4T = S4T[I]; // Pointer to the 4-th mantissa
293 lastB4T = B4T[I]; // Pointer to the 4-th slope
294 }
295 if(lastLP>lastPIN && lastLP<lPMax)
296 {
297 lastPIN=GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);// Can update upper logP-Limit in tabs
298 PIN[I]=lastPIN; // Remember the new P-Limit of the tables
299 }
300 }
301 else // This isotope wasn't initialized => CREATE
302 {
303 lastPAR = new G4double[nPoints]; // Allocate memory for parameters of CS function
304 lastPAR[nLast]=0; // Initialization for VALGRIND
305 lastCST = new G4double[nPoints]; // Allocate memory for Tabulated CS function
306 lastSST = new G4double[nPoints]; // Allocate memory for Tabulated first sqaredSlope
307 lastS1T = new G4double[nPoints]; // Allocate memory for Tabulated first mantissa
308 lastB1T = new G4double[nPoints]; // Allocate memory for Tabulated first slope
309 lastS2T = new G4double[nPoints]; // Allocate memory for Tabulated second mantissa
310 lastB2T = new G4double[nPoints]; // Allocate memory for Tabulated second slope
311 lastS3T = new G4double[nPoints]; // Allocate memory for Tabulated third mantissa
312 lastB3T = new G4double[nPoints]; // Allocate memory for Tabulated third slope
313 lastS4T = new G4double[nPoints]; // Allocate memory for Tabulated 4-th mantissa
314 lastB4T = new G4double[nPoints]; // Allocate memory for Tabulated 4-th slope
315 lastPIN = GetPTables(lastLP,lPMin,PDG,tgZ,tgN); // Returns the new P-limit for tables
316 PIN.push_back(lastPIN); // Fill parameters of CS function to AMDB
317 PAR.push_back(lastPAR); // Fill parameters of CS function to AMDB
318 CST.push_back(lastCST); // Fill Tabulated CS function to AMDB
319 SST.push_back(lastSST); // Fill Tabulated first sq.slope to AMDB
320 S1T.push_back(lastS1T); // Fill Tabulated first mantissa to AMDB
321 B1T.push_back(lastB1T); // Fill Tabulated first slope to AMDB
322 S2T.push_back(lastS2T); // Fill Tabulated second mantissa to AMDB
323 B2T.push_back(lastB2T); // Fill Tabulated second slope to AMDB
324 S3T.push_back(lastS3T); // Fill Tabulated third mantissa to AMDB
325 B3T.push_back(lastB3T); // Fill Tabulated third slope to AMDB
326 S4T.push_back(lastS4T); // Fill Tabulated 4-th mantissa to AMDB
327 B4T.push_back(lastB4T); // Fill Tabulated 4-th slope to AMDB
328 } // End of creation/update of the new set of parameters and tables
329 // =---------= NOW Update (if necessary) and Calculate the Cross Section =-----------=
330 if(lastLP>lastPIN && lastLP<lPMax)
331 {
332 lastPIN = GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);
333 }
334 if(!onlyCS) lastTM=GetQ2max(PDG, tgZ, tgN, pMom); // Calculate (-t)_max=Q2_max (GeV2)
335 if(lastLP>lPMin && lastLP<=lastPIN) // Linear fit is made using precalculated tables
336 {
337 if(lastLP==lastPIN)
338 {
339 G4double shift=(lastLP-lPMin)/dlnP+.000001; // Log distance from lPMin
340 G4int blast=static_cast<int>(shift); // this is a bin number of the lower edge (0)
341 if(blast<0 || blast>=nLast) G4cout<<"G4QaBarElCS::CCS:b="<<blast<<","<<nLast<<G4endl;
342 lastSIG = lastCST[blast];
343 if(!onlyCS) // Skip the differential cross-section parameters
344 {
345 theSS = lastSST[blast];
346 theS1 = lastS1T[blast];
347 theB1 = lastB1T[blast];
348 theS2 = lastS2T[blast];
349 theB2 = lastB2T[blast];
350 theS3 = lastS3T[blast];
351 theB3 = lastB3T[blast];
352 theS4 = lastS4T[blast];
353 theB4 = lastB4T[blast];
354 }
355 }
356 else
357 {
358 G4double shift=(lastLP-lPMin)/dlnP; // a shift from the beginning of the table
359 G4int blast=static_cast<int>(shift); // the lower bin number
360 if(blast<0) blast=0;
361 if(blast>=nLast) blast=nLast-1; // low edge of the last bin
362 shift-=blast; // step inside the unit bin
363 G4int lastL=blast+1; // the upper bin number
364 G4double SIGL=lastCST[blast]; // the basic value of the cross-section
365 lastSIG= SIGL+shift*(lastCST[lastL]-SIGL); // calculated total elastic cross-section
366 if(!onlyCS) // Skip the differential cross-section parameters
367 {
368 G4double SSTL=lastSST[blast]; // the low bin of the first squared slope
369 theSS=SSTL+shift*(lastSST[lastL]-SSTL); // the basic value of the first sq.slope
370 G4double S1TL=lastS1T[blast]; // the low bin of the first mantissa
371 theS1=S1TL+shift*(lastS1T[lastL]-S1TL); // the basic value of the first mantissa
372 G4double B1TL=lastB1T[blast]; // the low bin of the first slope
373 theB1=B1TL+shift*(lastB1T[lastL]-B1TL); // the basic value of the first slope
374 G4double S2TL=lastS2T[blast]; // the low bin of the second mantissa
375 theS2=S2TL+shift*(lastS2T[lastL]-S2TL); // the basic value of the second mantissa
376 G4double B2TL=lastB2T[blast]; // the low bin of the second slope
377 theB2=B2TL+shift*(lastB2T[lastL]-B2TL); // the basic value of the second slope
378 G4double S3TL=lastS3T[blast]; // the low bin of the third mantissa
379 theS3=S3TL+shift*(lastS3T[lastL]-S3TL); // the basic value of the third mantissa
380 G4double B3TL=lastB3T[blast]; // the low bin of the third slope
381 theB3=B3TL+shift*(lastB3T[lastL]-B3TL); // the basic value of the third slope
382 G4double S4TL=lastS4T[blast]; // the low bin of the 4-th mantissa
383 theS4=S4TL+shift*(lastS4T[lastL]-S4TL); // the basic value of the 4-th mantissa
384 G4double B4TL=lastB4T[blast]; // the low bin of the 4-th slope
385 theB4=B4TL+shift*(lastB4T[lastL]-B4TL); // the basic value of the 4-th slope
386 }
387 }
388 }
389 else lastSIG=GetTabValues(lastLP, PDG, tgZ, tgN); // Direct calculation beyond the table
390 if(lastSIG<0.) lastSIG = 0.; // @@ a Warning print can be added
391 return lastSIG;
392}
393
394// It has parameter sets for all tZ/tN/PDG, using them the tables can be created/updated
395G4double G4ChipsAntiBaryonElasticXS::GetPTables(G4double LP, G4double ILP, G4int PDG,
396 G4int tgZ, G4int tgN)
397{
398 // @@ At present all nA==pA ---------> Each neucleus can have not more than 51 parameters
399 static const G4double pwd=2727;
400 const G4int n_appel=30; // #of parameters for app-elastic (<nPoints=128)
401 // -0- -1- -2- -3- -4- -5- -6- -7- -8--9--10--11--12--13--14-
402 G4double app_el[n_appel]={1.25,3.5,80.,1.,.0557,6.72,5.,74.,3.,3.4,.2,.17,.001,8.,.055,
403 3.64,5.e-5,4000.,1500.,.46,1.2e6,3.5e6,5.e-5,1.e10,8.5e8,
404 1.e10,1.1,3.4e6,6.8e6,0.};
405 // -15- -16- -17- -18- -19- -20- -21- -22- -23- -24-
406 // -25- -26- -27- -28- -29-
407 if(PDG>-3334 && PDG<-1111)
408 {
409 // -- Total pp elastic cross section cs & s1/b1 (main), s2/b2 (tail1), s3/b3 (tail2) --
410 //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|
411 //CS=2.865/p2s/(1+.0022/p2s)+(18.9+.6461*dl1*dl1+9./p)/(1.+.425*lp)/(1.+.4276/p4);
412 // par(0) par(7) par(1) par(2) par(4) par(5) par(6)
413 //dl2=lp-5., s1=(74.+3.*dl2*dl2)/(1+3.4/p4/p)+(.2/p2+17.*p)/(p4+.001*sp),
414 // par(8) par(9) par(10) par(11) par(12)par(13) par(14)
415 // b1=8.*p**.055/(1.+3.64/p3); s2=5.e-5+4000./(p4+1500.*p); b2=.46+1.2e6/(p4+3.5e6/sp);
416 // par(15) par(16) par(17) par(18) par(19) par(20) par(21) par(22) par(23)
417 // s3=5.e-5+1.e10/(p4*p4+8.5e8*p2+1.e10); b3=1.1+3.4e6/(p4+6.8e6); ss=0.
418 // par(24) par(25) par(26) par(27) par(28) par(29) par(30) par(31)
419 //
420 if(lastPAR[nLast]!=pwd) // A unique flag to avoid the repeatable definition
421 {
422 if ( tgZ == 1 && tgN == 0 )
423 {
424 for (G4int ip=0; ip<n_appel; ip++) lastPAR[ip]=app_el[ip]; // PiMinus+P
425 }
426 else
427 {
428 G4double a=tgZ+tgN;
429 G4double sa=std::sqrt(a);
430 G4double ssa=std::sqrt(sa);
431 G4double asa=a*sa;
432 G4double a2=a*a;
433 G4double a3=a2*a;
434 G4double a4=a3*a;
435 G4double a5=a4*a;
436 G4double a6=a4*a2;
437 G4double a7=a6*a;
438 G4double a8=a7*a;
439 G4double a9=a8*a;
440 G4double a10=a5*a5;
441 G4double a12=a6*a6;
442 G4double a14=a7*a7;
443 G4double a16=a8*a8;
444 G4double a17=a16*a;
445 //G4double a20=a16*a4;
446 G4double a32=a16*a16;
447 // Reaction cross-section parameters (pel=peh_fit.f)
448 lastPAR[0]=.23*asa/(1.+a*.15); // p1
449 lastPAR[1]=2.8*asa/(1.+a*(.015+.05/ssa)); // p2
450 lastPAR[2]=15.*a/(1.+.005*a2); // p3
451 lastPAR[3]=.013*a2/(1.+a3*(.006+a*.00001)); // p4
452 lastPAR[4]=5.; // p5
453 lastPAR[5]=0.; // p6 not used
454 lastPAR[6]=0.; // p7 not used
455 lastPAR[7]=0.; // p8 not used
456 lastPAR[8]=0.; // p9 not used
457 // @@ the differential cross-section is parameterized separately for A>6 & A<7
458 if(a<6.5)
459 {
460 G4double a28=a16*a12;
461 // The main pre-exponent (pel_sg)
462 lastPAR[ 9]=4000*a; // p1
463 lastPAR[10]=1.2e7*a8+380*a17; // p2
464 lastPAR[11]=.7/(1.+4.e-12*a16); // p3
465 lastPAR[12]=2.5/a8/(a4+1.e-16*a32); // p4
466 lastPAR[13]=.28*a; // p5
467 lastPAR[14]=1.2*a2+2.3; // p6
468 lastPAR[15]=3.8/a; // p7
469 // The main slope (pel_sl)
470 lastPAR[16]=.01/(1.+.0024*a5); // p1
471 lastPAR[17]=.2*a; // p2
472 lastPAR[18]=9.e-7/(1.+.035*a5); // p3
473 lastPAR[19]=(42.+2.7e-11*a16)/(1.+.14*a); // p4
474 // The main quadratic (pel_sh)
475 lastPAR[20]=2.25*a3; // p1
476 lastPAR[21]=18.; // p2
477 lastPAR[22]=2.4e-3*a8/(1.+2.6e-4*a7); // p3
478 lastPAR[23]=3.5e-36*a32*a8/(1.+5.e-15*a32/a); // p4
479 // The 1st max pre-exponent (pel_qq)
480 lastPAR[24]=1.e5/(a8+2.5e12/a16); // p1
481 lastPAR[25]=8.e7/(a12+1.e-27*a28*a28); // p2
482 lastPAR[26]=.0006*a3; // p3
483 // The 1st max slope (pel_qs)
484 lastPAR[27]=10.+4.e-8*a12*a; // p1
485 lastPAR[28]=.114; // p2
486 lastPAR[29]=.003; // p3
487 lastPAR[30]=2.e-23; // p4
488 // The effective pre-exponent (pel_ss)
489 lastPAR[31]=1./(1.+.0001*a8); // p1
490 lastPAR[32]=1.5e-4/(1.+5.e-6*a12); // p2
491 lastPAR[33]=.03; // p3
492 // The effective slope (pel_sb)
493 lastPAR[34]=a/2; // p1
494 lastPAR[35]=2.e-7*a4; // p2
495 lastPAR[36]=4.; // p3
496 lastPAR[37]=64./a3; // p4
497 // The gloria pre-exponent (pel_us)
498 lastPAR[38]=1.e8*std::exp(.32*asa); // p1
499 lastPAR[39]=20.*std::exp(.45*asa); // p2
500 lastPAR[40]=7.e3+2.4e6/a5; // p3
501 lastPAR[41]=2.5e5*std::exp(.085*a3); // p4
502 lastPAR[42]=2.5*a; // p5
503 // The gloria slope (pel_ub)
504 lastPAR[43]=920.+.03*a8*a3; // p1
505 lastPAR[44]=93.+.0023*a12; // p2
506 }
507 else // A > Li6 (li7, ...)
508 {
509 G4double p1a10=2.2e-28*a10;
510 G4double r4a16=6.e14/a16;
511 G4double s4a16=r4a16*r4a16;
512 // a24
513 // a36
514 // The main pre-exponent (peh_sg)
515 lastPAR[ 9]=4.5*std::pow(a,1.15); // p1
516 lastPAR[10]=.06*std::pow(a,.6); // p2
517 lastPAR[11]=.6*a/(1.+2.e15/a16); // p3
518 lastPAR[12]=.17/(a+9.e5/a3+1.5e33/a32); // p4
519 lastPAR[13]=(.001+7.e-11*a5)/(1.+4.4e-11*a5); // p5
520 lastPAR[14]=(p1a10*p1a10+2.e-29)/(1.+2.e-22*a12); // p6
521 // The main slope (peh_sl)
522 lastPAR[15]=400./a12+2.e-22*a9; // p1
523 lastPAR[16]=1.e-32*a12/(1.+5.e22/a14); // p2
524 lastPAR[17]=1000./a2+9.5*sa*ssa; // p3
525 lastPAR[18]=4.e-6*a*asa+1.e11/a16; // p4
526 lastPAR[19]=(120./a+.002*a2)/(1.+2.e14/a16); // p5
527 lastPAR[20]=9.+100./a; // p6
528 // The main quadratic (peh_sh)
529 lastPAR[21]=.002*a3+3.e7/a6; // p1
530 lastPAR[22]=7.e-15*a4*asa; // p2
531 lastPAR[23]=9000./a4; // p3
532 // The 1st max pre-exponent (peh_qq)
533 lastPAR[24]=.0011*asa/(1.+3.e34/a32/a4); // p1
534 lastPAR[25]=1.e-5*a2+2.e14/a16; // p2
535 lastPAR[26]=1.2e-11*a2/(1.+1.5e19/a12); // p3
536 lastPAR[27]=.016*asa/(1.+5.e16/a16); // p4
537 // The 1st max slope (peh_qs)
538 lastPAR[28]=.002*a4/(1.+7.e7/std::pow(a-6.83,14)); // p1
539 lastPAR[29]=2.e6/a6+7.2/std::pow(a,.11); // p2
540 lastPAR[30]=11.*a3/(1.+7.e23/a16/a8); // p3
541 lastPAR[31]=100./asa; // p4
542 // The 2nd max pre-exponent (peh_ss)
543 lastPAR[32]=(.1+4.4e-5*a2)/(1.+5.e5/a4); // p1
544 lastPAR[33]=3.5e-4*a2/(1.+1.e8/a8); // p2
545 lastPAR[34]=1.3+3.e5/a4; // p3
546 lastPAR[35]=500./(a2+50.)+3; // p4
547 lastPAR[36]=1.e-9/a+s4a16*s4a16; // p5
548 // The 2nd max slope (peh_sb)
549 lastPAR[37]=.4*asa+3.e-9*a6; // p1
550 lastPAR[38]=.0005*a5; // p2
551 lastPAR[39]=.002*a5; // p3
552 lastPAR[40]=10.; // p4
553 // The effective pre-exponent (peh_us)
554 lastPAR[41]=.05+.005*a; // p1
555 lastPAR[42]=7.e-8/sa; // p2
556 lastPAR[43]=.8*sa; // p3
557 lastPAR[44]=.02*sa; // p4
558 lastPAR[45]=1.e8/a3; // p5
559 lastPAR[46]=3.e32/(a32+1.e32); // p6
560 // The effective slope (peh_ub)
561 lastPAR[47]=24.; // p1
562 lastPAR[48]=20./sa; // p2
563 lastPAR[49]=7.e3*a/(sa+1.); // p3
564 lastPAR[50]=900.*sa/(1.+500./a3); // p4
565 }
566 // Parameter for lowEnergyNeutrons
567 lastPAR[51]=1.e15+2.e27/a4/(1.+2.e-18*a16);
568 }
569 lastPAR[nLast]=pwd;
570 // and initialize the zero element of the table
571 G4double lp=lPMin; // ln(momentum)
572 G4bool memCS=onlyCS; // ??
573 onlyCS=false;
574 lastCST[0]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables
575 onlyCS=memCS;
576 lastSST[0]=theSS;
577 lastS1T[0]=theS1;
578 lastB1T[0]=theB1;
579 lastS2T[0]=theS2;
580 lastB2T[0]=theB2;
581 lastS3T[0]=theS3;
582 lastB3T[0]=theB3;
583 lastS4T[0]=theS4;
584 lastB4T[0]=theB4;
585 }
586 if(LP>ILP)
587 {
588 G4int ini = static_cast<int>((ILP-lPMin+.000001)/dlnP)+1; // already inited till this
589 if(ini<0) ini=0;
590 if(ini<nPoints)
591 {
592 G4int fin = static_cast<int>((LP-lPMin)/dlnP)+1; // final bin of initialization
593 if(fin>=nPoints) fin=nLast; // Limit of the tabular initialization
594 if(fin>=ini)
595 {
596 G4double lp=0.;
597 for(G4int ip=ini; ip<=fin; ip++) // Calculate tabular CS,S1,B1,S2,B2,S3,B3
598 {
599 lp=lPMin+ip*dlnP; // ln(momentum)
600 G4bool memCS=onlyCS;
601 onlyCS=false;
602 lastCST[ip]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables (ret CS)
603 onlyCS=memCS;
604 lastSST[ip]=theSS;
605 lastS1T[ip]=theS1;
606 lastB1T[ip]=theB1;
607 lastS2T[ip]=theS2;
608 lastB2T[ip]=theB2;
609 lastS3T[ip]=theS3;
610 lastB3T[ip]=theB3;
611 lastS4T[ip]=theS4;
612 lastB4T[ip]=theB4;
613 }
614 return lp;
615 }
616 else G4cout<<"*Warning*G4ChipsAntiBaryonElasticXS::GetPTables: PDG="<<PDG
617 <<", Z="<<tgZ<<", N="<<tgN<<", i="<<ini<<" > fin="<<fin<<", LP="<<LP
618 <<" > ILP="<<ILP<<" nothing is done!"<<G4endl;
619 }
620 else G4cout<<"*Warning*G4ChipsAntiBaryonElasticXS::GetPTables: PDG="<<PDG
621 <<", Z="<<tgZ<<", N="<<tgN<<", i="<<ini<<">= max="<<nPoints<<", LP="<<LP
622 <<" > ILP="<<ILP<<", lPMax="<<lPMax<<" nothing is done!"<<G4endl;
623 }
624 }
625 else
626 {
627 // G4cout<<"*Error*G4ChipsAntiBaryonElasticXS::GetPTables: PDG="<<PDG<<", Z="<<tgZ
628 // <<", N="<<tgN<<", while it is defined only for Anti Baryons"<<G4endl;
629 // throw G4QException("G4ChipsAntiBaryonElasticXS::GetPTables:onlyaBA implemented");
631 ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
632 << ", while it is defined only for Anti Baryons" << G4endl;
633 G4Exception("G4ChipsAntiBaryonElasticXS::GetPTables()", "HAD_CHPS_0000",
634 FatalException, ed);
635 }
636 return ILP;
637}
638
639// Returns Q2=-t in independent units (MeV^2) (all internal calculations are in GeV)
641{
642 static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
643 static const G4double third=1./3.;
644 static const G4double fifth=1./5.;
645 static const G4double sevth=1./7.;
646
647 if(PDG<-3334 || PDG>-1111)G4cout<<"*Warning*G4QAntiBaryonElCS::GetExT:PDG="<<PDG<<G4endl;
648 if(onlyCS)G4cout<<"WarningG4ChipsAntiBaryonElasticXS::GetExchanT:onlyCS=1"<<G4endl;
649 if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
650 G4double q2=0.;
651 if(tgZ==1 && tgN==0) // ===> p+p=p+p
652 {
653 G4double E1=lastTM*theB1;
654 G4double R1=(1.-std::exp(-E1));
655 G4double E2=lastTM*theB2;
656 G4double R2=(1.-std::exp(-E2*E2*E2));
657 G4double E3=lastTM*theB3;
658 G4double R3=(1.-std::exp(-E3));
659 G4double I1=R1*theS1/theB1;
660 G4double I2=R2*theS2;
661 G4double I3=R3*theS3;
662 G4double I12=I1+I2;
663 G4double rand=(I12+I3)*G4UniformRand();
664 if (rand<I1 )
665 {
666 G4double ran=R1*G4UniformRand();
667 if(ran>1.) ran=1.;
668 q2=-std::log(1.-ran)/theB1;
669 }
670 else if(rand<I12)
671 {
672 G4double ran=R2*G4UniformRand();
673 if(ran>1.) ran=1.;
674 q2=-std::log(1.-ran);
675 if(q2<0.) q2=0.;
676 q2=std::pow(q2,third)/theB2;
677 }
678 else
679 {
680 G4double ran=R3*G4UniformRand();
681 if(ran>1.) ran=1.;
682 q2=-std::log(1.-ran)/theB3;
683 }
684 }
685 else
686 {
687 G4double a=tgZ+tgN;
688 G4double E1=lastTM*(theB1+lastTM*theSS);
689 G4double R1=(1.-std::exp(-E1));
690 G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
691 G4double tm2=lastTM*lastTM;
692 G4double E2=lastTM*tm2*theB2; // power 3 for lowA, 5 for HighA (1st)
693 if(a>6.5)E2*=tm2; // for heavy nuclei
694 G4double R2=(1.-std::exp(-E2));
695 G4double E3=lastTM*theB3;
696 if(a>6.5)E3*=tm2*tm2*tm2; // power 1 for lowA, 7 (2nd) for HighA
697 G4double R3=(1.-std::exp(-E3));
698 G4double E4=lastTM*theB4;
699 G4double R4=(1.-std::exp(-E4));
700 G4double I1=R1*theS1;
701 G4double I2=R2*theS2;
702 G4double I3=R3*theS3;
703 G4double I4=R4*theS4;
704 G4double I12=I1+I2;
705 G4double I13=I12+I3;
706 G4double rand=(I13+I4)*G4UniformRand();
707 if(rand<I1)
708 {
709 G4double ran=R1*G4UniformRand();
710 if(ran>1.) ran=1.;
711 q2=-std::log(1.-ran)/theB1;
712 if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
713 }
714 else if(rand<I12)
715 {
716 G4double ran=R2*G4UniformRand();
717 if(ran>1.) ran=1.;
718 q2=-std::log(1.-ran)/theB2;
719 if(q2<0.) q2=0.;
720 if(a<6.5) q2=std::pow(q2,third);
721 else q2=std::pow(q2,fifth);
722 }
723 else if(rand<I13)
724 {
725 G4double ran=R3*G4UniformRand();
726 if(ran>1.) ran=1.;
727 q2=-std::log(1.-ran)/theB3;
728 if(q2<0.) q2=0.;
729 if(a>6.5) q2=std::pow(q2,sevth);
730 }
731 else
732 {
733 G4double ran=R4*G4UniformRand();
734 if(ran>1.) ran=1.;
735 q2=-std::log(1.-ran)/theB4;
736 if(a<6.5) q2=lastTM-q2; // u reduced for lightA (starts from 0)
737 }
738 }
739 if(q2<0.) q2=0.;
740 if(!(q2>=-1.||q2<=1.))G4cout<<"*NAN*G4QaBElasticCrossSect::GetExchangeT:-t="<<q2<<G4endl;
741 if(q2>lastTM)
742 {
743 q2=lastTM;
744 }
745 return q2*GeVSQ;
746}
747
748// Returns B in independent units (MeV^-2) (all internal calculations are in GeV) see ExT
749G4double G4ChipsAntiBaryonElasticXS::GetSlope(G4int tgZ, G4int tgN, G4int PDG)
750{
751 static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
752 if(onlyCS)G4cout<<"WarningG4ChipsAntiBaryonElasticXS::GetSlope:onlCS=true"<<G4endl;
753 if(lastLP<-4.3) return 0.; // S-wave for p<14 MeV/c (kinE<.1MeV)
754 if(PDG<-3334 || PDG>-1111)
755 {
756 // G4cout<<"*Error*G4ChipsAntiBaryonElasticXS::GetSlope: PDG="<<PDG<<", Z="<<tgZ
757 // <<", N="<<tgN<<", while it is defined only for Anti Baryons"<<G4endl;
758 // throw G4QException("G4ChipsAntiBaryonElasticXS::GetSlope: AnBa are implemented");
760 ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
761 << ", while it is defined only for Anti Baryons" << G4endl;
762 G4Exception("G4ChipsAntiBaryonElasticXS::GetSlope()", "HAD_CHPS_0000",
763 FatalException, ed);
764 }
765 if(theB1<0.) theB1=0.;
766 if(!(theB1>=-1.||theB1<=1.))G4cout<<"*NAN*G4QaBaElasticCrossS::Getslope:"<<theB1<<G4endl;
767 return theB1/GeVSQ;
768}
769
770// Returns half max(Q2=-t) in independent units (MeV^2)
771G4double G4ChipsAntiBaryonElasticXS::GetHMaxT()
772{
773 static const G4double HGeVSQ=gigaelectronvolt*gigaelectronvolt/2.;
774 return lastTM*HGeVSQ;
775}
776
777// lastLP is used, so calculating tables, one need to remember and then recover lastLP
778G4double G4ChipsAntiBaryonElasticXS::GetTabValues(G4double lp, G4int PDG, G4int tgZ,
779 G4int tgN)
780{
781 if(PDG<-3334 || PDG>-1111) G4cout<<"*Warning*G4QAntiBaryElCS::GetTabV:PDG="<<PDG<<G4endl;
782 if(tgZ<0 || tgZ>92)
783 {
784 G4cout<<"*Warning*G4QAntiBaryonElCS::GetTabValue:(1-92) NoIsotopesFor Z="<<tgZ<<G4endl;
785 return 0.;
786 }
787 G4int iZ=tgZ-1; // Z index
788 if(iZ<0)
789 {
790 iZ=0; // conversion of the neutron target to the proton target
791 tgZ=1;
792 tgN=0;
793 }
794 G4double p=std::exp(lp); // momentum
795 G4double sp=std::sqrt(p); // sqrt(p)
796 G4double p2=p*p;
797 G4double p3=p2*p;
798 G4double p4=p3*p;
799 if ( tgZ == 1 && tgN == 0 ) // PiMin+P
800 {
801 G4double dl2=lp-lastPAR[6]; // ld ?
802 theSS=lastPAR[29];
803 theS1=(lastPAR[7]+lastPAR[8]*dl2*dl2)/(1.+lastPAR[9]/p4/p)+
804 (lastPAR[10]/p2+lastPAR[11]*p)/(p4+lastPAR[12]*sp);
805 theB1=lastPAR[13]*std::pow(p,lastPAR[14])/(1.+lastPAR[15]/p3);
806 theS2=lastPAR[16]+lastPAR[17]/(p4+lastPAR[18]*p);
807 theB2=lastPAR[19]+lastPAR[20]/(p4+lastPAR[21]/sp);
808 theS3=lastPAR[22]+lastPAR[23]/(p4*p4+lastPAR[24]*p2+lastPAR[25]);
809 theB3=lastPAR[26]+lastPAR[27]/(p4+lastPAR[28]);
810 theS4=0.;
811 theB4=0.;
812 // Returns the total elastic pim-p cross-section (to avoid spoiling lastSIG)
813 G4double ye=std::exp(lp*lastPAR[0]);
814 G4double dp=lp-lastPAR[1];
815 return lastPAR[2]/(ye+lastPAR[3])+lastPAR[4]*dp*dp+lastPAR[5];
816 }
817 else
818 {
819 G4double p5=p4*p;
820 G4double p6=p5*p;
821 G4double p8=p6*p2;
822 G4double p10=p8*p2;
823 G4double p12=p10*p2;
824 G4double p16=p8*p8;
825 //G4double p24=p16*p8;
826 G4double dl=lp-5.;
827 G4double a=tgZ+tgN;
828 G4double pah=std::pow(p,a/2);
829 G4double pa=pah*pah;
830 G4double pa2=pa*pa;
831 if(a<6.5)
832 {
833 theS1=lastPAR[9]/(1.+lastPAR[10]*p4*pa)+lastPAR[11]/(p4+lastPAR[12]*p4/pa2)+
834 (lastPAR[13]*dl*dl+lastPAR[14])/(1.+lastPAR[15]/p2);
835 theB1=(lastPAR[16]+lastPAR[17]*p2)/(p4+lastPAR[18]/pah)+lastPAR[19];
836 theSS=lastPAR[20]/(1.+lastPAR[21]/p2)+lastPAR[22]/(p6/pa+lastPAR[23]/p16);
837 theS2=lastPAR[24]/(pa/p2+lastPAR[25]/p4)+lastPAR[26];
838 theB2=lastPAR[27]*std::pow(p,lastPAR[28])+lastPAR[29]/(p8+lastPAR[30]/p16);
839 theS3=lastPAR[31]/(pa*p+lastPAR[32]/pa)+lastPAR[33];
840 theB3=lastPAR[34]/(p3+lastPAR[35]/p6)+lastPAR[36]/(1.+lastPAR[37]/p2);
841 theS4=p2*(pah*lastPAR[38]*std::exp(-pah*lastPAR[39])+
842 lastPAR[40]/(1.+lastPAR[41]*std::pow(p,lastPAR[42])));
843 theB4=lastPAR[43]*pa/p2/(1.+pa*lastPAR[44]);
844 }
845 else
846 {
847 theS1=lastPAR[9]/(1.+lastPAR[10]/p4)+lastPAR[11]/(p4+lastPAR[12]/p2)+
848 lastPAR[13]/(p5+lastPAR[14]/p16);
849 theB1=(lastPAR[15]/p8+lastPAR[19])/(p+lastPAR[16]/std::pow(p,lastPAR[20]))+
850 lastPAR[17]/(1.+lastPAR[18]/p4);
851 theSS=lastPAR[21]/(p4/std::pow(p,lastPAR[23])+lastPAR[22]/p4);
852 theS2=lastPAR[24]/p4/(std::pow(p,lastPAR[25])+lastPAR[26]/p12)+lastPAR[27];
853 theB2=lastPAR[28]/std::pow(p,lastPAR[29])+lastPAR[30]/std::pow(p,lastPAR[31]);
854 theS3=lastPAR[32]/std::pow(p,lastPAR[35])/(1.+lastPAR[36]/p12)+
855 lastPAR[33]/(1.+lastPAR[34]/p6);
856 theB3=lastPAR[37]/p8+lastPAR[38]/p2+lastPAR[39]/(1.+lastPAR[40]/p8);
857 theS4=(lastPAR[41]/p4+lastPAR[46]/p)/(1.+lastPAR[42]/p10)+
858 (lastPAR[43]+lastPAR[44]*dl*dl)/(1.+lastPAR[45]/p12);
859 theB4=lastPAR[47]/(1.+lastPAR[48]/p)+lastPAR[49]*p4/(1.+lastPAR[50]*p5);
860 }
861 // Returns the total elastic (n/p)A cross-section (to avoid spoiling lastSIG)
862 G4double dlp=lp-lastPAR[4]; // ax
863 // p1 p2 p3 p4
864 return (lastPAR[0]*dlp*dlp+lastPAR[1]+lastPAR[2]/p)/(1.+lastPAR[3]/p);
865 }
866 return 0.;
867} // End of GetTableValues
868
869// Returns max -t=Q2 (GeV^2) for the momentum pP(GeV) and the target nucleus (tgN,tgZ)
870G4double G4ChipsAntiBaryonElasticXS::GetQ2max(G4int PDG, G4int tgZ, G4int tgN,
871 G4double pP)
872{
873 static const G4double mNeut= G4Neutron::Neutron()->GetPDGMass()*.001; // MeV to GeV
874 static const G4double mProt= G4Proton::Proton()->GetPDGMass()*.001; // MeV to GeV
875 static const G4double mNuc2= sqr((mProt+mNeut)/2);
876 G4double pP2=pP*pP; // squared momentum of the projectile
877 if(tgZ || tgN>-1) // ---> pipA
878 {
879 G4double mt=G4ParticleTable::GetParticleTable()->FindIon(tgZ,tgZ+tgN,0,tgZ)->GetPDGMass()*.001; // Target mass in GeV
880 G4double dmt=mt+mt;
881 G4double mds=dmt*std::sqrt(pP2+mNuc2)+mNuc2+mt*mt; // Mondelstam mds (@@ other AntiBar?)
882 return dmt*dmt*pP2/mds;
883 }
884 else
885 {
886 // G4cout<<"*Error*G4ChipsAntiBaryonElasticXS::GetQ2ma:PDG="<<PDG<<",Z="<<tgZ<<",N="
887 // <<tgN<<", while it is defined only for p projectiles & Z_target>0"<<G4endl;
888 // throw G4QException("G4ChipsAntiBaryonElasticXS::GetQ2max: only aBA implemented");
890 ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
891 << ", while it is defined only for p projectiles & Z_target>0" << G4endl;
892 G4Exception("G4ChipsAntiBaryonElasticXS::GetQ2max()", "HAD_CHPS_0000",
893 FatalException, ed);
894 return 0;
895 }
896}
#define G4_DECLARE_XS_FACTORY(cross_section)
@ FatalException
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 G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
static G4AntiLambda * AntiLambda()
static G4AntiNeutron * AntiNeutron()
static G4AntiOmegaMinus * AntiOmegaMinus()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4AntiSigmaZero * AntiSigmaZero()
static G4AntiXiMinus * AntiXiMinus()
static G4AntiXiZero * AntiXiZero()
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)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
G4ParticleDefinition * GetDefinition() const
G4double GetTotalMomentum() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4ParticleDefinition * FindIon(G4int atomicNumber, G4int atomicMass, G4double excitationEnergy)
static G4ParticleTable * GetParticleTable()
static G4Proton * Proton()
Definition: G4Proton.cc:93
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
T sqr(const T &x)
Definition: templates.hh:145