Geant4 11.1.1
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 ()
 
virtual ~G4MuonVDNuclearModel ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
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 G4MuonVDNuclearModel.hh.

Constructor & Destructor Documentation

◆ G4MuonVDNuclearModel()

G4MuonVDNuclearModel::G4MuonVDNuclearModel ( )
explicit

Definition at line 72 of file G4MuonVDNuclearModel.cc.

73 : G4HadronicInteraction("G4MuonVDNuclearModel"),isMaster(false)
74{
76 GetCrossSectionDataSet(G4KokoulinMuonNuclearXS::Default_Name());
77
78 SetMinEnergy(0.0);
79 SetMaxEnergy(1*CLHEP::PeV);
80 CutFixed = 0.2*CLHEP::GeV;
81
82 if(!fElementData && G4Threading::IsMasterThread()) {
83 fElementData = new G4ElementData();
84 MakeSamplingTable();
85 isMaster = true;
86 }
87
88 // reuse existing pre-compound model
89 G4GeneratorPrecompoundInterface* precoInterface
93 G4VPreCompoundModel* pre = static_cast<G4VPreCompoundModel*>(p);
94 if(!pre) { pre = new G4PreCompoundModel(); }
95 precoInterface->SetDeExcitation(pre);
96
97 // Build FTFP model
98 ftfp = new G4TheoFSGenerator();
99 ftfp->SetTransport(precoInterface);
100 theFragmentation = new G4LundStringFragmentation();
101 theStringDecay = new G4ExcitedStringDecay(theFragmentation);
102 G4FTFModel* theStringModel = new G4FTFModel;
103 theStringModel->SetFragmentationModel(theStringDecay);
104 ftfp->SetHighEnergyGenerator(theStringModel);
105
106 // Build Bertini cascade
107 bert = new G4CascadeInterface();
108
109 // Creator model ID
110 secID = G4PhysicsModelCatalog::GetModelID( "model_" + GetModelName() );
111}
static G4CrossSectionDataSetRegistry * Instance()
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
void SetMinEnergy(G4double anEnergy)
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static const char * Default_Name()
static G4int GetModelID(const G4int modelIndex)
void SetTransport(G4VIntraNuclearTransportModel *const value)
void SetHighEnergyGenerator(G4VHighEnergyGenerator *const value)
void SetDeExcitation(G4VPreCompoundModel *ptr)
void SetFragmentationModel(G4VStringFragmentation *aModel)
G4bool IsMasterThread()
Definition: G4Threading.cc:124

◆ ~G4MuonVDNuclearModel()

G4MuonVDNuclearModel::~G4MuonVDNuclearModel ( )
virtual

Definition at line 113 of file G4MuonVDNuclearModel.cc.

114{
115 delete theFragmentation;
116 delete theStringDecay;
117
118 if(isMaster) {
119 delete fElementData;
120 fElementData = nullptr;
121 }
122}

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 125 of file G4MuonVDNuclearModel.cc.

127{
129
130 // For very low energy, return initial track
131 G4double epmax = aTrack.GetTotalEnergy() - 0.5*proton_mass_c2;
132 if (epmax <= CutFixed) {
136 return &theParticleChange;
137 }
138
139 // Produce recoil muon and transferred photon
140 G4DynamicParticle* transferredPhoton = CalculateEMVertex(aTrack, targetNucleus);
141
142 // Interact the gamma with the nucleus
143 CalculateHadronicVertex(transferredPhoton, targetNucleus);
144 return &theParticleChange;
145}
@ isAlive
double G4double
Definition: G4Types.hh:83
Hep3Vector unit() const
Hep3Vector vect() const
void SetStatusChange(G4HadFinalStateStatus aS)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 397 of file G4MuonVDNuclearModel.cc.

398{
399 outFile << "G4MuonVDNuclearModel handles the inelastic scattering\n"
400 << "of mu- and mu+ from nuclei using the equivalent photon\n"
401 << "approximation in which the incoming lepton generates a\n"
402 << "virtual photon at the electromagnetic vertex, and the\n"
403 << "virtual photon is converted to a real photon. At low\n"
404 << "energies, the photon interacts directly with the nucleus\n"
405 << "using the Bertini cascade. At high energies the photon\n"
406 << "is converted to a pi0 which interacts using the FTFP\n"
407 << "model. The muon-nuclear cross sections of R. Kokoulin \n"
408 << "are used to generate the virtual photon spectrum\n";
409}

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