Geant4  10.01.p02
G4MuonVDNuclearModel.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: $
27 //
28 // Author: D.H. Wright
29 // Date: 2 February 2011
30 //
31 // Description: model of muon nuclear interaction in which a gamma from
32 // the virtual photon spectrum interacts in the nucleus as
33 // a real gamma at low energies and as a pi0 at high energies.
34 // Kokoulin's muon cross section and equivalent gamma spectrum
35 // are used.
36 //
37 
38 #include "G4MuonVDNuclearModel.hh"
39 
40 #include "Randomize.hh"
41 #include "G4PhysicalConstants.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "G4CascadeInterface.hh"
44 #include "G4TheoFSGenerator.hh"
46 #include "G4ExcitationHandler.hh"
47 #include "G4PreCompoundModel.hh"
49 #include "G4ExcitedStringDecay.hh"
50 #include "G4FTFModel.hh"
52 
54  : G4HadronicInteraction("G4MuonVDNuclearModel")
55 {
56  SetMinEnergy(0.0);
57  SetMaxEnergy(1*PeV);
58  CutFixed = 0.2*GeV;
59  NBIN = 1000;
60 
61  for (G4int k = 0; k < 5; k++) {
62  for (G4int j = 0; j < 8; j++) {
63  for (G4int i = 0; i < 1001; i++) {
64  proba[k][j][i] = 0.0;
65  ya[i] = 0.0;
66  }
67  }
68  }
69 
71 
72  // reuse existing pre-compound model
73  G4GeneratorPrecompoundInterface* precoInterface
77  G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
78  if(!pre) { pre = new G4PreCompoundModel(); }
79  precoInterface->SetDeExcitation(pre);
80 
81  // Build FTFP model
82  ftfp = new G4TheoFSGenerator();
83  ftfp->SetTransport(precoInterface);
86  G4FTFModel* theStringModel = new G4FTFModel;
87  theStringModel->SetFragmentationModel(theStringDecay);
88  ftfp->SetHighEnergyGenerator(theStringModel);
89 
90  // Build Bertini cascade
91  bert = new G4CascadeInterface();
92 }
93 
94 
96 {
97  delete theFragmentation;
98  delete theStringDecay;
99 }
100 
101 
104  G4Nucleus& targetNucleus)
105 {
107 
108  // For very low energy, return initial track
109  G4double epmax = aTrack.GetTotalEnergy() - 0.5*proton_mass_c2;
110  if (epmax <= CutFixed) {
113  theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
114  return &theParticleChange;
115  }
116 
117  // Produce recoil muon and transferred photon
118  G4DynamicParticle* transferredPhoton = CalculateEMVertex(aTrack, targetNucleus);
119 
120  // Interact the gamma with the nucleus
121  CalculateHadronicVertex(transferredPhoton, targetNucleus);
122  return &theParticleChange;
123 }
124 
125 
128  G4Nucleus& targetNucleus)
129 {
130  // Select sampling table
131  G4double KineticEnergy = aTrack.GetKineticEnergy();
132  G4double TotalEnergy = aTrack.GetTotalEnergy();
134  G4double lnZ = std::log(G4double(targetNucleus.GetZ_asInt() ) );
135 
136  G4double epmin = CutFixed;
137  G4double epmax = TotalEnergy - 0.5*proton_mass_c2;
138  G4double m0 = 0.2*GeV;
139 
140  G4double delmin = 1.e10;
141  G4double del;
142  G4int izz = 0;
143  G4int itt = 0;
144  G4int NBINminus1 = NBIN - 1;
145 
146  G4int nzdat = 5;
147  G4double zdat[] = {1.,4.,13.,29.,92.};
148  for (G4int iz = 0; iz < nzdat; iz++) {
149  del = std::abs(lnZ-std::log(zdat[iz]));
150  if (del < delmin) {
151  delmin = del;
152  izz = iz;
153  }
154  }
155 
156  G4int ntdat = 8;
157  G4double tdat[] = {1.e3,1.e4,1.e5,1.e6,1.e7,1.e8,1.e9,1.e10};
158  delmin = 1.e10;
159  for (G4int it = 0; it < ntdat; it++) {
160  del = std::abs(std::log(KineticEnergy)-std::log(tdat[it]) );
161  if (del < delmin) {
162  delmin = del;
163  itt = it;
164  }
165  }
166 
167  // Sample the energy transfer according to the probability table
168  G4double r = G4UniformRand();
169 
170  G4int iy = -1;
171  do {
172  iy += 1 ;
173  } while (((proba[izz][itt][iy]) < r)&&(iy < NBINminus1)) ;
174 
175  // Sampling is done uniformly in y in the bin
176 
177  G4double y;
178  if (iy < NBIN)
179  y = ya[iy] + G4UniformRand() * (ya[iy+1] - ya[iy]);
180  else
181  y = ya[iy];
182 
183  G4double x = std::exp(y);
184  G4double ep = epmin*std::exp(x*std::log(epmax/epmin) );
185 
186  // Sample scattering angle of mu, but first t should be sampled.
187  G4double yy = ep/TotalEnergy;
188  G4double tmin = Mass*Mass*yy*yy/(1.-yy);
189  G4double tmax = 2.*proton_mass_c2*ep;
190  G4double t1;
191  G4double t2;
192  if (m0 < ep) {
193  t1 = m0*m0;
194  t2 = ep*ep;
195  } else {
196  t1 = ep*ep;
197  t2 = m0*m0;
198  }
199 
200  G4double w1 = tmax*t1;
201  G4double w2 = tmax+t1;
202  G4double w3 = tmax*(tmin+t1)/(tmin*w2);
203  G4double y1 = 1.-yy;
204  G4double y2 = 0.5*yy*yy;
205  G4double y3 = y1+y2;
206 
207  G4double t;
208  G4double rej;
209 
210  // Now sample t
211  G4int ntry = 0;
212  do
213  {
214  ntry += 1;
215  t = w1/(w2*std::exp(G4UniformRand()*std::log(w3))-tmax);
216  rej = (1.-t/tmax)*(y1*(1.-tmin/t)+y2)/(y3*(1.-t/t2));
217  } while (G4UniformRand() > rej) ;
218 
219  // compute angle from t
220  G4double sinth2 =
221  0.5*(t-tmin)/(2.*(TotalEnergy*(TotalEnergy-ep)-Mass*Mass)-tmin);
222  G4double theta = std::acos(1. - 2.*sinth2);
223 
224  G4double phi = twopi*G4UniformRand();
225  G4double sinth = std::sin(theta);
226  G4double dirx = sinth*std::cos(phi);
227  G4double diry = sinth*std::sin(phi);
228  G4double dirz = std::cos(theta);
229  G4ThreeVector finalDirection(dirx,diry,dirz);
230  G4ThreeVector ParticleDirection(aTrack.Get4Momentum().vect().unit() );
231  finalDirection.rotateUz(ParticleDirection);
232 
233  G4double NewKinEnergy = KineticEnergy - ep;
234  G4double finalMomentum = std::sqrt(NewKinEnergy*(NewKinEnergy+2.*Mass) );
235  G4double Ef = NewKinEnergy + Mass;
236  G4double initMomentum = std::sqrt(KineticEnergy*(TotalEnergy+Mass) );
237 
238  // Set energy and direction of scattered primary in theParticleChange
240  theParticleChange.SetEnergyChange(NewKinEnergy);
241  theParticleChange.SetMomentumChange(finalDirection);
242 
243  // Now create the emitted gamma
244  G4LorentzVector primaryMomentum(initMomentum*ParticleDirection, TotalEnergy);
245  G4LorentzVector fsMomentum(finalMomentum*finalDirection, Ef);
246  G4LorentzVector momentumTransfer = primaryMomentum - fsMomentum;
247 
248  G4DynamicParticle* gamma =
249  new G4DynamicParticle(G4Gamma::Gamma(), momentumTransfer);
250 
251  return gamma;
252 }
253 
254 
255 void
257  G4Nucleus& target)
258 {
259  G4HadFinalState* hfs = 0;
260  G4double gammaE = incident->GetTotalEnergy();
261 
262  if (gammaE < 10*GeV) {
263  G4HadProjectile projectile(*incident);
264  hfs = bert->ApplyYourself(projectile, target);
265  } else {
266  // convert incident gamma to a pi0
268  G4double piKE = incident->GetTotalEnergy() - piMass;
269  G4double piMom = std::sqrt(piKE*(piKE + 2*piMass) );
270  G4ThreeVector piMomentum(incident->GetMomentumDirection() );
271  piMomentum *= piMom;
272  G4DynamicParticle theHadron(G4PionZero::PionZero(), piMomentum);
273  G4HadProjectile projectile(theHadron);
274  hfs = ftfp->ApplyYourself(projectile, target);
275  }
276 
277  delete incident;
278 
279  // Copy secondaries from sub-model to model
281 }
282 
283 
285 {
286  G4double adat[] = {1.01,9.01,26.98,63.55,238.03};
287  G4double zdat[] = {1.,4.,13.,29.,92.};
288  G4int nzdat = 5;
289 
290  G4double tdat[] = {1.e3,1.e4,1.e5,1.e6,1.e7,1.e8,1.e9,1.e10};
291  G4int ntdat = 8;
292 
293  G4int nbin;
294  G4double KineticEnergy;
295  G4double TotalEnergy;
296  G4double Maxep;
297  G4double CrossSection;
298 
299  G4double c;
300  G4double y;
301  G4double ymin,ymax;
302  G4double dy,yy;
303  G4double dx,x;
304  G4double ep;
305 
306  G4double AtomicNumber;
307  G4double AtomicWeight;
308 
309  for (G4int iz = 0; iz < nzdat; iz++) {
310  AtomicNumber = zdat[iz];
311  AtomicWeight = adat[iz]*(g/mole);
312 
313  for (G4int it = 0; it < ntdat; it++) {
314  KineticEnergy = tdat[it];
315  TotalEnergy = KineticEnergy + G4MuonMinus::MuonMinus()->GetPDGMass();
316  Maxep = TotalEnergy - 0.5*proton_mass_c2;
317 
318  CrossSection = 0.0;
319 
320  // Calculate the differential cross section
321  // numerical integration in log .........
322  c = std::log(Maxep/CutFixed);
323  ymin = -5.0;
324  ymax = 0.0;
325  dy = (ymax-ymin)/NBIN;
326 
327  nbin=-1;
328 
329  y = ymin - 0.5*dy;
330  yy = ymin - dy;
331  for (G4int i = 0; i < NBIN; i++) {
332  y += dy;
333  x = std::exp(y);
334  yy += dy;
335  dx = std::exp(yy+dy)-std::exp(yy);
336 
337  ep = CutFixed*std::exp(c*x);
338 
339  CrossSection +=
340  ep*dx*muNucXS.ComputeDDMicroscopicCrossSection(KineticEnergy,
341  AtomicNumber,
342  AtomicWeight, ep);
343  if (nbin < NBIN) {
344  nbin += 1;
345  ya[nbin] = y;
346  proba[iz][it][nbin] = CrossSection;
347  }
348  }
349  ya[NBIN] = 0.;
350 
351  if (CrossSection > 0.0) {
352  for (G4int ib = 0; ib <= nbin; ib++) proba[iz][it][ib] /= CrossSection;
353  }
354  } // loop on it
355  } // loop on iz
356 
357  // G4cout << " Kokoulin XS = "
358  // << muNucXS.ComputeDDMicroscopicCrossSection(1*GeV, 20.0, 40.0*g/mole, 0.3*GeV)/millibarn
359  // << G4endl;
360 }
361 
G4double proba[5][8][1001]
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
void AddSecondaries(const std::vector< G4HadSecondary > &addSecs)
CLHEP::Hep3Vector G4ThreeVector
G4double GetTotalEnergy() const
void CalculateHadronicVertex(G4DynamicParticle *incident, G4Nucleus &target)
void SetFragmentationModel(G4VStringFragmentation *aModel)
G4double ComputeDDMicroscopicCrossSection(G4double incidentKE, G4double Z, G4double A, G4double epsilon)
int G4int
Definition: G4Types.hh:78
void SetHighEnergyGenerator(G4VHighEnergyGenerator *const value)
G4TheoFSGenerator * ftfp
void SetStatusChange(G4HadFinalStateStatus aS)
void SetMinEnergy(G4double anEnergy)
#define G4UniformRand()
Definition: Randomize.hh:93
const G4ThreeVector & GetMomentumDirection() const
G4double iz
Definition: TRTMaterials.hh:39
G4double GetKineticEnergy() const
static const double GeV
Definition: G4SIunits.hh:196
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4LorentzVector & Get4Momentum() const
static G4PionZero * PionZero()
Definition: G4PionZero.cc:108
G4DynamicParticle * CalculateEMVertex(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4HadronicInteraction * FindModel(const G4String &name)
void SetEnergyChange(G4double anEnergy)
static const double PeV
Definition: G4SIunits.hh:198
G4double GetPDGMass() const
G4CascadeInterface * bert
static const double g
Definition: G4SIunits.hh:162
static G4HadronicInteractionRegistry * Instance()
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void SetMaxEnergy(const G4double anEnergy)
void SetDeExcitation(G4VPreCompoundModel *ptr)
static const double mole
Definition: G4SIunits.hh:265
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
void SetTransport(G4VIntraNuclearTransportModel *const value)
G4LundStringFragmentation * theFragmentation
double G4double
Definition: G4Types.hh:76
void SetMomentumChange(const G4ThreeVector &aV)
G4KokoulinMuonNuclearXS muNucXS
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4ExcitedStringDecay * theStringDecay
G4double GetTotalEnergy() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)
CLHEP::HepLorentzVector G4LorentzVector