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

#include <G4DNADiracRMatrixExcitationModel.hh>

+ Inheritance diagram for G4DNADiracRMatrixExcitationModel:

Public Member Functions

 G4DNADiracRMatrixExcitationModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNADiracRMatrixExcitationModel")
 
 ~G4DNADiracRMatrixExcitationModel () override
 
G4DNADiracRMatrixExcitationModeloperator= (const G4DNADiracRMatrixExcitationModel &right)=delete
 
 G4DNADiracRMatrixExcitationModel (const G4DNADiracRMatrixExcitationModel &)=delete
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector())) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 
virtual G4double GetExtendedTotalCrossSection (const G4Material *material, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double GetExtendedPartialCrossSection (const G4Material *material, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SelectStationary (G4bool input)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 50 of file G4DNADiracRMatrixExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNADiracRMatrixExcitationModel() [1/2]

G4DNADiracRMatrixExcitationModel::G4DNADiracRMatrixExcitationModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNADiracRMatrixExcitationModel" )

Definition at line 47 of file G4DNADiracRMatrixExcitationModel.cc.

48 :
49 G4VEmModel(nam)
50{
51 fpMaterialDensity = nullptr;
52 fHighEnergyLimit = 0;
53 fExperimentalEnergyLimit= 0;
54 fLowEnergyLimit = 0;
55 fParticleDefinition = nullptr;
56
57 verboseLevel = 0;
58
59 if (verboseLevel > 0)
60 {
61 G4cout << "Dirac R-matrix excitation model is constructed " << G4endl;
62 }
63
65 statCode = false;
66}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4DNADiracRMatrixExcitationModel()

G4DNADiracRMatrixExcitationModel::~G4DNADiracRMatrixExcitationModel ( )
override

Definition at line 70 of file G4DNADiracRMatrixExcitationModel.cc.

71{
72 delete fTableData;
73}

◆ G4DNADiracRMatrixExcitationModel() [2/2]

G4DNADiracRMatrixExcitationModel::G4DNADiracRMatrixExcitationModel ( const G4DNADiracRMatrixExcitationModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNADiracRMatrixExcitationModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double ekin,
G4double emin,
G4double emax )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 125 of file G4DNADiracRMatrixExcitationModel.cc.

131{
132 if (verboseLevel > 3)
133 {
134 G4cout <<
135 "Calling CrossSectionPerVolume() of G4DNADiracRMatrixExcitationModel"
136 << G4endl;
137 }
138
139 G4double atomicNDensity = material->GetAtomicNumDensityVector()[0];
140
141 // Protection: for single element
142 if(material->GetNumberOfElements()>1) return 0.;
143
144 G4double z = material->GetZ();
145
146 // Protection: for Gold
147 if(z!=79){return 0.;}
148
149 G4double sigma=0.;
150
151 if(atomicNDensity!= 0.0)
152 {
153 if (ekin >= fLowEnergyLimit && ekin < fExperimentalEnergyLimit)
154 {
155 sigma = fTableData->FindValue(ekin);
156 }
157 else if ((fExperimentalEnergyLimit <= ekin) && (ekin < fHighEnergyLimit))
158 {
159 sigma = GetExtendedTotalCrossSection(material,particleDefinition,ekin);
160 }
161
162 if (verboseLevel > 2)
163 {
164 G4cout<<"__________________________________" << G4endl;
165 G4cout<<"=== G4DNADiracRMatrixExcitationModel - XS INFO START"<<G4endl;
166 G4cout<<"=== Kinetic energy (eV)=" << ekin/eV << " particle : "
167 <<particleDefinition->GetParticleName() << G4endl;
168 G4cout<<"=== Cross section per atom for Z="<<z<<" is (cm^2)"
169 <<sigma/cm/cm << G4endl;
170 G4cout<<"=== Cross section per atom for Z="<<z<<" is (cm^-1)="
171 <<sigma*atomicNDensity/(1./cm) << G4endl;
172 G4cout<<"=== G4DNADiracRMatrixExcitationModel - XS INFO END"<<G4endl;
173 }
174 }
175
176 return sigma*atomicNDensity;
177}
double G4double
Definition G4Types.hh:83
G4double FindValue(G4double e, G4int componentId=0) const override
virtual G4double GetExtendedTotalCrossSection(const G4Material *material, const G4ParticleDefinition *, G4double kineticEnergy)
G4double GetZ() const
const G4double * GetAtomicNumDensityVector() const
std::size_t GetNumberOfElements() const

◆ GetExtendedPartialCrossSection()

G4double G4DNADiracRMatrixExcitationModel::GetExtendedPartialCrossSection ( const G4Material * material,
G4int level,
const G4ParticleDefinition * particle,
G4double kineticEnergy )
virtual

Definition at line 235 of file G4DNADiracRMatrixExcitationModel.cc.

240{
241 G4double value=0;
242
243 if(particle->GetParticleName()=="e-"){
244
245 if(level==0){
246 // y = [0]+[1]/pow(x-2,2)
247 value = paramFuncTCS_5dto6s1[0]+paramFuncTCS_5dto6s1[1]
248 /std::pow(kineticEnergy/eV-paramFuncTCS_5dto6s1[2],2);
249 }
250 else if(level==1){
251 // y = [0]+[1]/pow(x-2,2)
252 value = paramFuncTCS_5dto6s2[0]+paramFuncTCS_5dto6s2[1]
253 /std::pow(kineticEnergy/eV-paramFuncTCS_5dto6s2[2],2);
254 }
255 else if(level==2){
256 // y = [0]+[1]*log(x-2)/(x-[2])
257 value = paramFuncTCS_6sto6p1[0]+paramFuncTCS_6sto6p1[1]
258 *G4Log(kineticEnergy/eV-paramFuncTCS_6sto6p1[2])
259 /(kineticEnergy/eV-paramFuncTCS_6sto6p1[2]);
260 }
261 else if(level==3){
262 // y = [0]+[1]*log(x-2)/(x-[2])
263 value = paramFuncTCS_6sto6p2[0]+paramFuncTCS_6sto6p2[1]
264 *G4Log(kineticEnergy/eV-paramFuncTCS_6sto6p2[2])
265 /(kineticEnergy/eV-paramFuncTCS_6sto6p2[2]);
266 }
267 }
268
269 return value*cm*cm;
270}
G4double G4Log(G4double x)
Definition G4Log.hh:227
const G4String & GetParticleName() const

Referenced by GetExtendedTotalCrossSection().

◆ GetExtendedTotalCrossSection()

G4double G4DNADiracRMatrixExcitationModel::GetExtendedTotalCrossSection ( const G4Material * material,
const G4ParticleDefinition * particle,
G4double kineticEnergy )
virtual

Definition at line 216 of file G4DNADiracRMatrixExcitationModel.cc.

220{
221 G4double value=0;
222
223 size_t N=fTableData->NumberOfComponents();
224
225 for(int i=0;i<(int)N;i++){
226 value = value+GetExtendedPartialCrossSection(material,i,particle,
227 kineticEnergy);
228 }
229
230 return value;
231}
size_t NumberOfComponents() const override
virtual G4double GetExtendedPartialCrossSection(const G4Material *material, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
#define N
Definition crc32.c:57

Referenced by CrossSectionPerVolume().

◆ Initialise()

void G4DNADiracRMatrixExcitationModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector & = *(new G4DataVector()) )
overridevirtual

Implements G4VEmModel.

Definition at line 77 of file G4DNADiracRMatrixExcitationModel.cc.

79{
80
81 if (verboseLevel > 3)
82 {
83 G4cout <<
84 "Calling G4DNADiracRMatrixExcitationModel::Initialise()"
85 << G4endl;
86 }
87
88 fParticleDefinition = particle;
89
90 if(particle->GetParticleName() == "e-")
91 {
92 fTableFile = "dna/sigma_excitation_e_diracrmatrix_Z79";
93 fLowEnergyLimit = 10 * eV;
94 fExperimentalEnergyLimit = 577.* eV;
95 fHighEnergyLimit = 1.0 * GeV;
96 }
97 else
98 {
99 G4Exception("G4DNADiracRMatrixExcitationModel::Initialise","em0001",
100 FatalException,"Not defined for other particles than electrons.");
101 return;
102 }
103
104 G4double scaleFactor = 1. * cm * cm;
105 fTableData = new G4DNACrossSectionDataSet
106 (new G4LogLogInterpolation,eV,scaleFactor );
107 fTableData->LoadData(fTableFile);
108
109 if( verboseLevel>0 )
110 {
111 G4cout << "Dirac R-matrix excitation model is initialized " << G4endl
112 << "Energy range: "
113 << LowEnergyLimit() / eV << " eV - "<< HighEnergyLimit() / keV << " keV "
114 << " for "<< particle->GetParticleName()
115 << G4endl;
116 }
117
118 if (isInitialised){return;}
120 isInitialised = true;
121}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4bool LoadData(const G4String &argFileName) override
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const

◆ operator=()

G4DNADiracRMatrixExcitationModel & G4DNADiracRMatrixExcitationModel::operator= ( const G4DNADiracRMatrixExcitationModel & right)
delete

◆ SampleSecondaries()

void G4DNADiracRMatrixExcitationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * ,
const G4MaterialCutsCouple * couple,
const G4DynamicParticle * aDynamicParticle,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 181 of file G4DNADiracRMatrixExcitationModel.cc.

186{
187
188 if (verboseLevel > 3)
189 {
190 G4cout <<
191 "Calling SampleSecondaries() of G4DNADiracRMatrixExcitationModel"
192 << G4endl;
193 }
194
195 G4ParticleDefinition* particle = aDynamicParticle->GetDefinition();
196 G4double k = aDynamicParticle->GetKineticEnergy();
197
198 G4int level = RandomSelect(couple->GetMaterial(),particle,
199 k);
200 G4double excitationEnergy = ExcitationEnergyAu[level]*eV;
201 G4double newEnergy = k - excitationEnergy;
202
203 if (newEnergy > 0)
204 {
205 //Energy Loss
207 (aDynamicParticle->GetMomentumDirection());
209 if(!statCode) fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
211 }
212}
int G4int
Definition G4Types.hh:85
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNADiracRMatrixExcitationModel::SelectStationary ( G4bool input)
inline

Definition at line 124 of file G4DNADiracRMatrixExcitationModel.hh.

125{
126 statCode = input;
127}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNADiracRMatrixExcitationModel::fParticleChangeForGamma
protected

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