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

#include <G4NeutronHPorLCapture.hh>

+ Inheritance diagram for G4NeutronHPorLCapture:

Public Member Functions

 G4NeutronHPorLCapture ()
 
 ~G4NeutronHPorLCapture ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
G4int GetNiso ()
 
G4bool IsThisElementOK (G4String)
 
G4VCrossSectionDataSetGiveXSectionDataSet ()
 
- 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 57 of file G4NeutronHPorLCapture.hh.

Constructor & Destructor Documentation

◆ G4NeutronHPorLCapture()

G4NeutronHPorLCapture::G4NeutronHPorLCapture ( )

Definition at line 43 of file G4NeutronHPorLCapture.cc.

44 :G4HadronicInteraction("NeutronHPorLCapture")
45{
47 if(!getenv("G4NEUTRONHPDATA"))
48 throw G4HadronicException(__FILE__, __LINE__, "Please setenv G4NEUTRONHPDATA to point to the neutron cross-section files.");
49 dirName = getenv("G4NEUTRONHPDATA");
50 G4String tString = "/Capture/";
51 dirName = dirName + tString;
52// G4cout <<"G4NeutronHPorLCapture::G4NeutronHPorLCapture testit "<<dirName<<G4endl;
54 theCapture = new G4NeutronHPChannel[numEle];
55 unavailable_elements.clear();
56 for (G4int i=0; i<numEle; i++)
57 {
58 theCapture[i].Init((*(G4Element::GetElementTable()))[i], dirName);
59 //G4cout << (*(G4Element::GetElementTable()))[i] -> GetName() << G4endl;
60 //while(!theCapture[i].Register(theFS));
61 try { while(!theCapture[i].Register(theFS)) ; }
62 catch ( G4HadronicException )
63 {
64 unavailable_elements.insert ( (*(G4Element::GetElementTable()))[i]->GetName() );
65 }
66 }
67 delete theFS;
68 SetMinEnergy(0.*eV);
69 SetMaxEnergy(20.*MeV);
70 if ( unavailable_elements.size() > 0 )
71 {
72 std::set< G4String>::iterator it;
73 G4cout << "HP Capture data are not available for thess elements "<< G4endl;
74 for ( it = unavailable_elements.begin() ; it != unavailable_elements.end() ; it++ )
75 G4cout << *it << G4endl;
76 G4cout << "Low Energy Parameterization Models will be used."<< G4endl;
77 }
78
79 createXSectionDataSet();
80}
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
void Init(G4Element *theElement, const G4String dirName)

◆ ~G4NeutronHPorLCapture()

G4NeutronHPorLCapture::~G4NeutronHPorLCapture ( )

Definition at line 82 of file G4NeutronHPorLCapture.cc.

83{
84 delete [] theCapture;
85 delete theDataSet;
86}

Member Function Documentation

◆ ApplyYourself()

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

Implements G4HadronicInteraction.

Definition at line 90 of file G4NeutronHPorLCapture.cc.

91{
92 const G4Material * theMaterial = aTrack.GetMaterial();
93 G4int n = theMaterial->GetNumberOfElements();
94 G4int index = theMaterial->GetElement(0)->GetIndex();
95 if(n!=1)
96 {
97 G4int i;
98 xSec = new G4double[n];
99 G4double sum=0;
100 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
101 G4double rWeight;
102 G4NeutronHPThermalBoost aThermalE;
103 for (i=0; i<n; i++)
104 {
105 index = theMaterial->GetElement(i)->GetIndex();
106 rWeight = NumAtomsPerVolume[i];
107 G4double x = aThermalE.GetThermalEnergy(aTrack, theMaterial->GetElement(i), theMaterial->GetTemperature());
108
109 //xSec[i] = theCapture[index].GetXsec(aThermalE.GetThermalEnergy(aTrack,
110 // theMaterial->GetElement(i),
111 // theMaterial->GetTemperature()));
112 xSec[i] = theCapture[index].GetXsec(x);
113
114 xSec[i] *= rWeight;
115 sum+=xSec[i];
116 }
117 G4double random = G4UniformRand();
118 G4double running = 0;
119 for (i=0; i<n; i++)
120 {
121 running += xSec[i];
122 index = theMaterial->GetElement(i)->GetIndex();
123 if(random<=running/sum) break;
124 }
125 delete [] xSec;
126 // it is element-wise initialised.
127 }
128 return theCapture[index].ApplyYourself(aTrack);
129}
double G4double
Definition: G4Types.hh:64
#define G4UniformRand()
Definition: Randomize.hh:53
size_t GetIndex() const
Definition: G4Element.hh:182
const G4Material * GetMaterial() const
G4double GetTemperature() const
Definition: G4Material.hh:181
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:201
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)
G4double GetXsec(G4double energy)
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)

Referenced by G4NeutronHPorLCaptureModel::ApplyYourself().

◆ GetFatalEnergyCheckLevels()

const std::pair< G4double, G4double > G4NeutronHPorLCapture::GetFatalEnergyCheckLevels ( ) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 147 of file G4NeutronHPorLCapture.cc.

148{
149 //return std::pair<G4double, G4double>(10*perCent,10*GeV);
150 return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
151}
#define DBL_MAX
Definition: templates.hh:83

◆ GetNiso()

G4int G4NeutronHPorLCapture::GetNiso ( )
inline

Definition at line 69 of file G4NeutronHPorLCapture.hh.

69{return theCapture[0].GetNiso();}

◆ GiveXSectionDataSet()

G4VCrossSectionDataSet * G4NeutronHPorLCapture::GiveXSectionDataSet ( )
inline

Definition at line 83 of file G4NeutronHPorLCapture.hh.

83{ return theDataSet; };

Referenced by G4NeutronHPorLCaptureModel::GiveHPXSectionDataSet().

◆ IsThisElementOK()

G4bool G4NeutronHPorLCapture::IsThisElementOK ( G4String  name)

Definition at line 133 of file G4NeutronHPorLCapture.cc.

134{
135 if ( unavailable_elements.find( name ) == unavailable_elements.end() )
136 return TRUE;
137 else
138 return FALSE;
139}
#define TRUE
Definition: globals.hh:55
#define FALSE
Definition: globals.hh:52

Referenced by G4NeutronHPorLCaptureModel::ApplyYourself().


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