Geant4 11.1.1
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=0, const G4String &nam="DNADiracRMatrixExcitationModel")
 
virtual ~G4DNADiracRMatrixExcitationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual G4double GetExtendedTotalCrossSection (const G4Material *material, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double GetExtendedPartialCrossSection (const G4Material *material, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SelectStationary (G4bool input)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=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 LPMFlag () 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 SetLPMFlag (G4bool val)
 
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)
 
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
 
size_t currentCoupleIndex = 0
 
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()

G4DNADiracRMatrixExcitationModel::G4DNADiracRMatrixExcitationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNADiracRMatrixExcitationModel" 
)

Definition at line 49 of file G4DNADiracRMatrixExcitationModel.cc.

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

◆ ~G4DNADiracRMatrixExcitationModel()

G4DNADiracRMatrixExcitationModel::~G4DNADiracRMatrixExcitationModel ( )
virtual

Definition at line 72 of file G4DNADiracRMatrixExcitationModel.cc.

73{
74 if (fTableData) delete fTableData;
75}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 127 of file G4DNADiracRMatrixExcitationModel.cc.

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

◆ GetExtendedPartialCrossSection()

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

Definition at line 237 of file G4DNADiracRMatrixExcitationModel.cc.

242{
243 G4double value=0;
244
245 if(particle->GetParticleName()=="e-"){
246
247 if(level==0){
248 // y = [0]+[1]/pow(x-2,2)
249 value = paramFuncTCS_5dto6s1[0]+paramFuncTCS_5dto6s1[1]
250 /std::pow(kineticEnergy/eV-paramFuncTCS_5dto6s1[2],2);
251 }
252 else if(level==1){
253 // y = [0]+[1]/pow(x-2,2)
254 value = paramFuncTCS_5dto6s2[0]+paramFuncTCS_5dto6s2[1]
255 /std::pow(kineticEnergy/eV-paramFuncTCS_5dto6s2[2],2);
256 }
257 else if(level==2){
258 // y = [0]+[1]*log(x-2)/(x-[2])
259 value = paramFuncTCS_6sto6p1[0]+paramFuncTCS_6sto6p1[1]
260 *G4Log(kineticEnergy/eV-paramFuncTCS_6sto6p1[2])
261 /(kineticEnergy/eV-paramFuncTCS_6sto6p1[2]);
262 }
263 else if(level==3){
264 // y = [0]+[1]*log(x-2)/(x-[2])
265 value = paramFuncTCS_6sto6p2[0]+paramFuncTCS_6sto6p2[1]
266 *G4Log(kineticEnergy/eV-paramFuncTCS_6sto6p2[2])
267 /(kineticEnergy/eV-paramFuncTCS_6sto6p2[2]);
268 }
269 }
270
271 return value*cm*cm;
272}
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 218 of file G4DNADiracRMatrixExcitationModel.cc.

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

Referenced by CrossSectionPerVolume().

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 79 of file G4DNADiracRMatrixExcitationModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 183 of file G4DNADiracRMatrixExcitationModel.cc.

188{
189
190 if (verboseLevel > 3)
191 {
192 G4cout <<
193 "Calling SampleSecondaries() of G4DNADiracRMatrixExcitationModel"
194 << G4endl;
195 }
196
197 G4ParticleDefinition* particle = aDynamicParticle->GetDefinition();
198 G4double k = aDynamicParticle->GetKineticEnergy();
199
200 G4int level = RandomSelect(couple->GetMaterial(),particle,
201 k);
202 G4double excitationEnergy = ExcitationEnergyAu[level]*eV;
203 G4double newEnergy = k - excitationEnergy;
204
205 if (newEnergy > 0)
206 {
207 //Energy Loss
209 (aDynamicParticle->GetMomentumDirection());
211 if(!statCode) fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
213 }
214}
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 128 of file G4DNADiracRMatrixExcitationModel.hh.

129{
130 statCode = input;
131}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNADiracRMatrixExcitationModel::fParticleChangeForGamma
protected

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