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

#include <G4DNARelativisticIonisationModel.hh>

+ Inheritance diagram for G4DNARelativisticIonisationModel:

Public Member Functions

 G4DNARelativisticIonisationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNARelativisticIonisationModel")
 
virtual ~G4DNARelativisticIonisationModel ()
 
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 GetTotalCrossSection (const G4Material *material, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double GetPartialCrossSection (const G4Material *material, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double GetDifferentialCrossSection (const G4Material *material, const G4ParticleDefinition *particle, G4double kineticEnergy, G4double secondaryEnergy, G4int level)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual void LoadAtomicStates (G4int z, const char *path)
 
void SelectStationary (G4bool input)
 
void SelectFasterComputation (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 46 of file G4DNARelativisticIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4DNARelativisticIonisationModel()

G4DNARelativisticIonisationModel::G4DNARelativisticIonisationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNARelativisticIonisationModel" 
)

Definition at line 61 of file G4DNARelativisticIonisationModel.cc.

63 :
64G4VEmModel(nam), isInitialised(false),statCode(false),fasterCode(true)
65{
66 fHighEnergyLimit = 0;
67 fLowEnergyLimit = 0;
68
69 verboseLevel = 0;
70
72 fAtomDeexcitation = 0;
73 fMaterialDensity = 0;
74 fParticleDefinition = 0;
76
77 if (verboseLevel > 0)
78 {
79 G4cout << "Relativistic Ionisation Model is constructed " << G4endl;
80 }
81}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802

◆ ~G4DNARelativisticIonisationModel()

G4DNARelativisticIonisationModel::~G4DNARelativisticIonisationModel ( )
virtual

Definition at line 85 of file G4DNARelativisticIonisationModel.cc.

86{
87 // Cross section
88}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 232 of file G4DNARelativisticIonisationModel.cc.

238{
239 if (verboseLevel > 3)
240 {
241 G4cout <<
242 "Calling CrossSectionPerVolume() of G4DNARelativisticIonisationModel"
243 << G4endl;
244 }
245
246 if(particleDefinition != fParticleDefinition) return 0;
247
248 // Calculate total cross section for model
249 G4double sigma=0;
250
251 if(material->GetNumberOfElements()>1) return 0.; // Protection for Molecules
252 G4double atomicNDensity = material->GetAtomicNumDensityVector()[0];
253 G4double z = material->GetZ();
254
255 if(atomicNDensity!= 0.0)
256 {
257 if (ekin >= fLowEnergyLimit && ekin < fHighEnergyLimit)
258 {
259 sigma = GetTotalCrossSection(material,particleDefinition,ekin);
260 }
261
262 if (verboseLevel > 2)
263 {
264 G4cout << "__________________________________" << G4endl;
265 G4cout << "=== G4DNARelativisticIonisationModel - XS INFO START" <<G4endl;
266 G4cout << "=== Kinetic energy (eV)=" << ekin/eV << " particle : "
267 << particleDefinition->GetParticleName() << G4endl;
268 G4cout << "=== Cross section per atom for Z="<<z<<" is (cm^2)"
269 << sigma/cm/cm << G4endl;
270 G4cout << "=== Cross section per atom for Z="<<z<<" is (cm^-1)="
271 << sigma*atomicNDensity/(1./cm) << G4endl;
272 G4cout << "=== G4DNARelativisticIonisationModel - XS INFO END" << G4endl;
273 }
274 }
275 return sigma*atomicNDensity;
276}
double G4double
Definition: G4Types.hh:83
virtual G4double GetTotalCrossSection(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

◆ GetDifferentialCrossSection()

G4double G4DNARelativisticIonisationModel::GetDifferentialCrossSection ( const G4Material material,
const G4ParticleDefinition particle,
G4double  kineticEnergy,
G4double  secondaryEnergy,
G4int  level 
)
virtual

Definition at line 523 of file G4DNARelativisticIonisationModel.cc.

529{
530 G4double value=0.;
531 G4double constRy =13.6057E-6;//MeV
532
533 G4int z = material->GetZ();
534
536 if(particle==electronDef){
537 G4double w = secondaryEnergy /Ebinding[z].at(level);
538 G4double t = kineticEnergy /Ebinding[z].at(level);
539 G4double tdash = kineticEnergy /CLHEP::electron_mass_c2;
540 G4double udash = Ekinetic[z].at(level)/CLHEP::electron_mass_c2;
541 G4double bdash = Ebinding[z].at(level)/CLHEP::electron_mass_c2;
542 G4double beta_t2 = 1.-1./std::pow(1.+tdash,2);
543 G4double beta_u2 = 1.-1./std::pow(1.+udash,2);
544 G4double beta_b2 = 1.-1./std::pow(1.+bdash,2);
545 G4double alpha = std::sqrt(2*constRy/CLHEP::electron_mass_c2);
546 G4double phi = std::cos(std::sqrt(std::pow(alpha,2)/(beta_t2+beta_b2))
547 *G4Log(beta_t2/beta_b2));
548 G4double constS = 4*CLHEP::pi*std::pow(CLHEP::Bohr_radius,2)
549 *Nelectrons[z].at(level)*std::pow(alpha,4);
550
551 if(secondaryEnergy<=((kineticEnergy-Ebinding[z].at(level))/2.))
552 {
553 value = constS/((beta_t2+(beta_u2+beta_b2)/iShell[z].at(level))*2.*bdash)
554 *(-phi/(t+1.)*(1./std::pow(w+1.,1.)+1./std::pow(t-w,1.))
555 *(1.+2*tdash)/std::pow(1.+tdash/2.,2.)
556 +1./std::pow(w+1.,2.)+1./std::pow(t-w,2.)
557 +std::pow(bdash,2)/std::pow(1+tdash/2.,2)
558 +(1./std::pow(w+1.,3.)+1./std::pow(t-w,3.))
559 *(G4Log(beta_t2/(1.-beta_t2))-beta_t2-G4Log(2*bdash)));
560 }
561 }
562 return value;
563}
G4double G4Log(G4double x)
Definition: G4Log.hh:227
int G4int
Definition: G4Types.hh:85
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88

◆ GetPartialCrossSection()

G4double G4DNARelativisticIonisationModel::GetPartialCrossSection ( const G4Material material,
G4int  level,
const G4ParticleDefinition particle,
G4double  kineticEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 482 of file G4DNARelativisticIonisationModel.cc.

487{
488 G4double value = 0;
489 G4double constRy =13.6057E-6;//MeV
490
492 G4int z = material->GetZ();
493 if(particle==electronDef){
494
495 G4double t = kineticEnergy /Ebinding[z].at(level);
496 G4double tdash = kineticEnergy /CLHEP::electron_mass_c2;
497 G4double udash = Ekinetic[z].at(level)/CLHEP::electron_mass_c2;
498 G4double bdash = Ebinding[z].at(level)/CLHEP::electron_mass_c2;
499 G4double beta_t2 = 1.-1./std::pow(1.+tdash,2);
500 G4double beta_u2 = 1.-1./std::pow(1.+udash,2);
501 G4double beta_b2 = 1.-1./std::pow(1.+bdash,2);
502 G4double alpha = std::sqrt(2*constRy/CLHEP::electron_mass_c2);
503 G4double phi = std::cos(std::sqrt(std::pow(alpha,2)
504 /(beta_t2+beta_b2))*G4Log(beta_t2/beta_b2));
505 G4double constS = 4*CLHEP::pi*std::pow(CLHEP::Bohr_radius,2)
506 *Nelectrons[z].at(level)*std::pow(alpha,4);
507
508 if(Ebinding[z].at(level)<=kineticEnergy)
509 {
510 value =constS/((beta_t2+(beta_u2+beta_b2)/iShell[z].at(level))*2.*bdash)
511 *(1./2.*(G4Log(beta_t2/(1.-beta_t2))-beta_t2-G4Log(2.*bdash))
512 *(1.-1./std::pow(t,2.))
513 +1.-1./t-G4Log(t)/(t+1.)*(1.+2.*tdash)/(std::pow(1.+tdash/2.,2.))
514 *phi+std::pow(bdash,2)/(std::pow(1+tdash/2.,2))*(t-1)/2.);
515 }
516
517 }
518 return value;
519}

Referenced by GetTotalCrossSection().

◆ GetTotalCrossSection()

G4double G4DNARelativisticIonisationModel::GetTotalCrossSection ( const G4Material material,
const G4ParticleDefinition particle,
G4double  kineticEnergy 
)
virtual

Definition at line 463 of file G4DNARelativisticIonisationModel.cc.

467{
468 G4double value=0;
469 G4int z = material->GetZ();
470 if(z!=79){ return 0.;}
471 else {
472 std::size_t N=iState[z].size();
473 for(G4int i=0; i<(G4int)N; ++i){
474 value = value+GetPartialCrossSection(material,i,particle,kineticEnergy);
475 }
476 return value;
477 }
478}
virtual G4double GetPartialCrossSection(const G4Material *material, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
#define N
Definition: crc32.c:56

Referenced by CrossSectionPerVolume().

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 92 of file G4DNARelativisticIonisationModel.cc.

94{
95
96 if (verboseLevel > 3)
97 {
98 G4cout <<
99 "Calling G4DNARelativisticIonisationModel::Initialise()"
100 << G4endl;
101 }
102
103
104 if(fParticleDefinition != 0 && fParticleDefinition != particle)
105 {
106 G4Exception("G4DNARelativisticIonisationModel::Initialise","em0001",
107 FatalException,"Model already initialized for another particle type.");
108 }
109
110 fParticleDefinition = particle;
112 if(particle == electronDef)
113 {
114 fLowEnergyLimit = 10 * eV;
115 fHighEnergyLimit = 1.0 * GeV;
116
117 std::ostringstream eFullFileNameZ;
118
119 const char *path = G4FindDataDir("G4LEDATA");
120 if (!path)
121 {
122 G4Exception("G4DNARelativisticIonisationModel::Initialise","em0006",
123 FatalException,"G4LEDATA environment variable not set.");
124 return;
125 }
126
127
128 G4ProductionCutsTable *coupletable
130 G4int Ncouple = (G4int)coupletable ->GetTableSize();
131 for(G4int i=0; i<Ncouple; ++i)
132 {
133 const G4MaterialCutsCouple* couple
134 = coupletable->GetMaterialCutsCouple(i);
135 const G4Material * material = couple ->GetMaterial();
136 {
137 // Protection: only for single element
138 if(material->GetNumberOfElements()>1) continue;
139
140 G4int Z = material->GetZ();
141 // Protection: only for GOLD
142 if(Z!=79) continue;
143
144 iState [Z].clear();
145 iShell [Z].clear();
146 iSubShell [Z].clear();
147 Nelectrons[Z].clear();
148 Ebinding [Z].clear();
149 Ekinetic [Z].clear();
150 LoadAtomicStates(Z,path);
151
152 /////////////Load cumulated DCS////////////////
153 eVecEZ.clear();
154 eVecEjeEZ.clear();
155 eProbaShellMapZ.clear();
156 eDiffCrossSectionDataZ.clear();
157
158 eFullFileNameZ.str("");
159 eFullFileNameZ.clear(stringstream::goodbit);
160
161 eFullFileNameZ
162 << path
163 << "/dna/sigmadiff_cumulated_ionisation_e_RBEBV_Z"
164 << Z << ".dat";
165 std::ifstream eDiffCrossSectionZ(eFullFileNameZ.str().c_str());
166 if (!eDiffCrossSectionZ)
167 G4Exception("G4DNARelativisticIonisationModel::Initialise","em0003",
169 "Missing data file for cumulated DCS");
170
171 eVecEZ[Z].push_back(0.);
172 while(!eDiffCrossSectionZ.eof())
173 {
174 G4double tDummy;
175 G4double eDummy;
176 eDiffCrossSectionZ>>tDummy>>eDummy;
177 if (tDummy != eVecEZ[Z].back())
178 {
179 eVecEZ[Z].push_back(tDummy);
180 eVecEjeEZ[Z][tDummy].push_back(0.);
181 }
182
183 for(G4int istate=0;istate<(G4int)iState[Z].size();istate++)
184 {
185 eDiffCrossSectionZ>>
186 eDiffCrossSectionDataZ[Z][istate][tDummy][eDummy];
187 eEjectedEnergyDataZ[Z][istate][tDummy]
188 [eDiffCrossSectionDataZ[Z][istate][tDummy][eDummy]]
189 = eDummy;
190 eProbaShellMapZ[Z][istate][tDummy].push_back(
191 eDiffCrossSectionDataZ[Z][istate][tDummy][eDummy]);
192 }
193
194 if (eDummy != eVecEjeEZ[Z][tDummy].back()){
195 eVecEjeEZ[Z][tDummy].push_back(eDummy);
196 }
197 }
198 }
199 }
200 }
201 else
202 {
203 G4cout<<
204 "Error : No particle Definition is found in G4DNARelativisticIonisationModel"
205 <<G4endl;
206 return;
207 }
208
209 if( verboseLevel>0 )
210 {
211 G4cout << "Relativistic Ionisation model is initialized " << G4endl
212 << "Energy range: "
213 << LowEnergyLimit() / eV << " eV - "
214 << HighEnergyLimit() / keV << " keV for "
215 << particle->GetParticleName()
216 << G4endl;
217 }
218
219 // Initialise gold density pointer
220 fMaterialDensity = G4DNAMolecularMaterial::Instance()
222
223 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
225
226 if (isInitialised){return;}
227 isInitialised = true;
228}
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
const G4int Z[17]
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
virtual void LoadAtomicStates(G4int z, const char *path)
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:691
const G4String & GetParticleName() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634

◆ LoadAtomicStates()

void G4DNARelativisticIonisationModel::LoadAtomicStates ( G4int  z,
const char *  path 
)
virtual

Definition at line 394 of file G4DNARelativisticIonisationModel.cc.

396{
397
398 if (verboseLevel > 3)
399 {
400 G4cout <<
401 "Calling LoadAtomicStates() of G4DNARelativisticIonisationModel"
402 << G4endl;
403 }
404 const char *datadir = path;
405 if(!datadir)
406 {
407 datadir = G4FindDataDir("G4LEDATA");
408 if(!datadir)
409 {
410 G4Exception("G4DNARelativisticIonisationModel::LoadAtomicStates()",
411 "em0002",FatalException,"Enviroment variable G4LEDATA not defined");
412
413 return;
414 }
415 }
416 std::ostringstream targetfile;
417 targetfile << datadir <<"/dna/atomicstate_Z"<< z <<".dat";
418 std::ifstream fin(targetfile.str().c_str());
419 if(!fin)
420 {
421 G4cout<< " Error : "<< targetfile.str() <<" is not found "<<G4endl;
422 G4Exception("G4DNARelativisticIonisationModel::LoadAtomicStates()","em0002",
423 FatalException,"There is no target file");
424 return;
425 }
426
427 G4String buff0,buff1,buff2,buff3,buff4,buff5,buff6;
428 fin >> buff0 >>buff1>>buff2>>buff3>>buff4>>buff5>>buff6;
429 G4int iline=0;
430 while(true){
431 fin >> buff0 >>buff1>>buff2>>buff3>>buff4>>buff5>>buff6;
432 if(!fin.eof())
433 {
434 iState [z].push_back(stoi(buff0));
435 iShell [z].push_back(stoi(buff1));
436 iSubShell [z].push_back(stoi(buff2));
437 Nelectrons[z].push_back(stoi(buff3));
438 Ebinding [z].push_back(stod(buff4));
439 if(stod(buff5)==0.)
440 {// if there is no kinetic energy in the file, kinetic energy
441 // for Bhor atomic model will be calculated: !!! I's not realistic!!!
442 G4double radius = std::pow(iShell[z].at(iline),2)
443 *std::pow(CLHEP::hbar_Planck,2)*(4*CLHEP::pi*CLHEP::epsilon0)
444 /CLHEP::electron_mass_c2;
445 G4double momentum = iShell[z].at(iline)*CLHEP::hbar_Planck/radius;
446 Ekinetic[z].push_back(std::pow(momentum,2)/(2*CLHEP::electron_mass_c2));
447 }
448 else
449 {
450 Ekinetic [z].push_back(stod(buff5));
451 }
452 iline++;
453 }
454 else
455 {
456 break;
457 }
458 }
459}

Referenced by Initialise().

◆ SampleSecondaries()

void G4DNARelativisticIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle particle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 280 of file G4DNARelativisticIonisationModel.cc.

285{
286 if (verboseLevel > 3)
287 {
288 G4cout <<
289 "Calling SampleSecondaries() of G4DNARelativisticIonisationModel"
290 << G4endl;
291 }
292
293
294 G4ParticleDefinition* particleDef = particle->GetDefinition();
295 G4double k = particle->GetKineticEnergy();
296 G4double ejectedE = 0.*eV;
297
298 if(fLowEnergyLimit <= k && k<fHighEnergyLimit)
299 {
300 G4ThreeVector primaryDir = particle ->GetMomentumDirection();
301
302 G4double particleMass = particleDef->GetPDGMass();
303 G4double totalEnergy = k+particleMass;
304 G4double pSquare = k*(totalEnergy+particleMass);
305 G4double totalMomentum = std::sqrt(pSquare);
306
307 const G4Material *material = couple->GetMaterial();
308 G4int z = material->GetZ();
309 G4int level = RandomSelect(material,particleDef,k);
310
311 if(k<Ebinding[z].at(level)) return;
312
313 G4int NumSecParticlesInit =0;
314 G4int NumSecParticlesFinal=0;
315
316 if(fAtomDeexcitation){
318 const G4AtomicShell *shell = fAtomDeexcitation->GetAtomicShell(z,as);
319 NumSecParticlesInit = (G4int)fvect->size();
320 fAtomDeexcitation->GenerateParticles(fvect,shell,z,0,0);
321 NumSecParticlesFinal = (G4int)fvect->size();
322 }
323
324 ejectedE
325 = GetEjectedElectronEnergy (material,particleDef,k,level);
326 G4ThreeVector ejectedDir
327 = GetEjectedElectronDirection(particleDef,k,ejectedE);
328 ejectedDir.rotateUz(primaryDir);
329
330 G4double scatteredE = k - Ebinding[z].at(level) - ejectedE;
331
332 if(particleDef == G4Electron::ElectronDefinition()){
333 G4double secondaryTotMomentum
334 = std::sqrt(ejectedE*(ejectedE+2*CLHEP::electron_mass_c2));
335 G4double finalMomentumX
336 = totalMomentum*primaryDir.x()- secondaryTotMomentum*ejectedDir.x();
337 G4double finalMomentumY
338 = totalMomentum*primaryDir.y()- secondaryTotMomentum*ejectedDir.y();
339 G4double finalMomentumZ
340 = totalMomentum*primaryDir.z()- secondaryTotMomentum*ejectedDir.z();
341
342 G4ThreeVector scatteredDir(finalMomentumX,finalMomentumY,finalMomentumZ);
344
345 }
346 else
347 {
349 }
350
351 //G4double deexSecEnergy=0.;
352 G4double restEproduction = Ebinding[z].at(level);
353 for(G4int iparticle=NumSecParticlesInit;
354 iparticle<NumSecParticlesFinal;iparticle++)
355 {
356 //deexSecEnergy = deexSecEnergy + (*fvect)[iparticle]->GetKineticEnergy();
357 G4double Edeex = (*fvect)[iparticle]->GetKineticEnergy();
358 if(restEproduction>=Edeex){
359 restEproduction -= Edeex;
360 }
361 else{
362 delete (*fvect)[iparticle];
363 (*fvect)[iparticle]=0;
364 }
365 }
366 if(restEproduction < 0.0){
367 G4Exception("G4DNARelativisticIonisationModel::SampleSecondaries()",
368 "em0008",FatalException,"Negative local energy deposit");
369 }
370
371 if(!statCode)
372 {
373 if(scatteredE>0){
376 //fParticleChangeForGamma
377 //->ProposeLocalEnergyDeposit(k-scatteredE-ejectedE-deexSecEnergy);
378 }
379 }
380 else
381 {
384 }
385
386 if(ejectedE>0){
387 G4DynamicParticle* ejectedelectron
388 = new G4DynamicParticle(G4Electron::Electron(),ejectedDir,ejectedE);
389 fvect->push_back(ejectedelectron);
390 }
391 }
392}
G4AtomicShellEnumerator
double z() const
double x() const
double y() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectFasterComputation()

void G4DNARelativisticIonisationModel::SelectFasterComputation ( G4bool  input)
inline

Definition at line 83 of file G4DNARelativisticIonisationModel.hh.

83{fasterCode = input;};

◆ SelectStationary()

void G4DNARelativisticIonisationModel::SelectStationary ( G4bool  input)
inline

Definition at line 82 of file G4DNARelativisticIonisationModel.hh.

82{statCode = input;};

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNARelativisticIonisationModel::fParticleChangeForGamma
protected

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