Geant4-11
G4ElectroVDNuclearModel.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: 1 May 2012
29//
30// Description: model for electron and positron interaction with nuclei
31// using the equivalent photon spectrum. A real gamma is
32// produced from the virtual photon spectrum and is then
33// interacted hadronically by the Bertini cascade at low
34// energies. At high energies the gamma is treated as a
35// pi0 and interacted with the nucleus using the FTFP model.
36// The electro- and photo-nuclear cross sections of
37// M. Kossov are used to generate the virtual photon
38// spectrum.
39//
40
42
44#include "G4SystemOfUnits.hh"
45
49
50#include "G4CascadeInterface.hh"
51#include "G4TheoFSGenerator.hh"
54#include "G4PreCompoundModel.hh"
57#include "G4FTFModel.hh"
58
59#include "G4HadFinalState.hh"
62
64 : G4HadronicInteraction("G4ElectroVDNuclearModel"),
65 leptonKE(0.0), photonEnergy(0.0), photonQ2(0.0), secID(-1)
66{
67 SetMinEnergy(0.0);
69 electroXS =
71 GetCrossSectionDataSet(G4ElectroNuclearCrossSection::Default_Name());
72 gammaXS =
74 GetCrossSectionDataSet(G4PhotoNuclearCrossSection::Default_Name());
75
76 // reuse existing pre-compound model
77 G4GeneratorPrecompoundInterface* precoInterface
81 G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
82 if(!pre) { pre = new G4PreCompoundModel(); }
83 precoInterface->SetDeExcitation(pre);
84
85 // string model
86 ftfp = new G4TheoFSGenerator();
87 ftfp->SetTransport(precoInterface);
90 G4FTFModel* theStringModel = new G4FTFModel();
92 ftfp->SetHighEnergyGenerator(theStringModel);
93
94 // Build Bertini model
96
97 // Creator model ID
99}
100
102{
103 delete theFragmentation;
104 delete theStringDecay;
105}
106
107void G4ElectroVDNuclearModel::ModelDescription(std::ostream& outFile) const
108{
109 outFile << "G4ElectroVDNuclearModel handles the inelastic scattering\n"
110 << "of e- and e+ from nuclei using the equivalent photon\n"
111 << "approximation in which the incoming lepton generates a\n"
112 << "virtual photon at the electromagnetic vertex, and the\n"
113 << "virtual photon is converted to a real photon. At low\n"
114 << "energies, the photon interacts directly with the nucleus\n"
115 << "using the Bertini cascade. At high energies the photon\n"
116 << "is converted to a pi0 which interacts using the FTFP\n"
117 << "model. The electro- and gamma-nuclear cross sections of\n"
118 << "M. Kossov are used to generate the virtual photon spectrum\n";
119}
120
121
124 G4Nucleus& targetNucleus)
125{
126 // Set up default particle change (just returns initial state)
129 leptonKE = aTrack.GetKineticEnergy();
132
133 // Set up sanity checks for real photon production
134 G4DynamicParticle lepton(aTrack.GetDefinition(), aTrack.Get4Momentum() );
135
136 // Need to call GetElementCrossSection before calling GetEquivalentPhotonEnergy.
137 G4Material* mat = 0;
138 G4int targZ = targetNucleus.GetZ_asInt();
139 electroXS->GetElementCrossSection(&lepton, targZ, mat);
140
142 // Photon energy cannot exceed lepton energy
143 if (photonEnergy < leptonKE) {
146 // Photon
147 if (photonEnergy > photonQ2/dM) {
148 // Produce recoil lepton and transferred photon
149 G4DynamicParticle* transferredPhoton = CalculateEMVertex(aTrack, targetNucleus);
150 // Interact gamma with nucleus
151 if (transferredPhoton) CalculateHadronicVertex(transferredPhoton, targetNucleus);
152 }
153 }
154 return &theParticleChange;
155}
156
157
160 G4Nucleus& targetNucleus)
161{
163 G4ThreeVector(0.,0.,1.) );
164
165 // Get gamma cross section at Q**2 = 0 (real gamma)
166 G4int targZ = targetNucleus.GetZ_asInt();
167 G4Material* mat = 0;
168 G4double sigNu =
170
171 // Change real gamma energy to equivalent energy and get cross section at that energy
173 photon.SetKineticEnergy(photonEnergy - photonQ2/dM);
174 G4double sigK =
177
178 // No gamma produced, return null ptr
179 if (sigNu*G4UniformRand() > sigK*rndFraction) return 0;
180
181 // Scatter the lepton
182 G4double mProj = aTrack.GetDefinition()->GetPDGMass();
183 G4double mProj2 = mProj*mProj;
184 G4double iniE = leptonKE + mProj; // Total energy of incident lepton
185 G4double finE = iniE - photonEnergy; // Total energy of scattered lepton
187 G4double iniP = std::sqrt(iniE*iniE-mProj2); // Incident lepton momentum
188 G4double finP = std::sqrt(finE*finE-mProj2); // Scattered lepton momentum
189 G4double cost = (iniE*finE - mProj2 - photonQ2/2.)/iniP/finP; // cos(theta) from Q**2
190 if (cost > 1.) cost= 1.;
191 if (cost < -1.) cost=-1.;
192 G4double sint = std::sqrt(1.-cost*cost);
193
194 G4ThreeVector dir = aTrack.Get4Momentum().vect().unit();
195 G4ThreeVector ortx = dir.orthogonal().unit(); // Ortho-normal to scattering plane
196 G4ThreeVector orty = dir.cross(ortx); // Third unit vector
198 G4double sinx = sint*std::sin(phi);
199 G4double siny = sint*std::cos(phi);
200 G4ThreeVector findir = cost*dir+sinx*ortx+siny*orty;
201 theParticleChange.SetMomentumChange(findir); // change lepton direction
202
203 // Create a gamma with momentum equal to momentum transfer
204 G4ThreeVector photonMomentum = iniP*dir - finP*findir;
206 photonEnergy, photonMomentum);
207 return gamma;
208}
209
210
211void
213 G4Nucleus& target)
214{
215 G4HadFinalState* hfs = 0;
216 G4double gammaE = incident->GetTotalEnergy();
217
218 if (gammaE < 10*GeV) {
219 G4HadProjectile projectile(*incident);
220 hfs = bert->ApplyYourself(projectile, target);
221 } else {
222 // At high energies convert incident gamma to a pion
224 G4double piMom = std::sqrt(gammaE*gammaE - piMass*piMass);
225 G4ThreeVector piMomentum(incident->GetMomentumDirection() );
226 piMomentum *= piMom;
227 G4DynamicParticle theHadron(G4PionZero::PionZero(), piMomentum);
228 G4HadProjectile projectile(theHadron);
229 hfs = ftfp->ApplyYourself(projectile, target);
230 }
231
232 delete incident;
233
234 // Assign the creator model ID to the secondaries
235 for ( size_t i = 0; i < hfs->GetNumberOfSecondaries(); ++i ) {
237 }
238
239 // Copy secondaries from sub-model to model
241}
242
static const G4double dM
@ isAlive
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double PeV
Definition: G4SIunits.hh:205
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
Hep3Vector vect() const
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &theNucleus)
static G4CrossSectionDataSetRegistry * Instance()
const G4ThreeVector & GetMomentumDirection() const
G4double GetTotalEnergy() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat)
G4double GetVirtualFactor(G4double nu, G4double Q2)
void CalculateHadronicVertex(G4DynamicParticle *incident, G4Nucleus &target)
virtual void ModelDescription(std::ostream &outFile) const
G4ExcitedStringDecay * theStringDecay
G4LundStringFragmentation * theFragmentation
G4PhotoNuclearCrossSection * gammaXS
G4DynamicParticle * CalculateEMVertex(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
G4ElectroNuclearCrossSection * electroXS
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
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)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() 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 G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *)
static G4int GetModelID(const G4int modelIndex)
static G4PionZero * PionZero()
Definition: G4PionZero.cc:107
static G4Proton * Proton()
Definition: G4Proton.cc:92
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)