Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/coherent_elastic/src/G4LMsdGenerator.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/hadronic/models/coherent_elastic/src/G4LMsdGenerator.cc (Version 11.3.0) and /processes/hadronic/models/coherent_elastic/src/G4LMsdGenerator.cc (Version 9.0.p2)


  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 //                                                
 27 // G4LMsdGenerator                                
 28 //                                                
 29 //                                                
 30                                                   
 31 #include "G4DynamicParticle.hh"                   
 32 #include "G4LMsdGenerator.hh"                     
 33 #include "G4ReactionProductVector.hh"             
 34 #include "G4ReactionProduct.hh"                   
 35 #include "G4IonTable.hh"                          
 36 #include "G4NucleiProperties.hh"                  
 37 #include "G4ParticleDefinition.hh"                
 38 #include "G4HadFinalState.hh"                     
 39 #include "G4KineticTrack.hh"                      
 40 #include "G4DecayKineticTracks.hh"                
 41 #include "G4KineticTrackVector.hh"                
 42 #include "G4Log.hh"                               
 43 #include "G4PhysicsModelCatalog.hh"               
 44                                                   
 45                                                   
 46 G4LMsdGenerator::G4LMsdGenerator(const G4Strin    
 47   : G4HadronicInteraction(name), secID(-1)        
 48                                                   
 49 {                                                 
 50   fPDGencoding = 0;                               
 51   secID = G4PhysicsModelCatalog::GetModelID( "    
 52                                                   
 53   // theParticleChange = new G4HadFinalState;     
 54 }                                                 
 55                                                   
 56 G4LMsdGenerator::~G4LMsdGenerator()               
 57 {                                                 
 58   // delete theParticleChange;                    
 59 }                                                 
 60                                                   
 61 void G4LMsdGenerator::ModelDescription(std::os    
 62 {                                                 
 63   outFile << GetModelName() <<" consists of a     
 64     << " string model and a stage to de-excite    
 65     << "\n<p>"                                    
 66     << "The string model simulates the interac    
 67           << "an incident hadron with a nucleu    
 68           << "excited strings, decays these st    
 69           << "and leaves an excited nucleus. \    
 70           << "<p>The string model:\n";            
 71 }                                                 
 72                                                   
 73 //////////////////////////////////////////////    
 74 //                                                
 75 // Particle and kinematical limitation od diff    
 76                                                   
 77 G4bool                                            
 78 G4LMsdGenerator::IsApplicable( const G4HadProj    
 79                                       G4Nucleu    
 80 {                                                 
 81   G4bool applied = false;                         
 82                                                   
 83   if( ( aTrack.GetDefinition() == G4Proton::Pr    
 84   aTrack.GetDefinition() == G4Neutron::Neutron    
 85         targetNucleus.GetA_asInt() >= 1 &&        
 86       // aTrack.GetKineticEnergy() > 1800*CLHE    
 87         aTrack.GetKineticEnergy() > 300*CLHEP:    
 88     ) //  750*CLHEP::MeV )                        
 89   {                                               
 90     applied = true;                               
 91   }                                               
 92   else if( ( aTrack.GetDefinition() == G4PionP    
 93        aTrack.GetDefinition() == G4PionMinus::    
 94              targetNucleus.GetA_asInt() >= 1 &    
 95              aTrack.GetKineticEnergy() > 2340*    
 96   {                                               
 97     applied = true;                               
 98   }                                               
 99   else if( ( aTrack.GetDefinition() == G4KaonP    
100        aTrack.GetDefinition() == G4KaonMinus::    
101              targetNucleus.GetA_asInt() >= 1 &    
102              aTrack.GetKineticEnergy() > 1980*    
103   {                                               
104     applied = true;                               
105   }                                               
106   return applied;                                 
107 }                                                 
108                                                   
109 //////////////////////////////////////////////    
110 //                                                
111 // Return dissociated particle products and re    
112                                                   
113 G4HadFinalState*                                  
114 G4LMsdGenerator::ApplyYourself( const G4HadPro    
115                                       G4Nucleu    
116 {                                                 
117   theParticleChange.Clear();                      
118                                                   
119   const G4HadProjectile* aParticle = &aTrack;     
120   G4double eTkin = aParticle->GetKineticEnergy    
121                                                   
122   if( eTkin <= 1.*CLHEP::GeV && aTrack.GetDefi    
123   {                                               
124     theParticleChange.SetEnergyChange(eTkin);     
125     theParticleChange.SetMomentumChange(aTrack    
126     return &theParticleChange;                    
127   }                                               
128                                                   
129   G4int A = targetNucleus.GetA_asInt();           
130   G4int Z = targetNucleus.GetZ_asInt();           
131                                                   
132   G4double plab = aParticle->GetTotalMomentum(    
133   G4double plab2 = plab*plab;                     
134                                                   
135   const G4ParticleDefinition* theParticle = aP    
136   G4double partMass = theParticle->GetPDGMass(    
137                                                   
138   G4double oldE = partMass + eTkin;               
139                                                   
140   G4double targMass   = G4NucleiProperties::Ge    
141   G4double targMass2   = targMass*targMass;       
142                                                   
143   G4LorentzVector partLV = aParticle->Get4Mome    
144                                                   
145   G4double sumE = oldE + targMass;                
146   G4double sumE2 = sumE*sumE;                     
147                                                   
148   G4ThreeVector p1     = partLV.vect();           
149   // G4cout<<"p1 = "<<p1<<G4endl;                 
150   G4ParticleMomentum p1unit = p1.unit();          
151                                                   
152   G4double Mx = SampleMx(aParticle); // in GeV    
153   G4double t  = SampleT( aParticle, Mx); // in    
154                                                   
155   Mx *= CLHEP::GeV;                               
156                                                   
157   G4double Mx2 = Mx*Mx;                           
158                                                   
159   // equation for q|| based on sum-E-P and new    
160                                                   
161   G4double B = sumE2 + targMass2 - Mx2 - plab2    
162                                                   
163   G4double a = 4*(plab2 - sumE2);                 
164   G4double b = 4*plab*B;                          
165   G4double c = B*B - 4*sumE2*targMass2;           
166   G4double det2 = b*b - 4*a*c;                    
167   G4double qLong,  det, eRetard; // , x2, x3,     
168                                                   
169   if( det2 >= 0.)                                 
170   {                                               
171     det = std::sqrt(det2);                        
172     qLong = (-b - det)/2./a;                      
173     eRetard = std::sqrt((plab-qLong)*(plab-qLo    
174   }                                               
175   else                                            
176   {                                               
177     theParticleChange.SetEnergyChange(eTkin);     
178     theParticleChange.SetMomentumChange(aTrack    
179     return &theParticleChange;                    
180   }                                               
181   theParticleChange.SetStatusChange(stopAndKil    
182                                                   
183   plab -= qLong;                                  
184                                                   
185   G4ThreeVector pRetard = plab*p1unit;            
186                                                   
187   G4ThreeVector pTarg = p1 - pRetard;             
188                                                   
189   G4double eTarg = std::sqrt( targMass2 + pTar    
190                                                   
191   G4LorentzVector lvRetard(pRetard, eRetard);     
192   G4LorentzVector lvTarg(pTarg, eTarg);           
193                                                   
194   lvTarg += lvRetard; // sum LV                   
195                                                   
196   G4ThreeVector bst = lvTarg.boostVector();       
197                                                   
198   lvRetard.boost(-bst); // to CNS                 
199                                                   
200   G4ThreeVector pCMS   = lvRetard.vect();         
201   G4double momentumCMS = pCMS.mag();              
202   G4double tMax        = 4.0*momentumCMS*momen    
203                                                   
204   if( t > tMax ) t = tMax*G4UniformRand();        
205                                                   
206   G4double cost = 1. - 2.0*t/tMax;                
207                                                   
208                                                   
209   G4double phi  = G4UniformRand()*CLHEP::twopi    
210   G4double sint;                                  
211                                                   
212   if( cost > 1.0 || cost < -1.0 ) //              
213   {                                               
214     cost = 1.0;                                   
215     sint = 0.0;                                   
216   }                                               
217   else  // normal situation                       
218   {                                               
219     sint = std::sqrt( (1.0-cost)*(1.0+cost) );    
220   }                                               
221   G4ThreeVector v1( sint*std::cos(phi), sint*s    
222                                                   
223   v1 *= momentumCMS;                              
224                                                   
225   G4LorentzVector lvRes( v1.x(),v1.y(),v1.z(),    
226                                                   
227   lvRes.boost(bst); // to LS                      
228                                                   
229   lvTarg -= lvRes;                                
230                                                   
231   G4double eRecoil =  lvTarg.e() - targMass;      
232                                                   
233   if( eRecoil > 100.*CLHEP::MeV ) // add recoi    
234   {                                               
235     G4ParticleDefinition * recoilDef = 0;         
236                                                   
237     if      ( Z == 1 && A == 1 ) { recoilDef =    
238     else if ( Z == 1 && A == 2 ) { recoilDef =    
239     else if ( Z == 1 && A == 3 ) { recoilDef =    
240     else if ( Z == 2 && A == 3 ) { recoilDef =    
241     else if ( Z == 2 && A == 4 ) { recoilDef =    
242     else                                          
243     {                                             
244       recoilDef =                                 
245   G4ParticleTable::GetParticleTable()->GetIonT    
246     }                                             
247     G4DynamicParticle * aSec = new G4DynamicPa    
248     theParticleChange.AddSecondary(aSec, secID    
249   }                                               
250   else if( eRecoil > 0.0 )                        
251   {                                               
252     theParticleChange.SetLocalEnergyDeposit( e    
253   }                                               
254                                                   
255   G4ParticleDefinition* ddPart = G4ParticleTab    
256                                  FindParticle(    
257                                                   
258   // G4cout<<fPDGencoding<<", "<<ddPart->GetPa    
259                                                   
260   // G4DynamicParticle * aRes = new G4DynamicP    
261   //  theParticleChange.AddSecondary(aRes); //    
262                                                   
263                                                   
264                                                   
265   // Recursive decay using methods of G4Kineti    
266                                                   
267   G4KineticTrack ddkt( ddPart, 0., G4ThreeVect    
268   G4KineticTrackVector* ddktv = ddkt.Decay();     
269                                                   
270   G4DecayKineticTracks decay( ddktv );            
271                                                   
272   for( unsigned int i = 0; i < ddktv->size();     
273   {                                               
274     G4DynamicParticle * aNew =                    
275       new G4DynamicParticle( ddktv->operator[]    
276                              ddktv->operator[]    
277                                                   
278     // G4cout<<"       "<<i<<", "<<aNew->GetDe    
279                                                   
280     theParticleChange.AddSecondary(aNew, secID    
281     delete ddktv->operator[](i);                  
282   }                                               
283   delete ddktv;                                   
284                                                   
285   return &theParticleChange;                      
286 }                                                 
287                                                   
288 //////////////////////////////////////            
289 //                                                
290 // Sample Mx as Roper resonances, set PDG enco    
291                                                   
292 G4double G4LMsdGenerator::SampleMx(const G4Had    
293 {                                                 
294   G4double Mx = 0.;                               
295   G4int i;                                        
296   G4double rand = G4UniformRand();                
297                                                   
298   for( i = 0; i < 60; i++)                        
299   {                                               
300     if( rand >= fProbMx[i][1] ) break;            
301   }                                               
302   if(i <= 0)       Mx = fProbMx[0][0];            
303   else if(i >= 59) Mx = fProbMx[59][0];           
304   else             Mx = fProbMx[i][0];            
305                                                   
306   fPDGencoding = 0;                               
307                                                   
308   if ( Mx <=  1.45 )                              
309   {                                               
310     if( aParticle->GetDefinition() == G4Proton    
311     {                                             
312       Mx = 1.44;                                  
313       // fPDGencoding = 12212;                    
314       fPDGencoding = 2214;                        
315     }                                             
316     else if( aParticle->GetDefinition() == G4N    
317     {                                             
318       Mx = 1.44;                                  
319       fPDGencoding = 12112;                       
320     }                                             
321     else if( aParticle->GetDefinition() == G4P    
322     {                                             
323       // Mx = 1.3;                                
324       // fPDGencoding = 100211;                   
325       Mx = 1.26;                                  
326       fPDGencoding = 20213; // a1(1260)+          
327     }                                             
328     else if( aParticle->GetDefinition() == G4P    
329     {                                             
330       // Mx = 1.3;                                
331       // fPDGencoding = -100211;                  
332       Mx = 1.26;                                  
333       fPDGencoding = -20213; // a1(1260)-         
334     }                                             
335     else if( aParticle->GetDefinition() == G4K    
336     {                                             
337       Mx = 1.27;                                  
338       fPDGencoding = 10323;                       
339     }                                             
340     else if( aParticle->GetDefinition() == G4K    
341     {                                             
342       Mx = 1.27;                                  
343       fPDGencoding = -10323;                      
344     }                                             
345   }                                               
346   else if ( Mx <=  1.55 )                         
347   {                                               
348     if( aParticle->GetDefinition() == G4Proton    
349     {                                             
350       Mx = 1.52;                                  
351       // fPDGencoding = 2124;                     
352       fPDGencoding = 2214;                        
353     }                                             
354     else if( aParticle->GetDefinition() == G4N    
355     {                                             
356       Mx = 1.52;                                  
357       fPDGencoding = 1214;                        
358     }                                             
359     else if( aParticle->GetDefinition() == G4P    
360     {                                             
361       // Mx = 1.45;                               
362       // fPDGencoding = 10211;                    
363       Mx = 1.32;                                  
364       fPDGencoding = 215; // a2(1320)+            
365     }                                             
366     else if( aParticle->GetDefinition() == G4P    
367     {                                             
368       // Mx = 1.45;                               
369       // fPDGencoding = -10211;                   
370       Mx = 1.32;                                  
371       fPDGencoding = -215; // a2(1320)-           
372     }                                             
373     else if( aParticle->GetDefinition() == G4K    
374     {                                             
375       Mx = 1.46;                                  
376       fPDGencoding = 100321;                      
377     }                                             
378     else if( aParticle->GetDefinition() == G4K    
379     {                                             
380       Mx = 1.46;                                  
381       fPDGencoding = -100321;                     
382     }                                             
383   }                                               
384   else                                            
385   {                                               
386     if( aParticle->GetDefinition() == G4Proton    
387     {                                             
388       Mx = 1.68;                                  
389       // fPDGencoding = 12216;                    
390       fPDGencoding = 2214;                        
391     }                                             
392     else if( aParticle->GetDefinition() == G4N    
393     {                                             
394       Mx = 1.68;                                  
395       fPDGencoding = 12116;                       
396     }                                             
397     else if( aParticle->GetDefinition() == G4P    
398     {                                             
399       Mx = 1.67;                                  
400       fPDGencoding = 10215;  // pi2(1670)+        
401       // Mx = 1.45;                               
402       // fPDGencoding = 10211;                    
403     }                                             
404     else if( aParticle->GetDefinition() == G4P    
405     {                                             
406       Mx = 1.67;     // f0 problems->4pi vmg 2    
407       fPDGencoding = -10215;  // pi2(1670)-       
408       // Mx = 1.45;                               
409       // fPDGencoding = -10211;                   
410     }                                             
411     else if( aParticle->GetDefinition() == G4K    
412     {                                             
413       Mx = 1.68;                                  
414       fPDGencoding = 30323;                       
415     }                                             
416     else if( aParticle->GetDefinition() == G4K    
417     {                                             
418       Mx = 1.68;                                  
419       fPDGencoding = -30323;                      
420     }                                             
421   }                                               
422   if(fPDGencoding == 0)                           
423   {                                               
424       Mx = 1.44;                                  
425       // fPDGencoding = 12212;                    
426       fPDGencoding = 2214;                        
427   }                                               
428   G4ParticleDefinition* myResonance =             
429   G4ParticleTable::GetParticleTable()->FindPar    
430                                                   
431   if ( myResonance ) Mx = myResonance->GetPDGM    
432                                                   
433   // G4cout<<"PDG-ID = "<<fPDGencoding<<"; wit    
434                                                   
435   return Mx/CLHEP::GeV;                           
436 }                                                 
437                                                   
438 //////////////////////////////////////            
439 //                                                
440 // Sample t with kinematic limitations of Mx a    
441                                                   
442 G4double G4LMsdGenerator::SampleT(  const G4Ha    
443 G4double Mx)                                      
444 {                                                 
445   G4double t=0., b=0., rTkin = 50.*CLHEP::GeV,    
446   G4int i;                                        
447                                                   
448   for( i = 0; i < 23; ++i)                        
449   {                                               
450     if( Mx <= fMxBdata[i][0] ) break;             
451   }                                               
452   if( i <= 0 )       b = fMxBdata[0][1];          
453   else if( i >= 22 ) b = fMxBdata[22][1];         
454   else               b = fMxBdata[i][1];          
455                                                   
456   if( eTkin > rTkin ) b *= 1. + G4Log(eTkin/rT    
457                                                   
458   G4double rand = G4UniformRand();                
459                                                   
460   t = -G4Log(rand)/b;                             
461                                                   
462   t *= (CLHEP::GeV*CLHEP::GeV); // in G4 inter    
463                                                   
464   return t;                                       
465 }                                                 
466                                                   
467                                                   
468 //////////////////////////////////////////////    
469 //                                                
470 // Integral spectrum of Mx (GeV)                  
471                                                   
472 const G4double G4LMsdGenerator::fProbMx[60][2]    
473 {                                                 
474   {1.000000e+00,  1.000000e+00},                  
475   {1.025000e+00,  1.000000e+00},                  
476   {1.050000e+00,  1.000000e+00},                  
477   {1.075000e+00,  1.000000e+00},                  
478   {1.100000e+00,  9.975067e-01},                  
479   {1.125000e+00,  9.934020e-01},                  
480   {1.150000e+00,  9.878333e-01},                  
481   {1.175000e+00,  9.805002e-01},                  
482   {1.200000e+00,  9.716846e-01},                  
483   {1.225000e+00,  9.604761e-01},                  
484   {1.250000e+00,  9.452960e-01},                  
485   {1.275000e+00,  9.265278e-01},                  
486   {1.300000e+00,  9.053632e-01},                  
487   {1.325000e+00,  8.775566e-01},                  
488   {1.350000e+00,  8.441969e-01},                  
489   {1.375000e+00,  8.076336e-01},                  
490   {1.400000e+00,  7.682520e-01},                  
491   {1.425000e+00,  7.238306e-01},                  
492   {1.450000e+00,  6.769306e-01},                  
493   {1.475000e+00,  6.303898e-01},                  
494   {1.500000e+00,  5.824632e-01},                  
495   {1.525000e+00,  5.340696e-01},                  
496   {1.550000e+00,  4.873736e-01},                  
497   {1.575000e+00,  4.422901e-01},                  
498   {1.600000e+00,  3.988443e-01},                  
499   {1.625000e+00,  3.583727e-01},                  
500   {1.650000e+00,  3.205405e-01},                  
501   {1.675000e+00,  2.856655e-01},                  
502   {1.700000e+00,  2.537508e-01},                  
503   {1.725000e+00,  2.247863e-01},                  
504   {1.750000e+00,  1.985798e-01},                  
505   {1.775000e+00,  1.750252e-01},                  
506   {1.800000e+00,  1.539777e-01},                  
507   {1.825000e+00,  1.352741e-01},                  
508   {1.850000e+00,  1.187157e-01},                  
509   {1.875000e+00,  1.040918e-01},                  
510   {1.900000e+00,  9.118422e-02},                  
511   {1.925000e+00,  7.980909e-02},                  
512   {1.950000e+00,  6.979378e-02},                  
513   {1.975000e+00,  6.097771e-02},                  
514   {2.000000e+00,  5.322122e-02},                  
515   {2.025000e+00,  4.639628e-02},                  
516   {2.050000e+00,  4.039012e-02},                  
517   {2.075000e+00,  3.510275e-02},                  
518   {2.100000e+00,  3.044533e-02},                  
519   {2.125000e+00,  2.633929e-02},                  
520   {2.150000e+00,  2.271542e-02},                  
521   {2.175000e+00,  1.951295e-02},                  
522   {2.200000e+00,  1.667873e-02},                  
523   {2.225000e+00,  1.416633e-02},                  
524   {2.250000e+00,  1.193533e-02},                  
525   {2.275000e+00,  9.950570e-03},                  
526   {2.300000e+00,  8.181515e-03},                  
527   {2.325000e+00,  6.601664e-03},                  
528   {2.350000e+00,  5.188025e-03},                  
529   {2.375000e+00,  3.920655e-03},                  
530   {2.400000e+00,  2.782246e-03},                  
531   {2.425000e+00,  1.757765e-03},                  
532   {2.450000e+00,  8.341435e-04},                  
533   {2.475000e+00,  0.000000e+00}                   
534 };                                                
535                                                   
536 //////////////////////////////////////////////    
537 //                                                
538 // Slope b (1/GeV/GeV) vs Mx (GeV) for t-sampl    
539                                                   
540 const G4double G4LMsdGenerator::fMxBdata[23][2    
541 {                                                 
542   {1.09014,      17.8620},                        
543   {1.12590,      19.2831},                        
544   {1.18549,      17.6907},                        
545   {1.21693,      16.4760},                        
546   {1.25194,      15.3867},                        
547   {1.26932,      14.4236},                        
548   {1.29019,      13.2931},                        
549   {1.30755,      12.2882},                        
550   {1.31790,      11.4509},                        
551   {1.33888,      10.6969},                        
552   {1.34911,      9.44130},                        
553   {1.37711,      8.56148},                        
554   {1.39101,      7.76593},                        
555   {1.42608,      6.88582},                        
556   {1.48593,      6.13019},                        
557   {1.53179,      5.87723},                        
558   {1.58111,      5.37308},                        
559   {1.64105,      4.95217},                        
560   {1.69037,      4.44803},                        
561   {1.81742,      3.89879},                        
562   {1.88096,      3.68693},                        
563   {1.95509,      3.43278},                        
564   {2.02219,      3.30445}                         
565 };                                                
566                                                   
567                                                   
568                                                   
569 //                                                
570 //                                                
571 /////////////////////////////////////////////     
572