Geant4 9.6.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 ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
 
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
 
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
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) 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
 

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 54 of file G4MuonVDNuclearModel.hh.

Constructor & Destructor Documentation

◆ G4MuonVDNuclearModel()

G4MuonVDNuclearModel::G4MuonVDNuclearModel ( )

Definition at line 52 of file G4MuonVDNuclearModel.cc.

53 : G4HadronicInteraction("G4MuonVDNuclearModel")
54{
55 SetMinEnergy(0.0);
56 SetMaxEnergy(1*PeV);
57 CutFixed = 0.2*GeV;
58 NBIN = 1000;
59
60 for (G4int k = 0; k < 5; k++) {
61 for (G4int j = 0; j < 8; j++) {
62 for (G4int i = 0; i < 1001; i++) {
63 proba[k][j][i] = 0.0;
64 ya[i] = 0.0;
65 }
66 }
67 }
68
69 MakeSamplingTable();
70
71 // Build FTFP model
72 ftfp = new G4TheoFSGenerator();
73 precoInterface = new G4GeneratorPrecompoundInterface();
74 theHandler = new G4ExcitationHandler();
75 preEquilib = new G4PreCompoundModel(theHandler);
76 precoInterface->SetDeExcitation(preEquilib);
77 ftfp->SetTransport(precoInterface);
78 theFragmentation = new G4LundStringFragmentation();
79 theStringDecay = new G4ExcitedStringDecay(theFragmentation);
80 theStringModel = new G4FTFModel;
81 theStringModel->SetFragmentationModel(theStringDecay);
82 ftfp->SetHighEnergyGenerator(theStringModel);
83
84 // Build Bertini cascade
85 bert = new G4CascadeInterface();
86}
int G4int
Definition: G4Types.hh:66
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)

◆ ~G4MuonVDNuclearModel()

G4MuonVDNuclearModel::~G4MuonVDNuclearModel ( )

Definition at line 89 of file G4MuonVDNuclearModel.cc.

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

Member Function Documentation

◆ ApplyYourself()

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

Implements G4HadronicInteraction.

Definition at line 101 of file G4MuonVDNuclearModel.cc.

103{
105
106 // For very low energy, return initial track
107 G4double epmax = aTrack.GetTotalEnergy() - 0.5*proton_mass_c2;
108 if (epmax <= CutFixed) {
112 return &theParticleChange;
113 }
114
115 // Produce recoil muon and transferred photon
116 G4DynamicParticle* transferredPhoton = CalculateEMVertex(aTrack, targetNucleus);
117
118 // Interact the gamma with the nucleus
119 CalculateHadronicVertex(transferredPhoton, targetNucleus);
120 return &theParticleChange;
121}
@ isAlive
double G4double
Definition: G4Types.hh:64
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

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