Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/src/G4AtimaEnergyLossModel.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

Diff markup

Differences between /processes/electromagnetic/standard/src/G4AtimaEnergyLossModel.cc (Version 11.3.0) and /processes/electromagnetic/standard/src/G4AtimaEnergyLossModel.cc (Version 10.0.p3)


  1 //                                                  1 
  2 // *******************************************    
  3 // * License and Disclaimer                       
  4 // *                                              
  5 // * The  Geant4 software  is  copyright of th    
  6 // * the Geant4 Collaboration.  It is provided    
  7 // * conditions of the Geant4 Software License    
  8 // * LICENSE and available at  http://cern.ch/    
  9 // * include a list of copyright holders.         
 10 // *                                              
 11 // * Neither the authors of this software syst    
 12 // * institutes,nor the agencies providing fin    
 13 // * work  make  any representation or  warran    
 14 // * regarding  this  software system or assum    
 15 // * use.  Please see the license in the file     
 16 // * for the full disclaimer and the limitatio    
 17 // *                                              
 18 // * This  code  implementation is the result     
 19 // * technical work of the GEANT4 collaboratio    
 20 // * By using,  copying,  modifying or  distri    
 21 // * any work based  on the software)  you  ag    
 22 // * use  in  resulting  scientific  publicati    
 23 // * acceptance of all terms of the Geant4 Sof    
 24 // *******************************************    
 25 //                                                
 26 // GEANT4 Class header file                       
 27 //                                                
 28 //                                                
 29 // File name:     G4AtimaEnergyLossModel          
 30 //                                                
 31 // Author:        Jose Luis Rodriguez Sanchez     
 32 //                                                
 33 // Creation date: 16.01.2018                      
 34 //                                                
 35 // Modifications:                                 
 36 // 09/10/2018: Solved bug in the determination    
 37 //                                                
 38 // Class Description:                             
 39 //                                                
 40 // This model calculates the stopping power of    
 41 // to the ATIMA code developed at GSI, Darmsta    
 42 // For details: http://web-docs.gsi.de/~weick/    
 43 //                                                
 44 // Helmut Weick, GSI (responsible for fortran     
 45 // Andrej Prochazka, GSI (responsible for C ve    
 46 // Christoph Scheidenberger, GSI (project coor    
 47 //                                                
 48 // -------------------------------------------    
 49 //                                                
 50 //....oooOO0OOooo........oooOO0OOooo........oo    
 51 //....oooOO0OOooo........oooOO0OOooo........oo    
 52                                                   
 53 #include "G4AtimaEnergyLossModel.hh"              
 54 #include "Randomize.hh"                           
 55 #include "G4PhysicalConstants.hh"                 
 56 #include "G4SystemOfUnits.hh"                     
 57 #include "G4Electron.hh"                          
 58 #include "G4GenericIon.hh"                        
 59 #include "G4LossTableManager.hh"                  
 60 #include "G4EmCorrections.hh"                     
 61 #include "G4ParticleChangeForLoss.hh"             
 62 #include "G4DeltaAngle.hh"                        
 63 #include "G4Log.hh"                               
 64 #include "G4Exp.hh"                               
 65 #include "G4Pow.hh"                               
 66                                                   
 67 //....oooOO0OOooo........oooOO0OOooo........oo    
 68                                                   
 69 using namespace std;                              
 70                                                   
 71 G4double G4AtimaEnergyLossModel::stepE = 0.0;     
 72 G4double G4AtimaEnergyLossModel::tableE[] = {0    
 73                                                   
 74 G4AtimaEnergyLossModel::G4AtimaEnergyLossModel    
 75                                      const G4S    
 76   : G4VEmModel(nam),                              
 77     particle(nullptr)                             
 78 {                                                 
 79   g4calc = G4Pow::GetInstance();                  
 80   fParticleChange = nullptr;                      
 81   theElectron = G4Electron::Electron();           
 82   corr = G4LossTableManager::Instance()->EmCor    
 83   nist = G4NistManager::Instance();               
 84   SetLowEnergyLimit(2.0*MeV);                     
 85   MLN10 = 2.30258509299;                          
 86   atomic_mass_unit = 931.4940954; // MeV/c^2      
 87   dedx_constant = 0.3070749187;  //4*pi*Na*me*    
 88   electron_mass = 0.510998928;   // MeV/c^2       
 89   fine_structure = 1.0/137.035999139;             
 90   domega2dx_constant = dedx_constant*electron_    
 91   if(tableE[0] == 0.0) {                          
 92     G4double logmin = 0.;                         
 93     G4double logmax = 5.;                         
 94     stepE = (logmax-logmin)/199.;                 
 95     for(G4int i=0; i<200; ++i){                   
 96       tableE[i] = G4Exp(MLN10*(logmin + ((G4do    
 97     }                                             
 98   }                                               
 99 }                                                 
100                                                   
101 //....oooOO0OOooo........oooOO0OOooo........oo    
102                                                   
103 G4AtimaEnergyLossModel::~G4AtimaEnergyLossMode    
104                                                   
105 //....oooOO0OOooo........oooOO0OOooo........oo    
106                                                   
107 void G4AtimaEnergyLossModel::Initialise(const     
108                                    const G4Dat    
109 {                                                 
110   SetGenericIon(p);                               
111   SetParticle(p);                                 
112                                                   
113   //G4cout << "G4AtimaEnergyLossModel::Initial    
114   //         << "  isIon= " << isIon              
115   //         << G4endl;                           
116                                                   
117   // always false before the run                  
118   SetDeexcitationFlag(false);                     
119                                                   
120   if(nullptr == fParticleChange) {                
121     fParticleChange = GetParticleChangeForLoss    
122     if(UseAngularGeneratorFlag() && !GetAngula    
123       SetAngularDistribution(new G4DeltaAngle(    
124     }                                             
125   }                                               
126 }                                                 
127                                                   
128 //....oooOO0OOooo........oooOO0OOooo........oo    
129                                                   
130 G4double G4AtimaEnergyLossModel::GetChargeSqua    
131                                                   
132                                                   
133 {                                                 
134   // this method is called only for ions          
135   G4double q2 = corr->EffectiveChargeSquareRat    
136   corrFactor = q2*corr->EffectiveChargeCorrect    
137   return corrFactor;                              
138 }                                                 
139                                                   
140 //....oooOO0OOooo........oooOO0OOooo........oo    
141                                                   
142 G4double G4AtimaEnergyLossModel::GetParticleCh    
143                                                   
144                                                   
145 {                                                 
146   //G4cout<<"G4AtimaEnergyLossModel::GetPartic    
147   //  " q= " <<  corr->GetParticleCharge(p,mat    
148   // this method is called only for ions, so n    
149   return corr->GetParticleCharge(p,mat,kinetic    
150 }                                                 
151                                                   
152 //....oooOO0OOooo........oooOO0OOooo........oo    
153                                                   
154 void G4AtimaEnergyLossModel::SetupParameters()    
155 {                                                 
156   mass = particle->GetPDGMass();                  
157   spin = particle->GetPDGSpin();                  
158   G4double q = particle->GetPDGCharge()/eplus;    
159   chargeSquare = q*q;                             
160   corrFactor = chargeSquare;                      
161   ratio = electron_mass_c2/mass;                  
162   static const G4double aMag = 1./(0.5*eplus*h    
163   G4double magmom = particle->GetPDGMagneticMo    
164   magMoment2 = magmom*magmom - 1.0;               
165   formfact = 0.0;                                 
166   tlimit = DBL_MAX;                               
167   if(particle->GetLeptonNumber() == 0) {          
168     G4int iz = G4lrint(q);                        
169     if(iz <= 1) {                                 
170       formfact = (spin == 0.0 && mass < GeV) ?    
171     } else {                                      
172       G4double x = nist->GetA27(iz);              
173       formfact = 3.969e-6*x*x;                    
174     }                                             
175     tlimit = std::sqrt(0.414/formfact +           
176                        electron_mass_c2*electr    
177   }                                               
178 }                                                 
179                                                   
180 //....oooOO0OOooo........oooOO0OOooo........oo    
181                                                   
182 G4double G4AtimaEnergyLossModel::MinEnergyCut(    
183                                          const    
184 {                                                 
185   return couple->GetMaterial()->GetIonisation(    
186 }                                                 
187                                                   
188 //....oooOO0OOooo........oooOO0OOooo........oo    
189                                                   
190 G4double                                          
191 G4AtimaEnergyLossModel::ComputeCrossSectionPer    
192                                                   
193                                                   
194                                                   
195 {                                                 
196   G4double cross = 0.0;                           
197   G4double tmax = MaxSecondaryEnergy(p, kineti    
198   G4double maxEnergy = std::min(tmax,maxKinEne    
199   if(cutEnergy < maxEnergy) {                     
200                                                   
201     G4double totEnergy = kineticEnergy + mass;    
202     G4double energy2   = totEnergy*totEnergy;     
203     G4double beta2     = kineticEnergy*(kineti    
204                                                   
205     cross = (maxEnergy - cutEnergy)/(cutEnergy    
206       - beta2*G4Log(maxEnergy/cutEnergy)/tmax;    
207                                                   
208     // +term for spin=1/2 particle                
209     if( 0.0 < spin ) { cross += 0.5*(maxEnergy    
210                                                   
211     cross *= twopi_mc2_rcl2*chargeSquare/beta2    
212   }                                               
213                                                   
214    // G4cout << "BB: e= " << kineticEnergy <<     
215    //        << " tmax= " << tmax << " cross=     
216                                                   
217   return cross;                                   
218 }                                                 
219                                                   
220 //....oooOO0OOooo........oooOO0OOooo........oo    
221                                                   
222 G4double G4AtimaEnergyLossModel::ComputeCrossS    
223                                            con    
224                                                   
225                                                   
226                                                   
227                                                   
228 {                                                 
229   G4double cross = Z*ComputeCrossSectionPerEle    
230                                          (p,ki    
231   return cross;                                   
232 }                                                 
233                                                   
234 //....oooOO0OOooo........oooOO0OOooo........oo    
235                                                   
236 G4double G4AtimaEnergyLossModel::CrossSectionP    
237                                            con    
238                                            con    
239                                                   
240                                                   
241                                                   
242 {                                                 
243   G4double eDensity = material->GetElectronDen    
244   G4double cross = eDensity*ComputeCrossSectio    
245                                          (p,ki    
246   return cross;                                   
247 }                                                 
248                                                   
249 //....oooOO0OOooo........oooOO0OOooo........oo    
250                                                   
251 G4double G4AtimaEnergyLossModel::ComputeDEDXPe    
252                                                   
253                                                   
254                                                   
255 {                                                 
256   //Call to ATIMA Stopping Power function         
257   G4double zt = material->GetIonisation()->Get    
258   zt = std::min(zt,93.);                          
259   G4double at = nist->GetAtomicMassAmu(G4lrint    
260   G4double dedx = StoppingPower(p->GetPDGMass(    
261     kineticEnergy/(MeV), at, zt) *material->Ge    
262   dedx = std::max(dedx, 0.0);                     
263                                                   
264   //  G4cout << "E(MeV)= " << kineticEnergy/Me    
265   //          << "  " << material->GetName() <    
266   return dedx;                                    
267 }                                                 
268                                                   
269 //....oooOO0OOooo........oooOO0OOooo........oo    
270                                                   
271 void G4AtimaEnergyLossModel::CorrectionsAlongS    
272                                                   
273                                                   
274                                                   
275 {                                                 
276   if(isIon) {                                     
277     const G4ParticleDefinition* p = dp->GetDef    
278     const G4Material* mat = couple->GetMateria    
279     G4double cutEnergy = DBL_MAX;                 
280     G4double kineticEnergy = dp->GetKineticEne    
281     eloss = ComputeDEDXPerVolume(mat, p, kinet    
282   }                                               
283 }                                                 
284                                                   
285 //....oooOO0OOooo........oooOO0OOooo........oo    
286                                                   
287 void G4AtimaEnergyLossModel::SampleSecondaries    
288                                           cons    
289                                           cons    
290                                           G4do    
291                                           G4do    
292 {                                                 
293   G4double kineticEnergy = dp->GetKineticEnerg    
294   G4double tmax = MaxSecondaryEnergy(dp->GetDe    
295                                                   
296   G4double maxKinEnergy = std::min(maxEnergy,t    
297   if(minKinEnergy >= maxKinEnergy) { return; }    
298                                                   
299   //G4cout << "G4AtimaEnergyLossModel::SampleS    
300   //         << " Emax= " << maxKinEnergy << G    
301                                                   
302   G4double totEnergy     = kineticEnergy + mas    
303   G4double etot2         = totEnergy*totEnergy    
304   G4double beta2         = kineticEnergy*(kine    
305                                                   
306   G4double deltaKinEnergy, f;                     
307   G4double f1 = 0.0;                              
308   G4double fmax = 1.0;                            
309   if( 0.0 < spin ) { fmax += 0.5*maxKinEnergy*    
310                                                   
311   CLHEP::HepRandomEngine* rndmEngineMod = G4Ra    
312   G4double rndm[2];                               
313                                                   
314   // sampling without nuclear size effect         
315   do {                                            
316     rndmEngineMod->flatArray(2, rndm);            
317     deltaKinEnergy = minKinEnergy*maxKinEnergy    
318                     /(minKinEnergy*(1.0 - rndm    
319                                                   
320     f = 1.0 - beta2*deltaKinEnergy/tmax;          
321     if( 0.0 < spin ) {                            
322       f1 = 0.5*deltaKinEnergy*deltaKinEnergy/e    
323       f += f1;                                    
324     }                                             
325                                                   
326     // Loop checking, 03-Aug-2015, Vladimir Iv    
327   } while( fmax*rndm[1] > f);                     
328                                                   
329   // projectile formfactor - suppresion of hig    
330   // delta-electron production at high energy     
331                                                   
332   G4double x = formfact*deltaKinEnergy*(deltaK    
333   if(x > 1.e-6) {                                 
334                                                   
335     G4double x10 = 1.0 + x;                       
336     G4double grej  = 1.0/(x10*x10);               
337     if( 0.0 < spin ) {                            
338       G4double x2 = 0.5*electron_mass_c2*delta    
339       grej *= (1.0 + magMoment2*(x2 - f1/f)/(1    
340     }                                             
341     if(grej > 1.1) {                              
342       G4cout << "### G4AtimaEnergyLossModel WA    
343              << "  " << dp->GetDefinition()->G    
344              << " Ekin(MeV)= " <<  kineticEner    
345              << " delEkin(MeV)= " << deltaKinE    
346              << G4endl;                           
347     }                                             
348     if(rndmEngineMod->flat() > grej) { return;    
349   }                                               
350                                                   
351   G4ThreeVector deltaDirection;                   
352                                                   
353   if(UseAngularGeneratorFlag()) {                 
354                                                   
355     const G4Material* mat =  couple->GetMateri    
356     G4int Z = SelectRandomAtomNumber(mat);        
357                                                   
358     deltaDirection =                              
359       GetAngularDistribution()->SampleDirectio    
360                                                   
361   } else {                                        
362                                                   
363     G4double deltaMomentum =                      
364       std::sqrt(deltaKinEnergy * (deltaKinEner    
365     G4double cost = deltaKinEnergy * (totEnerg    
366       (deltaMomentum * dp->GetTotalMomentum())    
367     if(cost > 1.0) { cost = 1.0; }                
368     G4double sint = std::sqrt((1.0 - cost)*(1.    
369                                                   
370     G4double phi = twopi*rndmEngineMod->flat()    
371                                                   
372     deltaDirection.set(sint*cos(phi),sint*sin(    
373     deltaDirection.rotateUz(dp->GetMomentumDir    
374   }                                               
375   /*                                              
376     G4cout << "### G4AtimaEnergyLossModel "       
377            << dp->GetDefinition()->GetParticle    
378            << " Ekin(MeV)= " <<  kineticEnergy    
379            << " delEkin(MeV)= " << deltaKinEne    
380            << " tmin(MeV)= " << minKinEnergy      
381            << " tmax(MeV)= " << maxKinEnergy      
382            << " dir= " << dp->GetMomentumDirec    
383            << " dirDelta= " << deltaDirection     
384            << G4endl;                             
385   */                                              
386   // create G4DynamicParticle object for delta    
387   auto delta = new G4DynamicParticle(theElectr    
388                                                   
389   vdp->push_back(delta);                          
390                                                   
391   // Change kinematics of primary particle        
392   kineticEnergy -= deltaKinEnergy;                
393   G4ThreeVector finalP = dp->GetMomentum() - d    
394   finalP               = finalP.unit();           
395                                                   
396   fParticleChange->SetProposedKineticEnergy(ki    
397   fParticleChange->SetProposedMomentumDirectio    
398 }                                                 
399                                                   
400 //....oooOO0OOooo........oooOO0OOooo........oo    
401                                                   
402 G4double G4AtimaEnergyLossModel::MaxSecondaryE    
403                                                   
404 {                                                 
405   // here particle type is checked for any met    
406   SetParticle(pd);                                
407   G4double tau  = kinEnergy/mass;                 
408   G4double tmax = 2.0*electron_mass_c2*tau*(ta    
409                   (1. + 2.0*(tau + 1.)*ratio +    
410   return std::min(tmax,tlimit);                   
411 }                                                 
412                                                   
413 //....oooOO0OOooo........oooOO0OOooo........oo    
414                                                   
415 G4double G4AtimaEnergyLossModel::StoppingPower    
416          G4double ep, G4double at, G4double zt    
417   //                                              
418   if(ep==0)return 0.0;                            
419   ap=ap/atomic_mass_unit;                         
420   ep=ep/ap;                                       
421   G4double se = 0.0;                              
422   // ep in MeV                                    
423   if(ep<=10.){                                    
424         se = sezi_dedx_e(zp, ep, at, zt);         
425   }                                               
426   else if(ep>10. && ep<30.){                      
427     G4double factor = 0.05 * ( ep - 10.0 );       
428     se = (1.0-factor)*sezi_dedx_e(zp, ep, at,     
429   }                                               
430   else {                                          
431     se = Bethek_dedx_e(ap, zp, ep, at, zt);       
432   }                                               
433   return se + dedx_n(ap, zp, ep, at, zt);         
434 }                                                 
435                                                   
436 //....oooOO0OOooo........oooOO0OOooo........oo    
437                                                   
438 G4double G4AtimaEnergyLossModel::dedx_n(const     
439          const G4double ep, const G4double at,    
440                                                   
441   G4double zpowers = g4calc->powA(zp,0.23)+g4c    
442   G4double asum = ap + at;                        
443   G4double epsilon = 32.53*at*1000.*ep*ap/(zp*    
444   G4double sn=0.;                                 
445                                                   
446   if(epsilon<=30.0){                              
447     sn = G4Log(1.+(1.1383*epsilon))/ (2.*(epsi    
448   }                                               
449   else{                                           
450     sn = G4Log(epsilon)/(2.*epsilon);             
451   }                                               
452   sn = 100.*8.4621*zp*zt*ap*sn*Avogadro/1.e+23    
453   return sn;                                      
454 }                                                 
455                                                   
456 //....oooOO0OOooo........oooOO0OOooo........oo    
457                                                   
458 G4double G4AtimaEnergyLossModel::sezi_dedx_e(c    
459   G4double e = ep*1000.; // e in keV/u            
460   G4double se = 0.;                               
461                                                   
462   // heavy ion Z > 2                              
463   G4double h1,h4;                                 
464   G4double a,q,b;                                 
465   G4double l1,l0,l;                               
466   G4double YRmin = 0.130; //  YRmin = VR / ZP*    
467   G4double VRmin = 1.0;                           
468   G4double vfermi = atima_vfermi[(G4int)zt-1];    
469   G4double yr=0;                                  
470   G4double zeta = 0;                              
471                                                   
472   G4double v = std::sqrt(e/25.0)/vfermi;          
473   G4double v2= v*v;                               
474                                                   
475   G4double vr = (v >= 1)? v*vfermi*(1.+ 1./(5.    
476                                                   
477   h1= 1./g4calc->powA(zp,0.6667);                 
478   yr = std::max(YRmin,vr*h1);                     
479   yr = std::max(yr, VRmin*h1);                    
480                                                   
481   //--  CALCULATE ZEFF                            
482   a = -0.803*g4calc->powA(yr,0.3) + 1.3167*g4c    
483   q = std::min(1.0, std::max(0.0 , (1.0 - G4Ex    
484                                                   
485   //-- IONIZATION LEVEL TO EFFECTIVE CHARGE       
486   h1 = 1./g4calc->powA(zp,0.3333);                
487                                                   
488   b = (std::min(0.43, std::max(0.32,0.12 + 0.0    
489   l0 = (0.8 - q * std::min(1.2,0.6 +zp/30.0))*    
490         if(q < 0.2){                              
491             l1 = 0;                               
492             }                                     
493         else{                                     
494             if (q < std::max(0.0,0.9-0.025*zp)    
495                 l1 = b*(q-0.2)/std::abs(std::m    
496                 }                                 
497          else{                                    
498               if(q < std::max(0.0,1.0 - 0.025*    
499                 else l1 = b*(1.0 - q)/(0.025*s    
500                 }                                 
501         }                                         
502         // calculate screening                    
503         l = std::max(l1,l0*atima_lambda_screen    
504         h1 =4.0*l*vfermi/1.919;                   
505         zeta = q + (1./(2.*(vfermi*vfermi)))*(    
506                                                   
507          // ZP**3 EFFECT AS IN REF. 779?          
508         a = 7.6 - std::max(0.0, G4Log(e));        
509         zeta = zeta*(1. + (1./(zp*zp))*(0.18 +    
510                                                   
511         h1= 1./g4calc->powA(zp,0.6667);           
512         if (yr <= ( std::max(YRmin, VRmin*h1))    
513             VRmin=std::max(VRmin, YRmin/h1);      
514             //--C        ..CALCULATE VELOCITY     
515             G4double vmin =.5*(VRmin + std::sq    
516             G4double eee = 25.0*vmin*vmin;        
517             G4double eval = 1.0;                  
518             if((zt == 6) || (((zt == 14) || (z    
519             else eval = 0.5;                      
520                                                   
521             h1 = zeta *zp;                        
522             h4 = g4calc->powA(e / eee,eval);      
523             se = sezi_p_se(eee*0.001,at,zt) *     
524             return se;                            
525         }                                         
526         else {                                    
527             se = sezi_p_se(ep,at,zt)*g4calc->p    
528             return se;                            
529         }                                         
530                                                   
531    return se;                                     
532 }                                                 
533                                                   
534 //....oooOO0OOooo........oooOO0OOooo........oo    
535                                                   
536 G4double G4AtimaEnergyLossModel::sezi_p_se(con    
537   G4double sp = 0.;                               
538   G4double e = 1000*energy; //e in keV/u          
539   G4int i = zt - 1;                               
540                                                   
541   if(e<=25)e=25;                                  
542    G4double sl = (proton_stopping_coef[i][0]*g    
543    G4double sh = proton_stopping_coef[i][4]/g4    
544    sp = sl*sh/(sl+sh);                            
545    e=1000*energy;                                 
546    if(e<=25){                                     
547        sp *=(zt>6)?g4calc->powA(e/25.,0.45):g4    
548    }                                              
549   return 100.*sp*Avogadro/1.e+23/at;              
550 }                                                 
551                                                   
552                                                   
553 //....oooOO0OOooo........oooOO0OOooo........oo    
554                                                   
555 G4double G4AtimaEnergyLossModel::Bethek_dedx_e    
556 //                                                
557 // According to C. Scheidenberger et al., Phys    
558 //                                                
559     G4double gamma=1.0 + ep/atomic_mass_unit;     
560     G4double beta2=1.0-1.0/(gamma*gamma);         
561     //                                            
562     G4double beta = std::sqrt(beta2);             
563     G4double zp_eff = zp*(1.0-G4Exp(-0.95/fine    
564     G4int zi = std::min(std::max((G4int)zt, 1)    
565     G4double Ipot = ionisation_potentials_z[zi    
566     G4double f1 = dedx_constant*g4calc->powA(z    
567                                                   
568     //                                            
569     G4double f2 = G4Log(2.0*electron_mass*1000    
570     G4double eta = beta*gamma;                    
571                                                   
572     if(!(eta>=0.13)){ //shell corrections         
573         G4double cor = (+0.422377*g4calc->powA    
574                     +0.0304043*g4calc->powA(et    
575                     -0.00038106*g4calc->powA(e    
576                   +(+3.858019*g4calc->powA(eta    
577                     -0.1667989*(g4calc->powA(e    
578                     +0.00157955*(g4calc->powA(    
579         f2 = f2 -cor/zt;                          
580     }                                             
581                                                   
582     f2+=2*G4Log(gamma) -beta2;                    
583                                                   
584     //Barkas correction                           
585     G4double barkas=1.0;                          
586     G4double V2FVA[4]={0.33,0.30,0.26,0.23};      
587     G4double VA[4]={1.,2.,3.,4.};                 
588     G4double v1 = eta/(fine_structure*std::sqr    
589     G4double v2fv;                                
590     if(v1 >= 4.){                                 
591         v2fv = 0.45 / std::sqrt(v1);              
592         }                                         
593     else if(v1 > 1. && v1 < 4.){//VALUES FROM     
594         G4int i=1;                                
595         for(; i<4; i++){                          
596     if(VA[i] >= v1) break;                        
597   }                                               
598         i = std::min(i, 3);                       
599   v2fv = V2FVA[i-1]+(v1-VA[i-1])*(V2FVA[i]-V2F    
600     }                                             
601     else{                                         
602         v2fv=0;                                   
603     }                                             
604     barkas= 1.0+2.0 * zp_eff * v2fv /(v1*v1*st    
605                                                   
606     //Fermi-density effect for relativistic ve    
607     gamma = 1./std::sqrt(1-(beta*beta));          
608     G4double x = G4Log(beta * gamma) / 2.30258    
609     G4int i = std::min(std::max(zi-1,0), 91);     
610     G4double del = 0.;                            
611                                                   
612     if (x < x0[i] ){                              
613         if(del_0[i] > 0.)del = del_0[i] * g4ca    
614         }                                         
615     else {                                        
616         del = 4.6052 * x - c[i];                  
617         if ( x0[i]<= x &&  x <= x1[i] ) del +=    
618         }                                         
619                                                   
620     //Precalculated lindhard correction           
621     G4double LS = 0.;                             
622     G4int z = G4lrint(zp);                        
623     if(z>109)z=109;                               
624     if(ep<tableE[0])ep=tableE[0];                 
625                                                   
626     G4double da = (ap - element_atomic_weights    
627     z = z-1;                                      
628                                                   
629     G4double v3 = EnergyTable_interpolate(ep,l    
630     G4double v4 = EnergyTable_interpolate(ep,l    
631     G4double dif = v4 - v3;                       
632     LS = v3+(dif*da/0.05);                        
633                                                   
634     //Final stopping power                        
635     G4double result  = (f2)*barkas + LS - del/    
636     result *=f1;                                  
637                                                   
638     return result;                                
639 }                                                 
640                                                   
641 //....oooOO0OOooo........oooOO0OOooo........oo    
642                                                   
643 G4double G4AtimaEnergyLossModel::EnergyTable_i    
644      G4double r;                                  
645      G4int num=200;                               
646      G4double lxval = G4Log(xval)/MLN10;          
647      if(xval<tableE[0] || xval>tableE[num-1])r    
648      if(xval==tableE[num-1])return y[num-1];      
649      G4int i = (G4int)(lxval/stepE);              
650      i = std::min(std::max(i, 0), num-2);         
651      G4double linstep = tableE[i+1] - tableE[i    
652      G4double x = 1.0 - ((xval - tableE[i])/li    
653      r = (x*y[i]) + ((1-x)*y[i+1]);               
654      return r;                                    
655 }                                                 
656                                                   
657 //....oooOO0OOooo........oooOO0OOooo........oo    
658                                                   
659 const G4double G4AtimaEnergyLossModel::proton_    
660 {  .0091827, .0053496, .69741, .48493, 316.07,    
661 {  .11393, .0051984, 1.0822, .39252,  1081.0,     
662 {  .85837, .0050147, 1.6044, .38844,  1337.3,     
663 {  .8781,  .0051049, 5.4232, .2032,   1200.6,     
664 { 1.4608,  .0048836, 2.338,  .44249, 1801.3, 1    
665 { 3.2579,  .0049148, 2.7156, .36473, 2092.2, 1    
666 {  .59674, .0050837, 4.2073, .30612,  2394.2,     
667 {  .75253, .0050314, 4.0824, .30067, 2455.8, 1    
668 { 1.226,   .0051385, 3.2246, .32703, 2525.0, 1    
669 { 1.0332,  .0051645, 3.004,  .33889, 2338.6, .    
670 { 6.0972,  .0044292, 3.1929, .45763, 1363.3, .    
671 {14.013,   .0043646, 2.2641, .36326, 2187.4, .    
672 { .039001, .0045415, 5.5463, .39562, 1589.2, .    
673 { 2.072,   .0044516, 3.5585, .53933, 1515.2, .    
674 {17.575,   .0038346, .078694, 1.2388,2806.0, .    
675 {16.126,   .0038315, .054164, 1.3104,2813.3, .    
676 { 3.217,   .0044579, 3.6696, .5091,  2734.6, .    
677 { 2.0379,  .0044775, 3.0743, .54773, 3505.0, .    
678 { .74171, .0043051, 1.1515, .95083,   917.21,     
679 { 9.1316,  .0043809, 5.4611, .31327, 3891.8, .    
680 { 7.2247,  .0043718, 6.1017, .37511, 2829.2, .    
681 {  .147,   .0048456, 6.3485, .41057, 2164.1, .    
682 { 5.0611,  .0039867, 2.6174, .57957, 2218.9, .    
683 {  .53267, .0042968, .39005, 1.2725, 1872.7, .    
684 {  .47697, .0043038, .31452, 1.3289, 1920.5, .    
685 {  .027426, .0035443, .031563, 2.1755,  1919.5    
686 {  .16383, .0043042, .073454, 1.8592, 1918.4,     
687 { 4.2562,  .0043737, 1.5606, .72067, 1546.8, .    
688 { 2.3508,  .0043237, 2.882,  .50113, 1837.7, .    
689 { 3.1095,  .0038455, .11477, 1.5037, 2184.7, .    
690 {15.322,   .0040306, .65391, .67668, 3001.7, .    
691 { 3.6932,  .0044813, 8.608,  .27638,  2982.7,     
692 { 7.1373,  .0043134, 9.4247, .27937,  2725.8,     
693 { 4.8979,  .0042937, 3.7793, .50004,  2824.5,     
694 { 1.3683,  .0043024, 2.5679, .60822,  6907.8,     
695 { 1.8301,  .0042983, 2.9057, .6038,   4744.6,     
696 {  .42056, .0041169, .01695, 2.3616, 2252.7, .    
697 {30.78,    .0037736, .55813, .76816, 7113.2, .    
698 {11.576,   .0042119, 7.0244, .37764, 4713.5, .    
699 { 6.2406,  .0041916, 5.2701, .49453, 4234.6, .    
700 {  .33073, .0041243, 1.7246, 1.1062, 1930.2, .    
701 {  .017747, .0041715, .14586, 1.7305,1803.6, .    
702 { 3.7229,  .0041768, 4.6286, .56769, 1678.0, .    
703 {  .13998, .0041329, .25573, 1.4241, 1919.3, .    
704 {  .2859,  .0041386, .31301, 1.3424, 1954.8, .    
705 { .76,    .0042179, 3.386,  .76285,  1867.4, .    
706 { 6.3957,  .0041935, 5.4689, .41378, 1712.6, .    
707 { 3.4717,  .0041344, 3.2337, .63788, 1116.4, .    
708 { 2.5265,  .0042282, 4.532,  .53562, 1030.8, .    
709 { 7.3683,  .0041007, 4.6791, .51428, 1160.0, .    
710 { 7.7197,  .004388, 3.242,   .68434, 1428.1, .    
711 {16.78,    .0041918, 9.3198, .29568, 3370.9, .    
712 { 4.2132,  .0042098, 4.6753, .57945, 3503.9, .    
713 { 4.0818,  .004214, 4.4425,  .58393, 3945.3, .    
714 {.18517, .0036215, .00058788,3.5315, 2931.3, .    
715 { 4.8248,  .0041458, 6.0934, .57026, 2300.1, .    
716 {  .49857, .0041054, 1.9775, .95877, 786.55, .    
717 { 3.2754,  .0042177, 5.768,  .54054, 6631.3, .    
718 { 2.9978,  .0040901, 4.5299, .62025, 2161.2, .    
719 { 2.8701,  .004096, 4.2568,  .6138,  2130.4, .    
720 {10.853,   .0041149, 5.8907, .46834, 2857.2, .    
721 { 3.6407,  .0041782, 4.8742, .57861, 1267.7, .    
722 {17.645,   .0040992, 6.5855, .32734, 3931.3, .    
723 { 7.5309,  .0040814, 4.9389, .50679, 2519.7, .    
724 { 5.4742,  .0040829, 4.897,  .51113, 2340.1, .    
725 { 4.2661,  .0040667, 4.5032, .55257, 2076.4, .    
726 { 6.8313,  .0040486, 4.3987, .51675, 2003.0, .    
727 { 1.2707,  .0040553, 4.6295, .57428, 1626.3, .    
728 { 5.7561,  .0040491, 4.357,  .52496, 2207.3, .    
729 {14.127,   .0040596, 5.8304, .37755, 3645.9, .    
730 { 6.6948,  .0040603, 4.9361, .47961, 2719.0, .    
731 { 3.0619,  .0040511, 3.5803, .59082, 2346.1, .    
732 {10.811,   .0033008, 1.3767, .76512, 2003.7, .    
733 { 2.7101,  .0040961, 1.2289, .98598, 1232.4, .    
734 {  .52345, .0040244, 1.4038, .8551,  1461.4, .    
735 {  .4616,  .0040203, 1.3014, .87043, 1473.5, .    
736 {  .97814, .0040374, 2.0127, .7225,  1890.8, .    
737 { 3.2086,  .004051, 3.6658,  .53618, 3091.2, .    
738 { 2.0035,  .0040431, 7.4882, .3561,  4464.3, .    
739 {15.43,    .0039432, 1.1237, .70703, 4595.7, .    
740 { 3.1512,  .0040524, 4.0996, .5425,  3246.3, .    
741 { 7.1896,  .0040588, 8.6927, .35842, 4760.6, .    
742 { 9.3209,  .004054, 11.543,  .32027, 4866.2, .    
743 {29.242,   .0036195, .16864, 1.1226, 5688.0, .    
744 { 1.8522,  .0039973, 3.1556, .65096, 3755.0, .    
745 { 3.222,   .0040041, 5.9024, .52678, 4040.2, .    
746 { 9.3412,  .0039661, 7.921,  .42977, 5180.9, .    
747 {36.183,   .0036003, .58341, .86747, 6990.2, .    
748 { 5.9284,  .0039695, 6.4082, .52122, 4619.5, .    
749 { 5.2454,  .0039744, 6.7969, .48542, 4586.3, .    
750 {33.702,   .0036901, .47257, .89235, 5295.7, .    
751 {2.7589,  .0039806, 3.2092, .66122,  2505.4, .    
752 };                                                
753                                                   
754 //....oooOO0OOooo........oooOO0OOooo........oo    
755                                                   
756 const G4double G4AtimaEnergyLossModel::atima_v    
757 1.0309,                                           
758 0.15976,                                          
759 0.59782,                                          
760 1.0781,                                           
761 1.0486,                                           
762 1.00,                                             
763 1.058,                                            
764 0.93942,                                          
765 0.74562,                                          
766 0.3424,                                           
767 0.45259,                                          
768 0.71074,                                          
769 0.90519,                                          
770 0.97411,                                          
771 0.97184,                                          
772 0.89852,                                          
773 0.70827,                                          
774 0.39816,                                          
775 0.36552,                                          
776 0.62712,                                          
777 0.81707,                                          
778 0.9943,                                           
779 1.1423,                                           
780 1.2381,                                           
781 1.1222,                                           
782 0.92705,                                          
783 1.0047,                                           
784 1.2,                                              
785 1.0661,                                           
786 0.97411,                                          
787 0.84912,                                          
788 0.95,                                             
789 1.0903,                                           
790 1.0429,                                           
791 0.49715,                                          
792 0.37755,                                          
793 0.35211,                                          
794 0.57801,                                          
795 0.77773,                                          
796 1.0207,                                           
797 1.029,                                            
798 1.2542,                                           
799 1.122,                                            
800 1.1241,                                           
801 1.0882,                                           
802 1.2709,                                           
803 1.2542,                                           
804 0.90094,                                          
805 0.74093,                                          
806 0.86054,                                          
807 0.93155,                                          
808 1.0047,                                           
809 0.55379,                                          
810 0.43289,                                          
811 0.32636,                                          
812 0.5131,                                           
813 0.6950,                                           
814 0.72591,                                          
815 0.71202,                                          
816 0.67413,                                          
817 0.71418,                                          
818 0.71453,                                          
819 0.5911,                                           
820 0.70263,                                          
821 0.68049,                                          
822 0.68203,                                          
823 0.68121,                                          
824 0.68532,                                          
825 0.68715,                                          
826 0.61884,                                          
827 0.71801,                                          
828 0.83048,                                          
829 1.1222,                                           
830 1.2381,                                           
831 1.045,                                            
832 1.0733,                                           
833 1.0953,                                           
834 1.2381,                                           
835 1.2879,                                           
836 0.78654,                                          
837 0.66401,                                          
838 0.84912,                                          
839 0.88433,                                          
840 0.80746,                                          
841 0.43357,                                          
842 0.41923,                                          
843 0.43638,                                          
844 0.51464,                                          
845 0.73087,                                          
846 0.81065,                                          
847 1.9578,                                           
848 1.0257                                            
849 };                                                
850                                                   
851 //....oooOO0OOooo........oooOO0OOooo........oo    
852                                                   
853 const G4double G4AtimaEnergyLossModel::atima_l    
854  1.00,                                            
855  1.00,                                            
856  1.10,                                            
857  1.06,                                            
858  1.01,                                            
859  1.03,                                            
860  1.04,                                            
861  0.99,                                            
862  0.95,                                            
863  0.90,                                            
864  0.82,                                            
865  0.81,                                            
866  0.83,                                            
867  0.88,                                            
868  1.00,                                            
869  0.95,                                            
870  0.97,                                            
871  0.99,                                            
872  0.98,                                            
873  0.97,                                            
874  0.98,                                            
875  0.97,                                            
876  0.96,                                            
877  0.93,                                            
878  0.91,                                            
879  0.9,                                             
880  0.88,                                            
881  0.9,                                             
882  0.9,                                             
883  0.9,                                             
884  0.9,                                             
885  0.85,                                            
886  0.9,                                             
887  0.9,                                             
888  0.91,                                            
889  0.92,                                            
890  0.9,                                             
891  0.9,                                             
892  0.9,                                             
893  0.9,                                             
894  0.9,                                             
895  0.88,                                            
896  0.9,                                             
897  0.88,                                            
898  0.88,                                            
899  0.9,                                             
900  0.9,                                             
901  0.88,                                            
902  0.9,                                             
903  0.9,                                             
904  0.9,                                             
905  0.9,                                             
906  0.96,                                            
907  1.2,                                             
908  0.9,                                             
909  0.88,                                            
910  0.88,                                            
911  0.85,                                            
912  0.90,                                            
913  0.90,                                            
914  0.92,                                            
915  0.95,                                            
916  0.99,                                            
917  1.03,                                            
918  1.05,                                            
919  1.07,                                            
920  1.08,                                            
921  1.10,                                            
922  1.08,                                            
923  1.08,                                            
924  1.08,                                            
925  1.08,                                            
926  1.09,                                            
927  1.09,                                            
928  1.10,                                            
929  1.11,                                            
930  1.12,                                            
931  1.13,                                            
932  1.14,                                            
933  1.15,                                            
934  1.17,                                            
935  1.20,                                            
936  1.18,                                            
937  1.17,                                            
938  1.17,                                            
939  1.16,                                            
940  1.16,                                            
941  1.16,                                            
942  1.16,                                            
943  1.16,                                            
944  1.16,                                            
945  1.16                                             
946 };                                                
947                                                   
948 //....oooOO0OOooo........oooOO0OOooo........oo    
949                                                   
950 const G4double G4AtimaEnergyLossModel::element    
951  0.0,                                             
952  1.00794, //H                                     
953  4.0026, //He                                     
954  6.941, //Li                                      
955  9.01218, //Be                                    
956  10.811, //B                                      
957  12.0107, //C                                     
958  14.0067, //N                                     
959  15.9994, //O                                     
960  18.9984, //F                                     
961  20.1797, //Ne                                    
962  22.9898, //Na                                    
963  24.305, //Mg                                     
964  26.9815, //Al                                    
965  28.0855, //Si                                    
966  30.9738, //P                                     
967  32.065, //S                                      
968  35.453, //Cl                                     
969  39.948, //Ar                                     
970  39.0983, //K                                     
971  40.078, //Ca                                     
972  44.9559, //Sc                                    
973  47.867, //Ti                                     
974  50.9415, //V                                     
975  51.9961, //Cr                                    
976  54.938, //Mn                                     
977  55.845, //Fe                                     
978  58.9332, //Co                                    
979  58.6934, //Ni                                    
980  63.546, //Cu                                     
981  65.409, //Zn                                     
982  69.723, //Ga                                     
983  72.64, //Ge                                      
984  74.9216, //As                                    
985  78.96, //Se                                      
986  79.904, //Br                                     
987  83.798, //Kr                                     
988  85.4678, //Rb                                    
989  87.62, //Sr                                      
990  88.9059, //Y                                     
991  91.224, //Zr                                     
992  92.9064, //Nb                                    
993  95.94, //Mo                                      
994  97.9072, //Tc                                    
995  101.07, //Ru                                     
996  102.906, //Rh                                    
997  106.42, //Pd                                     
998  107.868, //Ag                                    
999  112.411, //Cd                                    
1000  114.818, //In                                   
1001  118.71, //Sn                                    
1002  121.76, //Sb                                    
1003  127.6, //Te                                     
1004  126.904, //I                                    
1005  131.293, //Xe                                   
1006  132.905, //Cs                                   
1007  137.327, //Ba                                   
1008  138.905, //La                                   
1009  140.116, //Ce                                   
1010  140.908, //Pr                                   
1011  144.24, //Nd                                    
1012  144.913, //Pm                                   
1013  150.36, //Sm                                    
1014  151.964, //Eu                                   
1015  157.25, //Gd                                    
1016  158.925, //Tb                                   
1017  162.5, //Dy                                     
1018  164.93, //Ho                                    
1019  167.259, //Er                                   
1020  168.934, //Tm                                   
1021  173.04, //Yb                                    
1022  174.967, //Lu                                   
1023  178.49, //Hf                                    
1024  180.948, //Ta                                   
1025  183.84, //W                                     
1026  186.207, //Re                                   
1027  190.23, //Os                                    
1028  192.217, //Ir                                   
1029  195.078, //Pt                                   
1030  196.967, //Au                                   
1031  200.59, //Hg                                    
1032  204.383, //Tl                                   
1033  207.2, //Pb                                     
1034  208.98, //Bi                                    
1035  208.982, //Po                                   
1036  209.987, //At                                   
1037  222.018, //Rn                                   
1038  223.02, //Fr                                    
1039  226.025, //Ra                                   
1040  227.028, //Ac                                   
1041  232.038, //Th                                   
1042  231.036, //Pa                                   
1043  238.029, //U                                    
1044  237.048, //Np                                   
1045  244.064, //Pu                                   
1046  243.061, //Am                                   
1047  247.07, //Cm                                    
1048  247.07, //Bk                                    
1049  251.08, //Cf                                    
1050  252.083, //Es                                   
1051  257.095, //Fm                                   
1052  258.098, //Md                                   
1053  259.101, //No                                   
1054  262.11, //Lr                                    
1055  261.109, //Rf                                   
1056  262.114, //Db                                   
1057  266.122, //Sg                                   
1058  264.125, //Bh                                   
1059  269.134, //Hs                                   
1060  268.139 //Mt                                    
1061 };                                               
1062                                                  
1063 //....oooOO0OOooo........oooOO0OOooo........o    
1064                                                  
1065 //LS coefficient for A=atomic weight             
1066 const G4double G4AtimaEnergyLossModel::ls_coe    
1067 { {-0.0286005,-0.0269762,-0.0254384,-0.023982    
1068   {-0.108683,-0.102922,-0.097437,-0.0922162,-    
1069   {-0.223884,-0.213009,-0.202568,-0.192553,-0    
1070   {-0.356328,-0.340686,-0.325539,-0.310886,-0    
1071   {-0.492346,-0.472912,-0.453947,-0.435458,-0    
1072   {-0.623826,-0.601642,-0.579854,-0.558476,-0    
1073   {-0.747011,-0.722922,-0.699148,-0.675705,-0    
1074   {-0.860779,-0.835389,-0.810244,-0.785358,-0    
1075   {-0.965337,-0.939056,-0.912961,-0.887067,-0    
1076   {-1.06145,-1.03454,-1.00778,-0.981172,-0.95    
1077   {-1.15003,-1.12268,-1.09543,-1.06831,-1.041    
1078   {-1.23199,-1.2043,-1.1767,-1.14919,-1.12179    
1079   {-1.30813,-1.28018,-1.25231,-1.22451,-1.196    
1080   {-1.37915,-1.351,-1.32292,-1.29489,-1.26693    
1081   {-1.44565,-1.41735,-1.38909,-1.36088,-1.332    
1082   {-1.50815,-1.47972,-1.45132,-1.42297,-1.394    
1083   {-1.56708,-1.53854,-1.51003,-1.48156,-1.453    
1084   {-1.62282,-1.59419,-1.56558,-1.537,-1.50845    
1085   {-1.67568,-1.64697,-1.61829,-1.58962,-1.560    
1086   {-1.72594,-1.69717,-1.66841,-1.63967,-1.610    
1087   {-1.77384,-1.74501,-1.71619,-1.68738,-1.658    
1088   {-1.81958,-1.7907,-1.76183,-1.73296,-1.7041    
1089   {-1.86335,-1.83442,-1.80551,-1.7766,-1.7476    
1090   {-1.90531,-1.87635,-1.84739,-1.81843,-1.789    
1091   {-1.9456,-1.9166,-1.88761,-1.85861,-1.82962    
1092   {-1.98434,-1.95532,-1.92629,-1.89726,-1.868    
1093   {-2.02166,-1.99261,-1.96355,-1.93449,-1.905    
1094   {-2.05765,-2.02857,-1.99949,-1.9704,-1.9413    
1095   {-2.0924,-2.0633,-2.03419,-2.00508,-1.97596    
1096   {-2.12599,-2.09687,-2.06774,-2.03861,-2.009    
1097   {-2.1585,-2.12936,-2.10021,-2.07106,-2.0418    
1098   {-2.18999,-2.16084,-2.13167,-2.1025,-2.0733    
1099   {-2.22053,-2.19136,-2.16218,-2.13299,-2.103    
1100   {-2.25017,-2.22099,-2.19179,-2.16258,-2.133    
1101   {-2.27896,-2.24977,-2.22056,-2.19133,-2.162    
1102   {-2.30695,-2.27775,-2.24852,-2.21929,-2.190    
1103   {-2.33419,-2.30497,-2.27573,-2.24649,-2.217    
1104   {-2.36071,-2.33148,-2.30223,-2.27297,-2.243    
1105   {-2.38654,-2.3573,-2.32805,-2.29877,-2.2694    
1106   {-2.41173,-2.38248,-2.35322,-2.32393,-2.294    
1107   {-2.4363,-2.40705,-2.37777,-2.34848,-2.3191    
1108   {-2.46029,-2.43103,-2.40174,-2.37244,-2.343    
1109   {-2.48372,-2.45445,-2.42516,-2.39585,-2.366    
1110   {-2.50661,-2.47733,-2.44804,-2.41872,-2.389    
1111   {-2.529,-2.49971,-2.47041,-2.44109,-2.41175    
1112   {-2.55089,-2.5216,-2.4923,-2.46297,-2.43362    
1113   {-2.57232,-2.54303,-2.51371,-2.48438,-2.455    
1114   {-2.59331,-2.56401,-2.53469,-2.50535,-2.475    
1115   {-2.61386,-2.58456,-2.55523,-2.52589,-2.496    
1116   {-2.634,-2.60469,-2.57536,-2.54601,-2.51664    
1117   {-2.65375,-2.62444,-2.5951,-2.56575,-2.5363    
1118   {-2.67311,-2.6438,-2.61446,-2.5851,-2.55572    
1119   {-2.69211,-2.66279,-2.63345,-2.60409,-2.574    
1120   {-2.71076,-2.68144,-2.65209,-2.62273,-2.593    
1121   {-2.72907,-2.69974,-2.6704,-2.64102,-2.6116    
1122   {-2.74705,-2.71772,-2.68837,-2.659,-2.6296,    
1123   {-2.76471,-2.73538,-2.70603,-2.67665,-2.647    
1124   {-2.78207,-2.75274,-2.72338,-2.694,-2.6646,    
1125   {-2.79914,-2.7698,-2.74044,-2.71106,-2.6816    
1126   {-2.81592,-2.78658,-2.75722,-2.72783,-2.698    
1127   {-2.83242,-2.80308,-2.77372,-2.74433,-2.714    
1128   {-2.84866,-2.81932,-2.78996,-2.76057,-2.731    
1129   {-2.86464,-2.8353,-2.80593,-2.77654,-2.7471    
1130   {-2.88037,-2.85103,-2.82166,-2.79227,-2.762    
1131   {-2.89586,-2.86652,-2.83715,-2.80775,-2.778    
1132   {-2.91112,-2.88177,-2.8524,-2.823,-2.79358,    
1133   {-2.92614,-2.8968,-2.86742,-2.83803,-2.8086    
1134   {-2.94095,-2.9116,-2.88223,-2.85283,-2.8234    
1135   {-2.95554,-2.92619,-2.89682,-2.86742,-2.837    
1136   {-2.96992,-2.94057,-2.9112,-2.8818,-2.85237    
1137   {-2.98411,-2.95475,-2.92538,-2.89598,-2.866    
1138   {-2.99809,-2.96874,-2.93936,-2.90996,-2.880    
1139   {-3.01188,-2.98253,-2.95316,-2.92375,-2.894    
1140   {-3.02549,-2.99614,-2.96676,-2.93736,-2.907    
1141   {-3.03892,-3.00956,-2.98019,-2.95079,-2.921    
1142   {-3.05217,-3.02281,-2.99344,-2.96404,-2.934    
1143   {-3.06524,-3.03589,-3.00652,-2.97711,-2.947    
1144   {-3.07816,-3.0488,-3.01943,-2.99003,-2.9606    
1145   {-3.0909,-3.06155,-3.03218,-3.00277,-2.9733    
1146   {-3.10349,-3.07414,-3.04477,-3.01537,-2.985    
1147   {-3.11593,-3.08658,-3.0572,-3.0278,-2.99837    
1148   {-3.12821,-3.09886,-3.06949,-3.04009,-3.010    
1149   {-3.14035,-3.111,-3.08163,-3.05223,-3.0228,    
1150   {-3.15234,-3.12299,-3.09362,-3.06422,-3.034    
1151   {-3.16419,-3.13485,-3.10547,-3.07608,-3.046    
1152   {-3.17591,-3.14656,-3.11719,-3.0878,-3.0583    
1153   {-3.18749,-3.15815,-3.12878,-3.09938,-3.069    
1154   {-3.19894,-3.1696,-3.14023,-3.11084,-3.0814    
1155   {-3.21027,-3.18092,-3.15156,-3.12216,-3.092    
1156   {-3.22147,-3.19212,-3.16276,-3.13337,-3.103    
1157   {-3.23254,-3.2032,-3.17384,-3.14445,-3.1150    
1158   {-3.2435,-3.21416,-3.1848,-3.15541,-3.12599    
1159   {-3.25434,-3.225,-3.19564,-3.16626,-3.13684    
1160   {-3.26507,-3.23573,-3.20637,-3.17699,-3.147    
1161   {-3.27568,-3.24635,-3.21699,-3.18761,-3.158    
1162   {-3.28619,-3.25686,-3.2275,-3.19812,-3.1687    
1163   {-3.29659,-3.26726,-3.2379,-3.20853,-3.1791    
1164   {-3.30688,-3.27755,-3.2482,-3.21883,-3.1894    
1165   {-3.31707,-3.28775,-3.2584,-3.22902,-3.1996    
1166   {-3.32716,-3.29784,-3.26849,-3.23912,-3.209    
1167   {-3.33716,-3.30784,-3.27849,-3.24912,-3.219    
1168   {-3.34705,-3.31773,-3.28839,-3.25903,-3.229    
1169   {-3.35685,-3.32754,-3.2982,-3.26884,-3.2394    
1170   {-3.36656,-3.33725,-3.30791,-3.27855,-3.249    
1171   {-3.37618,-3.34687,-3.31754,-3.28818,-3.258    
1172   {-3.38571,-3.3564,-3.32707,-3.29772,-3.2683    
1173   {-3.39515,-3.36584,-3.33652,-3.30717,-3.277    
1174   {-3.4045,-3.3752,-3.34588,-3.31653,-3.28716    
1175 };                                               
1176                                                  
1177 //....oooOO0OOooo........oooOO0OOooo........o    
1178                                                  
1179 //LS coefficient for A=atomic weight * 1.05      
1180 const G4double G4AtimaEnergyLossModel::ls_coe    
1181 {                                                
1182   {-0.0286005,-0.0269762,-0.0254384,-0.023982    
1183   {-0.108683,-0.102922,-0.097437,-0.0922162,-    
1184   {-0.223884,-0.213009,-0.202568,-0.192553,-0    
1185   {-0.356328,-0.340686,-0.325539,-0.310886,-0    
1186   {-0.492346,-0.472912,-0.453947,-0.435458,-0    
1187   {-0.623826,-0.601642,-0.579854,-0.558476,-0    
1188   {-0.747011,-0.722922,-0.699148,-0.675705,-0    
1189   {-0.860779,-0.835389,-0.810244,-0.785358,-0    
1190   {-0.965337,-0.939056,-0.912961,-0.887067,-0    
1191   {-1.06145,-1.03454,-1.00778,-0.981172,-0.95    
1192   {-1.15003,-1.12268,-1.09543,-1.06831,-1.041    
1193   {-1.23199,-1.2043,-1.1767,-1.14919,-1.12179    
1194   {-1.30813,-1.28018,-1.25231,-1.22451,-1.196    
1195   {-1.37915,-1.351,-1.32292,-1.29489,-1.26693    
1196   {-1.44565,-1.41735,-1.38909,-1.36088,-1.332    
1197   {-1.50815,-1.47972,-1.45132,-1.42297,-1.394    
1198   {-1.56708,-1.53854,-1.51003,-1.48156,-1.453    
1199   {-1.62282,-1.59419,-1.56558,-1.537,-1.50845    
1200   {-1.67568,-1.64697,-1.61829,-1.58962,-1.560    
1201   {-1.72594,-1.69717,-1.66841,-1.63967,-1.610    
1202   {-1.77384,-1.74501,-1.71619,-1.68738,-1.658    
1203   {-1.81958,-1.7907,-1.76183,-1.73296,-1.7041    
1204   {-1.86335,-1.83442,-1.80551,-1.7766,-1.7476    
1205   {-1.90531,-1.87635,-1.84739,-1.81843,-1.789    
1206   {-1.9456,-1.9166,-1.88761,-1.85861,-1.82962    
1207   {-1.98434,-1.95532,-1.92629,-1.89726,-1.868    
1208   {-2.02166,-1.99261,-1.96355,-1.93449,-1.905    
1209   {-2.05765,-2.02857,-1.99949,-1.9704,-1.9413    
1210   {-2.0924,-2.0633,-2.03419,-2.00508,-1.97596    
1211   {-2.12599,-2.09687,-2.06774,-2.03861,-2.009    
1212   {-2.1585,-2.12936,-2.10021,-2.07106,-2.0418    
1213   {-2.18999,-2.16084,-2.13167,-2.1025,-2.0733    
1214   {-2.22053,-2.19136,-2.16218,-2.13299,-2.103    
1215   {-2.25017,-2.22099,-2.19179,-2.16258,-2.133    
1216   {-2.27896,-2.24977,-2.22056,-2.19133,-2.162    
1217   {-2.30695,-2.27775,-2.24852,-2.21929,-2.190    
1218   {-2.33419,-2.30497,-2.27573,-2.24649,-2.217    
1219   {-2.36071,-2.33148,-2.30223,-2.27297,-2.243    
1220   {-2.38654,-2.3573,-2.32805,-2.29877,-2.2694    
1221   {-2.41173,-2.38248,-2.35322,-2.32393,-2.294    
1222   {-2.4363,-2.40705,-2.37777,-2.34848,-2.3191    
1223   {-2.46029,-2.43103,-2.40174,-2.37244,-2.343    
1224   {-2.48372,-2.45445,-2.42516,-2.39585,-2.366    
1225   {-2.50661,-2.47733,-2.44804,-2.41872,-2.389    
1226   {-2.529,-2.49971,-2.47041,-2.44109,-2.41175    
1227   {-2.55089,-2.5216,-2.4923,-2.46297,-2.43362    
1228   {-2.57232,-2.54303,-2.51371,-2.48438,-2.455    
1229   {-2.59331,-2.56401,-2.53469,-2.50535,-2.475    
1230   {-2.61386,-2.58456,-2.55523,-2.52589,-2.496    
1231   {-2.634,-2.60469,-2.57536,-2.54601,-2.51664    
1232   {-2.65375,-2.62444,-2.5951,-2.56575,-2.5363    
1233   {-2.67311,-2.6438,-2.61446,-2.5851,-2.55572    
1234   {-2.69211,-2.66279,-2.63345,-2.60409,-2.574    
1235   {-2.71076,-2.68144,-2.65209,-2.62273,-2.593    
1236   {-2.72907,-2.69974,-2.6704,-2.64102,-2.6116    
1237   {-2.74705,-2.71772,-2.68837,-2.659,-2.6296,    
1238   {-2.76471,-2.73538,-2.70603,-2.67665,-2.647    
1239   {-2.78207,-2.75274,-2.72338,-2.694,-2.6646,    
1240   {-2.79914,-2.7698,-2.74044,-2.71106,-2.6816    
1241   {-2.81592,-2.78658,-2.75722,-2.72783,-2.698    
1242   {-2.83242,-2.80308,-2.77372,-2.74433,-2.714    
1243   {-2.84866,-2.81932,-2.78996,-2.76057,-2.731    
1244   {-2.86464,-2.8353,-2.80593,-2.77654,-2.7471    
1245   {-2.88037,-2.85103,-2.82166,-2.79227,-2.762    
1246   {-2.89586,-2.86652,-2.83715,-2.80775,-2.778    
1247   {-2.91112,-2.88177,-2.8524,-2.823,-2.79358,    
1248   {-2.92614,-2.8968,-2.86742,-2.83803,-2.8086    
1249   {-2.94095,-2.9116,-2.88223,-2.85283,-2.8234    
1250   {-2.95554,-2.92619,-2.89682,-2.86742,-2.837    
1251   {-2.96992,-2.94057,-2.9112,-2.8818,-2.85237    
1252   {-2.98411,-2.95475,-2.92538,-2.89598,-2.866    
1253   {-2.99809,-2.96874,-2.93936,-2.90996,-2.880    
1254   {-3.01188,-2.98253,-2.95316,-2.92375,-2.894    
1255   {-3.02549,-2.99614,-2.96676,-2.93736,-2.907    
1256   {-3.03892,-3.00956,-2.98019,-2.95079,-2.921    
1257   {-3.05217,-3.02281,-2.99344,-2.96404,-2.934    
1258   {-3.06524,-3.03589,-3.00652,-2.97711,-2.947    
1259   {-3.07816,-3.0488,-3.01943,-2.99003,-2.9606    
1260   {-3.0909,-3.06155,-3.03218,-3.00277,-2.9733    
1261   {-3.10349,-3.07414,-3.04477,-3.01536,-2.985    
1262   {-3.11593,-3.08658,-3.0572,-3.0278,-2.99837    
1263   {-3.12821,-3.09886,-3.06949,-3.04009,-3.010    
1264   {-3.14035,-3.111,-3.08163,-3.05223,-3.0228,    
1265   {-3.15234,-3.12299,-3.09362,-3.06422,-3.034    
1266   {-3.16419,-3.13485,-3.10547,-3.07608,-3.046    
1267   {-3.17591,-3.14656,-3.11719,-3.08779,-3.058    
1268   {-3.18749,-3.15815,-3.12878,-3.09938,-3.069    
1269   {-3.19894,-3.1696,-3.14023,-3.11084,-3.0814    
1270   {-3.21027,-3.18092,-3.15156,-3.12216,-3.092    
1271   {-3.22147,-3.19212,-3.16276,-3.13337,-3.103    
1272   {-3.23254,-3.2032,-3.17384,-3.14445,-3.1150    
1273   {-3.2435,-3.21416,-3.1848,-3.15541,-3.12599    
1274   {-3.25434,-3.225,-3.19564,-3.16626,-3.13684    
1275   {-3.26507,-3.23573,-3.20637,-3.17699,-3.147    
1276   {-3.27568,-3.24635,-3.21699,-3.18761,-3.158    
1277   {-3.28619,-3.25686,-3.2275,-3.19812,-3.1687    
1278   {-3.29659,-3.26726,-3.2379,-3.20853,-3.1791    
1279   {-3.30688,-3.27755,-3.2482,-3.21883,-3.1894    
1280   {-3.31707,-3.28775,-3.2584,-3.22902,-3.1996    
1281   {-3.32716,-3.29784,-3.26849,-3.23912,-3.209    
1282   {-3.33716,-3.30784,-3.27849,-3.24912,-3.219    
1283   {-3.34705,-3.31773,-3.28839,-3.25903,-3.229    
1284   {-3.35685,-3.32754,-3.2982,-3.26884,-3.2394    
1285   {-3.36656,-3.33725,-3.30791,-3.27855,-3.249    
1286   {-3.37618,-3.34687,-3.31754,-3.28818,-3.258    
1287   {-3.38571,-3.3564,-3.32707,-3.29772,-3.2683    
1288   {-3.39515,-3.36584,-3.33652,-3.30717,-3.277    
1289   {-3.4045,-3.3752,-3.34588,-3.31653,-3.28716    
1290 };                                               
1291                                                  
1292 //....oooOO0OOooo........oooOO0OOooo........o    
1293                                                  
1294 const G4double G4AtimaEnergyLossModel::ionisa    
1295     9999999,  // 0                               
1296     19.2, // H                                   
1297     41.8, // He                                  
1298     40.0,                                        
1299     63.7,                                        
1300     76.0,                                        
1301     78.0, // C                                   
1302     82.0,                                        
1303     95.0,                                        
1304     115.0,                                       
1305     137.0,                                       
1306     149.0,                                       
1307     156.0,                                       
1308     166.0,                                       
1309     173.0,                                       
1310     173.0,                                       
1311     180.0,                                       
1312     174.0,                                       
1313     188.0,                                       
1314     190.0,                                       
1315     191.0, //Ca (20)                             
1316     216.0,                                       
1317     233.0,                                       
1318     245.0,                                       
1319     257.0,                                       
1320     272.0,                                       
1321     286.0, //Fe 26                               
1322     297.0,                                       
1323     311.0,                                       
1324     322.0,                                       
1325     330.0,                                       
1326     334.0,                                       
1327     350.0,                                       
1328     347.0,                                       
1329     348.0,                                       
1330     343.0,                                       
1331     352.0,                                       
1332     363.0,                                       
1333     366.0,                                       
1334     379.0,                                       
1335     393.0,                                       
1336     417.0,                                       
1337     424.0,                                       
1338     428.0,                                       
1339     441.0,                                       
1340     449.0,                                       
1341     470.0,                                       
1342     470.0,                                       
1343     469.0,                                       
1344     488.0,                                       
1345     488.0, //Sn 50                               
1346     487.0,                                       
1347     485.0,                                       
1348     491.0,                                       
1349     482.0, // Xe 54                              
1350     488.0,                                       
1351     491.0,                                       
1352     501.0,                                       
1353     523.0,                                       
1354     535.0,                                       
1355     546.0,                                       
1356     560.0,                                       
1357     574.0,                                       
1358     580.0,                                       
1359     591.0,                                       
1360     614.0,                                       
1361     628.0,                                       
1362     650.0,                                       
1363     658.0,                                       
1364     674.0,                                       
1365     684.0,                                       
1366     694.0,                                       
1367     705.0,                                       
1368     718.0,                                       
1369     727.0,                                       
1370     736.0,                                       
1371     746.0,                                       
1372     757.0,                                       
1373     790.0,                                       
1374     790.0, // Au 79                              
1375     800.0,                                       
1376     810.0,                                       
1377     823.0, //Pb 82                               
1378     823.0,                                       
1379     830.0,                                       
1380     825.0,                                       
1381     794.0,                                       
1382     827.0,                                       
1383     826.0,                                       
1384     841.0,                                       
1385     847.0,                                       
1386     878.0,                                       
1387     890.0, // U 92                               
1388       900.0,                                     
1389       910.0,                                     
1390       920.0,                                     
1391       930.0,                                     
1392       940.0,                                     
1393       950.0,                                     
1394       960.0,                                     
1395      970.0,                                      
1396      980.0,                                      
1397      990.0,                                      
1398     1000.0,                                      
1399     1010.0,                                      
1400     1020.0,                                      
1401     1030.0,                                      
1402     1040.0,                                      
1403     1050.0,                                      
1404     1060.0,                                      
1405     1070.0,                                      
1406     1080.0,                                      
1407     1090.0,                                      
1408     1100.0,                                      
1409     1110.0,                                      
1410     1120.0,                                      
1411     1130.0,                                      
1412     1140.0,                                      
1413     1150.0,                                      
1414     1160.0,                                      
1415     1170.0,                                      
1416     };                                           
1417                                                  
1418 //....oooOO0OOooo........oooOO0OOooo........o    
1419                                                  
1420 const G4double G4AtimaEnergyLossModel::x0[92]    
1421 1.8639, 2.2017,  0.1304, 0.0592, 0.0305, -.01    
1422 0.2880, 0.1499,  0.1708, 0.2014, 0.1696, 0.15    
1423 0.1640, 0.0957,  0.0691, 0.0340, 0.0447, -.00    
1424 0.2267, 0.3376,  0.1767, 0.2258, 1.5262, 1.71    
1425 0.1785, 0.2267,  0.0949, 0.0599, 0.0576, 0.05    
1426 0.3189, 0.3296,  0.0549, 1.5630, 0.5473, 0.41    
1427 0.1627, 0.1520,  0.1888, 0.1058, 0.0947, 0.08    
1428 0.1560, 0.1965,  0.2117, 0.2167, 0.0559, 0.08    
1429 0.3491, 0.3776,  0.4152, 0.4267, 0.4300, 1.53    
1430 0.3144, 0.2260 };                                
1431                                                  
1432 //....oooOO0OOooo........oooOO0OOooo........o    
1433                                                  
1434 const G4double G4AtimaEnergyLossModel::x1[92]    
1435 3.2718, 3.6122,  1.6397, 1.6922, 1.9688, 2.34    
1436 3.1962, 3.0668,  3.0127, 2.8715, 2.7815, 2.71    
1437 3.0593, 3.0386,  3.0322, 3.0451, 3.1074, 3.15    
1438 3.5434, 3.6096,  3.5702, 3.6264, 4.9899, 5.07    
1439 3.2201, 3.2784,  3.1253, 3.0834, 3.1069, 3.05    
1440 3.3489, 3.4418,  3.2596, 4.7371, 3.5914, 3.45    
1441 3.3199, 3.3460,  3.4633, 3.3932, 3.4224, 3.44    
1442 3.5218, 3.4337,  3.4805, 3.4960, 3.4845, 3.54    
1443 3.8044, 3.8073,  3.8248, 3.8293, 4.0000, 4.98    
1444 3.5079, 3.3721 };                                
1445                                                  
1446 //....oooOO0OOooo........oooOO0OOooo........o    
1447                                                  
1448 const G4double G4AtimaEnergyLossModel::afermi    
1449  0.14092, 0.13443,  0.95136, 0.80392, 0.56224    
1450  0.07772, 0.08163,  0.08024, 0.14921, 0.23610    
1451  0.15754, 0.15662,  0.15436, 0.15419, 0.14973    
1452  0.09440, 0.07188,  0.06633, 0.06568, 0.06335    
1453  0.13883, 0.10525,  0.16572, 0.19342, 0.19205    
1454  0.16652, 0.13815,  0.23766, 0.23314, 0.18233    
1455  0.24280, 0.24698,  0.24448, 0.25109, 0.24453    
1456  0.24033, 0.22918,  0.17798, 0.15509, 0.15184    
1457  0.09455, 0.09359,  0.09410, 0.09282, 0.09000    
1458  0.14770,  0.19677 };                            
1459                                                  
1460 //....oooOO0OOooo........oooOO0OOooo........o    
1461                                                  
1462 const G4double G4AtimaEnergyLossModel::c[92]=    
1463 9.5835, 11.1393, 3.1221, 2.7847, 2.8477, 2.86    
1464 5.0526, 4.5297,  4.2395, 4.4351, 4.5214, 4.66    
1465 4.6949, 4.4450,  4.2659, 4.1781, 4.2702, 4.29    
1466 4.9353, 5.1411,  5.0510, 5.3210, 11.7307,12.5    
1467 5.0141, 4.8793,  4.7769, 4.7694, 4.8008, 4.93    
1468 5.6241, 5.7131,  5.9488, 12.7281, 6.9135, 6.3    
1469 5.8224, 5.8597,  6.2278, 5.8738, 5.9045, 5.91    
1470 5.9785, 5.7139,  5.5262, 5.4059, 5.3445, 5.30    
1471 6.1365, 6.2018,  6.3505, 6.4003, 6.4, 13.2839    
1472 6.0327, 5.8694 };                                
1473                                                  
1474 //....oooOO0OOooo........oooOO0OOooo........o    
1475                                                  
1476 const G4double G4AtimaEnergyLossModel::m0[92]    
1477 5.7273, 5.8347,  2.4993, 2.4339, 2.4512, 2.86    
1478 3.6452, 3.6166,  3.6345, 3.2546, 2.9158, 2.64    
1479 3.0517, 3.0302,  3.0163, 2.9896, 2.9796, 2.96    
1480 3.1314, 3.3306,  3.4176, 3.4317, 3.4670, 3.40    
1481 3.0930, 3.2549,  2.9738, 2.8707, 2.8633, 2.72    
1482 2.9319, 3.0354,  2.7276, 2.7414, 2.8866, 2.89    
1483 2.6674, 2.6403,  2.6245, 2.5977, 2.6056, 2.58    
1484 2.5643, 2.6155,  2.7623, 2.8447, 2.8627, 2.96    
1485 3.1450, 3.1608,  3.1671, 3.1830, 1.1111, 2.74    
1486 2.9845, 2.8171 };                                
1487                                                  
1488 //....oooOO0OOooo........oooOO0OOooo........o    
1489                                                  
1490 const G4double G4AtimaEnergyLossModel::del_0[    
1491 0.0,  0.00,  0.14, 0.14, 0.14, 0.12, 0.00, 0.    
1492 0.08, 0.08,  0.12, 0.14, 0.14, 0.14, 0.00, 0.    
1493 0.10, 0.12,  0.14, 0.14, 0.14, 0.12, 0.12, 0.    
1494 0.14, 0.14,  0.08, 0.10, 0.00, 0.00, 0.14, 0.    
1495 0.14, 0.14,  0.14, 0.14, 0.14, 0.14, 0.14, 0.    
1496 0.14, 0.14,  0.00, 0.00, 0.14, 0.14, 0.14, 0.    
1497 0.14, 0.14,  0.14, 0.14, 0.14, 0.14, 0.14, 0.    
1498 0.14, 0.14,  0.14, 0.14, 0.08, 0.10, 0.10, 0.    
1499 0.14, 0.14,  0.14, 0.14, 0.00, 0.00, 0.00, 0.    
1500 0.14, 0.14 };                                    
1501                                                  
1502 //....oooOO0OOooo........oooOO0OOooo........o    
1503