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

#include <G4LENDCapture.hh>

+ Inheritance diagram for G4LENDCapture:

Public Member Functions

 G4LENDCapture (G4ParticleDefinition *pd)
 
 ~G4LENDCapture ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
- Public Member Functions inherited from G4LENDModel
 G4LENDModel (G4String name="LENDModel")
 
 ~G4LENDModel ()
 
void ChangeDefaultEvaluation (G4String name)
 
void AllowNaturalAbundanceTarget ()
 
void AllowAnyCandidateTarget ()
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DumpLENDTargetInfo (G4bool force=false)
 
- 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 ModelDescription (std::ostream &outFile) const
 
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 G4LENDModel
void create_used_target_map ()
 
void recreate_used_target_map ()
 
G4GIDI_targetget_target_from_map (G4int nuclear_code)
 
G4HadFinalStatereturnUnchanged (const G4HadProjectile &aTrack, G4HadFinalState *theResult)
 
- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4LENDModel
G4ParticleDefinitionproj
 
G4LENDManagerlend_manager
 
std::map< G4int, G4LENDUsedTarget * > usedTarget_map
 
G4int secID
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 44 of file G4LENDCapture.hh.

Constructor & Destructor Documentation

◆ G4LENDCapture()

G4LENDCapture::G4LENDCapture ( G4ParticleDefinition * pd)
inline

Definition at line 49 of file G4LENDCapture.hh.

50 :G4LENDModel( "LENDCapture" )
51 {
52 proj = pd;
53
54// theModelName = "LENDCapture for ";
55// theModelName += proj->GetParticleName();
57 };
G4LENDModel(G4String name="LENDModel")
void create_used_target_map()
G4ParticleDefinition * proj

◆ ~G4LENDCapture()

G4LENDCapture::~G4LENDCapture ( )
inline

Definition at line 59 of file G4LENDCapture.hh.

59{;};

Member Function Documentation

◆ ApplyYourself()

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

Reimplemented from G4LENDModel.

Definition at line 35 of file G4LENDCapture.cc.

