Geant4 11.1.1
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//
27// Author: D.H. Wright
28// Date: 2 February 2011
29//
30// Description: model of muon nuclear interaction in which a gamma from
31// the virtual photon spectrum interacts in the nucleus as
32// a real gamma at low energies and as a pi0 at high energies.
33// Kokoulin's muon cross section and equivalent gamma spectrum
34// are used.
35//
36
38
39#include "Randomize.hh"
40#include "G4Log.hh"
42#include "G4SystemOfUnits.hh"
43#include "G4CascadeInterface.hh"
44#include "G4TheoFSGenerator.hh"
46#include "G4PreCompoundModel.hh"
49#include "G4FTFModel.hh"
53#include "G4ElementData.hh"
54#include "G4Physics2DVector.hh"
55#include "G4Pow.hh"
57
58const G4int G4MuonVDNuclearModel::zdat[] = {1, 4, 13, 29, 92};
59const G4double G4MuonVDNuclearModel::adat[] = {1.01,9.01,26.98,63.55,238.03};
60const G4double G4MuonVDNuclearModel::tdat[] = {
61 1.e3,2.e3,3.e3,4.e3,5.e3,6.e3,7.e3,8.e3,9.e3,
62 1.e4,2.e4,3.e4,4.e4,5.e4,6.e4,7.e4,8.e4,9.e4,
63 1.e5,2.e5,3.e5,4.e5,5.e5,6.e5,7.e5,8.e5,9.e5,
64 1.e6,2.e6,3.e6,4.e6,5.e6,6.e6,7.e6,8.e6,9.e6,
65 1.e7,2.e7,3.e7,4.e7,5.e7,6.e7,7.e7,8.e7,9.e7,
66 1.e8,2.e8,3.e8,4.e8,5.e8,6.e8,7.e8,8.e8,9.e8,
67 1.e9,2.e9,3.e9,4.e9,5.e9,6.e9,7.e9,8.e9,9.e9,
68 1.e10,2.e10,3.e10,4.e10,5.e10,6.e10,7.e10,8.e10,9.e10,1.e11};
69
70G4ElementData* G4MuonVDNuclearModel::fElementData = nullptr;
71
73 : G4HadronicInteraction("G4MuonVDNuclearModel"),isMaster(false)
74{
76 GetCrossSectionDataSet(G4KokoulinMuonNuclearXS::Default_Name());
77
78 SetMinEnergy(0.0);
79 SetMaxEnergy(1*CLHEP::PeV);
80 CutFixed = 0.2*CLHEP::GeV;
81
82 if(!fElementData && G4Threading::IsMasterThread()) {
83 fElementData = new G4ElementData();
84 MakeSamplingTable();
85 isMaster = true;
86 }
87
88 // reuse existing pre-compound model
89 G4GeneratorPrecompoundInterface* precoInterface
93 G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
94 if(!pre) { pre = new G4PreCompoundModel(); }
95 precoInterface->SetDeExcitation(pre);
96
97 // Build FTFP model
98 ftfp = new G4TheoFSGenerator();
99 ftfp->SetTransport(precoInterface);
100 theFragmentation = new G4LundStringFragmentation();
101 theStringDecay = new G4ExcitedStringDecay(theFragmentation);
102 G4FTFModel* theStringModel = new G4FTFModel;
103 theStringModel->SetFragmentationModel(theStringDecay);
104 ftfp->SetHighEnergyGenerator(theStringModel);
105
106 // Build Bertini cascade
107 bert = new G4CascadeInterface();
108
109 // Creator model ID
110 secID = G4PhysicsModelCatalog::GetModelID( "model_" + GetModelName() );
111}
112
114{
115 delete theFragmentation;
116 delete theStringDecay;
117
118 if(isMaster) {
119 delete fElementData;
120 fElementData = nullptr;
121 }
122}
123
126 G4Nucleus& targetNucleus)
127{
129
130 // For very low energy, return initial track
131 G4double epmax = aTrack.GetTotalEnergy() - 0.5*proton_mass_c2;
132 if (epmax <= CutFixed) {
136 return &theParticleChange;
137 }
138
139 // Produce recoil muon and transferred photon
140 G4DynamicParticle* transferredPhoton = CalculateEMVertex(aTrack, targetNucleus);
141
142 // Interact the gamma with the nucleus
143 CalculateHadronicVertex(transferredPhoton, targetNucleus);
144 return &theParticleChange;
145}
146
148G4MuonVDNuclearModel::CalculateEMVertex(const G4HadProjectile& aTrack,
149 G4Nucleus& targetNucleus)
150{
151 // Select sampling table
152 G4double KineticEnergy = aTrack.GetKineticEnergy();
153 G4double TotalEnergy = aTrack.GetTotalEnergy();
155 G4Pow* g4calc = G4Pow::GetInstance();
156 G4double lnZ = g4calc->logZ(targetNucleus.GetZ_asInt());
157
158 G4double epmin = CutFixed;
159 G4double epmax = TotalEnergy - 0.5*proton_mass_c2;
160 G4double m0 = CutFixed;
161
162 G4double delmin = 1.e10;
163 G4double del;
164 G4int izz = 0;
165 G4int itt = 0;
166
167 for (G4int iz = 0; iz < nzdat; ++iz) {
168 del = std::abs(lnZ - g4calc->logZ(zdat[iz]));
169 if (del < delmin) {
170 delmin = del;
171 izz = iz;
172 }
173 }
174
175 delmin = 1.e10;
176 for (G4int it = 0; it < ntdat; ++it) {
177 del = std::abs(G4Log(KineticEnergy)-G4Log(tdat[it]) );
178 if (del < delmin) {
179 delmin = del;
180 itt = it;
181 }
182 }
183
184 // Sample the energy transfer according to the probability table
186
187 G4int iy;
188
189 G4int Z = zdat[izz];
190
191 for(iy = 0; iy<NBIN; ++iy) {
192
193 G4double pvv = fElementData->GetElement2DData(Z)->GetValue(iy, itt);
194 if(pvv >= r) { break; }
195 }
196
197 // Sampling is done uniformly in y in the bin
198 G4double pvx = fElementData->GetElement2DData(Z)->GetX(iy);
199 G4double pvx1 = fElementData->GetElement2DData(Z)->GetX(iy+1);
200 G4double y = pvx + G4UniformRand() * (pvx1 - pvx);
201
202 G4double x = G4Exp(y);
203 G4double ep = epmin*G4Exp(x*G4Log(epmax/epmin) );
204
205 // Sample scattering angle of mu, but first t should be sampled.
206 G4double yy = ep/TotalEnergy;
207 G4double tmin = Mass*Mass*yy*yy/(1.-yy);
208 G4double tmax = 2.*proton_mass_c2*ep;
209 G4double t1;
210 G4double t2;
211 if (m0 < ep) {
212 t1 = m0*m0;
213 t2 = ep*ep;
214 } else {
215 t1 = ep*ep;
216 t2 = m0*m0;
217 }
218
219 G4double w1 = tmax*t1;
220 G4double w2 = tmax+t1;
221 G4double w3 = tmax*(tmin+t1)/(tmin*w2);
222 G4double y1 = 1.-yy;
223 G4double y2 = 0.5*yy*yy;
224 G4double y3 = y1+y2;
225
226 G4double t = 0.0;
227 G4double rej = 0.0;
228
229 // Now sample t
230 G4int ntry = 0;
231 do
232 {
233 ntry += 1;
234 if (ntry > 10000) {
236 eda << " While count exceeded " << G4endl;
237 G4Exception("G4MuonVDNuclearModel::CalculateEMVertex()", "HAD_RPG_100", JustWarning, eda);
238 break;
239 }
240
241 t = w1/(w2*G4Exp(G4UniformRand()*G4Log(w3))-tmax);
242 rej = (1.-t/tmax)*(y1*(1.-tmin/t)+y2)/(y3*(1.-t/t2));
243 } while (G4UniformRand() > rej) ; /* Loop checking, 01.09.2015, D.Wright */
244
245 // compute angle from t
246 G4double sinth2 =
247 0.5*(t-tmin)/(2.*(TotalEnergy*(TotalEnergy-ep)-Mass*Mass)-tmin);
248 G4double theta = std::acos(1. - 2.*sinth2);
249
250 G4double phi = twopi*G4UniformRand();
251 G4double sinth = std::sin(theta);
252 G4double dirx = sinth*std::cos(phi);
253 G4double diry = sinth*std::sin(phi);
254 G4double dirz = std::cos(theta);
255 G4ThreeVector finalDirection(dirx,diry,dirz);
256 G4ThreeVector ParticleDirection(aTrack.Get4Momentum().vect().unit() );
257 finalDirection.rotateUz(ParticleDirection);
258
259 G4double NewKinEnergy = KineticEnergy - ep;
260 G4double finalMomentum = std::sqrt(NewKinEnergy*(NewKinEnergy+2.*Mass) );
261 G4double Ef = NewKinEnergy + Mass;
262 G4double initMomentum = std::sqrt(KineticEnergy*(TotalEnergy+Mass) );
263
264 // Set energy and direction of scattered primary in theParticleChange
267 theParticleChange.SetMomentumChange(finalDirection);
268
269 // Now create the emitted gamma
270 G4LorentzVector primaryMomentum(initMomentum*ParticleDirection, TotalEnergy);
271 G4LorentzVector fsMomentum(finalMomentum*finalDirection, Ef);
272 G4LorentzVector momentumTransfer = primaryMomentum - fsMomentum;
273
274 G4DynamicParticle* gamma =
275 new G4DynamicParticle(G4Gamma::Gamma(), momentumTransfer);
276
277 return gamma;
278}
279
280void
281G4MuonVDNuclearModel::CalculateHadronicVertex(G4DynamicParticle* incident,
282 G4Nucleus& target)
283{
284 G4HadFinalState* hfs = 0;
285 G4double gammaE = incident->GetTotalEnergy();
286
287 if (gammaE < 10*GeV) {
288 G4HadProjectile projectile(*incident);
289 hfs = bert->ApplyYourself(projectile, target);
290 } else {
291 // convert incident gamma to a pi0
293 G4double piKE = incident->GetTotalEnergy() - piMass;
294 G4double piMom = std::sqrt(piKE*(piKE + 2*piMass) );
295 G4ThreeVector piMomentum(incident->GetMomentumDirection() );
296 piMomentum *= piMom;
297 G4DynamicParticle theHadron(G4PionZero::PionZero(), piMomentum);
298 G4HadProjectile projectile(theHadron);
299 hfs = ftfp->ApplyYourself(projectile, target);
300 }
301
302 delete incident;
303
304 // Assign the creator model ID to the secondaries
305 for ( size_t i = 0; i < hfs->GetNumberOfSecondaries(); ++i ) {
306 hfs->GetSecondary( i )->SetCreatorModelID( secID );
307 }
308
309 // Copy secondaries from sub-model to model
311}
312
313
314void G4MuonVDNuclearModel::MakeSamplingTable()
315{
316 G4int nbin;
318 G4double TotalEnergy;
319 G4double Maxep;
320 G4double CrossSection;
321
322 G4double c;
323 G4double y;
324 G4double ymin,ymax;
325 G4double dy,yy;
326 G4double dx,x;
327 G4double ep;
328
329 G4double AtomicNumber;
330 G4double AtomicWeight;
331
333
334 for (G4int iz = 0; iz < nzdat; ++iz) {
335 AtomicNumber = zdat[iz];
336 AtomicWeight = adat[iz]*(g/mole);
337
338 G4Physics2DVector* pv = new G4Physics2DVector(NBIN+1,ntdat+1);
339 G4double pvv;
340
341 for (G4int it = 0; it < ntdat; ++it) {
342 KineticEnergy = tdat[it];
343 TotalEnergy = KineticEnergy + mumass;
344 Maxep = TotalEnergy - 0.5*proton_mass_c2;
345
346 CrossSection = 0.0;
347
348 // Calculate the differential cross section
349 // numerical integration in log .........
350 c = G4Log(Maxep/CutFixed);
351 ymin = -5.0;
352 ymax = 0.0;
353 dy = (ymax-ymin)/NBIN;
354
355 nbin=-1;
356
357 y = ymin - 0.5*dy;
358 yy = ymin - dy;
359 for (G4int i = 0; i < NBIN; ++i) {
360 y += dy;
361 x = G4Exp(y);
362 yy += dy;
363 dx = G4Exp(yy+dy)-G4Exp(yy);
364
365 ep = CutFixed*G4Exp(c*x);
366
367 CrossSection +=
368 ep*dx*muNucXS->ComputeDDMicroscopicCrossSection(KineticEnergy,
369 AtomicNumber,
370 AtomicWeight, ep);
371 if (nbin < NBIN) {
372 ++nbin;
373 pv->PutValue(nbin, it, CrossSection);
374 pv->PutX(nbin, y);
375 }
376 }
377 pv->PutX(NBIN, 0.);
378
379 if (CrossSection > 0.0) {
380 for (G4int ib = 0; ib <= nbin; ++ib) {
381 pvv = pv->GetValue(ib, it);
382 pvv = pvv/CrossSection;
383 pv->PutValue(ib, it, pvv);
384 }
385 }
386 } // loop on it
387
388 fElementData->InitialiseForElement(zdat[iz], pv);
389 } // loop on iz
390
391 // G4cout << " Kokoulin XS = "
392 // << muNucXS->ComputeDDMicroscopicCrossSection(1*GeV, 20.0,
393 // 40.0*g/mole, 0.3*GeV)/millibarn
394 // << G4endl;
395}
396
397void G4MuonVDNuclearModel::ModelDescription(std::ostream& outFile) const
398{
399 outFile << "G4MuonVDNuclearModel handles the inelastic scattering\n"
400 << "of mu- and mu+ from nuclei using the equivalent photon\n"
401 << "approximation in which the incoming lepton generates a\n"
402 << "virtual photon at the electromagnetic vertex, and the\n"
403 << "virtual photon is converted to a real photon. At low\n"
404 << "energies, the photon interacts directly with the nucleus\n"
405 << "using the Bertini cascade. At high energies the photon\n"
406 << "is converted to a pi0 which interacts using the FTFP\n"
407 << "model. The muon-nuclear cross sections of R. Kokoulin \n"
408 << "are used to generate the virtual photon spectrum\n";
409}
410
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
@ isAlive
G4double G4Log(G4double x)
Definition: G4Log.hh:227
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector vect() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
static G4CrossSectionDataSetRegistry * Instance()
const G4ThreeVector & GetMomentumDirection() const
G4double GetTotalEnergy() const
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
G4Physics2DVector * GetElement2DData(G4int Z)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondaries(const std::vector< G4HadSecondary > &addSecs)
std::size_t GetNumberOfSecondaries() const
void SetEnergyChange(G4double anEnergy)
G4HadSecondary * GetSecondary(size_t i)
void SetMomentumChange(const G4ThreeVector &aV)
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
void SetCreatorModelID(G4int id)
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
void SetMinEnergy(G4double anEnergy)
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static const char * Default_Name()
G4double ComputeDDMicroscopicCrossSection(G4double incidentKE, G4double Z, G4double A, G4double epsilon)
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
virtual void ModelDescription(std::ostream &outFile) const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
G4double GetValue(std::size_t idx, std::size_t idy) const
G4double GetX(std::size_t index) const
static G4int GetModelID(const G4int modelIndex)
static G4PionZero * PionZero()
Definition: G4PionZero.cc:107
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double logZ(G4int Z) const
Definition: G4Pow.hh:137
void SetTransport(G4VIntraNuclearTransportModel *const value)
void SetHighEnergyGenerator(G4VHighEnergyGenerator *const value)
G4HadFinalState * ApplyYourself(const G4HadProjectile &thePrimary, G4Nucleus &theNucleus) override
void SetDeExcitation(G4VPreCompoundModel *ptr)
void SetFragmentationModel(G4VStringFragmentation *aModel)
G4bool IsMasterThread()
Definition: G4Threading.cc:124