Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChipsNeutronElasticXS.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: G4ChipsNeutronElasticXS for nA elastic cross sections
31// Created: M.V. Kossov, CERN/ITEP(Moscow), 10-OCT-01
32// The last update: M.V. Kossov, CERN/ITEP (Moscow) 12-Jan-10 (from G4QElCrSect)
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 "G4Neutron.hh"
45#include "G4Nucleus.hh"
46#include "G4ParticleTable.hh"
47#include "G4NucleiProperties.hh"
48
49// factory
51//
53
54G4ChipsNeutronElasticXS::G4ChipsNeutronElasticXS():G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
55{
56 lPMin=-8.; // Min tabulated log Momentum (D)
57 lPMax= 8.; // Max tabulated log Momentum (D)
58 dlnP=(lPMax-lPMin)/nLast;// Log step in table (D)
59 onlyCS=true;// Flag to calc only CS (not Si/Bi)(L)
60 lastSIG=0.; // Last calculated cross section (L)
61 lastLP=-10.;// Last log(momOfIncidentHadron) (L)
62 lastTM=0.; // Last t_maximum (L)
63 theSS=0.; // The Last sq.slope of 1st difMax (L)
64 theS1=0.; // The Last mantissa of 1st difMax (L)
65 theB1=0.; // The Last slope of 1st difr. Max (L)
66 theS2=0.; // The Last mantissa of 2nd difMax (L)
67 theB2=0.; // The Last slope of 2nd difr. Max (L)
68 theS3=0.; // The Last mantissa of 3d difrMax (L)
69 theB3=0.; // The Last slope of 3d difructMax (L)
70 theS4=0.; // The Last mantissa of 4th difMax (L)
71 theB4=0.; // The Last slope of 4th difr. Max (L)
72 lastTZ=0; // Last atomic number of the target
73 lastTN=0; // Last # of neutrons in the target
74 lastPIN=0.; // Last initialized max momentum
75 lastCST=0; // Elastic cross-section table
76 lastPAR=0; // Parameters of FunctionalCalculation
77 lastSST=0; // E-dep of sq.slope of the 1st difMax
78 lastS1T=0; // E-dep of mantissa of the 1st difMax
79 lastB1T=0; // E-dep of theSlope of the 1st difMax
80 lastS2T=0; // E-dep of mantissa of the 2nd difMax
81 lastB2T=0; // E-dep of theSlope of the 2nd difMax
82 lastS3T=0; // E-dep of mantissa of the 3d difrMax
83 lastB3T=0; // E-dep of the slope of the 3d difMax
84 lastS4T=0; // E-dep of mantissa of the 4th difMax
85 lastB4T=0; // E-dep of theSlope of the 4th difMax
86 lastN=0; // The last N of calculated nucleus
87 lastZ=0; // The last Z of calculated nucleus
88 lastP=0.; // Last used in cross section Momentum
89 lastTH=0.; // Last threshold momentum
90 lastCS=0.; // Last value of the Cross Section
91 lastI=0; // The last position in the DAMDB
92}
93
95{
96 std::vector<G4double*>::iterator pos;
97 for (pos=CST.begin(); pos<CST.end(); pos++)
98 { delete [] *pos; }
99 CST.clear();
100 for (pos=PAR.begin(); pos<PAR.end(); pos++)
101 { delete [] *pos; }
102 PAR.clear();
103 for (pos=SST.begin(); pos<SST.end(); pos++)
104 { delete [] *pos; }
105 SST.clear();
106 for (pos=S1T.begin(); pos<S1T.end(); pos++)
107 { delete [] *pos; }
108 S1T.clear();
109 for (pos=B1T.begin(); pos<B1T.end(); pos++)
110 { delete [] *pos; }
111 B1T.clear();
112 for (pos=S2T.begin(); pos<S2T.end(); pos++)
113 { delete [] *pos; }
114 S2T.clear();
115 for (pos=B2T.begin(); pos<B2T.end(); pos++)
116 { delete [] *pos; }
117 B2T.clear();
118 for (pos=S3T.begin(); pos<S3T.end(); pos++)
119 { delete [] *pos; }
120 S3T.clear();
121 for (pos=B3T.begin(); pos<B3T.end(); pos++)
122 { delete [] *pos; }
123 B3T.clear();
124 for (pos=S4T.begin(); pos<S4T.end(); pos++)
125 { delete [] *pos; }
126 S4T.clear();
127 for (pos=B4T.begin(); pos<B4T.end(); pos++)
128 { delete [] *pos; }
129 B4T.clear();
130}
131
133 const G4Element*,
134 const G4Material*)
135{
136 G4ParticleDefinition* particle = Pt->GetDefinition();
137 if (particle == G4Proton::Proton() ) return true;
138 return false;
139}
140
142 const G4Isotope*,
143 const G4Element*,
144 const G4Material*)
145{
146 G4double pMom=Pt->GetTotalMomentum();
147 G4int tgN = A - tgZ;
148
149 return GetChipsCrossSection(pMom, tgZ, tgN, 2112);
150}
151
152// The main member function giving the collision cross section (P is in IU, CS is in mb)
153// Make pMom in independent units ! (Now it is MeV)
155{
156 static std::vector <G4int> colN; // Vector of N for calculated nuclei (isotops)
157 static std::vector <G4int> colZ; // Vector of Z for calculated nuclei (isotops)
158 static std::vector <G4double> colP; // Vector of last momenta for the reaction
159 static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
160 static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
161 // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
162
163 G4double pEn=pMom;
164 onlyCS=false;
165
166 G4bool in=false; // By default the isotope must be found in the AMDB
167 lastP = 0.; // New momentum history (nothing to compare with)
168 lastN = tgN; // The last N of the calculated nucleus
169 lastZ = tgZ; // The last Z of the calculated nucleus
170 lastI = colN.size(); // Size of the Associative Memory DB in the heap
171 if(lastI) for(G4int i=0; i<lastI; i++) // Loop over proj/tgZ/tgN lines of DB
172 { // The nucleus with projPDG is found in AMDB
173 if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
174 {
175 lastI=i;
176 lastTH =colTH[i]; // Last THreshold (A-dependent)
177 if(pEn<=lastTH)
178 {
179 return 0.; // Energy is below the Threshold value
180 }
181 lastP =colP [i]; // Last Momentum (A-dependent)
182 lastCS =colCS[i]; // Last CrossSect (A-dependent)
183 // if(std::fabs(lastP/pMom-1.)<tolerance) //VI (do not use tolerance)
184 if(lastP == pMom) // Do not recalculate
185 {
186 CalculateCrossSection(false,-1,i,2112,lastZ,lastN,pMom); // Update param's only
187 return lastCS*millibarn; // Use theLastCS
188 }
189 in = true; // This is the case when the isotop is found in DB
190 // Momentum pMom is in IU ! @@ Units
191 lastCS=CalculateCrossSection(false,-1,i,2112,lastZ,lastN,pMom); // read & update
192 if(lastCS<=0. && pEn>lastTH) // Correct the threshold
193 {
194 lastTH=pEn;
195 }
196 break; // Go out of the LOOP with found lastI
197 }
198 }
199 if(!in) // This nucleus has not been calculated previously
200 {
201 //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
202 lastCS=CalculateCrossSection(false,0,lastI,2112,lastZ,lastN,pMom);//calculate&create
203 if(lastCS<=0.)
204 {
205 lastTH = 0; // ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
206 if(pEn>lastTH)
207 {
208 lastTH=pEn;
209 }
210 }
211 colN.push_back(tgN);
212 colZ.push_back(tgZ);
213 colP.push_back(pMom);
214 colTH.push_back(lastTH);
215 colCS.push_back(lastCS);
216 return lastCS*millibarn;
217 } // End of creation of the new set of parameters
218 else
219 {
220 colP[lastI]=pMom;
221 colCS[lastI]=lastCS;
222 }
223 return lastCS*millibarn;
224}
225
226// Calculation of total elastic cross section (p in IU, CS in mb) @@ Units (?)
227// F=0 - create AMDB, F=-1 - read&update AMDB, F=1 - update AMDB (sinchro with higher AMDB)
228G4double G4ChipsNeutronElasticXS::CalculateCrossSection(G4bool CS, G4int F,G4int I,
229 G4int PDG, G4int tgZ, G4int tgN, G4double pIU)
230{
231 // *** Begin of Associative Memory DB for acceleration of the cross section calculations
232 static std::vector <G4double> PIN; // Vector of max initialized log(P) in the table
233 // *** End of Static Definitions (Associative Memory Data Base) ***
234 G4double pMom=pIU/GeV; // All calculations are in GeV
235 onlyCS=CS; // Flag to calculate only CS (not Si/Bi)
236 lastLP=std::log(pMom); // Make a logarithm of the momentum for calculation
237 if(F) // This isotope was found in AMDB =>RETRIEVE/UPDATE
238 {
239 if(F<0) // the AMDB must be loded
240 {
241 lastPIN = PIN[I]; // Max log(P) initialised for this table set
242 lastPAR = PAR[I]; // Pointer to the parameter set
243 lastCST = CST[I]; // Pointer to the total sross-section table
244 lastSST = SST[I]; // Pointer to the first squared slope
245 lastS1T = S1T[I]; // Pointer to the first mantissa
246 lastB1T = B1T[I]; // Pointer to the first slope
247 lastS2T = S2T[I]; // Pointer to the second mantissa
248 lastB2T = B2T[I]; // Pointer to the second slope
249 lastS3T = S3T[I]; // Pointer to the third mantissa
250 lastB3T = B3T[I]; // Pointer to the rhird slope
251 lastS4T = S4T[I]; // Pointer to the 4-th mantissa
252 lastB4T = B4T[I]; // Pointer to the 4-th slope
253 }
254 if(lastLP>lastPIN && lastLP<lPMax)
255 {
256 lastPIN=GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);// Can update upper logP-Limit in tabs
257 PIN[I]=lastPIN; // Remember the new P-Limit of the tables
258 }
259 }
260 else // This isotope wasn't initialized => CREATE
261 {
262 lastPAR = new G4double[nPoints]; // Allocate memory for parameters of CS function
263 lastPAR[nLast]=0; // Initialization for VALGRIND
264 lastCST = new G4double[nPoints]; // Allocate memory for Tabulated CS function
265 lastSST = new G4double[nPoints]; // Allocate memory for Tabulated first sqaredSlope
266 lastS1T = new G4double[nPoints]; // Allocate memory for Tabulated first mantissa
267 lastB1T = new G4double[nPoints]; // Allocate memory for Tabulated first slope
268 lastS2T = new G4double[nPoints]; // Allocate memory for Tabulated second mantissa
269 lastB2T = new G4double[nPoints]; // Allocate memory for Tabulated second slope
270 lastS3T = new G4double[nPoints]; // Allocate memory for Tabulated third mantissa
271 lastB3T = new G4double[nPoints]; // Allocate memory for Tabulated third slope
272 lastS4T = new G4double[nPoints]; // Allocate memory for Tabulated 4-th mantissa
273 lastB4T = new G4double[nPoints]; // Allocate memory for Tabulated 4-th slope
274 lastPIN = GetPTables(lastLP,lPMin,PDG,tgZ,tgN); // Returns the new P-limit for tables
275 PIN.push_back(lastPIN); // Fill parameters of CS function to AMDB
276 PAR.push_back(lastPAR); // Fill parameters of CS function to AMDB
277 CST.push_back(lastCST); // Fill Tabulated CS function to AMDB
278 SST.push_back(lastSST); // Fill Tabulated first sq.slope to AMDB
279 S1T.push_back(lastS1T); // Fill Tabulated first mantissa to AMDB
280 B1T.push_back(lastB1T); // Fill Tabulated first slope to AMDB
281 S2T.push_back(lastS2T); // Fill Tabulated second mantissa to AMDB
282 B2T.push_back(lastB2T); // Fill Tabulated second slope to AMDB
283 S3T.push_back(lastS3T); // Fill Tabulated third mantissa to AMDB
284 B3T.push_back(lastB3T); // Fill Tabulated third slope to AMDB
285 S4T.push_back(lastS4T); // Fill Tabulated 4-th mantissa to AMDB
286 B4T.push_back(lastB4T); // Fill Tabulated 4-th slope to AMDB
287 } // End of creation/update of the new set of parameters and tables
288 // =-------= NOW Update (if necessary) and Calculate the Cross Section =---------=
289 if(lastLP>lastPIN && lastLP<lPMax)
290 {
291 lastPIN = GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);
292 }
293 if(!onlyCS) lastTM=GetQ2max(PDG, tgZ, tgN, pMom); // Calculate (-t)_max=Q2_max (GeV2)
294 if(lastLP>lPMin && lastLP<=lastPIN) // Linear fit is made using precalculated tables
295 {
296 if(lastLP==lastPIN)
297 {
298 G4double shift=(lastLP-lPMin)/dlnP+.000001; // Log distance from lPMin
299 G4int blast=static_cast<int>(shift); // this is a bin number of the lower edge (0)
300 if(blast<0 || blast>=nLast) G4cout<<"G4QNeutElCS::CCS:b="<<blast<<","<<nLast<<G4endl;
301 lastSIG = lastCST[blast];
302 if(!onlyCS) // Skip the differential cross-section parameters
303 {
304 theSS = lastSST[blast];
305 theS1 = lastS1T[blast];
306 theB1 = lastB1T[blast];
307 theS2 = lastS2T[blast];
308 theB2 = lastB2T[blast];
309 theS3 = lastS3T[blast];
310 theB3 = lastB3T[blast];
311 theS4 = lastS4T[blast];
312 theB4 = lastB4T[blast];
313 }
314 }
315 else
316 {
317 G4double shift=(lastLP-lPMin)/dlnP; // a shift from the beginning of the table
318 G4int blast=static_cast<int>(shift); // the lower bin number
319 if(blast<0) blast=0;
320 if(blast>=nLast) blast=nLast-1; // low edge of the last bin
321 shift-=blast; // step inside the unit bin
322 G4int lastL=blast+1; // the upper bin number
323 G4double SIGL=lastCST[blast]; // the basic value of the cross-section
324 lastSIG= SIGL+shift*(lastCST[lastL]-SIGL); // calculated total elastic cross-section
325 if(!onlyCS) // Skip the differential cross-section parameters
326 {
327 G4double SSTL=lastSST[blast]; // the low bin of the first squared slope
328 theSS=SSTL+shift*(lastSST[lastL]-SSTL); // the basic value of the first sq.slope
329 G4double S1TL=lastS1T[blast]; // the low bin of the first mantissa
330 theS1=S1TL+shift*(lastS1T[lastL]-S1TL); // the basic value of the first mantissa
331 G4double B1TL=lastB1T[blast]; // the low bin of the first slope
332 theB1=B1TL+shift*(lastB1T[lastL]-B1TL); // the basic value of the first slope
333 G4double S2TL=lastS2T[blast]; // the low bin of the second mantissa
334 theS2=S2TL+shift*(lastS2T[lastL]-S2TL); // the basic value of the second mantissa
335 G4double B2TL=lastB2T[blast]; // the low bin of the second slope
336 theB2=B2TL+shift*(lastB2T[lastL]-B2TL); // the basic value of the second slope
337 G4double S3TL=lastS3T[blast]; // the low bin of the third mantissa
338 theS3=S3TL+shift*(lastS3T[lastL]-S3TL); // the basic value of the third mantissa
339 G4double B3TL=lastB3T[blast]; // the low bin of the third slope
340 theB3=B3TL+shift*(lastB3T[lastL]-B3TL); // the basic value of the third slope
341 G4double S4TL=lastS4T[blast]; // the low bin of the 4-th mantissa
342 theS4=S4TL+shift*(lastS4T[lastL]-S4TL); // the basic value of the 4-th mantissa
343 G4double B4TL=lastB4T[blast]; // the low bin of the 4-th slope
344 theB4=B4TL+shift*(lastB4T[lastL]-B4TL); // the basic value of the 4-th slope
345 }
346 }
347 }
348 else lastSIG=GetTabValues(lastLP, PDG, tgZ, tgN); // Direct calculation beyond the table
349 if(lastSIG<0.) lastSIG = 0.; // @@ a Warning print can be added
350 return lastSIG;
351}
352
353// It has parameter sets for all tZ/tN/PDG, using them the tables can be created/updated
354G4double G4ChipsNeutronElasticXS::GetPTables(G4double LP, G4double ILP, G4int PDG,
355 G4int tgZ, G4int tgN)
356{
357 // @@ At present all nA==pA ---------> Each neucleus can have not more than 51 parameters
358 static const G4double pwd=2727;
359 const G4int n_npel=24; // #of parameters for np-elastic (<nPoints=128)
360 const G4int n_ppel=32; // #of parameters for pp-elastic (<nPoints=128)
361 // -0- -1- -2- -3- -4- -5- -6- -7- -8- -9--10--11--12--13- -14-
362 G4double np_el[n_npel]={12.,.05,.0001,5.,.35,6.75,.14,19.,.6,6.75,.14,13.,.14,.6,.00013,
363 75.,.001,7.2,4.32,.012,2.5,0.0,12.,.34};
364 // -15--16--17- -18- -19--20--21--22--23-
365 // -0- -1- -2- -3- -4- -5- -6- -7- -8--9--10--11--12--13-
366 G4double pp_el[n_ppel]={2.865,18.9,.6461,3.,9.,.425,.4276,.0022,5.,74.,3.,3.4,.2,.17,
367 .001,8.,.055,3.64,5.e-5,4000.,1500.,.46,1.2e6,3.5e6,5.e-5,1.e10,
368 8.5e8,1.e10,1.1,3.4e6,6.8e6,0.};
369 // -14--15- -16- -17- -18- -19- -20- -21- -22- -23- -24- -25-
370 // -26- -27- -28- -29- -30- -31-
371 //==> n (Z=0)
372 static const G4int N0=1; // *** Not used (fake)***
373 static const G4double pZ0N1[7]={0., 0., 0., 0., 0., 0., 0.}; // Not used (fake)
374 static const std::pair<G4int, const G4double*> Z0N1(1,pZ0N1);
375 static const std::pair<G4int, const G4double*> Z0[N0]={Z0N1};
376 //==> H (Z=1) *** protons are treated separately ***
377 static const G4int N1=3;
378 static const G4double pZ1N0[7]={0., 0., 0., 0., 0., 0., 0.}; // Not used (fake)
379 static const std::pair<G4int, const G4double*> Z1N0(0,pZ1N0);
380 static const G4double pZ1N1[7]={6.E-5, 4., .055, 1.1E-8, .008, 1.2E-8, .019};
381 static const std::pair<G4int, const G4double*> Z1N1(1,pZ1N1);
382 static const G4double pZ1N2[7]={6.E-5, 2.2, .051, 1.E-8, .04, 9.E-8, .0075};
383 static const std::pair<G4int, const G4double*> Z1N2(2,pZ1N2);
384 static const std::pair<G4int, const G4double*> Z1[N1]={Z1N0, Z1N1, Z1N2};
385 //==> He(Z=2)
386 static const G4int N2=2;
387 static const G4double pZ2N1[7]={6.E-5, 3., .06, 4.E-9, .03, 7.E-8, .015};
388 static const std::pair<G4int, const G4double*> Z2N1(1,pZ2N1);
389 static const G4double pZ2N2[7]={3.E-4, .23, 1., 1.5E-9, 2.E-02, 1.E-8, .003};
390 static const std::pair<G4int, const G4double*> Z2N2(2,pZ2N2);
391 static const std::pair<G4int, const G4double*> Z2[N2]={Z2N1, Z2N2};
392 //==> Li(Z=3)
393 static const G4int N3=2;
394 static const G4double pZ3N3[7]={3.1E-7, 1.7, 1.3E-4, 1.E-8, .02, 1.1E-7, .0023};
395 static const std::pair<G4int, const G4double*> Z3N1(3,pZ3N3);
396 static const G4double pZ3N4[7]={1.3E-6, 1.8, 7.6E-4, 9.E-9, .03, 1.E-7, .0029};
397 static const std::pair<G4int, const G4double*> Z3N2(4,pZ3N4);
398 static const std::pair<G4int, const G4double*> Z3[N3]={Z3N1, Z3N2};
399 //==> Be(Z=4)
400 static const G4int N4=2;
401 static const G4double pZ4N3[7]={2.E-4, 1.4, 2.7, 0., .02, 5.E-8, 0.};
402 static const std::pair<G4int, const G4double*> Z4N3(3,pZ4N3);
403 static const G4double pZ4N5[7]={1.E-6, 5.7, .0011, 3.E-9, .007, 2.E-8, .016};
404 static const std::pair<G4int, const G4double*> Z4N5(5,pZ4N5);
405 static const std::pair<G4int, const G4double*> Z4[N4]={Z4N3,Z4N5};
406 //==> B (Z=5)
407 static const G4int N5=2;
408 static const G4double pZ5N5[7]={8.E-7, 5., 3.4E-4, 7.E-9, 1.E-02, 1.E-07, .0053};
409 static const std::pair<G4int, const G4double*> Z5N5(5,pZ5N5);
410 static const G4double pZ5N6[7]={4.8E-6, 6.6, .0035, 4.E-9, .003, 1.E-8, .012};
411 static const std::pair<G4int, const G4double*> Z5N6(6,pZ5N6);
412 static const std::pair<G4int, const G4double*> Z5[N5]={Z5N5, Z5N6};
413 //==> C (Z=6) *** Only nat (C13=C12=C_nat) ***
414 static const G4int N6=2;
415 static const G4double pZ6N6[7]={4.9E-6, 6.6, .0035, 4.E-9, .002, 6.E-9, .011};
416 static const std::pair<G4int, const G4double*> Z6N6(6,pZ6N6);
417 static const G4double pZ6N7[7]={4.9E-6, 6.6, .0035, 4.E-9, .002, 6.E-9, .011};
418 static const std::pair<G4int, const G4double*> Z6N7(7,pZ6N7);
419 static const std::pair<G4int, const G4double*> Z6[N6]={Z6N6, Z6N7};
420 //==> N (Z=7)
421 static const G4int N7=2;
422 static const G4double pZ7N7[7]={4.9E-6, 1.6, .03, .4E-9, .02, 6.E-8, .021};
423 static const std::pair<G4int, const G4double*> Z7N7(7,pZ7N7);
424 static const G4double pZ7N8[7]={2.5E-6, 5., .0021, 2.5E-9, .015, 5.E-8, .009};
425 static const std::pair<G4int, const G4double*> Z7N8(8,pZ7N8);
426 static const std::pair<G4int, const G4double*> Z7[N7]={Z7N7, Z7N8};
427 //==> O (Z=8) (O18=O17, No data)
428 static const G4int N8=3;
429 static const G4double pZ8N8[7]={2.5E-6, 5.3, .0018, 3.E-9, .01, 1.5E-8, .0075};
430 static const std::pair<G4int, const G4double*> Z8N8(8,pZ8N8);
431 static const G4double pZ8N9[7]={1.4E-6, 2.1, .0025, 1.3E-9, .02, 1.3E-7, .0072};
432 static const std::pair<G4int, const G4double*> Z8N9(9,pZ8N9);
433 static const G4double pZ8N10[7]={1.4E-6, 2.1, .0025, 1.3E-9, .02, 1.3E-7, .0072};
434 static const std::pair<G4int, const G4double*> Z8N10(10,pZ8N10);
435 static const std::pair<G4int, const G4double*> Z8[N8]={Z8N8, Z8N9, Z8N10};
436 //==> F (Z=9)
437 static const G4int N9=1;
438 static const G4double pZ9N10[7]={1.4E-6, 7.5, 6.7E-4, 4.E-9, 2.E-5, 7.E-12, .0066};
439 static const std::pair<G4int, const G4double*> Z9N10(10,pZ9N10);
440 static const std::pair<G4int, const G4double*> Z9[N9]={Z9N10};
441 //==> Ne(Z=10) *** No data *** (Ne20=Na22, Ne21=F19, Ne22=Na22)
442 static const G4int N10=3;
443 static const G4double pZ10N10[7]={1.4E-5, 7., .01, 0., 3.E-11, 7.E-24, .12};
444 static const std::pair<G4int, const G4double*> Z10N10(10,pZ10N10);
445 static const G4double pZ10N11[7]={1.4E-6, 7.5, 6.7E-4, 4.E-9, 2.E-5, 7.E-12, .0066};
446 static const std::pair<G4int, const G4double*> Z10N11(11,pZ10N11);
447 static const G4double pZ10N12[7]={1.4E-5, 7., .01, 0., 3.E-11, 7.E-24, .12};
448 static const std::pair<G4int, const G4double*> Z10N12(12,pZ10N12);
449 static const std::pair<G4int, const G4double*> Z10[N10]={Z10N10, Z10N11, Z10N12};
450 //==> Na(Z=11)
451 static const G4int N11=2;
452 static const G4double pZ11N11[7]={1.4E-5, 7., .01, 0., 3.E-11, 7.E-24, .12};
453 static const std::pair<G4int, const G4double*> Z11N11(11,pZ11N11);
454 static const G4double pZ11N12[7]={1.4E-6, 7.6, 6.E-4, 5.E-9, 7.E-9, 3.E-18, .0056};
455 static const std::pair<G4int, const G4double*> Z11N12(12,pZ11N12);
456 static const std::pair<G4int, const G4double*> Z11[N11]={Z11N11, Z11N12};
457 //==> Mg(Z=12)
458 static const G4int N12=3;
459 static const G4double pZ12N12[7]={8.E-7, 3., .001, 1.8E-9, .0015, .2E-9, .006};
460 static const std::pair<G4int, const G4double*> Z12N12(12,pZ12N12);
461 static const G4double pZ12N13[7]={8.E-7, 7., 3.E-4, 6.E-9, .006, 4.E-8, .0042};
462 static const std::pair<G4int, const G4double*> Z12N13(13,pZ12N13);
463 static const G4double pZ12N14[7]={1.2E-6, 6.8, 5.E-4, 5.E-9, .007, 2.E-8, .0044};
464 static const std::pair<G4int, const G4double*> Z12N14(14,pZ12N14);
465 static const std::pair<G4int, const G4double*> Z12[N12]={Z12N12, Z12N13, Z12N14};
466 //==> Al(Z=13)
467 static const G4int N13=1;
468 static const G4double pZ13N14[7]={3.E-7, 5., 8.4E-5, 7.E-9, .008, 2.E-8, .0022};
469 static const std::pair<G4int, const G4double*> Z13N14(14,pZ13N14);
470 static const std::pair<G4int, const G4double*> Z13[N13]={Z13N14};
471 //==> Si(Z=14)
472 static const G4int N14=3;
473 static const G4double pZ14N14[7]={1.2E-6, 6., 4.E-4, 6.E-9, .012, 8.E-8, .0029};
474 static const std::pair<G4int, const G4double*> Z14N14(14,pZ14N14);
475 static const G4double pZ14N15[7]={2.4E-6, 4., .0016, 3.E-9, .018, 6.E-8, .0037};
476 static const std::pair<G4int, const G4double*> Z14N15(15,pZ14N15);
477 static const G4double pZ14N16[7]={6.E-7, 4., 3.7E-4, 3.E-9, .018, 6.E-8, .0036};
478 static const std::pair<G4int, const G4double*> Z14N16(16,pZ14N16);
479 static const std::pair<G4int, const G4double*> Z14[N14]={Z14N14, Z14N15, Z14N16};
480 //==> P (Z=15)
481 static const G4int N15=1;
482 static const G4double pZ15N16[7]={6.E-7, 3., 8.2E-4, 1.4E-9, .03, 8.E-8, .0059};
483 static const std::pair<G4int, const G4double*> Z15N16(16,pZ15N16);
484 static const std::pair<G4int, const G4double*> Z15[N15]={Z15N16};
485 //==> S (Z=16)
486 static const G4int N16=4;
487 static const G4double pZ16N16[7]={6.E-7, 3., 1.9E-4, 5.E-9, .03, 6.E-8, .0013};
488 static const std::pair<G4int, const G4double*> Z16N16(16,pZ16N16);
489 static const G4double pZ16N17[7]={2.4E-6, 3., .0023, 2.E-9, .03, 6.5E-8, .004};
490 static const std::pair<G4int, const G4double*> Z16N17(17,pZ16N17);
491 static const G4double pZ16N18[7]={2.4E-6, 1.6, .0031, 1.4E-9, .03, 4.E-08, .0028};
492 static const std::pair<G4int, const G4double*> Z16N18(18,pZ16N18);
493 static const G4double pZ16N20[7]={2.4E-6, 3.1, .0017, 2.5E-9, .03, 5.E-08, .0029};
494 static const std::pair<G4int, const G4double*> Z16N20(20,pZ16N20);
495 static const std::pair<G4int, const G4double*> Z16[N16]={Z16N16, Z16N17, Z16N18, Z16N20};
496 //==> Cl(Z=17)
497 static const G4int N17=2;
498 static const G4double pZ17N18[7]={1.2E-7, .04, .062, 3.E-12, 3.E-02, 3.E-08, .027};
499 static const std::pair<G4int, const G4double*> Z17N18(18,pZ17N18);
500 static const G4double pZ17N20[7]={1.2E-7, 2., 6.8E-5, 2.7E-9, .03, 4.E-8, .0015};
501 static const std::pair<G4int, const G4double*> Z17N20(20,pZ17N20);
502 static const std::pair<G4int, const G4double*> Z17[N17]={Z17N18, Z17N20};
503 //==> Ar(Z=18)
504 static const G4int N18=3;
505 static const G4double pZ18N18[7]={1.2E-7, .52, .017, 1.1E-11, .03, 3.E-8, .095};
506 static const std::pair<G4int, const G4double*> Z18N18(18,pZ18N18);
507 static const G4double pZ18N20[7]={1.2E-07, .09, .012, 1.8E-11, .03, 3.E-8, .011};
508 static const std::pair<G4int, const G4double*> Z18N20(20,pZ18N20);
509 static const G4double pZ18N22[7]={1.2E-7, .65, 1.2E-4, 1.5E-9, .03, 5.E-8, 8.E-4};
510 static const std::pair<G4int, const G4double*> Z18N22(22,pZ18N22);
511 static const std::pair<G4int, const G4double*> Z18[N18]={Z18N18, Z18N20, Z18N22};
512 //==> K (Z=19)
513 static const G4int N19=3;
514 static const G4double pZ19N20[7]={1.2E-7, 1.3, 1.9E-4, .9E-9, .04, 5.5E-8, .0026};
515 static const std::pair<G4int, const G4double*> Z19N20(20,pZ19N20);
516 static const G4double pZ19N21[7]={1.6E-7, 1.2, 3.7E-4, .8E-9, .04, 6.5E-8, .0034};
517 static const std::pair<G4int, const G4double*> Z19N21(21,pZ19N21);
518 static const G4double pZ19N22[7]={6.E-8, 1.3, 1.2E-4, .9E-9, .04, 6.E-8, .0031};
519 static const std::pair<G4int, const G4double*> Z19N22(22,pZ19N22);
520 static const std::pair<G4int, const G4double*> Z19[N19]={Z19N20, Z19N21, Z19N22};
521 //==> Ca(Z=20)
522 static const G4int N20=6;
523 static const G4double pZ20N20[7]={2.4E-7, 3.4, 2.1E-4, 1.5E-9, .035, 6.E-8, .0037};
524 static const std::pair<G4int, const G4double*> Z20N20(20,pZ20N20);
525 static const G4double pZ20N22[7]={6.E-8, 2.7, 2.7E-5, 3.E-9, .035, 6.E-8, .0014};
526 static const std::pair<G4int, const G4double*> Z20N22(22,pZ20N22);
527 static const G4double pZ20N23[7]={1.5E-8, 1.8, 3.4E-5, .6E-9, .04, 6.E-8, .0049};
528 static const std::pair<G4int, const G4double*> Z20N23(23,pZ20N23);
529 static const G4double pZ20N24[7]={3.E-6, 5., .002, 2.E-9, .03, 7.E-8, .0038};
530 static const std::pair<G4int, const G4double*> Z20N24(24,pZ20N24);
531 static const G4double pZ20N26[7]={1.7E-5, 18., .0027, 1.E-8, 2.E-7, 7.E-17, .0047};
532 static const std::pair<G4int, const G4double*> Z20N26(26,pZ20N26);
533 static const G4double pZ20N28[7]={7.6E-6, .4, .07, .13E-9, .05, 4.E-8, .0042};
534 static const std::pair<G4int, const G4double*> Z20N28(28,pZ20N28);
535 static const std::pair<G4int, const G4double*> Z20[N20]={Z20N20, Z20N22, Z20N23,
536 Z20N24, Z20N26, Z20N28};
537 //==> Sc(Z=21)
538 static const G4int N21=1;
539 static const G4double pZ21N24[7]={3.6E-9, 1.5, 5.2E-5, .1E-9, .05, 1.E-7, .025};
540 static const std::pair<G4int, const G4double*> Z21N24(24,pZ21N24);
541 static const std::pair<G4int, const G4double*> Z21[N21]={Z21N24};
542 //==> Ti(Z=22)
543 static const G4int N22=5;
544 static const G4double pZ22N24[7]={2.8E-8, 1.8, 5.6E-5, .6E-9, .05, 8.E-8, .0042};
545 static const std::pair<G4int, const G4double*> Z22N24(24,pZ22N24);
546 static const G4double pZ22N25[7]={3.1E-9, 1.6, 6.E-6, .8E-9, .04, 8.E-8, .0036};
547 static const std::pair<G4int, const G4double*> Z22N25(25,pZ22N25);
548 static const G4double pZ22N26[7]={3.E-9, 4., 3.2E-6, 1.4E-9, .05, 2.E-7, .0048};
549 static const std::pair<G4int, const G4double*> Z22N26(26,pZ22N26);
550 static const G4double pZ22N27[7]={1.E-8, 2., 3.4E-6, 4.5E-9, .05, 8.E-8, 7.7E-4};
551 static const std::pair<G4int, const G4double*> Z22N27(27,pZ22N27);
552 static const G4double pZ22N28[7]={4.E-7, 4., 3.7E-4, 1.E-09, .05, 1.E-7, .0041};
553 static const std::pair<G4int, const G4double*> Z22N28(28,pZ22N28);
554 static const std::pair<G4int, const G4double*> Z22[N22]={Z22N24, Z22N25, Z22N26,
555 Z22N27, Z22N28};
556 //==> V (Z=23) *** Only nat *** (v50=v51=v_nat)
557 static const G4int N23=2;
558 static const G4double pZ23N27[7]={.3E-9, 2., 8.E-7, .55E-9, .07, 1.7E-7, .0055};
559 static const std::pair<G4int, const G4double*> Z23N27(27,pZ23N27);
560 static const G4double pZ23N28[7]={.3E-9, 2., 8.E-7, .55E-9, .07, 1.7E-7, .0055};
561 static const std::pair<G4int, const G4double*> Z23N28(28,pZ23N28);
562 static const std::pair<G4int, const G4double*> Z23[N23]={Z23N27, Z23N28};
563 //==> Cr(Z=24)
564 static const G4int N24=4;
565 static const G4double pZ24N26[7]={1.2E-9, 2.8, 1.E-6, 1.7E-9, .07, 1.7E-7, .0026};
566 static const std::pair<G4int, const G4double*> Z24N26(26,pZ24N26);
567 static const G4double pZ24N28[7]={4.4E-6, 11., .0012, 5.E-9, .04, 3.E-7, .0032};
568 static const std::pair<G4int, const G4double*> Z24N28(28,pZ24N28);
569 static const G4double pZ24N29[7]={1.8E-9, 2.4, 6.3E-6, .5E-9, .07, 2.E-7, .0085};
570 static const std::pair<G4int, const G4double*> Z24N29(29,pZ24N29);
571 static const G4double pZ24N30[7]={4.8E-8, 2.8, 4.4E-5, 1.4E-9, .07, 2.E-7, .0027};
572 static const std::pair<G4int, const G4double*> Z24N30(30,pZ24N30);
573 static const std::pair<G4int, const G4double*> Z24[N24]={Z24N26, Z24N28, Z24N29, Z24N30};
574 //==> Mn(Z=25)
575 static const G4int N25=1;
576 static const G4double pZ25N30[7]={6.5E-11, 1.4, 1.E-7, .8E-9, .07, 1.7E-7, .0022};
577 static const std::pair<G4int, const G4double*> Z25N30(30,pZ25N30);
578 static const std::pair<G4int, const G4double*> Z25[N25]={Z25N30};
579 //==> Fe(Z=26)
580 static const G4int N26=4;
581 static const G4double pZ26N28[7]={3.9E-8, 5., 1.7E-5, 3.E-9, .07, 3.E-7, .0023};
582 static const std::pair<G4int, const G4double*> Z26N28(28,pZ26N28);
583 static const G4double pZ26N30[7]={5.E-9, .4, 1.5E-4, 4.E-11, .1, 3.E-7, .012};
584 static const std::pair<G4int, const G4double*> Z26N30(30,pZ26N30);
585 static const G4double pZ26N31[7]={.5E-9, .5, 2.6E-6, .3E-9, .11, 5.E-7, .0027};
586 static const std::pair<G4int, const G4double*> Z26N31(31,pZ26N31);
587 static const G4double pZ26N32[7]={1.E-7, 3.1, 1.E-4, 1.3E-9, .11, 5.E-7, .0031};
588 static const std::pair<G4int, const G4double*> Z26N32(32,pZ26N32);
589 static const std::pair<G4int, const G4double*> Z26[N26]={Z26N28, Z26N30, Z26N31, Z26N32};
590 //==> Co(Z=27)
591 static const G4int N27=2;
592 static const G4double pZ27N31[7]={4.E-7, 3., .004, 0., .11, 4.5E-7, .07};
593 static const std::pair<G4int, const G4double*> Z27N31(31,pZ27N31);
594 static const G4double pZ27N32[7]={4.E-7, 5., 5.E-4, 1.2E-9, .13, 6.E-7, .006};
595 static const std::pair<G4int, const G4double*> Z27N32(32,pZ27N32);
596 static const std::pair<G4int, const G4double*> Z27[N27]={Z27N32, Z27N31};
597 //==> Ni(Z=28)
598 static const G4int N28=6;
599 static const G4double pZ28N30[7]={1.E-7, 2.5, .001, .14E-9, .13, 6.E-7, .025};
600 static const std::pair<G4int, const G4double*> Z28N30(30,pZ28N30);
601 static const G4double pZ28N31[7]={1.E-7, 19., 1.2E-5, 1.E-8, 4.E-12, 3.E-22, .0024};
602 static const std::pair<G4int, const G4double*> Z28N31(31,pZ28N31);
603 static const G4double pZ28N32[7]={1.E-8, 2.5, 3.9E-6, 3.5E-9, .13, 6.E-7, .001};
604 static const std::pair<G4int, const G4double*> Z28N32(32,pZ28N32);
605 static const G4double pZ28N33[7]={5.E-9, 2.6, 1.5E-5, .42E-9, .13, 7.E-7, .008};
606 static const std::pair<G4int, const G4double*> Z28N33(33,pZ28N33);
607 static const G4double pZ28N34[7]={.24E-9, 2., 1.2E-6, .25E-9, .13, 6.E-7, .0094};
608 static const std::pair<G4int, const G4double*> Z28N34(34,pZ28N34);
609 static const G4double pZ28N36[7]={1.E-8, 3., 5.5E-8, 2.8E-7, .12, 6.E-7, 1.6E-5};
610 static const std::pair<G4int, const G4double*> Z28N36(36,pZ28N36);
611 static const std::pair<G4int, const G4double*> Z28[N28]={Z28N30, Z28N31, Z28N32, Z28N33,
612 Z28N34, Z28N36};
613 //==> Cu(Z=29)
614 static const G4int N29=2;
615 static const G4double pZ29N34[7]={1.1E-7, 3.5, 1.6E-4, .9E-9, .13, 7.E-7, .005};
616 static const std::pair<G4int, const G4double*> Z29N34(34,pZ29N34);
617 static const G4double pZ29N36[7]={1.1E-7, 3.5, 4.3E-4, .3E-9, .13, 7.E-7, .013};
618 static const std::pair<G4int, const G4double*> Z29N36(36,pZ29N36);
619 static const std::pair<G4int, const G4double*> Z29[N29]={Z29N34, Z29N36};
620 //==> Zn(Z=30) *** Only nat *** (zn64=zn66=zn67=zn68=zn70=zn_nat)
621 static const G4int N30=5;
622 static const G4double pZ30N34[7]={1.1E-7, 4., 1.2E-4, 1.2E-9, .17, 7.E-7, .004};
623 static const std::pair<G4int, const G4double*> Z30N34(34,pZ30N34);
624 static const G4double pZ30N36[7]={1.1E-7, 4., 1.2E-4, 1.2E-9, .17, 7.E-7, .004};
625 static const std::pair<G4int, const G4double*> Z30N36(36,pZ30N36);
626 static const G4double pZ30N37[7]={1.1E-7, 4., 1.2E-4, 1.2E-9, .17, 7.E-7, .004};
627 static const std::pair<G4int, const G4double*> Z30N37(37,pZ30N37);
628 static const G4double pZ30N38[7]={1.1E-7, 4., 1.2E-4, 1.2E-9, .17, 7.E-7, .004};
629 static const std::pair<G4int, const G4double*> Z30N38(38,pZ30N38);
630 static const G4double pZ30N40[7]={1.1E-7, 4., 1.2E-4, 1.2E-9, .17, 7.E-7, .004};
631 static const std::pair<G4int, const G4double*> Z30N40(40,pZ30N40);
632 static const std::pair<G4int, const G4double*> Z30[N30]={Z30N34, Z30N36, Z30N37,
633 Z30N38, Z30N40};
634 //==> Ga(Z=31)
635 static const G4int N31=2;
636 static const G4double pZ31N38[7]={5.E-8, 3.7, 1.1E-4, .55E-9, .17, 8.4E-7, .0076};
637 static const std::pair<G4int, const G4double*> Z31N38(38,pZ31N38);
638 static const G4double pZ31N40[7]={1.E-8, 3.1, 1.7E-5, .7E-9, .17, 9.E-7, .0048};
639 static const std::pair<G4int, const G4double*> Z31N40(40,pZ31N40);
640 static const std::pair<G4int, const G4double*> Z31[N31]={Z31N38, Z31N40};
641 //==> Ge(Z=32)
642 static const G4int N32=5;
643 static const G4double pZ32N38[7]={5.E-5, 4., .17, .35E-9, .17, 9.E-7, .013};
644 static const std::pair<G4int, const G4double*> Z32N38(38,pZ32N38);
645 static const G4double pZ32N40[7]={5.E-7, 4.4, .001, .6E-9, .17, 9.E-7, .008};
646 static const std::pair<G4int, const G4double*> Z32N40(40,pZ32N40);
647 static const G4double pZ32N41[7]={5.E-9, 3., 8.E-6, .7E-9, .17, 1.E-6, .0043};
648 static const std::pair<G4int, const G4double*> Z32N41(41,pZ32N41);
649 static const G4double pZ32N42[7]={1.E-7, 4.2, 1.7E-4, .7E-9, .17, 1.E-6, .0065};
650 static const std::pair<G4int, const G4double*> Z32N42(42,pZ32N42);
651 static const G4double pZ32N44[7]={1.E-6, 4.6, .0018, .6E-9, .17, 1.E-6, .0073};
652 static const std::pair<G4int, const G4double*> Z32N44(44,pZ32N44);
653 static const std::pair<G4int, const G4double*> Z32[N32]={Z32N38, Z32N40, Z32N41,
654 Z32N42, Z32N44};
655 //==> As(Z=33)
656 static const G4int N33=2;
657 static const G4double pZ33N41[7]={1.E-8, 3.4, 1.5E-5, .72E-9, .17, 1.E-6, .0045};
658 static const std::pair<G4int, const G4double*> Z33N41(41,pZ33N41);
659 static const G4double pZ33N42[7]={1.E-8, 4.1, 1.3E-5, .75E-9, .2, 1.2E-6, .0048};
660 static const std::pair<G4int, const G4double*> Z33N42(42,pZ33N42);
661 static const std::pair<G4int, const G4double*> Z33[N33]={Z33N41, Z33N42};
662 //==> Se(Z=34)
663 static const G4int N34=7;
664 static const G4double pZ34N40[7]={6.E-8, 7.2, 6.E-5, 1.E-9, .32, 2.E-6, .0063};
665 static const std::pair<G4int, const G4double*> Z34N40(40,pZ34N40);
666 static const G4double pZ34N42[7]={4.E-5, 7.4, .1, .43E-9, .34, 2.1E-6, .016};
667 static const std::pair<G4int, const G4double*> Z34N42(42,pZ34N42);
668 static const G4double pZ34N43[7]={1.E-7, 6.2, 1.4E-4, .9E-9, .34, 2.1E-6, .0075};
669 static const std::pair<G4int, const G4double*> Z34N43(43,pZ34N43);
670 static const G4double pZ34N44[7]={1.E-7, 6.6, 1.3E-4, .9E-9, .34, 2.1E-6, .0075};
671 static const std::pair<G4int, const G4double*> Z34N44(44,pZ34N44);
672 static const G4double pZ34N45[7]={5.E-8, 6.6, 4.8E-5, 1.2E-9, .4, 2.6E-6, .0055};
673 static const std::pair<G4int, const G4double*> Z34N45(45,pZ34N45);
674 static const G4double pZ34N46[7]={2.E-7, 7.7, 1.3E-4, 1.7E-9, .34, 2.1E-6, .0043};
675 static const std::pair<G4int, const G4double*> Z34N46(46,pZ34N46);
676 static const G4double pZ34N48[7]={2.E-7, 8.3, 1.2E-4, 1.7E-9, .34, 2.1E-6, .0043};
677 static const std::pair<G4int, const G4double*> Z34N48(48,pZ34N48);
678 static const std::pair<G4int, const G4double*> Z34[N34]={Z34N40, Z34N42, Z34N43, Z34N44,
679 Z34N45, Z34N46, Z34N48};
680 //==> Br(Z=35)
681 static const G4int N35=2;
682 static const G4double pZ35N44[7]={5.E-8, 6., 2.8E-5, 2.E-9, .34, 2.1E-6, .0028};
683 static const std::pair<G4int, const G4double*> Z35N44(44,pZ35N44);
684 static const G4double pZ35N46[7]={4.E-8, 6.2, 3.7E-5, 1.1E-9, .34, 2.1E-6, .0049};
685 static const std::pair<G4int, const G4double*> Z35N46(46,pZ35N46);
686 static const std::pair<G4int, const G4double*> Z35[N35]={Z35N44, Z35N46};
687 //==> Kr(Z=36)
688 static const G4int N36=7;
689 static const G4double pZ36N42[7]={1.6E-7, 6.8, 2.E-4, .8E-9, .35, 2.1E-6, .0076};
690 static const std::pair<G4int, const G4double*> Z36N42(42,pZ36N42);
691 static const G4double pZ36N44[7]={1.6E-7, 7.3, 1.6E-4, 1.E-9, .35, 2.1E-6, .0062};
692 static const std::pair<G4int, const G4double*> Z36N44(44,pZ36N44);
693 static const G4double pZ36N46[7]={1.6E-7, 7.3, 3.3E-4, .7E-9, .35, 2.1E-6, .013};
694 static const std::pair<G4int, const G4double*> Z36N46(46,pZ36N46);
695 static const G4double pZ36N47[7]={1.6E-6, 7.3, .003, .6E-9, .35, 2.1E-6, .011};
696 static const std::pair<G4int, const G4double*> Z36N47(47,pZ36N47);
697 static const G4double pZ36N48[7]={1.6E-7, 7.8, 7.6E-5, 2.E-9, .35, 2.1E-6, .0031};
698 static const std::pair<G4int, const G4double*> Z36N48(48,pZ36N48);
699 static const G4double pZ36N49[7]={6.E-7, 8., 4.8E-4, 1.4E-9, .27, 2.1E-6, .0053};
700 static const std::pair<G4int, const G4double*> Z36N49(49,pZ36N49);
701 static const G4double pZ36N50[7]={4.E-7, 8.1, 2.7E-4, 1.6E-9, .35, 2.1E-6, .0045};
702 static const std::pair<G4int, const G4double*> Z36N50(50,pZ36N50);
703 static const std::pair<G4int, const G4double*> Z36[N36]={Z36N42, Z36N44, Z36N46,
704 Z36N47, Z36N48, Z36N49, Z36N50};
705 //==> Rb(Z=37)
706 static const G4int N37=3;
707 static const G4double pZ37N48[7]={1.6E-7, 7.2, 1.4E-4, 1.2E-9, .35, 2.1E-6, .0052};
708 static const std::pair<G4int, const G4double*> Z37N48(48,pZ37N48);
709 static const G4double pZ37N49[7]={8.E-8, 7.1, 4.7E-5, 1.6E-9, .27, 2.1E-6, .0034};
710 static const std::pair<G4int, const G4double*> Z37N49(49,pZ37N49);
711 static const G4double pZ37N50[7]={1.E-7, 8., 5.5E-5, 1.9E-9, .27, 1.5E-6, .0036};
712 static const std::pair<G4int, const G4double*> Z37N50(50,pZ37N50);
713 static const std::pair<G4int, const G4double*> Z37[N37]={Z37N48, Z37N49, Z37N50};
714 //==> Sr(Z=38)
715 static const G4int N38=6;
716 static const G4double pZ38N46[7]={8.E-8, 7.3, 6.E-5, 1.3E-9, .27, 2.E-6, .0045};
717 static const std::pair<G4int, const G4double*> Z38N46(46,pZ38N46);
718 static const G4double pZ38N48[7]={8.E-8, 9.7, 2.3E-5, 3.5E-9, .4, 2.7E-6, .0023};
719 static const std::pair<G4int, const G4double*> Z38N48(48,pZ38N48);
720 static const G4double pZ38N49[7]={2.6E-7, 9.5, 1.9E-4, 1.5E-9, .4, 2.7E-6, .0057};
721 static const std::pair<G4int, const G4double*> Z38N49(49,pZ38N49);
722 static const G4double pZ38N50[7]={2.6E-7, 9.5, 2.E-4, 1.4E-9, .37, 3.2E-6, .0059};
723 static const std::pair<G4int, const G4double*> Z38N50(50,pZ38N50);
724 static const G4double pZ38N51[7]={1.3E-7, 9.9, 7.5E-5, 1.7E-9, .37, 3.2E-6, .0046};
725 static const std::pair<G4int, const G4double*> Z38N51(51,pZ38N51);
726 static const G4double pZ38N52[7]={2.6E-7, 9.6, 1.6E-4, 1.8E-9, .37, 2.7E06, .0047};
727 static const std::pair<G4int, const G4double*> Z38N52(52,pZ38N52);
728 static const std::pair<G4int, const G4double*> Z38[N38]={Z38N46, Z38N48, Z38N49, Z38N50,
729 Z38N51, Z38N52};
730 //==> Y (Z=39)
731 static const G4int N39=3;
732 static const G4double pZ39N50[7]={2.6E-7, 9.9, 2.E-4, 1.1E-9, .37, 3.2E-6, .0062};
733 static const std::pair<G4int, const G4double*> Z39N50(50,pZ39N50);
734 static const G4double pZ39N51[7]={2.7E-5, 20., .013, 2.2E-9, .37, 1.9E-6, .0078};
735 static const std::pair<G4int, const G4double*> Z39N51(51,pZ39N51);
736 static const G4double pZ39N52[7]={2.E-7, 9.6, 1.2E-4, 1.7E-9, .37, 3.E-6, .0046};
737 static const std::pair<G4int, const G4double*> Z39N52(52,pZ39N52);
738 static const std::pair<G4int, const G4double*> Z39[N39]={Z39N50, Z39N51, Z39N52};
739 //==> Zr(Z=40)
740 static const G4int N40=7;
741 static const G4double pZ40N50[7]={1.E-7, 9., 6.2E-5, 1.7E-9, .3, 3.E-6, .0044};
742 static const std::pair<G4int, const G4double*> Z40N50(50,pZ40N50);
743 static const G4double pZ40N51[7]={5.E-7, 9.8, 5.E-4, 1.E-9, .25, 2.1E-6, .0079};
744 static const std::pair<G4int, const G4double*> Z40N51(51,pZ40N51);
745 static const G4double pZ40N52[7]={3.E-7, 9.6, 2.2E-4, 1.2E-9, .25, 2.2E-6, .0056};
746 static const std::pair<G4int, const G4double*> Z40N52(52,pZ40N52);
747 static const G4double pZ40N53[7]={2.E-7, 9.6, 1.2E-4, 1.6E-9, .38, 2.9E-6, .0046};
748 static const std::pair<G4int, const G4double*> Z40N53(53,pZ40N53);
749 static const G4double pZ40N54[7]={2.E-7, 9.6, 1.8E-4, 1.1E-9, .25, 2.1E-6, .0067};
750 static const std::pair<G4int, const G4double*> Z40N54(54,pZ40N54);
751 static const G4double pZ40N55[7]={1.8E-7, 9.4, 1.1E-4, 1.6E-9, .33, 2.3E-6, .0045};
752 static const std::pair<G4int, const G4double*> Z40N55(55,pZ40N55);
753 static const G4double pZ40N56[7]={1.8E-7, 9.4, 1.1E-4, 1.6E-9, .2, 1.5E-6, .0045};
754 static const std::pair<G4int, const G4double*> Z40N56(56,pZ40N56);
755 static const std::pair<G4int, const G4double*> Z40[N40]={Z40N50, Z40N51, Z40N52, Z40N53,
756 Z40N54, Z40N55, Z40N56};
757 //==> Nb(Z=41)
758 static const G4int N41=3;
759 static const G4double pZ41N52[7]={2.6E-7, 8.3, 2.E-4, 1.2E-9, .4, 4.E-6, .0051};
760 static const std::pair<G4int, const G4double*> Z41N52(52,pZ41N52);
761 static const G4double pZ41N53[7]={2.E-7, 8.3, 1.6E-4, 1.4E-9, .35, 2.5E-6, .0051};
762 static const std::pair<G4int, const G4double*> Z41N53(53,pZ41N53);
763 static const G4double pZ41N54[7]={1.5E-7, 8.6, 1.E-4, 1.5E-9, .35, 2.5E-6, .0045};
764 static const std::pair<G4int, const G4double*> Z41N54(54,pZ41N54);
765 static const std::pair<G4int, const G4double*> Z41[N41]={Z41N52, Z41N52, Z41N52};
766 //==> Mo(Z=42)
767 static const G4int N42=8;
768 static const G4double pZ42N50[7]={2.E-7, 10., 1.1E-4, 1.8E-9, .3, 2.7E-6, .0044};
769 static const std::pair<G4int, const G4double*> Z42N50(50,pZ42N50);
770 static const G4double pZ42N52[7]={2.1E-7, 10., 1.2E-4, 1.7E-9, .3, 2.8E-6, .0046};
771 static const std::pair<G4int, const G4double*> Z42N52(52,pZ42N52);
772 static const G4double pZ42N53[7]={3.E-7, 10., 1.9E-4, 1.5E-9, .29, 3.E-6, .005};
773 static const std::pair<G4int, const G4double*> Z42N53(53,pZ42N53);
774 static const G4double pZ42N54[7]={1.5E-7, 10., 7.1E-5, 2.1E-9, .29, 2.9E-6, .0037};
775 static const std::pair<G4int, const G4double*> Z42N54(54,pZ42N54);
776 static const G4double pZ42N55[7]={1.9E-7, 9., 1.4E-4, 1.3E-9, .29, 2.8E-6, .0052};
777 static const std::pair<G4int, const G4double*> Z42N55(55,pZ42N55);
778 static const G4double pZ42N56[7]={1.9E-7, 9.9, 1.1E-4, 1.8E-9, .29, 2.4E-6, .0044};
779 static const std::pair<G4int, const G4double*> Z42N56(56,pZ42N56);
780 static const G4double pZ42N57[7]={1.4E-7, 8., 1.E-4, 1.4E-9, .34, 2.5E-6, .0044};
781 static const std::pair<G4int, const G4double*> Z42N57(57,pZ42N57);
782 static const G4double pZ42N58[7]={1.8E-7, 9.5, 1.E-4, 1.7E-9, .27, 2.2E-6, .0041};
783 static const std::pair<G4int, const G4double*> Z42N58(58,pZ42N58);
784 static const std::pair<G4int, const G4double*> Z42[N42]={Z42N50, Z42N52, Z42N53, Z42N54,
785 Z42N55, Z42N56, Z42N57, Z42N58};
786 //==> Tc(Z=43)
787 static const G4int N43=1;
788 static const G4double pZ43N56[7]={1.E-7, 8., 7.2E-5, 1.4E-9, .24, 2.5E-6, .0044};
789 static const std::pair<G4int, const G4double*> Z43N56(56,pZ43N56);
790 static const std::pair<G4int, const G4double*> Z43[N43]={Z43N56};
791 //==> Ru(Z=44)
792 static const G4int N44=9;
793 static const G4double pZ44N52[7]={1.9E-7, 10., 1.E-4, 2.1E-9, .4, 4.E-6, .004};
794 static const std::pair<G4int, const G4double*> Z44N52(52,pZ44N52);
795 static const G4double pZ44N54[7]={1.5E-7, 10., 7.7E-5, 2.1E-9, .29, 2.9E-6, .004};
796 static const std::pair<G4int, const G4double*> Z44N54(54,pZ44N54);
797 static const G4double pZ44N55[7]={1.8E-7, 10., 6.6E-5, 2.6E-9, .47, 4.6E-6, .0028};
798 static const std::pair<G4int, const G4double*> Z44N55(55,pZ44N55);
799 static const G4double pZ44N56[7]={1.8E-6, 10., .0017, 1.1E-9, .47, 6.E-6, .0073};
800 static const std::pair<G4int, const G4double*> Z44N56(56,pZ44N56);
801 static const G4double pZ44N57[7]={1.8E-7, 7.8, 1.3E-4, 1.4E-9, .42, 5.E-6, .0043};
802 static const std::pair<G4int, const G4double*> Z44N57(57,pZ44N57);
803 static const G4double pZ44N58[7]={1.7E-6, 9.8, .0015, 1.2E-9, .32, 4.E-6, .0065};
804 static const std::pair<G4int, const G4double*> Z44N58(58,pZ44N58);
805 static const G4double pZ44N59[7]={3.3E-7, 8.7, 1.9E-4, 1.6E-9, .32, 3.8E-6, .0038};
806 static const std::pair<G4int, const G4double*> Z44N59(59,pZ44N59);
807 static const G4double pZ44N60[7]={3.E-7, 8.7, 1.8E-4, 1.6E-9, .3, 3.2E-6, .004};
808 static const std::pair<G4int, const G4double*> Z44N60(60,pZ44N60);
809 static const G4double pZ44N61[7]={3.E-7, 8.8, 1.4E-4, 1.9E-9, .3, 3.2E-6, .003};
810 static const std::pair<G4int, const G4double*> Z44N61(61,pZ44N61);
811 static const std::pair<G4int, const G4double*> Z44[N44]={Z44N52, Z44N54, Z44N55, Z44N56,
812 Z44N57, Z44N58, Z44N59, Z44N60,
813 Z44N61};
814 //==> Rh(Z=45)
815 static const G4int N45=2;
816 static const G4double pZ45N58[7]={8.E-8, 8.7, 4.E-5, 1.8E-9, .29, 2.9E-6, .0033};
817 static const std::pair<G4int, const G4double*> Z45N58(58,pZ45N58);
818 static const G4double pZ45N60[7]={8.E-8, 8.7, .09, 1.3E-12, .29, 2.9E-6, 7.};
819 static const std::pair<G4int, const G4double*> Z45N60(60,pZ45N60);
820 static const std::pair<G4int, const G4double*> Z45[N45]={Z45N58, Z45N60};
821 //==> Pd(Z=46)
822 static const G4int N46=7;
823 static const G4double pZ46N56[7]={2.E-7, 9.9, 1.2E-4, 1.5E-9, .35, 3.3E-6, .0045};
824 static const std::pair<G4int, const G4double*> Z46N56(56,pZ46N56);
825 static const G4double pZ46N58[7]={2.E-7, 9.9, 9.5E-5, 1.8E-9, .4, 4.E-6, .0036};
826 static const std::pair<G4int, const G4double*> Z46N58(58,pZ46N58);
827 static const G4double pZ46N59[7]={5.6E-7, 9., 4.6E-4, 1.2E-9, .5, 4.8E-6, .0056};
828 static const std::pair<G4int, const G4double*> Z46N59(59,pZ46N59);
829 static const G4double pZ46N60[7]={2.4E-7, 9.2, 1.2E-4, 1.8E-9, .47, 4.6E-6, .0035};
830 static const std::pair<G4int, const G4double*> Z46N60(60,pZ46N60);
831 static const G4double pZ46N61[7]={1.2E-7, 9.2, 4.4E-5, 2.8E-9, .5, 4.3E-6, .0025};
832 static const std::pair<G4int, const G4double*> Z46N61(61,pZ46N61);
833 static const G4double pZ46N62[7]={1.2E-7, 9.2, 3.2E-5, 3.4E-9, .48, 4.5E-6, .0018};
834 static const std::pair<G4int, const G4double*> Z46N62(62,pZ46N62);
835 static const G4double pZ46N64[7]={4.E-7, 9.1, 2.5E-4, 1.5E-9, .48, 4.7E-6, .0042};
836 static const std::pair<G4int, const G4double*> Z46N64(64,pZ46N64);
837 static const std::pair<G4int, const G4double*> Z46[N46]={Z46N56, Z46N58, Z46N59,
838 Z46N60, Z46N61, Z46N62, Z46N64};
839 //==> Ag(Z=47)
840 static const G4int N47=4;
841 static const G4double pZ47N60[7]={1.4E-6, 9.7, .0011, 1.4E-9, .55, 5.E-6, .0056};
842 static const std::pair<G4int, const G4double*> Z47N60(60,pZ47N60);
843 static const G4double pZ47N62[7]={3.E-8, 8.7, 8.5E-6, 3.5E-9, .6, 5.6E-6, .0018};
844 static const std::pair<G4int, const G4double*> Z47N62(62,pZ47N62);
845 static const G4double pZ47N63[7]={3.E-6, 9.5, .002, 1.5E-9, .58, 5.E-6, .0047};
846 static const std::pair<G4int, const G4double*> Z47N63(63,pZ47N63);
847 static const G4double pZ47N64[7]={1.5E-7, 9., 9.E-5, 1.7E-9, .58, 5.6E-6, .0039};
848 static const std::pair<G4int, const G4double*> Z47N64(64,pZ47N64);
849 static const std::pair<G4int, const G4double*> Z47[N47]={Z47N60, Z47N62, Z47N63, Z47N64};
850 //==> Cd(Z=48)
851 static const G4int N48=9;
852 static const G4double pZ48N58[7]={2.9E-7, 10., 1.3E-4, 1.9E-9, .4, 3.8E-6, .0034};
853 static const std::pair<G4int, const G4double*> Z48N58(58,pZ48N58);
854 static const G4double pZ48N60[7]={2.3E-7, 10., 8.2E-5, 2.5E-9, .5, 4.7E-6, .0026};
855 static const std::pair<G4int, const G4double*> Z48N60(60,pZ48N60);
856 static const G4double pZ48N62[7]={2.3E-7, 10., 9.9E-5, 2.5E-9, .5, 4.7E-6, .0031};
857 static const std::pair<G4int, const G4double*> Z48N62(62,pZ48N62);
858 static const G4double pZ48N63[7]={8.4E-7, 11., 4.3E-4, 1.8E-9, .5, 4.5E-6, .0042};
859 static const std::pair<G4int, const G4double*> Z48N63(63,pZ48N63);
860 static const G4double pZ48N64[7]={4.E-7, 11., 1.8E-4, 1.8E-9, .5, 4.6E-6, .0036};
861 static const std::pair<G4int, const G4double*> Z48N64(64,pZ48N64);
862 static const G4double pZ48N65[7]={1.6E-6, 12., .001, 0., .5, 4.6E-6, .013};
863 static const std::pair<G4int, const G4double*> Z48N65(65,pZ48N65);
864 static const G4double pZ48N66[7]={3.E-7, 11., 1.2E-4, 1.9E-9, .5, 4.6E-6, .0031};
865 static const std::pair<G4int, const G4double*> Z48N66(66,pZ48N66);
866 static const G4double pZ48N67[7]={3.8E-7, 11., 1.7E-4, 2.E-9, .6, 6.6E-6, .0035};
867 static const std::pair<G4int, const G4double*> Z48N67(67,pZ48N67);
868 static const G4double pZ48N68[7]={6.E-7, 11., 3.3E-4, 1.9E-9, .5, 4.6E-6, .0043};
869 static const std::pair<G4int, const G4double*> Z48N68(68,pZ48N68);
870 static const std::pair<G4int, const G4double*> Z48[N48]={Z48N58, Z48N60, Z48N62, Z48N63,
871 Z48N64, Z48N65, Z48N66, Z48N67,
872 Z48N68};
873 //==> In(Z=49)
874 static const G4int N49=2;
875 static const G4double pZ49N64[7]={2.7E-7, 12., 8.1E-5, 2.7E-9, .5, 5.E-6, .0026};
876 static const std::pair<G4int, const G4double*> Z49N64(64,pZ49N64);
877 static const G4double pZ49N66[7]={2.7E-7, 12., 5.5E-5, 4.E-9, .5, 5.E-6, .0018};
878 static const std::pair<G4int, const G4double*> Z49N66(66,pZ49N66);
879 static const std::pair<G4int, const G4double*> Z49[N49]={Z49N64, Z49N66};
880 //==> Sn(Z=50)
881 static const G4int N50=14;
882 static const G4double pZ50N62[7]={4.E-7, 11., 1.6E-4, 2.2E-9, .5, 4.5E-6, .0032};
883 static const std::pair<G4int, const G4double*> Z50N62(62,pZ50N62);
884 static const G4double pZ50N63[7]={4.1E-7, 11., 1.6E-4, 2.4E-9, .54, 6.E-6, .0031};
885 static const std::pair<G4int, const G4double*> Z50N63(63,pZ50N63);
886 static const G4double pZ50N64[7]={5.E-7, 12., 1.9E-4, 2.2E-9, .5, 4.4E-6, .0032};
887 static const std::pair<G4int, const G4double*> Z50N64(64,pZ50N64);
888 static const G4double pZ50N65[7]={1.E-5, 12., .0077, 1.4E-9, .5, 5.E-6, .0066};
889 static const std::pair<G4int, const G4double*> Z50N65(65,pZ50N65);
890 static const G4double pZ50N66[7]={5.E-7, 12., 1.8E-4, 2.4E-9, .5, 5.E-6, .0031};
891 static const std::pair<G4int, const G4double*> Z50N66(66,pZ50N66);
892 static const G4double pZ50N67[7]={1.E-6, 12., 4.4E-4, 1.8E-9, .5, 5.E-6, .0037};
893 static const std::pair<G4int, const G4double*> Z50N67(67,pZ50N67);
894 static const G4double pZ50N68[7]={5.E-7, 12., 2.E-4, 2.4E-9, .5, 5.E-6, .0033};
895 static const std::pair<G4int, const G4double*> Z50N68(68,pZ50N68);
896 static const G4double pZ50N69[7]={6.E-7, 12., 2.5E-4, 2.E-9, .5, 5.E-6, .0035};
897 static const std::pair<G4int, const G4double*> Z50N69(69,pZ50N69);
898 static const G4double pZ50N70[7]={1.E-6, 12., 4.7E-4, 2.E-9, .5, 5.E-6, .0039};
899 static const std::pair<G4int, const G4double*> Z50N70(70,pZ50N70);
900 static const G4double pZ50N72[7]={1.E-6, 12., 3.7E-4, 2.2E-9, .5, 5.E-6, .0031};
901 static const std::pair<G4int, const G4double*> Z50N72(72,pZ50N72);
902 static const G4double pZ50N73[7]={5.E-7, 12., 1.7E-4, 2.8E-9, .5, 5.E-6, .0028};
903 static const std::pair<G4int, const G4double*> Z50N73(73,pZ50N73);
904 static const G4double pZ50N74[7]={5.E-7, 12., 2.E-4, 2.E-9, .5, 5.E-6, .0033};
905 static const std::pair<G4int, const G4double*> Z50N74(74,pZ50N74);
906 static const G4double pZ50N75[7]={5.E-7, 12., 1.9E-4, 2.8E-9, .5, 5.E-6, .003};
907 static const std::pair<G4int, const G4double*> Z50N75(75,pZ50N75);
908 static const G4double pZ50N76[7]={5.E-7, 12., 1.7E-4, 2.8E-9, .5, 5.E-6, .0028};
909 static const std::pair<G4int, const G4double*> Z50N76(76,pZ50N76);
910 static const std::pair<G4int, const G4double*> Z50[N50]={Z50N62, Z50N63, Z50N64, Z50N65,
911 Z50N66, Z50N67, Z50N68, Z50N69,
912 Z50N70, Z50N72, Z50N73, Z50N74,
913 Z50N75, Z50N76};
914 //==> Sb(Z=51)
915 static const G4int N51=5;
916 static const G4double pZ51N70[7]={6.E-7, 12., 2.E-4, 2.8E-9, .5, 5.E-6, .0028};
917 static const std::pair<G4int, const G4double*> Z51N70(70,pZ51N70);
918 static const G4double pZ51N72[7]={6.E-7, 12., 1.9E-4, 3.E-9, .5, 5.E-6, .0025};
919 static const std::pair<G4int, const G4double*> Z51N72(72,pZ51N72);
920 static const G4double pZ51N73[7]={1.1E-6, 12., 3.5E-4, 2.9E-9, .5, 5.E-6, .0026};
921 static const std::pair<G4int, const G4double*> Z51N73(73,pZ51N73);
922 static const G4double pZ51N74[7]={5.5E-7, 12., 1.9E-4, 2.9E-9, .5, 5.E-6, .0027};
923 static const std::pair<G4int, const G4double*> Z51N74(74,pZ51N74);
924 static const G4double pZ51N75[7]={6.E-7, 12., 2.E-4, 2.9E-9, .5, 5.E-6, .0027};
925 static const std::pair<G4int, const G4double*> Z51N75(75,pZ51N75);
926 static const std::pair<G4int, const G4double*> Z51[N51]={Z51N70, Z51N72, Z51N73, Z51N74,
927 Z51N75};
928 //==> Te(Z=52)
929 static const G4int N52=11;
930 static const G4double pZ52N68[7]={2.7E-7, 12., 8.4E-5, 3.2E-9, 1., 8.E-6, .0026};
931 static const std::pair<G4int, const G4double*> Z52N68(68,pZ52N68);
932 static const G4double pZ52N70[7]={2.7E-7, 12., 3.8E-5, 6.E-9, 1., 8.E-6, .0012};
933 static const std::pair<G4int, const G4double*> Z52N70(70,pZ52N70);
934 static const G4double pZ52N71[7]={2.7E-8, 12., 1.8E-6, 2.E-8, 1., 8.E-6, 4.8E-4};
935 static const std::pair<G4int, const G4double*> Z52N71(71,pZ52N71);
936 static const G4double pZ52N72[7]={2.6E-6, 14., .0014, 2.E-9, 1., 9.E-6, .005};
937 static const std::pair<G4int, const G4double*> Z52N72(72,pZ52N72);
938 static const G4double pZ52N73[7]={1.E-6, 14., 2.4E-4, 3.9E-9, 1., 9.E-6, .0022};
939 static const std::pair<G4int, const G4double*> Z52N73(73,pZ52N73);
940 static const G4double pZ52N74[7]={8.E-7, 14., 2.3E-4, 3.6E-9, 1.4, 1.3E-5, .0028};
941 static const std::pair<G4int, const G4double*> Z52N74(74,pZ52N74);
942 static const G4double pZ52N75[7]={8.E-7, 14., 2.1E-4, 3.6E-9, 1.4, 1.3E-5, .0025};
943 static const std::pair<G4int, const G4double*> Z52N75(75,pZ52N75);
944 static const G4double pZ52N76[7]={8.E-7, 14., 2.5E-4, 3.E-9, 1.4, 1.3E-5, .003};
945 static const std::pair<G4int, const G4double*> Z52N76(76,pZ52N76);
946 static const G4double pZ52N77[7]={5.E-7, 15., 1.2E-4, 4.3E-9, 1.4, 1.4E-5, .0023};
947 static const std::pair<G4int, const G4double*> Z52N77(77,pZ52N77);
948 static const G4double pZ52N78[7]={8.E-7, 14., 2.7E-4, 2.7E-9, 1.4, 1.3E-5, .0031};
949 static const std::pair<G4int, const G4double*> Z52N78(78,pZ52N78);
950 static const G4double pZ52N80[7]={4.7E-7, 14., 1.8E-4, 2.2E-9, .83, 1.E-5, .0036};
951 static const std::pair<G4int, const G4double*> Z52N80(80,pZ52N80);
952 static const std::pair<G4int, const G4double*> Z52[N52]={Z52N68, Z52N70, Z52N71, Z52N72,
953 Z52N73, Z52N74, Z52N75, Z52N76,
954 Z52N77, Z52N78, Z52N80};
955 //==> I (Z=53)
956 static const G4int N53=5;
957 static const G4double pZ53N74[7]={9.4E-7, 14., 2.5E-4, 3.E-9, .7, 7.3E-6, .0025};
958 static const std::pair<G4int, const G4double*> Z53N74(74,pZ53N74);
959 static const G4double pZ53N76[7]={2.1E-5, 14., .015, 1.1E-9, 1.1, 1.E-5, .007};
960 static const std::pair<G4int, const G4double*> Z53N76(76,pZ53N76);
961 static const G4double pZ53N77[7]={1.1E-6, 14., 2.4E-4, 3.3E-9, .9, 1.E-5, .0021};
962 static const std::pair<G4int, const G4double*> Z53N77(77,pZ53N77);
963 static const G4double pZ53N78[7]={5.5E-7, 14., 1.5E-4, 3.7E-9, 1.2, 1.2E-5, .0024};
964 static const std::pair<G4int, const G4double*> Z53N78(78,pZ53N78);
965 static const G4double pZ53N82[7]={3.2E-6, 14., .0017, 1.8E-9, .8, 8.E-6, .0024};
966 static const std::pair<G4int, const G4double*> Z53N82(82,pZ53N82);
967 static const std::pair<G4int, const G4double*> Z53[N53]={Z53N74, Z53N76, Z53N77, Z53N78,
968 Z53N82};
969 //==> Xe(Z=54)
970 static const G4int N54=12;
971 static const G4double pZ54N69[7]={3.E-6, 14., 8.E-4, 3.7E-9, .9, 1.1E-5, .15};
972 static const std::pair<G4int, const G4double*> Z54N69(69,pZ54N69);
973 static const G4double pZ54N70[7]={1.5E-7, 14., 1.4E-6, 9.E-8, .7, 8.E-6, 9.5E-5};
974 static const std::pair<G4int, const G4double*> Z54N70(70,pZ54N70);
975 static const G4double pZ54N72[7]={1.5E-6, 14., 5.6E-4, 3.E-9, 1.2, 1.1E-5, .0036};
976 static const std::pair<G4int, const G4double*> Z54N72(72,pZ54N72);
977 static const G4double pZ54N74[7]={1.8E-6, 14., 8.8E-4, 2.E-9, 1.3, 1.2E-5, .0047};
978 static const std::pair<G4int, const G4double*> Z54N74(74,pZ54N74);
979 static const G4double pZ54N75[7]={1.E-6, 14., 2.6E-4, 3.7E-9, 1.5, 1.4E-5, .0024};
980 static const std::pair<G4int, const G4double*> Z54N75(75,pZ54N75);
981 static const G4double pZ54N76[7]={1.8E-6, 14., 8.E-4, 2.E-9, 1.2, 1.4E-5, .0042};
982 static const std::pair<G4int, const G4double*> Z54N76(76,pZ54N76);
983 static const G4double pZ54N77[7]={2.3E-7, 14., 1.9E-5, 9.E-9, 1.2, 1.4E-5, 7.7E-4};
984 static const std::pair<G4int, const G4double*> Z54N77(77,pZ54N77);
985 static const G4double pZ54N78[7]={6.E-7, 14., 1.6E-4, 3.E-9, 1.2, 1.4E-5, .0025};
986 static const std::pair<G4int, const G4double*> Z54N78(78,pZ54N78);
987 static const G4double pZ54N79[7]={6.E-7, 14., 1.6E-4, 3.3E-9, 1.6, 1.5E-5, .0024};
988 static const std::pair<G4int, const G4double*> Z54N79(79,pZ54N79);
989 static const G4double pZ54N80[7]={6.6E-7, 14., 2.1E-4, 2.5E-9, 1.2, 1.4E-5, .003};
990 static const std::pair<G4int, const G4double*> Z54N80(80,pZ54N80);
991 static const G4double pZ54N81[7]={.03, 40., 2.1, 2.5E-9, 1.E-16, 6.E-36, 140.};
992 static const std::pair<G4int, const G4double*> Z54N81(81,pZ54N81);
993 static const G4double pZ54N82[7]={3.1E-6, 14., .0019, 1.6E-9, 1., 1.3E-5, .0054};
994 static const std::pair<G4int, const G4double*> Z54N82(82,pZ54N82);
995 static const std::pair<G4int, const G4double*> Z54[N54]={Z54N69, Z54N70, Z54N72, Z54N74,
996 Z54N75, Z54N76, Z54N77, Z54N78,
997 Z54N79, Z54N80, Z54N81, Z54N82};
998 //==> Cs(Z=55)
999 static const G4int N55=5;
1000 static const G4double pZ55N78[7]={1.4E-6, 14., 4.E-4, 3.E-9, 1.2, 1.4E-5, .0026};
1001 static const std::pair<G4int, const G4double*> Z55N78(78,pZ55N78);
1002 static const G4double pZ55N79[7]={.028, 14., 44., .5E-9, 1.2, 1.3E-5, .015};
1003 static const std::pair<G4int, const G4double*> Z55N79(79,pZ55N79);
1004 static const G4double pZ55N80[7]={2.E-6, 14., 9.5E-4, 2.E-9, 1.2, 1.4E-5, .0042};
1005 static const std::pair<G4int, const G4double*> Z55N80(80,pZ55N80);
1006 static const G4double pZ55N81[7]={2.5E-7, 14., 6.5E-5, 3.8E-9, 1.2, 1.4E-5, .0023};
1007 static const std::pair<G4int, const G4double*> Z55N81(81,pZ55N81);
1008 static const G4double pZ55N82[7]={2.5E-7, 14., 6.5E-5, 3.8E-9, 1.2, 1.4E-5, .0023};
1009 static const std::pair<G4int, const G4double*> Z55N82(82,pZ55N82);
1010 static const std::pair<G4int, const G4double*> Z55[N55]={Z55N78, Z55N79, Z55N80, Z55N81,
1011 Z55N82};
1012 //==> Ba(Z=56)
1013 static const G4int N56=9;
1014 static const G4double pZ56N74[7]={4.E-7, 14., 2.8E-5, 1.2E-8, 1.2, 1.5E-5, 6.6E-4};
1015 static const std::pair<G4int, const G4double*> Z56N74(74,pZ56N74);
1016 static const G4double pZ56N76[7]={4.E-6, 14., .0022, 1.4E-9, 1.3, 1.6E-5, .0053};
1017 static const std::pair<G4int, const G4double*> Z56N76(76,pZ56N76);
1018 static const G4double pZ56N77[7]={2.E-7, 14., 3.7E-5, 5.E-9, 1.1, 1.5E-5, .0016};
1019 static const std::pair<G4int, const G4double*> Z56N77(77,pZ56N77);
1020 static const G4double pZ56N78[7]={1.6E-6, 14., 6.E-4, 2.E-9, 1.3, 1.6E-5, .0033};
1021 static const std::pair<G4int, const G4double*> Z56N78(78,pZ56N78);
1022 static const G4double pZ56N79[7]={5.E-7, 17., 8.E-5, 4.5E-9, 1.3, 1.6E-5, .0018};
1023 static const std::pair<G4int, const G4double*> Z56N79(79,pZ56N79);
1024 static const G4double pZ56N80[7]={2.E-6, 20., 3.E-4, 6.E-9, 1.3, 1.8E-5, .0019};
1025 static const std::pair<G4int, const G4double*> Z56N80(80,pZ56N80);
1026 static const G4double pZ56N81[7]={5.8E-6, 20., .0018, 3.E-9, 1.3, 1.7E-5, .0041};
1027 static const std::pair<G4int, const G4double*> Z56N81(81,pZ56N81);
1028 static const G4double pZ56N82[7]={2.7E-6, 20., 5.5E-4, 4.E-9, 1.4, 2.E-5, .0027};
1029 static const std::pair<G4int, const G4double*> Z56N82(82,pZ56N82);
1030 static const G4double pZ56N84[7]={1.1E-6, 21., 1.E-4, 9.E-9, 2., 2.7E-5, .0012};
1031 static const std::pair<G4int, const G4double*> Z56N84(84,pZ56N84);
1032 static const std::pair<G4int, const G4double*> Z56[N56]={Z56N74, Z56N76, Z56N77, Z56N78,
1033 Z56N79, Z56N80, Z56N81, Z56N82,
1034 Z56N84};
1035 //==> La(Z=57)
1036 static const G4int N57=3;
1037 static const G4double pZ57N81[7]={2.7E-6, 20., .0017, 1.1E-9, 1.4, 2.E-5, .0083};
1038 static const std::pair<G4int, const G4double*> Z57N81(81,pZ57N81);
1039 static const G4double pZ57N82[7]={5.4E-6, 20., .0027, 1.5E-9, 1., 1.5E-5, .0065};
1040 static const std::pair<G4int, const G4double*> Z57N82(82,pZ57N82);
1041 static const G4double pZ57N83[7]={2.7E-6, 20., 2.6E-4, 6.E-9, 1.4, 2.E-5, .0012};
1042 static const std::pair<G4int, const G4double*> Z57N83(83,pZ57N83);
1043 static const std::pair<G4int, const G4double*> Z57[N57]={Z57N81, Z57N82, Z57N83};
1044 //==> Ce(Z=58)
1045 static const G4int N58=8;
1046 static const G4double pZ58N78[7]={1.8E-6, 20., 3.7E-4, 4.E-9, 1.4, 2.E-5, .0028};
1047 static const std::pair<G4int, const G4double*> Z58N78(78,pZ58N78);
1048 static const G4double pZ58N80[7]={1.8E-6, 18., 2.6E-4, 6.E-9, 1.3, 2.1E-5, .0017};
1049 static const std::pair<G4int, const G4double*> Z58N80(80,pZ58N80);
1050 static const G4double pZ58N81[7]={.0018, 18., 2.9, .6E-9, 1.3, 2.E-5, .02};
1051 static const std::pair<G4int, const G4double*> Z58N81(81,pZ58N81);
1052 static const G4double pZ58N82[7]={1.8E-6, 18., 3.7E-4, 5.E-9, 1.3, 2.E-5, .0024};
1053 static const std::pair<G4int, const G4double*> Z58N82(82,pZ58N82);
1054 static const G4double pZ58N83[7]={7.2E-6, 20., .0025, 2.3E-9, 1.1, 1.7E-5, .0045};
1055 static const std::pair<G4int, const G4double*> Z58N83(83,pZ58N83);
1056 static const G4double pZ58N84[7]={1.2E-6, 18., 2.E-4, 6.E-9, 1.5, 1.9E-5, .0018};
1057 static const std::pair<G4int, const G4double*> Z58N84(84,pZ58N84);
1058 static const G4double pZ58N85[7]={1.2E-6, 16., 6.E-4, 1.7E-9, 1.4, 2.E-5, .0053};
1059 static const std::pair<G4int, const G4double*> Z58N85(85,pZ58N85);
1060 static const G4double pZ58N86[7]={6.E-7, 18., 1.E-4, 6.E-9, 1.5, 1.9E-5, .0018};
1061 static const std::pair<G4int, const G4double*> Z58N86(86,pZ58N86);
1062 static const std::pair<G4int, const G4double*> Z58[N58]={Z58N78, Z58N80, Z58N81, Z58N82,
1063 Z58N83, Z58N84, Z58N85, Z58N86};
1064 //==> Pr(Z=59)
1065 static const G4int N59=3;
1066 static const G4double pZ59N82[7]={9.5E-7, 16., 1.6E-4, 4.E-9, 1.4, 2.E-5, .0017};
1067 static const std::pair<G4int, const G4double*> Z59N82(82,pZ59N82);
1068 static const G4double pZ59N83[7]={9.5E-7, 16., 1.9E-4, 4.E-9, 1.4, 1.9E-5, .0021};
1069 static const std::pair<G4int, const G4double*> Z59N83(83,pZ59N83);
1070 static const G4double pZ59N84[7]={9.5E-6, 16., .019, .4E-9, 2., 2.4E-5, .021};
1071 static const std::pair<G4int, const G4double*> Z59N84(84,pZ59N84);
1072 static const std::pair<G4int, const G4double*> Z59[N59]={Z59N82, Z59N83, Z59N84};
1073 //==> Nd(Z=60)
1074 static const G4int N60=8;
1075 static const G4double pZ60N82[7]={9.6E-6, 21., .0036, 2.3E-9, 1.4, 2.E-5, .0052};
1076 static const std::pair<G4int, const G4double*> Z60N82(82,pZ60N82);
1077 static const G4double pZ60N83[7]={9.6E-4, 20., 4., .25E-9, 1.4, 2.E-5, .052};
1078 static const std::pair<G4int, const G4double*> Z60N83(83,pZ60N83);
1079 static const G4double pZ60N84[7]={4.8E-7, 21., 3.3E-5, 1.E-08, 1.3, 2.2E-5, 9.E-4};
1080 static const std::pair<G4int, const G4double*> Z60N84(84,pZ60N84);
1081 static const G4double pZ60N85[7]={.0048, 20., 4.5, .9E-9, 1.3, 2.E-5, .012};
1082 static const std::pair<G4int, const G4double*> Z60N85(85,pZ60N85);
1083 static const G4double pZ60N86[7]={1.2E-6, 16., 7.7E-4, 1.5E-9, 1.3, 1.8E-5, .0066};
1084 static const std::pair<G4int, const G4double*> Z60N86(86,pZ60N86);
1085 static const G4double pZ60N87[7]={.0012, 15., 8.4, .1E-9, 1.3, 1.5E-5, .071};
1086 static const std::pair<G4int, const G4double*> Z60N87(87,pZ60N87);
1087 static const G4double pZ60N88[7]={1.5E-7, 16., 4.1E-5, 2.5E-9, 1.3, 1.6E-5, .0027};
1088 static const std::pair<G4int, const G4double*> Z60N88(88,pZ60N88);
1089 static const G4double pZ60N90[7]={1.5E-7, 16., 4.3E-5, 2.5E-9, 1.3, 1.6E-5, .0029};
1090 static const std::pair<G4int, const G4double*> Z60N90(90,pZ60N90);
1091 static const std::pair<G4int, const G4double*> Z60[N60]={Z60N82, Z60N83, Z60N84, Z60N85,
1092 Z60N86, Z60N87, Z60N88, Z60N90};
1093 //==> Pm(Z=61)
1094 static const G4int N61=3;
1095 static const G4double pZ61N86[7]={6.E-7, 16., 8.E-4, .6E-9, 3., 2.8E-5, .014};
1096 static const std::pair<G4int, const G4double*> Z61N86(86,pZ61N86);
1097 static const G4double pZ61N87[7]={6.2E-8, 16., 1.2E-5, 4.E-9, 2.2, 2.5E-5, .0019};
1098 static const std::pair<G4int, const G4double*> Z61N87(87,pZ61N87);
1099 static const G4double pZ61N88[7]={3.2E-8, 16., 6.4E-6, 4.E-9, 2.2, 2.5E-5, .002};
1100 static const std::pair<G4int, const G4double*> Z61N88(88,pZ61N88);
1101 static const std::pair<G4int, const G4double*> Z61[N61]={Z61N86, Z61N87, Z61N88};
1102 //==> Sm(Z=62)
1103 static const G4int N62=9;
1104 static const G4double pZ62N82[7]={1.2E-7, 16., 2.1E-5, 5.E-9, 1.4, 2.E-5, .0017};
1105 static const std::pair<G4int, const G4double*> Z62N82(82,pZ62N82);
1106 static const G4double pZ62N85[7]={1.2E-7, 16., 5.3E-5, 1.5E-9, 1.3, 1.7E-5, .0045};
1107 static const std::pair<G4int, const G4double*> Z62N85(85,pZ62N85);
1108 static const G4double pZ62N86[7]={6.E-8, 16., 1.7E-5, 3.E-9, 1.3, 1.7E-5, .0028};
1109 static const std::pair<G4int, const G4double*> Z62N86(86,pZ62N86);
1110 static const G4double pZ62N87[7]={6.E-8, 15., 5.2E-4, .11E-9, 1.3, 1.5E-5, .074};
1111 static const std::pair<G4int, const G4double*> Z62N87(87,pZ62N87);
1112 static const G4double pZ62N88[7]={6.E-7, 16., 8.6E-4, .7E-9, 1.3, 1.7E-5, .015};
1113 static const std::pair<G4int, const G4double*> Z62N88(88,pZ62N88);
1114 static const G4double pZ62N89[7]={6.E-7, 16., 5.E-4, 0., 1.3, 1.6E-5, .053};
1115 static const std::pair<G4int, const G4double*> Z62N89(89,pZ62N89);
1116 static const G4double pZ62N90[7]={6.E-8, 15., 1.3E-5, 4.5E-9, 1.3, 1.7E-5, .0019};
1117 static const std::pair<G4int, const G4double*> Z62N90(90,pZ62N90);
1118 static const G4double pZ62N91[7]={6.E-8, 15., 1.5E-5, 2.E-9, 1.3, 1.6E-5, .0024};
1119 static const std::pair<G4int, const G4double*> Z62N91(91,pZ62N91);
1120 static const G4double pZ62N92[7]={1.2E-7, 15., 8.6E-5, 1.2E-9, 1.3, 1.6E-5, .007};
1121 static const std::pair<G4int, const G4double*> Z62N92(92,pZ62N92);
1122 static const std::pair<G4int, const G4double*> Z62[N62]={Z62N82, Z62N85, Z62N86, Z62N87,
1123 Z62N88, Z62N89, Z62N90, Z62N91,
1124 Z62N92};
1125 //==> Eu(Z=63)
1126 static const G4int N63=7;
1127 static const G4double pZ63N88[7]={6.E-8, 15., 2.8E-5, 2.E-9, 1.3, 1.5E-5, .0046};
1128 static const std::pair<G4int, const G4double*> Z63N88(88,pZ63N88);
1129 static const G4double pZ63N89[7]={6.E-7, 15., .0011, .5E-9, 2.4, 2.4E-5, .017};
1130 static const std::pair<G4int, const G4double*> Z63N89(89,pZ63N89);
1131 static const G4double pZ63N90[7]={3.E-7, 15., 1.8E-4, 1.1E-9, 1., 1.2E-5, .0054};
1132 static const std::pair<G4int, const G4double*> Z63N90(90,pZ63N90);
1133 static const G4double pZ63N91[7]={4.1E-7, 15., 1.4E-4, 1.9E-9, 1., 1.4E-5, .0032};
1134 static const std::pair<G4int, const G4double*> Z63N91(91,pZ63N91);
1135 static const G4double pZ63N92[7]={5.E-8, 15., 2.4E-5, 2.8E-9, 1., 1.3E-5, .0037};
1136 static const std::pair<G4int, const G4double*> Z63N92(92,pZ63N92);
1137 static const G4double pZ63N93[7]={4.1E-8, 17., 1.6E-5, 2.E-9, 3.3, 3.4E-5, .004};
1138 static const std::pair<G4int, const G4double*> Z63N93(93,pZ63N93);
1139 static const G4double pZ63N94[7]={4.2E-8, 17., 1.6E-5, 2.E-9, 1.2, 1.6E-5, .004};
1140 static const std::pair<G4int, const G4double*> Z63N94(94,pZ63N94);
1141 static const std::pair<G4int, const G4double*> Z63[N63]={Z63N88, Z63N89, Z63N90, Z63N91,
1142 Z63N92, Z63N93, Z63N94};
1143 //==> Gd(Z=64)
1144 static const G4int N64=8;
1145 static const G4double pZ64N88[7]={4.2E-8, 14., 2.E-4, 0., 1.2, 1.3E-5, .19};
1146 static const std::pair<G4int, const G4double*> Z64N88(88,pZ64N88);
1147 static const G4double pZ64N89[7]={2.E-6, 14., .0016, 1.4E-9, 1.6, 1.6E-5, .0057};
1148 static const std::pair<G4int, const G4double*> Z64N89(89,pZ64N89);
1149 static const G4double pZ64N90[7]={1.7E-7, 12., 8.4E-5, 2.E-9, 1.8, 2.2E-5, .0035};
1150 static const std::pair<G4int, const G4double*> Z64N90(90,pZ64N90);
1151 static const G4double pZ64N91[7]={1.7E-7, 13., 5.E-4, .3E-9, 1.7, 1.9E-5, .026};
1152 static const std::pair<G4int, const G4double*> Z64N91(91,pZ64N91);
1153 static const G4double pZ64N92[7]={1.7E-7, 13., 7.E-5, 2.5E-9, 1.8, 2.2E-5, .003};
1154 static const std::pair<G4int, const G4double*> Z64N92(92,pZ64N92);
1155 static const G4double pZ64N93[7]={1.7E-6, 12., .002, 0., 1.7, 1.8E-5, .47};
1156 static const std::pair<G4int, const G4double*> Z64N93(93,pZ64N93);
1157 static const G4double pZ64N94[7]={3.4E-7, 13., 1.5E-4, 2.E-9, 1.8, 2.3E-5, .0034};
1158 static const std::pair<G4int, const G4double*> Z64N94(94,pZ64N94);
1159 static const G4double pZ64N96[7]={2.6E-6, 13., .0019, 1.2E-9, 1., 1.2E-5, .0056};
1160 static const std::pair<G4int, const G4double*> Z64N96(96,pZ64N96);
1161 static const std::pair<G4int, const G4double*> Z64[N64]={Z64N88, Z64N89, Z64N90, Z64N91,
1162 Z64N92, Z64N93, Z64N94, Z64N96};
1163 //==> Tb(Z=65)
1164 static const G4int N65=2;
1165 static const G4double pZ65N94[7]={9.E-7, 16., 3.9E-4, 1.7E-9, 2., 2.2E-5, .0042};
1166 static const std::pair<G4int, const G4double*> Z65N94(94,pZ65N94);
1167 static const G4double pZ65N95[7]={4.5E-7, 16., 1.1E-4, 3.E-9, 1.7, 2.2E-5, .0024};
1168 static const std::pair<G4int, const G4double*> Z65N95(95,pZ65N95);
1169 static const std::pair<G4int, const G4double*> Z65[N65]={Z65N94, Z65N94};
1170 //==> Dy(Z=66)
1171 static const G4int N66=7;
1172 static const G4double pZ66N90[7]={1.2E-7, 13., 4.E-5, 3.E-9, 1.2, 1.4E-5, .0025};
1173 static const std::pair<G4int, const G4double*> Z66N90(90,pZ66N90);
1174 static const G4double pZ66N92[7]={1.2E-7, 13., 6.7E-5, 2.E-9, 1., 1.1E-5, .004};
1175 static const std::pair<G4int, const G4double*> Z66N92(92,pZ66N92);
1176 static const G4double pZ66N94[7]={1.2E-7, 13., 5.3E-5, 1.6E-9, 1., 1.1E-5, .0034};
1177 static const std::pair<G4int, const G4double*> Z66N94(94,pZ66N94);
1178 static const G4double pZ66N95[7]={1.2E-6, 13., .0017, .7E-9, 1.3, 1.3E-5, .011};
1179 static const std::pair<G4int, const G4double*> Z66N95(95,pZ66N95);
1180 static const G4double pZ66N96[7]={1.2E-7, 13., 8.E-6, 1.5E-7, 1., 1.1E-5, 1.E-4};
1181 static const std::pair<G4int, const G4double*> Z66N96(96,pZ66N96);
1182 static const G4double pZ66N97[7]={1.5E-7, 13., 4.E-5, 4.E-9, 1.3, 1.3E-5, .002};
1183 static const std::pair<G4int, const G4double*> Z66N97(97,pZ66N97);
1184 static const G4double pZ66N98[7]={3.E-7, 13., .001, 4.E-9, 1.3, 1.3E-5, .23};
1185 static const std::pair<G4int, const G4double*> Z66N98(98,pZ66N98);
1186 static const std::pair<G4int, const G4double*> Z66[N66]={Z66N90, Z66N92, Z66N94, Z66N95,
1187 Z66N96, Z66N97, Z66N98};
1188 //==> Ho(Z=67)
1189 static const G4int N67=2;
1190 static const G4double pZ67N98[7]={3.E-7, 13., 2.2E-4, 1.5E-9, 1., 1.E-5, .0054};
1191 static const std::pair<G4int, const G4double*> Z67N98(98,pZ67N98);
1192 static const G4double pZ67N99[7]={7.5E-8, 13., 2.6E-5, 4.5E-9, 1.5, 1.5E-5, .0021};
1193 static const std::pair<G4int, const G4double*> Z67N99(99,pZ67N99);
1194 static const std::pair<G4int, const G4double*> Z67[N67]={Z67N98, Z67N99};
1195 //==> Er(Z=68)
1196 static const G4int N68=6;
1197 static const G4double pZ68N94[7]={1.2E-7, 13., 7.8E-5, 1.6E-9, .9, 9.E-6, .005};
1198 static const std::pair<G4int, const G4double*> Z68N94(94,pZ68N94);
1199 static const G4double pZ68N96[7]={1.2E-7, 13., 8.5E-5, 1.2E-9, .9, 8.E-6, .0055};
1200 static const std::pair<G4int, const G4double*> Z68N96(96,pZ68N96);
1201 static const G4double pZ68N98[7]={1.E-6, 13., .0011, .8E-9, .9, 8.E-6, .0087};
1202 static const std::pair<G4int, const G4double*> Z68N98(98,pZ68N98);
1203 static const G4double pZ68N99[7]={1.2E-7, 13., 2.5E-5, 4.5E-9, .9, 9.E-6, .0015};
1204 static const std::pair<G4int, const G4double*> Z68N99(99,pZ68N99);
1205 static const G4double pZ68N100[7]={2.E-6, 13., .0015, 1.1E-9, .9, 8.E-6, .0058};
1206 static const std::pair<G4int, const G4double*> Z68N100(100,pZ68N100);
1207 static const G4double pZ68N102[7]={2.E-6, 13., .0018, 1.E-9, .9, 8.E-6, .007};
1208 static const std::pair<G4int, const G4double*> Z68N102(102,pZ68N102);
1209 static const std::pair<G4int, const G4double*> Z68[N68]={Z68N94, Z68N96, Z68N98,
1210 Z68N99, Z68N100, Z68N102};
1211 //==> Tm(Z=69) *** No data *** (Tm169=Er167)
1212 static const G4int N69=1;
1213 static const G4double pZ69N100[7]={1.2E-7, 13., 2.5E-5, 4.5E-9, .9, 9.E-6, .0015};
1214 static const std::pair<G4int, const G4double*> Z69N100(100,pZ69N100);
1215 static const std::pair<G4int, const G4double*> Z69[N69]={Z69N100};
1216 //==> Yb(Z=70) *** No data *** (Yb168=Er166, Yb170=Er168, Yb171=Er167, Yb172=Er170,
1217 // Yb173=Hf177, Yb174=Hf176, Yb176=Hf178)
1218 static const G4int N70=7;
1219 static const G4double pZ70N98[7]={1.E-6, 13., .0011, .8E-9, .9, 8.E-6, .0087};
1220 static const std::pair<G4int, const G4double*> Z70N98(98,pZ70N98);
1221 static const G4double pZ70N100[7]={2.E-6, 13., .0015, 1.1E-9, .9, 8.E-6, .0058};
1222 static const std::pair<G4int, const G4double*> Z70N100(100,pZ70N100);
1223 static const G4double pZ70N101[7]={1.2E-7, 13., 2.5E-5, 4.5E-9, .9, 9.E-6, .0015};
1224 static const std::pair<G4int, const G4double*> Z70N101(101,pZ70N101);
1225 static const G4double pZ70N102[7]={2.E-6, 13., .0018, 1.E-9, .9, 8.E-6, .007};
1226 static const std::pair<G4int, const G4double*> Z70N102(102,pZ70N102);
1227 static const G4double pZ70N103[7]={5.E-7, 18., .01, 1.7E-6, 1.2, 1.4E-5, 1.4E-5};
1228 static const std::pair<G4int, const G4double*> Z70N103(103,pZ70N103);
1229 static const G4double pZ70N104[7]={5.E-7, 18., 1.9E-4, 2.5E-9, 1.2, 1.4E-5, .004};
1230 static const std::pair<G4int, const G4double*> Z70N104(104,pZ70N104);
1231 static const G4double pZ70N106[7]={5.E-7, 18., 1.3E-4, 2.E-9, 1.2, 1.4E-5, .0027};
1232 static const std::pair<G4int, const G4double*> Z70N106(106,pZ70N106);
1233 static const std::pair<G4int, const G4double*> Z70[N70]={Z70N98, Z70N100, Z70N101,
1234 Z70N102, Z70N103, Z70N104,
1235 Z70N106};
1236 //==> Lu(Z=71)
1237 static const G4int N71=2;
1238 static const G4double pZ71N104[7]={5.E-7, 18., 1.8E-4, 2.E-9, .9, 9.E-6, .0036};
1239 static const std::pair<G4int, const G4double*> Z71N104(104,pZ71N104);
1240 static const G4double pZ71N105[7]={2.5E-7, 18., 9.E-5, 1.E-8, .9, 9.E-6, .0016};
1241 static const std::pair<G4int, const G4double*> Z71N105(105,pZ71N105);
1242 static const std::pair<G4int, const G4double*> Z71[N71]={Z71N104, Z71N105};
1243 //==> Hf(Z=72)
1244 static const G4int N72=6;
1245 static const G4double pZ72N102[7]={1.E-6, 18., 8.8E-4, .8E-9, 1., 1.1E-5, .0092};
1246 static const std::pair<G4int, const G4double*> Z72N102(102,pZ72N102);
1247 static const G4double pZ72N104[7]={5.E-7, 18., 1.9E-4, 2.5E-9, 1.2, 1.4E-5, .004};
1248 static const std::pair<G4int, const G4double*> Z72N104(104,pZ72N104);
1249 static const G4double pZ72N105[7]={5.E-7, 18., .01, 1.7E-6, 1.2, 1.4E-5, 1.4E-5};
1250 static const std::pair<G4int, const G4double*> Z72N105(105,pZ72N105);
1251 static const G4double pZ72N106[7]={5.E-7, 18., 1.3E-4, 2.E-9, 1.2, 1.4E-5, .0027};
1252 static const std::pair<G4int, const G4double*> Z72N106(106,pZ72N106);
1253 static const G4double pZ72N107[7]={2.5E-7, 18., 1.E-4, 2.E-9, 1.2, 1.5E-5, .0041};
1254 static const std::pair<G4int, const G4double*> Z72N107(107,pZ72N107);
1255 static const G4double pZ72N108[7]={1.E-6, 18., .0012, .6E-9, 1.2, 1.5E-5, .012};
1256 static const std::pair<G4int, const G4double*> Z72N108(108,pZ72N108);
1257 static const std::pair<G4int, const G4double*> Z72[N72]={Z72N102, Z72N104, Z72N105,
1258 Z72N106, Z72N107, Z72N108};
1259 //==> Ta(Z=73)
1260 static const G4int N73=2;
1261 static const G4double pZ73N108[7]={5.E-7, 18., 1.7E-4, 2.E-9, 1.2, 1.4E-5, .0035};
1262 static const std::pair<G4int, const G4double*> Z73N108(108,pZ73N108);
1263 static const G4double pZ73N109[7]={1.E-6, 14., .002, .3E-9, 1.3, 1.5E-5, .016};
1264 static const std::pair<G4int, const G4double*> Z73N109(109,pZ73N109);
1265 static const std::pair<G4int, const G4double*> Z73[N73]={Z73N108, Z73N108};
1266 //==> W (Z=74) *** W180 only bad TENDL-2008 *** (W180=Hf178)
1267 static const G4int N74=5;
1268 static const G4double pZ74N106[7]={5.E-7, 18., 1.3E-4, 2.E-9, 1.2, 1.4E-5, .0027};
1269 static const std::pair<G4int, const G4double*> Z74N106(106,pZ74N106);
1270 static const G4double pZ74N108[7]={4.E-6, 14., .0034, .9E-9, 1.3, 1.4E-5, .0067};
1271 static const std::pair<G4int, const G4double*> Z74N108(108,pZ74N108);
1272 static const G4double pZ74N109[7]={1.2E-7, 14., 3.E-5, 3.E-9, 1.3, 1.4E-5, .0019};
1273 static const std::pair<G4int, const G4double*> Z74N109(109,pZ74N109);
1274 static const G4double pZ74N110[7]={1.2E-7, 14., 3.6E-5, 2.E-9, 1.3, 1.3E-5, .0024};
1275 static const std::pair<G4int, const G4double*> Z74N110(110,pZ74N110);
1276 static const G4double pZ74N112[7]={1.2E-7, 14., 3.E-5, 1.3E-7, 1.3, 1.3E-5, 1.4E-4};
1277 static const std::pair<G4int, const G4double*> Z74N112(112,pZ74N112);
1278 static const std::pair<G4int, const G4double*> Z74[N74]={Z74N106, Z74N108, Z74N109,
1279 Z74N110, Z74N112};
1280 //==> Re(Z=75)
1281 static const G4int N75=2;
1282 static const G4double pZ75N110[7]={1.2E-7, 14., 8.E-5, 1.2E-9, 1.3, 1.5E-5, .005};
1283 static const std::pair<G4int, const G4double*> Z75N110(110,pZ75N110);
1284 static const G4double pZ75N112[7]={1.2E-7, 14., 8.8E-5, 1.1E-9, 1.3, 1.5E-5, .0055};
1285 static const std::pair<G4int, const G4double*> Z75N112(112,pZ75N112);
1286 static const std::pair<G4int, const G4double*> Z75[N75]={Z75N110, Z75N112};
1287 //==> Os(Z=76) *** No data *** (Os184=W182, Os186=W184, Os187=Re187, Os188=W186,
1288 // Os189=Re187, Os190=Os192=W186)
1289 static const G4int N76=7;
1290 static const G4double pZ76N108[7]={4.E-6, 14., .0034, .9E-9, 1.3, 1.4E-5, .0067};
1291 static const std::pair<G4int, const G4double*> Z76N108(108,pZ76N108);
1292 static const G4double pZ76N110[7]={1.2E-7, 14., 3.6E-5, 2.E-9, 1.3, 1.3E-5, .0024};
1293 static const std::pair<G4int, const G4double*> Z76N110(110,pZ76N110);
1294 static const G4double pZ76N111[7]={1.2E-7, 14., 8.8E-5, 1.1E-9, 1.3, 1.5E-5, .0055};
1295 static const std::pair<G4int, const G4double*> Z76N111(111,pZ76N111);
1296 static const G4double pZ76N112[7]={1.2E-7, 14., 3.E-5, 1.3E-7, 1.3, 1.3E-5, 1.4E-4};
1297 static const std::pair<G4int, const G4double*> Z76N112(112,pZ76N112);
1298 static const G4double pZ76N113[7]={1.2E-7, 14., 8.8E-5, 1.1E-9, 1.3, 1.5E-5, .0055};
1299 static const std::pair<G4int, const G4double*> Z76N113(113,pZ76N113);
1300 static const G4double pZ76N114[7]={1.2E-7, 14., 3.E-5, 1.3E-7, 1.3, 1.3E-5, 1.4E-4};
1301 static const std::pair<G4int, const G4double*> Z76N114(114,pZ76N114);
1302 static const G4double pZ76N116[7]={1.2E-7, 14., 3.E-5, 1.3E-7, 1.3, 1.3E-5, 1.4E-4};
1303 static const std::pair<G4int, const G4double*> Z76N116(116,pZ76N116);
1304 static const std::pair<G4int, const G4double*> Z76[N76]={Z76N108, Z76N110, Z76N111,
1305 Z76N112, Z76N113, Z76N114,
1306 Z76N116};
1307 //==> Ir(Z=77)
1308 static const G4int N77=2;
1309 static const G4double pZ77N114[7]={4.8E-7, 14., 5.2E-4, .7E-9, 1.5, 1.7E-5, .0082};
1310 static const std::pair<G4int, const G4double*> Z77N114(114,pZ77N114);
1311 static const G4double pZ77N116[7]={4.8E-7, 14., 4.5E-4, .8E-9, 1.8, 2.3E-5, .0073};
1312 static const std::pair<G4int, const G4double*> Z77N116(116,pZ77N116);
1313 static const std::pair<G4int, const G4double*> Z77[N77]={Z77N114, Z77N116};
1314 //==> Pt(Z=78) *** No data *** (Pt190=Pt192=Pt194=Hg196, Pt195=Hg199, Pt196=Hg198,
1315 // Pt198=Hg200)
1316 static const G4int N78=6;
1317 static const G4double pZ78N112[7]={6.E-8, 19., 3.3E-4, .1E-9, 1.6, 1.8E-5, .06};
1318 static const std::pair<G4int, const G4double*> Z78N112(112,pZ78N112);
1319 static const G4double pZ78N114[7]={6.E-8, 19., 3.3E-4, .1E-9, 1.6, 1.8E-5, .06};
1320 static const std::pair<G4int, const G4double*> Z78N114(114,pZ78N114);
1321 static const G4double pZ78N116[7]={6.E-8, 19., 3.3E-4, .1E-9, 1.6, 1.8E-5, .06};
1322 static const std::pair<G4int, const G4double*> Z78N116(116,pZ78N116);
1323 static const G4double pZ78N117[7]={9.6E-7, 20., .001, .2E-9, 1.6, 2.E-5, .037};
1324 static const std::pair<G4int, const G4double*> Z78N117(117,pZ78N117);
1325 static const G4double pZ78N118[7]={2.4E-7, 20., 1.6E-4, 1.3E-9, 1.6, 1.8E-5, .007};
1326 static const std::pair<G4int, const G4double*> Z78N118(118,pZ78N118);
1327 static const G4double pZ78N120[7]={2.E-6, 19., .0015, .9E-9, 1.6, 1.8E-5, .0078};
1328 static const std::pair<G4int, const G4double*> Z78N120(120,pZ78N120);
1329 static const std::pair<G4int, const G4double*> Z78[N78]={Z78N112, Z78N114, Z78N116,
1330 Z78N117, Z78N118, Z78N120};
1331 //==> Au(Z=79)
1332 static const G4int N79=1;
1333 static const G4double pZ79N118[7]={2.4E-7, 19., 1.E-4, 1.4E-9, 1.3, 1.7E-5, .0042};
1334 static const std::pair<G4int, const G4double*> Z79N118(118,pZ79N118);
1335 static const std::pair<G4int, const G4double*> Z79[N79]={Z79N118};
1336 //==> Hg(Z=80)
1337 static const G4int N80=7;
1338 static const G4double pZ80N116[7]={6.E-8, 19., 3.3E-4, .1E-9, 1.6, 1.8E-5, .06};
1339 static const std::pair<G4int, const G4double*> Z80N116(116,pZ80N116);
1340 static const G4double pZ80N118[7]={2.4E-7, 20., 1.6E-4, 1.3E-9, 1.6, 1.8E-5, .007};
1341 static const std::pair<G4int, const G4double*> Z80N118(118,pZ80N118);
1342 static const G4double pZ80N119[7]={9.6E-7, 20., .001, .2E-9, 1.6, 2.E-5, .037};
1343 static const std::pair<G4int, const G4double*> Z80N119(119,pZ80N119);
1344 static const G4double pZ80N120[7]={2.E-6, 19., .0015, .9E-9, 1.6, 1.8E-5, .0078};
1345 static const std::pair<G4int, const G4double*> Z80N120(120,pZ80N120);
1346 static const G4double pZ80N121[7]={1.E-6, 20., 7.E-4, 1.E-9, 1.6, 1.8E-5, .0076};
1347 static const std::pair<G4int, const G4double*> Z80N121(121,pZ80N121);
1348 static const G4double pZ80N122[7]={2.E-6, 18., .0016, .8E-9, 1.6, 1.8E-5, .0078};
1349 static const std::pair<G4int, const G4double*> Z80N122(122,pZ80N122);
1350 static const G4double pZ80N124[7]={2.0E-6, 18., .0032, .4E-9, 1.6, 1.8E-5, .016};
1351 static const std::pair<G4int, const G4double*> Z80N124(124,pZ80N124);
1352 static const std::pair<G4int, const G4double*> Z80[N80]={Z80N116, Z80N118, Z80N119,
1353 Z80N120, Z80N121, Z80N122,
1354 Z80N124};
1355 //==> Tl(Z=81) *** No data *** (Tl203=Au196, Tl198=Bi209)
1356 static const G4int N81=2;
1357 static const G4double pZ81N122[7]={2.4E-7, 19., 1.E-4, 1.4E-9, 1.3, 1.7E-5, .0042};
1358 static const std::pair<G4int, const G4double*> Z81N122(122,pZ81N122);
1359 static const G4double pZ81N124[7]={};
1360 static const std::pair<G4int, const G4double*> Z81N124(124,pZ81N124);
1361 static const std::pair<G4int, const G4double*> Z81[N81]={Z81N122, Z81N124};
1362 //==> Pb(Z=82)
1363 static const G4int N82=4;
1364 static const G4double pZ82N122[7]={4.E-6, 20., .0022, 1.E-9, 1.6, 1.8E-5, .0058};
1365 static const std::pair<G4int, const G4double*> Z82N122(122,pZ82N122);
1366 static const G4double pZ82N124[7]={4.E-6, 20., .0022, 1.E-9, 1.6, 1.8E-5, .0058};
1367 static const std::pair<G4int, const G4double*> Z82N124(124,pZ82N124);
1368 static const G4double pZ82N125[7]={2.E-6, 20., .0011, 1.2E-9, 1.6, 1.8E-5, .0056};
1369 static const std::pair<G4int, const G4double*> Z82N125(125,pZ82N125);
1370 static const G4double pZ82N126[7]={4.E-6, 20., .0023, 1.2E-9, 1.6, 1.8E-5, .0058};
1371 static const std::pair<G4int, const G4double*> Z82N126(126,pZ82N126);
1372 static const std::pair<G4int, const G4double*> Z82[N82]={Z82N122, Z82N124, Z82N125,
1373 Z82N126};
1374 //==> Bi(Z=83)
1375 static const G4int N83=1;
1376 static const G4double pZ83N126[7]={8.E-7, 23., 3.3E-4, 1.8E-9, 1.6, 1.8E-5, .005};
1377 static const std::pair<G4int, const G4double*> Z83N126(126,pZ83N126);
1378 static const std::pair<G4int, const G4double*> Z83[N83]={Z83N126};
1379 //==> Po(Z=84) *** No data *** (Po209=Pb207)
1380 static const G4int N84=1;
1381 static const G4double pZ84N125[7]={2.E-6, 20., .0011, 1.2E-9, 1.6, 1.8E-5, .0056};
1382 static const std::pair<G4int, const G4double*> Z84N125(125,pZ84N125);
1383 static const std::pair<G4int, const G4double*> Z84[N84]={Z84N125};
1384 //==> At(Z=85) *** No data *** (At210=Pb207)
1385 static const G4int N85=1;
1386 static const G4double pZ85N125[7]={2.E-6, 20., .0011, 1.2E-9, 1.6, 1.8E-5, .0056};
1387 static const std::pair<G4int, const G4double*> Z85N125(125,pZ85N125);
1388 static const std::pair<G4int, const G4double*> Z85[N85]={Z85N125};
1389 //==> Rn(Z=86) *** No data *** (Rn222=Ra224)
1390 static const G4int N86=1;
1391 static const G4double pZ86N136[7]={1.E-7, 23., 5.5E-5, 1.2E-9, 1.6, 1.8E-5, .0062};
1392 static const std::pair<G4int, const G4double*> Z86N136(136,pZ86N136);
1393 static const std::pair<G4int, const G4double*> Z86[N86]={Z86N136};
1394 //==> Fr(Z=87) *** No data *** (Fr223=Ac225)
1395 static const G4int N87=1;
1396 static const G4double pZ87N136[7]={2.E-7, 23., 1.1E-4, 1.2E-9, 1.6, 1.8E-5, .0062};
1397 static const std::pair<G4int, const G4double*> Z87N136(136,pZ87N136);
1398 static const std::pair<G4int, const G4double*> Z87[N87]={Z87N136};
1399 //==> Ra(Z=88)
1400 static const G4int N88=4;
1401 static const G4double pZ88N135[7]={1.E-7, 23., 5.5E-5, 1.2E-9, 1.6, 1.8E-5, .0062};
1402 static const std::pair<G4int, const G4double*> Z88N135(135,pZ88N135);
1403 static const G4double pZ88N136[7]={1.E-7, 23., 5.5E-5, 1.2E-9, 1.6, 1.8E-5, .0062};
1404 static const std::pair<G4int, const G4double*> Z88N136(136,pZ88N136);
1405 static const G4double pZ88N137[7]={1.E-7, 23., 5.5E-5, 1.2E-9, 1.6, 1.8E-5, .0062};
1406 static const std::pair<G4int, const G4double*> Z88N137(137,pZ88N137);
1407 static const G4double pZ88N138[7]={4.E-7, 23., 1.7E-4, 1.5E-9, 1.6, 1.8E-5, .005};
1408 static const std::pair<G4int, const G4double*> Z88N138(138,pZ88N138);
1409 static const std::pair<G4int, const G4double*> Z88[N88]={Z88N135, Z88N136, Z88N137,
1410 Z88N138};
1411 //==> Ac(Z=89)
1412 static const G4int N89=3;
1413 static const G4double pZ89N136[7]={2.E-7, 23., 1.1E-4, 1.2E-9, 1.6, 1.8E-5, .0062};
1414 static const std::pair<G4int, const G4double*> Z89N136(136,pZ89N136);
1415 static const G4double pZ89N137[7]={4.E-7, 23., 2.2E-4, 1.2E-9, 1.6, 1.8E-5, .0062};
1416 static const std::pair<G4int, const G4double*> Z89N137(137,pZ89N137);
1417 static const G4double pZ89N138[7]={1.E-7, 23., 5.5E-5, 1.2E-9, 1.6, 1.8E-5, .0062};
1418 static const std::pair<G4int, const G4double*> Z89N138(138,pZ89N138);
1419 static const std::pair<G4int, const G4double*> Z89[N89]={Z89N136, Z89N137, Z89N138};
1420 //==> Th(Z=90)
1421 static const G4int N90=7;
1422 static const G4double pZ90N137[7]={4.E-7, 23., 2.2E-4, 1.2E-9, 1.6, 1.8E-5, .0062};
1423 static const std::pair<G4int, const G4double*> Z90N137(137,pZ90N137);
1424 static const G4double pZ90N138[7]={1.E-6, 23., .0016, .4E-9, 3., 3.E-5, .019};
1425 static const std::pair<G4int, const G4double*> Z90N138(138,pZ90N138);
1426 static const G4double pZ90N139[7]={2.5E-7, 23., 1.1E-4, 1.4E-9, 2.4, 2.7E-5, .0049};
1427 static const std::pair<G4int, const G4double*> Z90N139(139,pZ90N139);
1428 static const G4double pZ90N140[7]={1.2E-7, 23., 3.E-5, 2.E-9, 3., 3.E-5, .003};
1429 static const std::pair<G4int, const G4double*> Z90N140(140,pZ90N140);
1430 static const G4double pZ90N142[7]={4.E-6, 23., .0023, 1.1E-9, 1.8, 2.3E-5, .0064};
1431 static const std::pair<G4int, const G4double*> Z90N142(142,pZ90N142);
1432 static const G4double pZ90N143[7]={9.4E-7, 23., 5.4E-4, 1.1E-9, 3., 3.E-5, .0066};
1433 static const std::pair<G4int, const G4double*> Z90N143(143,pZ90N143);
1434 static const G4double pZ90N144[7]={2.5E-7, 23., 1.4E-4, 1.1E-9, 3., 3.E-5, .0066};
1435 static const std::pair<G4int, const G4double*> Z90N144(144,pZ90N144);
1436 static const std::pair<G4int, const G4double*> Z90[N90]={Z90N137, Z90N138, Z90N139,
1437 Z90N140, Z90N142, Z90N143,
1438 Z90N144};
1439 //==> Pa(Z=91)
1440 static const G4int N91=3;
1441 static const G4double pZ91N140[7]={1.E-5, 23., .0052, 1.6E-9, 1.8, 2.3E-5, .0057};
1442 static const std::pair<G4int, const G4double*> Z91N140(140,pZ91N140);
1443 static const G4double pZ91N141[7]={8.E-6, 23., .006, 0., 3.5, 3.5E-5, .021};
1444 static const std::pair<G4int, const G4double*> Z91N141(141,pZ91N141);
1445 static const G4double pZ91N142[7]={8.E-6, 23., .0042, 1.E-9, 2., 2.5E-5, .006};
1446 static const std::pair<G4int, const G4double*> Z91N142(142,pZ91N142);
1447 static const std::pair<G4int, const G4double*> Z91[N91]={Z91N140, Z91N141, Z91N142};
1448 //==> U (Z=92)
1449 static const G4int N92=10;
1450 static const G4double pZ92N140[7]={1.4E-6, 20., 8.E-4, 1.5E-9, 2.5, 2.8E-5, .0055};
1451 static const std::pair<G4int, const G4double*> Z92N140(140,pZ92N140);
1452 static const G4double pZ92N141[7]={5.6E-6, 20., .0033, 1.E-9, 2.5, 2.8E-5, .006};
1453 static const std::pair<G4int, const G4double*> Z92N141(141,pZ92N141);
1454 static const G4double pZ92N142[7]={5.6E-6, 20., .0034, 0., 2.5, 2.8E-5, .0072};
1455 static const std::pair<G4int, const G4double*> Z92N142(142,pZ92N142);
1456 static const G4double pZ92N143[7]={5.6E-6, 20., .0032, 0., 2., 2.3E-5, .006};
1457 static const std::pair<G4int, const G4double*> Z92N143(143,pZ92N143);
1458 static const G4double pZ92N144[7]={3.6E-7, 20., 1.6E-4, 1.3E-9, 2.2, 2.7E-5, .0043};
1459 static const std::pair<G4int, const G4double*> Z92N144(144,pZ92N144);
1460 static const G4double pZ92N145[7]={3.6E-6, 20., .003, 0., 2.2, 2.7E-5, .045};
1461 static const std::pair<G4int, const G4double*> Z92N145(145,pZ92N145);
1462 static const G4double pZ92N146[7]={3.6E-7, 20., 1.6E-4, 1.3E-9, 2.2, 2.7E-5, .0043};
1463 static const std::pair<G4int, const G4double*> Z92N146(146,pZ92N146);
1464 static const G4double pZ92N147[7]={3.6E-6, 20., .0014, 1.3E-9, 2.2, 2.7E-5, 12.};
1465 static const std::pair<G4int, const G4double*> Z92N147(147,pZ92N147);
1466 static const G4double pZ92N148[7]={3.4E-7, 20., 1.3E-4, 1.3E-9, 2.2, 2.8E-5, .0036};
1467 static const std::pair<G4int, const G4double*> Z92N148(148,pZ92N148);
1468 static const G4double pZ92N149[7]={3.3E-7, 20., 1.5E-4, 1.2E-9, 3., 3.4E-5, .0044};
1469 static const std::pair<G4int, const G4double*> Z92N149(149,pZ92N149);
1470 static const std::pair<G4int, const G4double*> Z92[N92]={Z92N140, Z92N141, Z92N142,
1471 Z92N143, Z92N144, Z92N145,
1472 Z92N146, Z92N147, Z92N148,
1473 Z92N146};
1474 //==> Np(Z=93)
1475 static const G4int N93=5;
1476 static const G4double pZ93N142[7]={3.4E-6, 20., .002, 1.3E-9, 3., 3.3E-5, .0056};
1477 static const std::pair<G4int, const G4double*> Z93N142(142,pZ93N142);
1478 static const G4double pZ93N143[7]={3.4E-6, 20., .002, 1.6E-9, 3.5, 3.6E-5, .005};
1479 static const std::pair<G4int, const G4double*> Z93N143(143,pZ93N143);
1480 static const G4double pZ93N144[7]={6.8E-6, 18., .0052, .8E-9, 2.4, 3.E-5, .0072};
1481 static const std::pair<G4int, const G4double*> Z93N144(144,pZ93N144);
1482 static const G4double pZ93N145[7]={3.4E-6, 20., .002, 1.E-9, 3.5, 3.6E-5, .006};
1483 static const std::pair<G4int, const G4double*> Z93N145(145,pZ93N145);
1484 static const G4double pZ93N146[7]={3.4E-6, 20., .002, 1.5E-9, 3.5, 3.6E-5, .0053};
1485 static const std::pair<G4int, const G4double*> Z93N146(146,pZ93N146);
1486 static const std::pair<G4int, const G4double*> Z93[N93]={Z93N142, Z93N143, Z93N144,
1487 Z93N145, Z93N146};
1488 //==> Pu(Z=94)
1489 static const G4int N94=10;
1490 static const G4double pZ94N142[7]={6.8E-7, 16., 4.5E-4, 1.7E-9, 2.6, 3.E-5, .0047};
1491 static const std::pair<G4int, const G4double*> Z94N142(142,pZ94N142);
1492 static const G4double pZ94N143[7]={6.8E-6, 18., .0044, .9E-9, 3.3, 3.5E-5, .0058};
1493 static const std::pair<G4int, const G4double*> Z94N143(143,pZ94N143);
1494 static const G4double pZ94N144[7]={6.8E-7, 16., 6.E-4, 0., 2.7, 2.6E-5, .0082};
1495 static const std::pair<G4int, const G4double*> Z94N144(144,pZ94N144);
1496 static const G4double pZ94N145[7]={2.6E-6, 16., .0017, 1.8E-9, 1.8, 2.E-5, .004};
1497 static const std::pair<G4int, const G4double*> Z94N145(145,pZ94N145);
1498 static const G4double pZ94N146[7]={2.5E-7, 20., 9.E-5, 3.6E-8, 3.4, 3.8E-5, 5.4E-4};
1499 static const std::pair<G4int, const G4double*> Z94N146(146,pZ94N146);
1500 static const G4double pZ94N147[7]={1.4E-5, 16., .01, .8E-9, 2.7, 2.6E-5, .0055};
1501 static const std::pair<G4int, const G4double*> Z94N147(147,pZ94N147);
1502 static const G4double pZ94N148[7]={3.4E-7, 20., 1.3E-4, 1.2E-9, 3.2, 3.E-5, .0036};
1503 static const std::pair<G4int, const G4double*> Z94N148(148,pZ94N148);
1504 static const G4double pZ94N149[7]={5.2E-6, 20., .0035, .4E-9, 2.3, 3.E-5, .0095};
1505 static const std::pair<G4int, const G4double*> Z94N149(149,pZ94N149);
1506 static const G4double pZ94N150[7]={3.3E-7, 20., 1.6E-4, 1.2E-9, 3., 3.E-5, .0046};
1507 static const std::pair<G4int, const G4double*> Z94N150(150,pZ94N150);
1508 static const G4double pZ94N152[7]={2.5E-6, 16., .0018, 1.2E-9, 3., 3.1E-5, .0052};
1509 static const std::pair<G4int, const G4double*> Z94N152(152,pZ94N152);
1510 static const std::pair<G4int, const G4double*> Z94[N94]={Z94N142, Z94N143, Z94N144,
1511 Z94N145, Z94N146, Z94N147,
1512 Z94N148, Z94N149, Z94N150,
1513 Z94N152};
1514 //==> Am(Z=95)
1515 static const G4int N95=4;
1516 static const G4double pZ95N156[7]={2.5E-6, 18., .0016, .9E-9, 2., 2.3E-5, .0058};
1517 static const std::pair<G4int, const G4double*> Z95N156(156,pZ95N156);
1518 static const G4double pZ95N157[7]={5.E-6, 18., .003, 2.7E-9, 2., 2.3E-5, .0039};
1519 static const std::pair<G4int, const G4double*> Z95N157(157,pZ95N157);
1520 static const G4double pZ95N158[7]={5.E-6, 19., .0033, 2.6E-9, 2., 2.3E-5, .0044};
1521 static const std::pair<G4int, const G4double*> Z95N158(158,pZ95N158);
1522 static const G4double pZ95N159[7]={5.E-5, 20., .029, 1.1E-9, 2., 2.3E-5, .0057};
1523 static const std::pair<G4int, const G4double*> Z95N159(159,pZ95N159);
1524 static const std::pair<G4int, const G4double*> Z95[N95]={Z95N156, Z95N157, Z95N158,
1525 Z95N159};
1526 //==> Cm(Z=96)
1527 static const G4int N96=10;
1528 static const G4double pZ96N145[7]={5.E-5, 22., .027, 1.1E-9, 2.2, 2.2E-5, .006};
1529 static const std::pair<G4int, const G4double*> Z96N145(145,pZ96N145);
1530 static const G4double pZ96N146[7]={5.E-5, 24., .027, 2.E-9, 2.2, 2.2E-5, .0055};
1531 static const std::pair<G4int, const G4double*> Z96N146(146,pZ96N146);
1532 static const G4double pZ96N147[7]={5.E-5, 22., .025, 2.5E-9, 2.2, 2.4E-5, .0044};
1533 static const std::pair<G4int, const G4double*> Z96N147(147,pZ96N147);
1534 static const G4double pZ96N148[7]={5.E-5, 23., .028, 1.9E-9, 2.2, 3.E-5, .0055};
1535 static const std::pair<G4int, const G4double*> Z96N148(148,pZ96N148);
1536 static const G4double pZ96N149[7]={5.E-5, 23., .025, 1.6E-9, 3., 3.5E-5, .0054};
1537 static const std::pair<G4int, const G4double*> Z96N149(149,pZ96N149);
1538 static const G4double pZ96N150[7]={5.E-5, 24., .026, 2.E-9, 3., 3.6E-5, .0045};
1539 static const std::pair<G4int, const G4double*> Z96N150(150,pZ96N150);
1540 static const G4double pZ96N151[7]={5.E-5, 24., .022, 2.4E-9, 3., 3.6E-5, .0039};
1541 static const std::pair<G4int, const G4double*> Z96N151(151,pZ96N151);
1542 static const G4double pZ96N152[7]={6.5E-7, 25., 2.E-4, 3.4E-9, 3., 3.6E-5, .003};
1543 static const std::pair<G4int, const G4double*> Z96N152(152,pZ96N152);
1544 static const G4double pZ96N153[7]={1.6E-6, 21., 7.E-4, 1.4E-9, 3., 3.6E-5, .0045};
1545 static const std::pair<G4int, const G4double*> Z96N153(153,pZ96N153);
1546 static const G4double pZ96N154[7]={1.3E-5, 16., .016, 0., 3., 3.6E-5, .017};
1547 static const std::pair<G4int, const G4double*> Z96N154(154,pZ96N154);
1548 static const std::pair<G4int, const G4double*> Z96[N96]={Z96N145, Z96N146, Z96N147,
1549 Z96N148, Z96N149, Z96N150,
1550 Z96N151, Z96N152, Z96N153,
1551 Z96N154};
1552 //==> Bk(Z=97)
1553 static const G4int N97=2;
1554 static const G4double pZ97N152[7]={6.5E-7, 22., 3.5E-4, 2.7E-9, 3., 4.E-5, .004};
1555 static const std::pair<G4int, const G4double*> Z97N152(152,pZ97N152);
1556 static const G4double pZ97N153[7]={6.5E-6, 22., .0036, 1.E-9, 2.7, 4.E-5, .006};
1557 static const std::pair<G4int, const G4double*> Z97N153(153,pZ97N153);
1558 static const std::pair<G4int, const G4double*> Z97[N97]={Z97N152, Z97N153};
1559 //==> Cf(Z=98)
1560 static const G4int N98=6;
1561 static const G4double pZ98N151[7]={6.5E-6, 22., .0035, .9E-9, 3., 4.E-5, .0068};
1562 static const std::pair<G4int, const G4double*> Z98N151(151,pZ98N151);
1563 static const G4double pZ98N152[7]={1.3E-6, 22., 7.E-4, 2.E-9, 2.7, 4.E-5, .0045};
1564 static const std::pair<G4int, const G4double*> Z98N152(152,pZ98N152);
1565 static const G4double pZ98N153[7]={2.6E-6, 22., .0014, 2.1E-9, 2.7, 4.E-5, .0044};
1566 static const std::pair<G4int, const G4double*> Z98N153(153,pZ98N153);
1567 static const G4double pZ98N154[7]={2.6E-6, 22., .0014, 1.3E-9, 2.7, 4.E-5, .0054};
1568 static const std::pair<G4int, const G4double*> Z98N154(154,pZ98N154);
1569 static const G4double pZ98N155[7]={2.6E-5, 22., .03, 0., 2.7, 4.E-5, .03};
1570 static const std::pair<G4int, const G4double*> Z98N155(155,pZ98N155);
1571 static const G4double pZ98N156[7]={5.2E-7, 22., 2.6E-4, 1.3E-9, 2.7, 4.E-5, .005};
1572 static const std::pair<G4int, const G4double*> Z98N156(156,pZ98N156);
1573 static const std::pair<G4int, const G4double*> Z98[N98]={Z98N151, Z98N152, Z98N153,
1574 Z98N154, Z98N155, Z98N156};
1575
1576 static const G4int NZ=99; // #of Elements covered by CHIPS elastic
1577 static const std::pair<G4int, const G4double*>* Pars[NZ]={Z0,Z1,Z2,Z3,Z4,Z5,Z6,Z7,Z8,Z9,
1578 Z10,Z11,Z12,Z13,Z14,Z15,Z16,Z17,Z18,Z19,Z20,Z21,Z22,Z23,Z24,Z25,Z26,Z27,Z28,Z29,Z30,
1579 Z31,Z32,Z33,Z34,Z35,Z36,Z37,Z38,Z39,Z40,Z41,Z42,Z43,Z44,Z45,Z46,Z47,Z48,Z49,Z50,Z51,
1580 Z52,Z53,Z54,Z55,Z56,Z57,Z58,Z59,Z60,Z61,Z62,Z63,Z64,Z65,Z66,Z67,Z68,Z69,Z70,Z71,Z72,
1581 Z73,Z74,Z75,Z76,Z77,Z78,Z79,Z80,Z81,Z82,Z83,Z84,Z85,Z86,Z87,Z88,Z89,Z90,Z91,Z92,Z93,
1582 Z94,Z95,Z96,Z97,Z98};
1583 static const G4int NIso[NZ]={N0,N1,N2,N3,N4,N5,N6,N7,N8,N9,N10,N11,N12,N13,N14,N15,N16,
1584 N17,N18,N19,N20,N21,N22,N23,N24,N25,N26,N27,N28,N29,N30,N31,N32,N33,N34,N35,N36,N37,
1585 N38,N39,N40,N41,N42,N43,N44,N45,N46,N47,N48,N49,N50,N51,N52,N53,N54,N55,N56,N57,N58,
1586 N59,N60,N61,N62,N63,N64,N65,N66,N67,N68,N69,N70,N71,N72,N73,N74,N75,N76,N77,N78,N79,
1587 N80,N81,N82,N83,N84,N85,N86,N87,N88,N89,N90,N91,N92,N93,N94,N95,N96,N97,N98};
1588 if(PDG==2112)
1589 {
1590 // --- Total np elastic cross section cs & s1/b1 (t), s2/b2 (u) --- NotTuned for highE
1591 //p2=p*p;p3=p2*p;sp=sqrt(p);p2s=p2*sp;lp=log(p);dl1=lp-(5.=par(3));p4=p2*p2; p=|3-mom|
1592 //CS=12./(p2s+.05*p+.0001/sqrt(sp))+.35/p+(6.75+.14*dl1*dl1+19./p)/(1.+.6/p4);
1593 // par(0) par(1) par(2) par(4) par(5) par(6) par(7) par(8)
1594 //s1=(6.75+.14*dl2*dl2+13./p)/(1.+.14/p4)+.6/(p4+.00013), s2=(75.+.001/p4/p)/p3
1595 // par(9) par(10) par(11) par(12) par(13) par(14) par(15) par(16)
1596 //b1=(7.2+4.32/(p4*p4+.012*p3))/(1.+2.5/p4), ss=0., b2=12./(p*sp+.34)
1597 //par(17) par(18) par(19) par(20) par(21) par(22) par(23)
1598 //
1599 if(lastPAR[nLast]!=pwd) // A unique flag to avoid the repeatable definition
1600 {
1601 if ( tgZ == 1 && tgN == 0 )
1602 {
1603 for (G4int ip=0; ip<n_npel; ip++) lastPAR[ip]=np_el[ip]; // np
1604 }
1605 else if ( tgZ == 0 && tgN == 1 )
1606 {
1607 for (G4int ip=0; ip<n_ppel; ip++) lastPAR[ip]=pp_el[ip]; // nn
1608 }
1609 else
1610 {
1611 G4double a=tgZ+tgN;
1612 G4double ala=std::log(a); // for powers of a
1613 G4double sa=std::sqrt(a);
1614 G4double ssa=std::sqrt(sa);
1615 G4double asa=a*sa;
1616 G4double a2=a*a;
1617 G4double a3=a2*a;
1618 G4double a4=a3*a;
1619 G4double a5=a4*a;
1620 G4double a6=a4*a2;
1621 G4double a7=a6*a;
1622 G4double a8=a7*a;
1623 G4double a9=a8*a;
1624 G4double a10=a5*a5;
1625 G4double a12=a6*a6;
1626 G4double a14=a7*a7;
1627 G4double a16=a8*a8;
1628 G4double a17=a16*a;
1629 G4double a32=a16*a16;
1630 // Reaction cross-section parameters (na_el.f)
1631 lastPAR[ 0]=5./(1.+22./asa); // p1
1632 lastPAR[ 1]=4.8*std::exp(ala*1.14)/(1.+3.6/a3); // p2
1633 lastPAR[ 2]=1./(1.+.004*a4)+2.E-6*a3/(1.+1.3E-6*a3); // p3
1634 lastPAR[ 3]=.07*asa/(1.+.009*a2); // o4
1635 lastPAR[ 5]=1.7*a; // p5
1636 lastPAR[ 6]=5.5E-6*std::exp(ala*1.3); // p6
1637 lastPAR[13]=0.; // reserved
1638 lastPAR[14]=0.; // reserved
1639 G4int nn=NIso[tgZ];
1640 G4bool nfound=true;
1641 if(nn) for (G4int in=0; in<nn; in++)
1642 {
1643 std::pair<G4int, const G4double*> curIs=Pars[tgZ][in];
1644 G4int cn=curIs.first;
1645 if(cn == tgN)
1646 {
1647 const G4double* curT=curIs.second;
1648 lastPAR[ 4]=curT[0]; // p4
1649 lastPAR[ 7]=curT[1]; // p7
1650 lastPAR[ 8]=curT[2]; // p8
1651 lastPAR[ 9]=curT[3]; // p9
1652 lastPAR[10]=curT[4]; // p10
1653 lastPAR[11]=curT[5]; // p11
1654 lastPAR[12]=curT[6]; // p12
1655 nfound = false;
1656 break;
1657 }
1658 }
1659 if(nfound)
1660 {
1661 G4cout<<"-Warning-G4ChipsNeutronElasticXS::CalcCS: Z="<<tgZ<<", N="<<tgN
1662 <<" isotope is not implemented in CHIPS"<<G4endl; // Put default values:
1663 lastPAR[ 4]=5.2E-7; // p4
1664 lastPAR[ 7]=22.; // p7
1665 lastPAR[ 8]=.00026; // p8
1666 lastPAR[ 9]=1.3E-9; // p9
1667 lastPAR[10]=2.7; // p10
1668 lastPAR[11]=4.E-5; // p11
1669 lastPAR[12]=.005; // p12
1670 }
1671 // @@ the differential cross-section is parameterized separately for A>6 & A<7
1672 if(a<6.5)
1673 {
1674 G4double a28=a16*a12;
1675 // The main pre-exponent (pel_sg)
1676 lastPAR[15]=4000*a; // p1
1677 lastPAR[16]=1.2e7*a8+380*a17; // p2
1678 lastPAR[17]=.7/(1.+4.e-12*a16); // p3
1679 lastPAR[18]=2.5/a8/(a4+1.e-16*a32); // p4
1680 lastPAR[19]=.28*a; // p5
1681 lastPAR[20]=1.2*a2+2.3; // p6
1682 lastPAR[21]=3.8/a; // p7
1683 // The main slope (pel_sl)
1684 lastPAR[22]=.01/(1.+.0024*a5); // p1
1685 lastPAR[23]=.2*a; // p2
1686 lastPAR[24]=9.e-7/(1.+.035*a5); // p3
1687 lastPAR[25]=(42.+2.7e-11*a16)/(1.+.14*a); // p4
1688 // The main quadratic (pel_sh)
1689 lastPAR[26]=2.25*a3; // p1
1690 lastPAR[27]=18.; // p2
1691 lastPAR[28]=.0024*a8/(1.+2.6e-4*a7); // p3
1692 lastPAR[29]=3.5e-36*a32*a8/(1.+5.e-15*a32/a); // p4
1693 lastPAR[30]=1.e5/(a8+2.5e12/a16); // p1
1694 lastPAR[31]=8.e7/(a12+1.e-27*a28*a28); // p2
1695 lastPAR[32]=.0006*a3; // p3
1696 // The 1st max slope (pel_qs)
1697 lastPAR[33]=10.+4.e-8*a12*a; // p1
1698 lastPAR[34]=.114; // p2
1699 lastPAR[35]=.003; // p3
1700 lastPAR[36]=2.e-23; // p4
1701 // The effective pre-exponent (pel_ss)
1702 lastPAR[37]=1./(1.+.0001*a8); // p1
1703 lastPAR[38]=1.5e-4/(1.+5.e-6*a12); // p2
1704 lastPAR[39]=.03; // p3
1705 // The effective slope (pel_sb)
1706 lastPAR[40]=a/2; // p1
1707 lastPAR[41]=2.e-7*a4; // p2
1708 lastPAR[42]=4.; // p3
1709 lastPAR[43]=64./a3; // p4
1710 // The gloria pre-exponent (pel_us)
1711 lastPAR[44]=1.e8*std::exp(.32*asa); // p1
1712 lastPAR[45]=20.*std::exp(.45*asa); // p2
1713 lastPAR[46]=7.e3+2.4e6/a5; // p3
1714 lastPAR[47]=2.5e5*std::exp(.085*a3); // p4
1715 lastPAR[48]=2.5*a; // p5
1716 // The gloria slope (pel_ub)
1717 lastPAR[49]=920.+.03*a8*a3; // p1
1718 lastPAR[50]=93.+.0023*a12; // p2
1719 }
1720 else
1721 {
1722 G4double p1a10=2.2e-28*a10;
1723 G4double r4a16=6.e14/a16;
1724 G4double s4a16=r4a16*r4a16;
1725 // a24
1726 // a36
1727 // The main pre-exponent (peh_sg)
1728 lastPAR[15]=4.5*std::pow(a,1.15); // p1
1729 lastPAR[16]=.06*std::pow(a,.6); // p2
1730 lastPAR[17]=.6*a/(1.+2.e15/a16); // p3
1731 lastPAR[18]=.17/(a+9.e5/a3+1.5e33/a32); // p4
1732 lastPAR[19]=(.001+7.e-11*a5)/(1.+4.4e-11*a5); // p5
1733 lastPAR[20]=(p1a10*p1a10+2.e-29)/(1.+2.e-22*a12); // p6
1734 // The main slope (peh_sl)
1735 lastPAR[21]=400./a12+2.e-22*a9; // p1
1736 lastPAR[22]=1.e-32*a12/(1.+5.e22/a14); // p2
1737 lastPAR[23]=1000./a2+9.5*sa*ssa; // p3
1738 lastPAR[24]=4.e-6*a*asa+1.e11/a16; // p4
1739 lastPAR[25]=(120./a+.002*a2)/(1.+2.e14/a16); // p5
1740 lastPAR[26]=9.+100./a; // p6
1741 // The main quadratic (peh_sh)
1742 lastPAR[27]=.002*a3+3.e7/a6; // p1
1743 lastPAR[28]=7.e-15*a4*asa; // p2
1744 lastPAR[29]=9000./a4; // p3
1745 // The 1st max pre-exponent (peh_qq)
1746 lastPAR[30]=.0011*asa/(1.+3.e34/a32/a4); // p1
1747 lastPAR[31]=1.e-5*a2+2.e14/a16; // p2
1748 lastPAR[32]=1.2e-11*a2/(1.+1.5e19/a12); // p3
1749 lastPAR[33]=.016*asa/(1.+5.e16/a16); // p4
1750 // The 1st max slope (peh_qs)
1751 lastPAR[34]=.002*a4/(1.+7.e7/std::pow(a-6.83,14)); // p1
1752 lastPAR[35]=2.e6/a6+7.2/std::pow(a,.11); // p2
1753 lastPAR[36]=11.*a3/(1.+7.e23/a16/a8); // p3
1754 lastPAR[37]=100./asa; // p4
1755 // The 2nd max pre-exponent (peh_ss)
1756 lastPAR[38]=(.1+4.4e-5*a2)/(1.+5.e5/a4); // p1
1757 lastPAR[39]=3.5e-4*a2/(1.+1.e8/a8); // p2
1758 lastPAR[40]=1.3+3.e5/a4; // p3
1759 lastPAR[41]=500./(a2+50.)+3; // p4
1760 lastPAR[42]=1.e-9/a+s4a16*s4a16; // p5
1761 // The 2nd max slope (peh_sb)
1762 lastPAR[43]=.4*asa+3.e-9*a6; // p1
1763 lastPAR[44]=.0005*a5; // p2
1764 lastPAR[45]=.002*a5; // p3
1765 lastPAR[46]=10.; // p4
1766 // The effective pre-exponent (peh_us)
1767 lastPAR[47]=.05+.005*a; // p1
1768 lastPAR[48]=7.e-8/sa; // p2
1769 lastPAR[49]=.8*sa; // p3
1770 lastPAR[50]=.02*sa; // p4
1771 lastPAR[51]=1.e8/a3; // p5
1772 lastPAR[52]=3.e32/(a32+1.e32); // p6
1773 // The effective slope (peh_ub)
1774 lastPAR[53]=24.; // p1
1775 lastPAR[54]=20./sa; // p2
1776 lastPAR[55]=7.e3*a/(sa+1.); // p3
1777 lastPAR[56]=900.*sa/(1.+500./a3); // p4
1778 }
1779 // Parameter for lowEnergyNeutrons
1780 lastPAR[57]=1.e15+2.e27/a4/(1.+2.e-18*a16);
1781 }
1782 lastPAR[nLast]=pwd;
1783 // and initialize the zero element of the table
1784 G4double lp=lPMin; // ln(momentum)
1785 G4bool memCS=onlyCS; // ??
1786 onlyCS=false;
1787 lastCST[0]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables
1788 onlyCS=memCS;
1789 lastSST[0]=theSS;
1790 lastS1T[0]=theS1;
1791 lastB1T[0]=theB1;
1792 lastS2T[0]=theS2;
1793 lastB2T[0]=theB2;
1794 lastS3T[0]=theS3;
1795 lastB3T[0]=theB3;
1796 lastS4T[0]=theS4;
1797 lastB4T[0]=theB4;
1798 }
1799 if(LP>ILP)
1800 {
1801 G4int ini = static_cast<int>((ILP-lPMin+.000001)/dlnP)+1; // already inited till this
1802 if(ini<0) ini=0;
1803 if(ini<nPoints)
1804 {
1805 G4int fin = static_cast<int>((LP-lPMin)/dlnP)+1; // final bin of initialization
1806 if(fin>=nPoints) fin=nLast; // Limit of the tabular initialization
1807 if(fin>=ini)
1808 {
1809 G4double lp=0.;
1810 for(G4int ip=ini; ip<=fin; ip++) // Calculate tabular CS,S1,B1,S2,B2,S3,B3
1811 {
1812 lp=lPMin+ip*dlnP; // ln(momentum)
1813 G4bool memCS=onlyCS;
1814 onlyCS=false;
1815 lastCST[ip]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables (ret CS)
1816 onlyCS=memCS;
1817 lastSST[ip]=theSS;
1818 lastS1T[ip]=theS1;
1819 lastB1T[ip]=theB1;
1820 lastS2T[ip]=theS2;
1821 lastB2T[ip]=theB2;
1822 lastS3T[ip]=theS3;
1823 lastB3T[ip]=theB3;
1824 lastS4T[ip]=theS4;
1825 lastB4T[ip]=theB4;
1826 }
1827 return lp;
1828 }
1829 else G4cout<<"*Warning*G4ChipsNeutronElasticXS::GetPTables: PDG="<<PDG
1830 <<", Z="<<tgZ<<", N="<<tgN<<", i="<<ini<<" > fin="<<fin<<", LP="<<LP
1831 <<" > ILP="<<ILP<<" nothing is done!"<<G4endl;
1832 }
1833 else G4cout<<"*Warning*G4ChipsNeutronElasticXS::GetPTables: PDG="<<PDG<<", Z="
1834 <<tgZ<<", N="<<tgN<<", i="<<ini<<">= max="<<nPoints<<", LP="<<LP
1835 <<" > ILP="<<ILP<<", lPMax="<<lPMax<<" nothing is done!"<<G4endl;
1836 }
1837 }
1838 else
1839 {
1840 // G4cout<<"*Error*G4ChipsNeutronElasticXS::GetPTables: PDG="<<PDG<<", Z="<<tgZ
1841 // <<", N="<<tgN<<", while it is defined only for PDG=2112(n)"<<G4endl;
1842 // throw G4QException("G4ChipsNeutronElasticXS::GetPTables:only nA're implemented");
1844 ed << "PDG = " << PDG << ", Z = " << tgZ <<", N = " << tgN
1845 << ", while it is defined only for PDG=2112 (n)" << G4endl;
1846 G4Exception("G4ChipsNeutronElasticXS::GetPTables()", "HAD_CHPS_0000",
1847 FatalException, ed);
1848 }
1849 return ILP;
1850}
1851
1852// Returns Q2=-t in independent units (MeV^2) (all internal calculations are in GeV)
1854{
1855 static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
1856 static const G4double third=1./3.;
1857 static const G4double fifth=1./5.;
1858 static const G4double sevth=1./7.;
1859 if(PDG!=2112) G4cout<<"*Warning*G4ChipsNeutronElasticXS::GetExT:PDG="<<PDG<<G4endl;
1860 if(onlyCS) G4cout<<"*Warning*G4ChipsNeutronElasticXS::GetExchangeT:onCS=1"<<G4endl;
1861 if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
1862 G4double q2=0.;
1863 if(tgZ==1 && tgN==0) // ===> n+p=n+p
1864 {
1865 G4double E1=lastTM*theB1;
1866 G4double R1=(1.-std::exp(-E1));
1867 G4double E2=lastTM*theB2;
1868 G4double R2=(1.-std::exp(-E2));
1869 G4double I1=R1*theS1;
1870 G4double I2=R2*theS2/theB2;
1871 //G4double I3=R3*theS3/theB3;
1872 G4double I12=I1+I2;
1873 //G4double rand=(I12+I3)*G4UniformRand();
1874 G4double rand=I12*G4UniformRand();
1875 if (rand<I1 )
1876 {
1877 G4double ran=R1*G4UniformRand();
1878 if(ran>1.) ran=1.;
1879 q2=-std::log(1.-ran)/theB1; // t-chan
1880 }
1881 else
1882 {
1883 G4double ran=R2*G4UniformRand();
1884 if(ran>1.) ran=1.;
1885 q2=lastTM+std::log(1.-ran)/theB2; // u-chan (ChEx)
1886 }
1887 }
1888 else
1889 {
1890 G4double a=tgZ+tgN;
1891 G4double E1=lastTM*(theB1+lastTM*theSS);
1892 G4double R1=(1.-std::exp(-E1));
1893 G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
1894 G4double tm2=lastTM*lastTM;
1895 G4double E2=lastTM*tm2*theB2; // power 3 for lowA, 5 for HighA (1st)
1896 if(a>6.5)E2*=tm2; // for heavy nuclei
1897 G4double R2=(1.-std::exp(-E2));
1898 G4double E3=lastTM*theB3;
1899 if(a>6.5)E3*=tm2*tm2*tm2; // power 1 for lowA, 7 (2nd) for HighA
1900 G4double R3=(1.-std::exp(-E3));
1901 G4double E4=lastTM*theB4;
1902 G4double R4=(1.-std::exp(-E4));
1903 G4double I1=R1*theS1;
1904 G4double I2=R2*theS2;
1905 G4double I3=R3*theS3;
1906 G4double I4=R4*theS4;
1907 G4double I12=I1+I2;
1908 G4double I13=I12+I3;
1909 G4double rand=(I13+I4)*G4UniformRand();
1910 if(rand<I1)
1911 {
1912 G4double ran=R1*G4UniformRand();
1913 if(ran>1.) ran=1.;
1914 q2=-std::log(1.-ran)/theB1;
1915 if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
1916 }
1917 else if(rand<I12)
1918 {
1919 G4double ran=R2*G4UniformRand();
1920 if(ran>1.) ran=1.;
1921 q2=-std::log(1.-ran)/theB2;
1922 if(q2<0.) q2=0.;
1923 if(a<6.5) q2=std::pow(q2,third);
1924 else q2=std::pow(q2,fifth);
1925 }
1926 else if(rand<I13)
1927 {
1928 G4double ran=R3*G4UniformRand();
1929 if(ran>1.) ran=1.;
1930 q2=-std::log(1.-ran)/theB3;
1931 if(q2<0.) q2=0.;
1932 if(a>6.5) q2=std::pow(q2,sevth);
1933 }
1934 else
1935 {
1936 G4double ran=R4*G4UniformRand();
1937 if(ran>1.) ran=1.;
1938 q2=-std::log(1.-ran)/theB4;
1939 if(a<6.5) q2=lastTM-q2; // u reduced for lightA (starts from 0)
1940 }
1941 }
1942 if(q2<0.) q2=0.;
1943 if(!(q2>=-1.||q2<=1.)) G4cout<<"*NAN*G4QNeutronElCroSect::GetExchangeT: -t="<<q2<<G4endl;
1944 if(q2>lastTM)
1945 {
1946 q2=lastTM;
1947 }
1948 return q2*GeVSQ;
1949}
1950
1951// Returns B in independent units (MeV^-2) (all internal calculations are in GeV) see ExT
1952G4double G4ChipsNeutronElasticXS::GetSlope(G4int tgZ, G4int tgN, G4int PDG)
1953{
1954 static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
1955
1956 if(onlyCS) G4cout<<"Warning*G4ChipsNeutronElasticXS::GetSlope:onlyCS=true"<<G4endl;
1957 if(lastLP<-4.3) return 0.; // S-wave for p<14 MeV/c (kinE<.1MeV)
1958 if(PDG!=2112)
1959 {
1960 // G4cout<<"*Error*G4ChipsNeutronElasticXS::GetSlope: PDG="<<PDG<<", Z="<<tgZ
1961 // <<", N="<<tgN<<", while it is defined only for PDG=2112"<<G4endl;
1962 // throw G4QException("G4ChipsNeutronElasticXS::GetSlope: only nA are implemented");
1964 ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
1965 <<", while it is defined only for PDG=2112 (n) " << G4endl;
1966 G4Exception("G4ChipsNeutronElasticXS::GetSlope()", "HAD_CHPS_0000",
1967 FatalException, ed);
1968 }
1969 if(theB1<0.) theB1=0.;
1970 if(!(theB1>=-1.||theB1<=1.))G4cout<<"*NAN*G4QNeutElasticCrosS::Getslope:"<<theB1<<G4endl;
1971 return theB1/GeVSQ;
1972}
1973
1974// Returns half max(Q2=-t) in independent units (MeV^2)
1976{
1977 static const G4double HGeVSQ=gigaelectronvolt*gigaelectronvolt/2.;
1978 return lastTM*HGeVSQ;
1979}
1980
1981// lastLP is used, so calculating tables, one need to remember and then recover lastLP
1982G4double G4ChipsNeutronElasticXS::GetTabValues(G4double lp, G4int PDG, G4int tgZ,
1983 G4int tgN)
1984{
1985 if(PDG!=2112) G4cout<<"*Warning*G4ChipsNeutronElasticXS::GetTaV:PDG="<<PDG<<G4endl;
1986 if(tgZ<0 || tgZ>92)
1987 {
1988 G4cout<<"*Warning*G4QNElasticCrS::GetTabValue: (1-92) No isotopes for Z="<<tgZ<<G4endl;
1989 return 0.;
1990 }
1991 G4int iZ=tgZ-1; // Z index
1992 if(iZ<0)
1993 {
1994 iZ=0; // conversion of the neutron target to the proton target
1995 tgZ=1;
1996 tgN=0;
1997 }
1998 G4double p=std::exp(lp); // momentum
1999 G4double sp=std::sqrt(p); // sqrt(p)
2000 G4double p2=p*p;
2001 G4double p3=p2*p;
2002 G4double p4=p3*p;
2003 if ( tgZ == 1 && tgN == 0 ) // np
2004 {
2005 G4double ssp=std::sqrt(sp); // sqrt(sqrt(p))=p^.25
2006 G4double p2s=p2*sp;
2007 G4double dl1=lp-lastPAR[3];
2008 theSS=lastPAR[27];
2009 theS1=(lastPAR[9]+lastPAR[10]*dl1*dl1+lastPAR[11]/p)/(1.+lastPAR[12]/p4)
2010 +lastPAR[13]/(p4+lastPAR[14]);
2011 theB1=(lastPAR[17]+lastPAR[18]/(p4*p4+lastPAR[19]*p3))/(1.+lastPAR[20]/p4);
2012 theS2=(lastPAR[15]+lastPAR[16]/p4/p)/p3;
2013 theB2=lastPAR[22]/(p*sp+lastPAR[23]);
2014 theS3=0.;
2015 theB3=0.;
2016 theS4=0.;
2017 theB4=0.;
2018 // Returns the total elastic pp cross-section (to avoid spoiling lastSIG)
2019 return lastPAR[0]/(p2s+lastPAR[1]*p+lastPAR[2]/ssp)+lastPAR[4]/p
2020 +(lastPAR[5]+lastPAR[6]*dl1*dl1+lastPAR[7]/p)/(1.+lastPAR[8]/p4);
2021
2022 }
2023 else
2024 {
2025 G4double p5=p4*p;
2026 G4double p6=p5*p;
2027 G4double p8=p6*p2;
2028 G4double p10=p8*p2;
2029 G4double p12=p10*p2;
2030 G4double p16=p8*p8;
2031 G4double dl=lp-5.;
2032 G4double a=tgZ+tgN;
2033 G4double pah=std::pow(p,a/2);
2034 G4double pa=pah*pah;
2035 G4double pa2=pa*pa;
2036 if(a<6.5)
2037 {
2038 theS1=lastPAR[15]/(1.+lastPAR[16]*p4*pa)+lastPAR[17]/(p4+lastPAR[18]*p4/pa2)+
2039 (lastPAR[19]*dl*dl+lastPAR[20])/(1.+lastPAR[21]/p2);
2040 theB1=(lastPAR[22]+lastPAR[23]*p2)/(p4+lastPAR[24]/pah)+lastPAR[25];
2041 theSS=lastPAR[26]/(1.+lastPAR[27]/p2)+lastPAR[28]/(p6/pa+lastPAR[29]/p16);
2042 theS2=lastPAR[30]/(pa/p2+lastPAR[31]/p4)+lastPAR[32];
2043 theB2=lastPAR[33]*std::pow(p,lastPAR[34])+lastPAR[35]/(p8+lastPAR[36]/p16);
2044 theS3=lastPAR[37]/(pa*p+lastPAR[38]/pa)+lastPAR[39];
2045 theB3=lastPAR[40]/(p3+lastPAR[41]/p6)+lastPAR[42]/(1.+lastPAR[43]/p2);
2046 theS4=p2*(pah*lastPAR[44]*std::exp(-pah*lastPAR[45])+
2047 lastPAR[46]/(1.+lastPAR[47]*std::pow(p,lastPAR[48])));
2048 theB4=lastPAR[49]*pa/p2/(1.+pa*lastPAR[50]);
2049 }
2050 else
2051 {
2052 theS1=lastPAR[15]/(1.+lastPAR[16]/p4)+lastPAR[17]/(p4+lastPAR[18]/p2)+
2053 lastPAR[19]/(p5+lastPAR[20]/p16);
2054 theB1=(lastPAR[21]/p8+lastPAR[25])/(p+lastPAR[22]/std::pow(p,lastPAR[26]))+
2055 lastPAR[23]/(1.+lastPAR[24]/p4);
2056 theSS=lastPAR[27]/(p4/std::pow(p,lastPAR[29])+lastPAR[28]/p4);
2057 theS2=lastPAR[30]/p4/(std::pow(p,lastPAR[31])+lastPAR[32]/p12)+lastPAR[33];
2058 theB2=lastPAR[34]/std::pow(p,lastPAR[35])+lastPAR[36]/std::pow(p,lastPAR[37]);
2059 theS3=lastPAR[38]/std::pow(p,lastPAR[41])/(1.+lastPAR[42]/p12)+
2060 lastPAR[39]/(1.+lastPAR[40]/p6);
2061 theB3=lastPAR[43]/p8+lastPAR[44]/p2+lastPAR[45]/(1.+lastPAR[46]/p8);
2062 theS4=(lastPAR[47]/p4+lastPAR[52]/p)/(1.+lastPAR[48]/p10)+
2063 (lastPAR[49]+lastPAR[50]*dl*dl)/(1.+lastPAR[51]/p12);
2064 theB4=lastPAR[53]/(1.+lastPAR[54]/p)+lastPAR[55]*p4/(1.+lastPAR[56]*p5);
2065 }
2066 // Returns the total elastic (n/p)A cross-section (to avoid spoiling lastSIG)
2067 // p1(p6) p2(p7) p3(p4) o4(p8) (p9)p5
2068 return (lastPAR[0]*dl*dl+lastPAR[1])/(1.+lastPAR[2]/p+lastPAR[3]/p4)+lastPAR[5]/
2069 (p3+lastPAR[6]/p3)+lastPAR[7]/(p2+lastPAR[4]/(p2+lastPAR[8])+lastPAR[9]/p)+
2070 lastPAR[10]/(p5+lastPAR[11]/p2)+lastPAR[12]/p;
2071 // p10 p11 p12
2072 }
2073 return 0.;
2074} // End of GetTableValues
2075
2076// Returns max -t=Q2 (GeV^2) for the momentum pP(GeV) and the target nucleus (tgN,tgZ)
2077G4double G4ChipsNeutronElasticXS::GetQ2max(G4int PDG, G4int tgZ, G4int tgN,
2078 G4double pP)
2079{
2080 static const G4double mNeut= G4Neutron::Neutron()->GetPDGMass()*.001; // MeV to GeV
2081 static const G4double mProt= G4Proton::Proton()->GetPDGMass()*.001; // MeV to GeV
2082
2083 static const G4double mNeut2= mNeut*mNeut;
2084
2085 G4double pP2=pP*pP; // squared momentum of the projectile
2086 if(tgZ==0 && tgN==1)
2087 {
2088 G4double tMid=std::sqrt(pP2+mNeut2)*mNeut-mNeut2; // CMS 90deg value of -t=Q2 (GeV^2)
2089 return tMid+tMid;
2090 }
2091 else if(tgZ || tgN) // ---> nA
2092 {
2093 G4double mt=mProt; // Target mass in GeV
2094 if(tgN||tgZ>1) mt=G4ParticleTable::GetParticleTable()->FindIon(tgZ,tgZ+tgN,0,tgZ)->GetPDGMass()*.001; // Target mass in GeV
2095 G4double dmt=mt+mt;
2096 G4double mds=dmt*std::sqrt(pP2+mNeut2)+mNeut2+mt*mt; // Mondelstam mds (in GeV^2)
2097 return dmt*dmt*pP2/mds;
2098 }
2099 else
2100 {
2101 // G4cout<<"*Error*G4ChipsNeutronElasticXS::GetQ2max:PDG="<<PDG<<", Z="<<tgZ<<", N="
2102 // <<tgN<<", while it is defined only for n projectiles & Z_target>0"<<G4endl;
2103 // throw G4QException("G4ChipsNeutronElasticXS::GetQ2max: only nA implemented");
2105 ed << "PDG = " << PDG << ", Z = " << tgZ << ", N =" << tgN
2106 <<", while it is defined only for n projectiles & Z_target>0" << G4endl;
2107 G4Exception("G4ChipsNeutronElasticXS::GetQ2max()", "HAD_CHPS_0000",
2108 FatalException, ed);
2109 return 0;
2110 }
2111}
#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
G4double GetExchangeT(G4int tZ, G4int tN, G4int pPDG)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
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)
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