Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ElectroVDNuclearModel Class Reference

#include <G4ElectroVDNuclearModel.hh>

+ Inheritance diagram for G4ElectroVDNuclearModel:

Public Member Functions

 G4ElectroVDNuclearModel ()
 
 ~G4ElectroVDNuclearModel ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
virtual void ModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 49 of file G4ElectroVDNuclearModel.hh.

Constructor & Destructor Documentation

◆ G4ElectroVDNuclearModel()

G4ElectroVDNuclearModel::G4ElectroVDNuclearModel ( )

Definition at line 62 of file G4ElectroVDNuclearModel.cc.

63 : G4HadronicInteraction("G4ElectroVDNuclearModel"),
64 leptonKE(0.0), photonEnergy(0.0), photonQ2(0.0)
65{
66 SetMinEnergy(0.0);
67 SetMaxEnergy(1*PeV);
68 electroXS =
70 GetCrossSectionDataSet(G4ElectroNuclearCrossSection::Default_Name());
71 gammaXS =
73 GetCrossSectionDataSet(G4PhotoNuclearCrossSection::Default_Name());
74
75 // reuse existing pre-compound model
76 G4GeneratorPrecompoundInterface* precoInterface
80 G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
81 if(!pre) { pre = new G4PreCompoundModel(); }
82 precoInterface->SetDeExcitation(pre);
83
84 // string model
85 ftfp = new G4TheoFSGenerator();
86 ftfp->SetTransport(precoInterface);
87 theFragmentation = new G4LundStringFragmentation();
88 theStringDecay = new G4ExcitedStringDecay(theFragmentation);
89 G4FTFModel* theStringModel = new G4FTFModel();
90 theStringModel->SetFragmentationModel(theStringDecay);
91 ftfp->SetHighEnergyGenerator(theStringModel);
92
93 // Build Bertini model
94 bert = new G4CascadeInterface();
95}
static G4CrossSectionDataSetRegistry * Instance()
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
void SetTransport(G4VIntraNuclearTransportModel *const value)
void SetHighEnergyGenerator(G4VHighEnergyGenerator *const value)
void SetDeExcitation(G4VPreCompoundModel *ptr)
void SetFragmentationModel(G4VStringFragmentation *aModel)

◆ ~G4ElectroVDNuclearModel()

G4ElectroVDNuclearModel::~G4ElectroVDNuclearModel ( )

Definition at line 97 of file G4ElectroVDNuclearModel.cc.

98{
99 delete theFragmentation;
100 delete theStringDecay;
101}

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4ElectroVDNuclearModel::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus aTargetNucleus 
)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 119 of file G4ElectroVDNuclearModel.cc.

121{
122 // Set up default particle change (just returns initial state)
125 leptonKE = aTrack.GetKineticEnergy();
128
129 // Set up sanity checks for real photon production
130 G4DynamicParticle lepton(aTrack.GetDefinition(), aTrack.Get4Momentum() );
131
132 // Need to call GetElementCrossSection before calling GetEquivalentPhotonEnergy.
133 G4Material* mat = 0;
134 G4int targZ = targetNucleus.GetZ_asInt();
135 electroXS->GetElementCrossSection(&lepton, targZ, mat);
136
137 photonEnergy = electroXS->GetEquivalentPhotonEnergy();
138 // Photon energy cannot exceed lepton energy
139 if (photonEnergy < leptonKE) {
140 photonQ2 = electroXS->GetEquivalentPhotonQ2(photonEnergy);
142 // Photon
143 if (photonEnergy > photonQ2/dM) {
144 // Produce recoil lepton and transferred photon
145 G4DynamicParticle* transferredPhoton = CalculateEMVertex(aTrack, targetNucleus);
146 // Interact gamma with nucleus
147 if (transferredPhoton) CalculateHadronicVertex(transferredPhoton, targetNucleus);
148 }
149 }
150 return &theParticleChange;
151}
@ isAlive
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
Hep3Vector unit() const
Hep3Vector vect() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat)
void SetStatusChange(G4HadFinalStateStatus aS)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4Proton * Proton()
Definition: G4Proton.cc:92

◆ ModelDescription()

void G4ElectroVDNuclearModel::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 103 of file G4ElectroVDNuclearModel.cc.

104{
105 outFile << "G4ElectroVDNuclearModel handles the inelastic scattering\n"
106 << "of e- and e+ from nuclei using the equivalent photon\n"
107 << "approximation in which the incoming lepton generates a\n"
108 << "virtual photon at the electromagnetic vertex, and the\n"
109 << "virtual photon is converted to a real photon. At low\n"
110 << "energies, the photon interacts directly with the nucleus\n"
111 << "using the Bertini cascade. At high energies the photon\n"
112 << "is converted to a pi0 which interacts using the FTFP\n"
113 << "model. The electro- and gamma-nuclear cross sections of\n"
114 << "M. Kossov are used to generate the virtual photon spectrum\n";
115}

The documentation for this class was generated from the following files: