Geant4  10.02
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 "G4Log.hh"
42 #include "G4PhysicalConstants.hh"
43 #include "G4SystemOfUnits.hh"
44 #include "G4CascadeInterface.hh"
45 #include "G4TheoFSGenerator.hh"
47 #include "G4ExcitationHandler.hh"
48 #include "G4PreCompoundModel.hh"
50 #include "G4ExcitedStringDecay.hh"
51 #include "G4FTFModel.hh"
55 
57  : G4HadronicInteraction("G4MuonVDNuclearModel")
58 {
60  GetCrossSectionDataSet(G4KokoulinMuonNuclearXS::Default_Name());
61 
62  SetMinEnergy(0.0);
63  SetMaxEnergy(1*PeV);
64  CutFixed = 0.2*GeV;
65  NBIN = 1000;
66 
67  for (G4int k = 0; k < 5; k++) {
68  for (G4int j = 0; j < 8; j++) {
69  for (G4int i = 0; i < 1001; i++) {
70  proba[k][j][i] = 0.0;
71  ya[i] = 0.0;
72  }
73  }
74  }
75 
77 
78  // reuse existing pre-compound model
79  G4GeneratorPrecompoundInterface* precoInterface
83  G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
84  if(!pre) { pre = new G4PreCompoundModel(); }
85  precoInterface->SetDeExcitation(pre);
86 
87  // Build FTFP model
88  ftfp = new G4TheoFSGenerator();
89  ftfp->SetTransport(precoInterface);
92  G4FTFModel* theStringModel = new G4FTFModel;
93  theStringModel->SetFragmentationModel(theStringDecay);
94  ftfp->SetHighEnergyGenerator(theStringModel);
95 
96  // Build Bertini cascade
97  bert = new G4CascadeInterface();
98 }
99 
100 
102 {
103  delete theFragmentation;
104  delete theStringDecay;
105 }
106 
107 
110  G4Nucleus& targetNucleus)
111 {
113 
114  // For very low energy, return initial track
115  G4double epmax = aTrack.GetTotalEnergy() - 0.5*proton_mass_c2;
116  if (epmax <= CutFixed) {
119  theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
120  return &theParticleChange;
121  }
122 
123  // Produce recoil muon and transferred photon
124  G4DynamicParticle* transferredPhoton = CalculateEMVertex(aTrack, targetNucleus);
125 
126  // Interact the gamma with the nucleus
127  CalculateHadronicVertex(transferredPhoton, targetNucleus);
128  return &theParticleChange;
129 }
130 
131 
134  G4Nucleus& targetNucleus)
135 {
136  // Select sampling table
137  G4double KineticEnergy = aTrack.GetKineticEnergy();
138  G4double TotalEnergy = aTrack.GetTotalEnergy();
140  G4double lnZ = G4Log(G4double(targetNucleus.GetZ_asInt() ) );
141 
142  G4double epmin = CutFixed;
143  G4double epmax = TotalEnergy - 0.5*proton_mass_c2;
144  G4double m0 = 0.2*GeV;
145 
146  G4double delmin = 1.e10;
147  G4double del;
148  G4int izz = 0;
149  G4int itt = 0;
150  G4int NBINminus1 = NBIN - 1;
151 
152  G4int nzdat = 5;
153  G4double zdat[] = {1.,4.,13.,29.,92.};
154  for (G4int iz = 0; iz < nzdat; iz++) {
155  del = std::abs(lnZ-G4Log(zdat[iz]));
156  if (del < delmin) {
157  delmin = del;
158  izz = iz;
159  }
160  }
161 
162  G4int ntdat = 8;
163  G4double tdat[] = {1.e3,1.e4,1.e5,1.e6,1.e7,1.e8,1.e9,1.e10};
164  delmin = 1.e10;
165  for (G4int it = 0; it < ntdat; it++) {
166  del = std::abs(G4Log(KineticEnergy)-G4Log(tdat[it]) );
167  if (del < delmin) {
168  delmin = del;
169  itt = it;
170  }
171  }
172 
173  // Sample the energy transfer according to the probability table
174  G4double r = G4UniformRand();
175 
176  G4int iy = -1;
178  ed << " While count exceeded " << G4endl;
179 
180  do {
181  iy += 1 ;
182  if (iy > 10000) {
183  G4Exception("G4RPGAntiProtonInelastic::Cascade()", "HAD_RPG_100", JustWarning, ed);
184  break;
185  }
186  } while (((proba[izz][itt][iy]) < r)&&(iy < NBINminus1)) ; /* Loop checking, 01.09.2015, D.Wright */
187 
188  // Sampling is done uniformly in y in the bin
189 
190  G4double y;
191  if (iy < NBIN)
192  y = ya[iy] + G4UniformRand() * (ya[iy+1] - ya[iy]);
193  else
194  y = ya[iy];
195 
196  G4double x = G4Exp(y);
197  G4double ep = epmin*G4Exp(x*G4Log(epmax/epmin) );
198 
199  // Sample scattering angle of mu, but first t should be sampled.
200  G4double yy = ep/TotalEnergy;
201  G4double tmin = Mass*Mass*yy*yy/(1.-yy);
202  G4double tmax = 2.*proton_mass_c2*ep;
203  G4double t1;
204  G4double t2;
205  if (m0 < ep) {
206  t1 = m0*m0;
207  t2 = ep*ep;
208  } else {
209  t1 = ep*ep;
210  t2 = m0*m0;
211  }
212 
213  G4double w1 = tmax*t1;
214  G4double w2 = tmax+t1;
215  G4double w3 = tmax*(tmin+t1)/(tmin*w2);
216  G4double y1 = 1.-yy;
217  G4double y2 = 0.5*yy*yy;
218  G4double y3 = y1+y2;
219 
220  G4double t;
221  G4double rej;
222 
223  // Now sample t
224  G4int ntry = 0;
226  eda << " While count exceeded " << G4endl;
227 
228  do
229  {
230  ntry += 1;
231  if (ntry > 10000) {
232  G4Exception("G4MuonVDNuclearModel::CalculateEMVertex()", "HAD_RPG_100", JustWarning, eda);
233  break;
234  }
235 
236  t = w1/(w2*G4Exp(G4UniformRand()*G4Log(w3))-tmax);
237  rej = (1.-t/tmax)*(y1*(1.-tmin/t)+y2)/(y3*(1.-t/t2));
238  } while (G4UniformRand() > rej) ; /* Loop checking, 01.09.2015, D.Wright */
239 
240  // compute angle from t
241  G4double sinth2 =
242  0.5*(t-tmin)/(2.*(TotalEnergy*(TotalEnergy-ep)-Mass*Mass)-tmin);
243  G4double theta = std::acos(1. - 2.*sinth2);
244 
245  G4double phi = twopi*G4UniformRand();
246  G4double sinth = std::sin(theta);
247  G4double dirx = sinth*std::cos(phi);
248  G4double diry = sinth*std::sin(phi);
249  G4double dirz = std::cos(theta);
250  G4ThreeVector finalDirection(dirx,diry,dirz);
251  G4ThreeVector ParticleDirection(aTrack.Get4Momentum().vect().unit() );
252  finalDirection.rotateUz(ParticleDirection);
253 
254  G4double NewKinEnergy = KineticEnergy - ep;
255  G4double finalMomentum = std::sqrt(NewKinEnergy*(NewKinEnergy+2.*Mass) );
256  G4double Ef = NewKinEnergy + Mass;
257  G4double initMomentum = std::sqrt(KineticEnergy*(TotalEnergy+Mass) );
258 
259  // Set energy and direction of scattered primary in theParticleChange
261  theParticleChange.SetEnergyChange(NewKinEnergy);
262  theParticleChange.SetMomentumChange(finalDirection);
263 
264  // Now create the emitted gamma
265  G4LorentzVector primaryMomentum(initMomentum*ParticleDirection, TotalEnergy);
266  G4LorentzVector fsMomentum(finalMomentum*finalDirection, Ef);
267  G4LorentzVector momentumTransfer = primaryMomentum - fsMomentum;
268 
269  G4DynamicParticle* gamma =
270  new G4DynamicParticle(G4Gamma::Gamma(), momentumTransfer);
271 
272  return gamma;
273 }
274 
275 
276 void
278  G4Nucleus& target)
279 {
280  G4HadFinalState* hfs = 0;
281  G4double gammaE = incident->GetTotalEnergy();
282 
283  if (gammaE < 10*GeV) {
284  G4HadProjectile projectile(*incident);
285  hfs = bert->ApplyYourself(projectile, target);
286  } else {
287  // convert incident gamma to a pi0
289  G4double piKE = incident->GetTotalEnergy() - piMass;
290  G4double piMom = std::sqrt(piKE*(piKE + 2*piMass) );
291  G4ThreeVector piMomentum(incident->GetMomentumDirection() );
292  piMomentum *= piMom;
293  G4DynamicParticle theHadron(G4PionZero::PionZero(), piMomentum);
294  G4HadProjectile projectile(theHadron);
295  hfs = ftfp->ApplyYourself(projectile, target);
296  }
297 
298  delete incident;
299 
300  // Copy secondaries from sub-model to model
302 }
303 
304 
306 {
307  G4double adat[] = {1.01,9.01,26.98,63.55,238.03};
308  G4double zdat[] = {1.,4.,13.,29.,92.};
309  G4int nzdat = 5;
310 
311  G4double tdat[] = {1.e3,1.e4,1.e5,1.e6,1.e7,1.e8,1.e9,1.e10};
312  G4int ntdat = 8;
313 
314  G4int nbin;
315  G4double KineticEnergy;
316  G4double TotalEnergy;
317  G4double Maxep;
318  G4double CrossSection;
319 
320  G4double c;
321  G4double y;
322  G4double ymin,ymax;
323  G4double dy,yy;
324  G4double dx,x;
325  G4double ep;
326 
327  G4double AtomicNumber;
328  G4double AtomicWeight;
329 
330  for (G4int iz = 0; iz < nzdat; iz++) {
331  AtomicNumber = zdat[iz];
332  AtomicWeight = adat[iz]*(g/mole);
333 
334  for (G4int it = 0; it < ntdat; it++) {
335  KineticEnergy = tdat[it];
336  TotalEnergy = KineticEnergy + G4MuonMinus::MuonMinus()->GetPDGMass();
337  Maxep = TotalEnergy - 0.5*proton_mass_c2;
338 
339  CrossSection = 0.0;
340 
341  // Calculate the differential cross section
342  // numerical integration in log .........
343  c = G4Log(Maxep/CutFixed);
344  ymin = -5.0;
345  ymax = 0.0;
346  dy = (ymax-ymin)/NBIN;
347 
348  nbin=-1;
349 
350  y = ymin - 0.5*dy;
351  yy = ymin - dy;
352  for (G4int i = 0; i < NBIN; i++) {
353  y += dy;
354  x = G4Exp(y);
355  yy += dy;
356  dx = G4Exp(yy+dy)-G4Exp(yy);
357 
358  ep = CutFixed*G4Exp(c*x);
359 
360  CrossSection +=
361  ep*dx*muNucXS->ComputeDDMicroscopicCrossSection(KineticEnergy,
362  AtomicNumber,
363  AtomicWeight, ep);
364  if (nbin < NBIN) {
365  nbin += 1;
366  ya[nbin] = y;
367  proba[iz][it][nbin] = CrossSection;
368  }
369  }
370  ya[NBIN] = 0.;
371 
372  if (CrossSection > 0.0) {
373  for (G4int ib = 0; ib <= nbin; ib++) proba[iz][it][ib] /= CrossSection;
374  }
375  } // loop on it
376  } // loop on iz
377 
378  // G4cout << " Kokoulin XS = "
379  // << muNucXS->ComputeDDMicroscopicCrossSection(1*GeV, 20.0,
380  // 40.0*g/mole, 0.3*GeV)/millibarn
381  // << G4endl;
382 }
383 
G4double proba[5][8][1001]
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
void AddSecondaries(const std::vector< G4HadSecondary > &addSecs)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
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:97
static const char * Default_Name()
const G4ThreeVector & GetMomentumDirection() const
G4double iz
Definition: TRTMaterials.hh:39
G4double GetKineticEnergy() const
static const double twopi
Definition: G4SIunits.hh:75
static G4CrossSectionDataSetRegistry * Instance()
static const double GeV
Definition: G4SIunits.hh:214
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)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4HadronicInteraction * FindModel(const G4String &name)
void SetEnergyChange(G4double anEnergy)
static const double PeV
Definition: G4SIunits.hh:216
G4double GetPDGMass() const
G4CascadeInterface * bert
static const double g
Definition: G4SIunits.hh:180
static G4HadronicInteractionRegistry * Instance()
const G4double x[NPOINTSGL]
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4KokoulinMuonNuclearXS * muNucXS
void SetMaxEnergy(const G4double anEnergy)
void SetDeExcitation(G4VPreCompoundModel *ptr)
static const double mole
Definition: G4SIunits.hh:283
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
#define G4endl
Definition: G4ios.hh:61
void SetTransport(G4VIntraNuclearTransportModel *const value)
G4LundStringFragmentation * theFragmentation
double G4double
Definition: G4Types.hh:76
void SetMomentumChange(const G4ThreeVector &aV)
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4ExcitedStringDecay * theStringDecay
G4double GetTotalEnergy() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)
CLHEP::HepLorentzVector G4LorentzVector