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