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

#include <G4LowEGammaNuclearModel.hh>

+ Inheritance diagram for G4LowEGammaNuclearModel:

Public Member Functions

 G4LowEGammaNuclearModel ()
 
 ~G4LowEGammaNuclearModel () override
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) final
 
void InitialiseModel () final
 
- 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 46 of file G4LowEGammaNuclearModel.hh.

Constructor & Destructor Documentation

◆ G4LowEGammaNuclearModel()

G4LowEGammaNuclearModel::G4LowEGammaNuclearModel ( )
explicit

Definition at line 46 of file G4LowEGammaNuclearModel.cc.

47 : G4HadronicInteraction("GammaNPreco"),lab4mom(0.,0.,0.,0.)
48{
49 SetMinEnergy( 0.0*CLHEP::GeV );
51
52 // reuse existing pre-compound model
55 fPreco = static_cast<G4PreCompoundModel*>(p);
56 if(!fPreco) { fPreco = new G4PreCompoundModel(); }
57}
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()

◆ ~G4LowEGammaNuclearModel()

G4LowEGammaNuclearModel::~G4LowEGammaNuclearModel ( )
override

Definition at line 59 of file G4LowEGammaNuclearModel.cc.

60{}

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4LowEGammaNuclearModel::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
finalvirtual

Reimplemented from G4HadronicInteraction.

Definition at line 65 of file G4LowEGammaNuclearModel.cc.

67{
69
70 G4int A = theNucleus.GetA_asInt();
71 G4int Z = theNucleus.GetZ_asInt();
72
73 // Create initial state
74 lab4mom.set(0.,0.,0.,G4NucleiProperties::GetNuclearMass(A, Z));
75 lab4mom += aTrack.Get4Momentum();
76
77 G4Fragment frag(A, Z, lab4mom);
78
79 if (verboseLevel > 1) {
80 G4cout << "G4LowEGammaNuclearModel::ApplyYourself initial G4Fragmet:"
81 << G4endl;
82 G4cout << frag << G4endl;
83 }
84 G4ReactionProductVector* res = fPreco->DeExcite(frag);
85
86 // secondaries produced
87 if(res) {
88
90 G4int nsec = res->size();
91 if (verboseLevel > 1) {
92 G4cout << "G4LowEGammaNuclearModel: " << nsec << " secondaries" << G4endl;
93 }
94 for(G4int i=0; i<nsec; ++i) {
95 if((*res)[i]) {
96 G4double ekin = (*res)[i]->GetKineticEnergy();
97 G4ThreeVector dir(0.,0.,1.);
98 if(ekin > 0.0) { dir = (*res)[i]->GetMomentum().unit(); }
100 new G4DynamicParticle((*res)[i]->GetDefinition(), dir, ekin));
101 news->SetTime((*res)[i]->GetTOF());
102 news->SetCreatorModelType((*res)[i]->GetCreatorModel());
104 if (verboseLevel > 1) {
105 G4cout << i << ". " << (*res)[i]->GetDefinition()->GetParticleName()
106 << " Ekin(MeV)= " << ekin/MeV
107 << " dir: " << dir << G4endl;
108 }
109 delete (*res)[i];
110 delete news;
111 }
112 }
113 delete res;
114 }
115 return &theParticleChange;
116}
double A(double temperature)
@ stopAndKill
std::vector< G4ReactionProduct * > G4ReactionProductVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void set(double x, double y, double z, double t)
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
const G4LorentzVector & Get4Momentum() const
void SetTime(G4double aT)
void SetCreatorModelType(G4int idx)
static G4double GetNuclearMass(const G4double A, const G4double Z)
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment) final

◆ InitialiseModel()

void G4LowEGammaNuclearModel::InitialiseModel ( )
finalvirtual

Reimplemented from G4HadronicInteraction.

Definition at line 62 of file G4LowEGammaNuclearModel.cc.

63{}

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