Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4ElectroVDNuclearModel Class Reference

#include <G4ElectroVDNuclearModel.hh>

Inheritance diagram for G4ElectroVDNuclearModel:
G4HadronicInteraction

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 G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
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)
 
const G4HadronicInteractionGetMyPointer () const
 
virtual G4int GetVerboseLevel () const
 
virtual 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
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair
< G4double, G4double
GetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 

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 56 of file G4ElectroVDNuclearModel.hh.

Constructor & Destructor Documentation

G4ElectroVDNuclearModel::G4ElectroVDNuclearModel ( )

Definition at line 63 of file G4ElectroVDNuclearModel.cc.

References G4PhotoNuclearCrossSection::Default_Name(), G4ElectroNuclearCrossSection::Default_Name(), G4CrossSectionDataSetRegistry::Instance(), python.hepunit::PeV, G4VIntraNuclearTransportModel::SetDeExcitation(), G4VPartonStringModel::SetFragmentationModel(), G4TheoFSGenerator::SetHighEnergyGenerator(), G4HadronicInteraction::SetMaxEnergy(), G4HadronicInteraction::SetMinEnergy(), and G4TheoFSGenerator::SetTransport().

64  : G4HadronicInteraction("G4ElectroVDNuclearModel"),
65  leptonKE(0.0), photonEnergy(0.0), photonQ2(0.0)
66 {
67  SetMinEnergy(0.0);
68  SetMaxEnergy(1*PeV);
69 
71  GetCrossSectionDataSet(G4ElectroNuclearCrossSection::Default_Name());
73  GetCrossSectionDataSet(G4PhotoNuclearCrossSection::Default_Name());
74  ftfp = new G4TheoFSGenerator();
75  precoInterface = new G4GeneratorPrecompoundInterface();
76  theHandler = new G4ExcitationHandler();
77  preEquilib = new G4PreCompoundModel(theHandler);
78  precoInterface->SetDeExcitation(preEquilib);
79  ftfp->SetTransport(precoInterface);
80  theFragmentation = new G4LundStringFragmentation();
81  theStringDecay = new G4ExcitedStringDecay(theFragmentation);
82  theStringModel = new G4FTFModel();
83  theStringModel->SetFragmentationModel(theStringDecay);
84  ftfp->SetHighEnergyGenerator(theStringModel);
85 
86  // Build Bertini model
87  bert = new G4CascadeInterface();
88 }
void SetFragmentationModel(G4VStringFragmentation *aModel)
void SetHighEnergyGenerator(G4VHighEnergyGenerator *const value)
void SetMinEnergy(G4double anEnergy)
static G4CrossSectionDataSetRegistry * Instance()
G4HadronicInteraction(const G4String &modelName="HadronicModel")
void SetMaxEnergy(const G4double anEnergy)
void SetDeExcitation(G4VPreCompoundModel *ptr)
void SetTransport(G4VIntraNuclearTransportModel *const value)
G4ElectroVDNuclearModel::~G4ElectroVDNuclearModel ( )

Definition at line 91 of file G4ElectroVDNuclearModel.cc.

92 {
93  delete ftfp;
94  delete preEquilib;
95  delete theFragmentation;
96  delete theStringDecay;
97  delete theStringModel;
98  delete bert;
99 }

Member Function Documentation

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

Implements G4HadronicInteraction.

Definition at line 118 of file G4ElectroVDNuclearModel.cc.

References G4HadFinalState::Clear(), G4HadProjectile::Get4Momentum(), G4HadProjectile::GetDefinition(), G4ElectroNuclearCrossSection::GetElementCrossSection(), G4ElectroNuclearCrossSection::GetEquivalentPhotonEnergy(), G4ElectroNuclearCrossSection::GetEquivalentPhotonQ2(), G4HadProjectile::GetKineticEnergy(), G4ParticleDefinition::GetPDGMass(), G4Nucleus::GetZ_asInt(), isAlive, G4Neutron::Neutron(), G4Proton::Proton(), G4HadFinalState::SetEnergyChange(), G4HadFinalState::SetMomentumChange(), G4HadFinalState::SetStatusChange(), G4HadronicInteraction::theParticleChange, CLHEP::Hep3Vector::unit(), and CLHEP::HepLorentzVector::vect().

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

Reimplemented from G4HadronicInteraction.

Definition at line 102 of file G4ElectroVDNuclearModel.cc.

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

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