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

#include <G4MuonVDNuclearModel.hh>

Inheritance diagram for G4MuonVDNuclearModel:

Public Member Functions

 G4MuonVDNuclearModel ()
 ~G4MuonVDNuclearModel () override
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
void ModelDescription (std::ostream &outFile) const override
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 &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 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 G4MuonVDNuclearModel.hh.

Constructor & Destructor Documentation

◆ G4MuonVDNuclearModel()

G4MuonVDNuclearModel::G4MuonVDNuclearModel ( )

Definition at line 72 of file G4MuonVDNuclearModel.cc.

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();
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}
static G4CrossSectionDataSetRegistry * Instance()
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static const char * Default_Name()
static G4int GetModelID(const G4int modelIndex)
void SetDeExcitation(G4VPreCompoundModel *ptr)
void SetFragmentationModel(G4VStringFragmentation *aModel)

◆ ~G4MuonVDNuclearModel()

G4MuonVDNuclearModel::~G4MuonVDNuclearModel ( )
override

Definition at line 112 of file G4MuonVDNuclearModel.cc.

113{
114 delete theFragmentation;
115 delete theStringDecay;
116}

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4MuonVDNuclearModel::ApplyYourself ( const G4HadProjectile & aTrack,
G4Nucleus & targetNucleus )
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 119 of file G4MuonVDNuclearModel.cc.

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}
@ isAlive
double G4double
Definition G4Types.hh:83
Hep3Vector unit() const
Hep3Vector vect() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const

◆ ModelDescription()

void G4MuonVDNuclearModel::ModelDescription ( std::ostream & outFile) const
overridevirtual

Reimplemented from G4HadronicInteraction.

Definition at line 391 of file G4MuonVDNuclearModel.cc.

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}

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