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