36{
37
38 G4double temp = aTrack.GetMaterial()->GetTemperature();
39
40 //G4int iZ = int ( aTarg.GetZ() );
41 //G4int iA = int ( aTarg.GetN() );
42 //migrate to integer A and Z (GetN_asInt returns number of neutrons in the nucleus since this)
43 G4int iZ = aTarg.GetZ_asInt();
44 G4int iA = aTarg.GetA_asInt();
45 G4int iM = 0;
46 if ( aTarg.GetIsotope() != NULL ) {
47 iM = aTarg.GetIsotope()->Getm();
48 }
49
50 G4double ke = aTrack.GetKineticEnergy();
51
53 theResult->Clear();
54
56 if ( aTarget == NULL ) return returnUnchanged( aTrack , theResult );
57 std::vector<G4GIDI_Product>* products = aTarget->getCaptureFinalState( ke*MeV, temp, MyRNG, NULL );
58
59 G4int ipZ = aTrack.GetDefinition()->GetAtomicNumber();
60 G4int ipA = aTrack.GetDefinition()->GetAtomicMass();
61
62 G4bool needResidual=true;
63
64 G4ThreeVector p(0,0,0);
65 if ( products != NULL )
66 {
67
68 G4int totN = 0;
69
70 for ( G4int j = 0; j < int( products->size() ); j++ )
71 {
72 G4int jZ = (*products)[j].Z;
73 G4int jA = (*products)[j].A;
74
75 //G4cout << "ZA = " << 1000 * (*products)[j].Z + (*products)[j].A << " EK = "
76 // << (*products)[j].kineticEnergy
77 // << " px " << (*products)[j].px
78 // << " py " << (*products)[j].py
79 // << " pz " << (*products)[j].pz
80 // << G4endl;
81
82 if ( jZ == iZ + ipZ && jA == iA + ipA ) needResidual = false;
83
84 G4ThreeVector dp((*products)[j].px,(*products)[j].py,(*products)[j].pz);
85 p += dp;
86
88
89 if ( jA == 1 && jZ == 1 ) {
91 totN += 1;
92 }
93 else if ( jA == 1 && jZ == 0 )
94 {
96 totN += 1;
97 }
98 else if ( jZ > 0 ) {
99 if ( jA != 0 )
100 {
101 theSec->SetDefinition( G4IonTable::GetIonTable()->GetIon( jZ , jA , iM ) );
102 totN += jA;
103 }
104 else
105 {
106 theSec->SetDefinition( G4IonTable::GetIonTable()->GetIon( jZ , iA+1-totN , iM ) );
107 }
108 }
109 else {
110 theSec->SetDefinition( G4Gamma::Gamma() );
111 }
112
113 theSec->SetMomentum( G4ThreeVector( (*products)[j].px*MeV , (*products)[j].py*MeV , (*products)[j].pz*MeV ) );
114
115/*
116 if ( dp.mag() == 0 )
117 {
118 //theSec->SetMomentum( -p*MeV );
119 }
120*/
121
122 theResult->AddSecondary( theSec, secID );
123 }
124 }
125 else
126 {
127
128 //For the case data does not provide final states
129
130 //G4cout << "products != NULL; iZ = " << iZ << ", iA = " << iA << G4endl;
131
132 // TK comment
133 // aTarg->ReturnTargetParticle()->Get4Momentum has trouble, thus we use following
134 G4Fragment nucleus( iA + ipA , iZ + ipZ , aTrack.Get4Momentum() + G4LorentzVector( G4ThreeVector(0,0,0) , G4IonTable::GetIonTable()->GetIon( iZ + ipZ , iA )->GetPDGMass() ) );
135 G4PhotonEvaporation photonEvaporation;
136 photonEvaporation.SetICM( TRUE );
137 G4FragmentVector* products_from_PE = photonEvaporation.BreakItUp(nucleus);
138 G4FragmentVector::iterator it;
139
140 for (it = products_from_PE->begin(); it != products_from_PE->end(); it++)
141 {
142 if ( (*it)->GetZ_asInt() == iZ + ipZ && (*it)->GetA_asInt() == iA + ipA ) needResidual = false;
144 if ( (*it)->GetParticleDefinition() != NULL ) {
145 //G4cout << (*it)->GetParticleDefinition()->GetParticleName() << G4endl;
146 theSec->SetDefinition( (*it)->GetParticleDefinition() );
147 theSec->Set4Momentum( (*it)->GetMomentum() );
148 } else {
149 //G4cout << (*it)->GetZ_asInt() << " " << (*it)->GetA_asInt() << G4endl;
150 theSec->SetDefinition( G4IonTable::GetIonTable()->GetIon( (*it)->GetZ_asInt() , (*it)->GetA_asInt() ) );
151 theSec->Set4Momentum( (*it)->GetMomentum() );
152 }
153 theResult->AddSecondary( theSec, secID );
154 }
155 delete products_from_PE;
156 }
157
158 //if necessary, generate residual nucleus
159 if ( needResidual ) {
160 G4DynamicParticle* residual = new G4DynamicParticle;
161 residual->SetDefinition( G4IonTable::GetIonTable()->GetIon( iZ + ipZ , iA + ipA ) );
162 residual->SetMomentum( -p*MeV );
163 theResult->AddSecondary( residual, secID );
164 }
165
166 delete products;
167
168 theResult->SetStatusChange( stopAndKill );
169
170 return theResult;
171
172}
std::vector< G4Fragment * > G4FragmentVector
Definition G4Fragment.hh:65
@ stopAndKill
double MyRNG(void *)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void Set4Momentum(const G4LorentzVector &momentum)
void SetMomentum(const G4ThreeVector &momentum)
std::vector< G4GIDI_Product > * getCaptureFinalState(double e_in, double temperature, double(*rng)(void *), void *rngState)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
void SetStatusChange(G4HadFinalStateStatus aS)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
static G4IonTable * GetIonTable()
G4int GetNucleusEncoding(G4int iZ, G4int iA, G4int iM)
G4LENDManager * lend_manager
G4HadFinalState * returnUnchanged(const G4HadProjectile &aTrack, G4HadFinalState *theResult)
G4GIDI_target * get_target_from_map(G4int nuclear_code)
G4double GetTemperature() const
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
void SetICM(G4bool) override
G4FragmentVector * BreakItUp(const G4Fragment &theNucleus)
static G4Proton * Proton()
Definition G4Proton.cc:90
#define TRUE
Definition globals.hh:41

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