Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/lepto_nuclear/src/G4MuonVDNuclearModel.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 ]

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