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

#include <G4DNAMeltonAttachmentModel.hh>

+ Inheritance diagram for G4DNAMeltonAttachmentModel:

Public Member Functions

 G4DNAMeltonAttachmentModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNAMeltonAttachmentModel")
 
 ~G4DNAMeltonAttachmentModel () override
 
G4DNAMeltonAttachmentModeloperator= (const G4DNAMeltonAttachmentModel &right)=delete
 
 G4DNAMeltonAttachmentModel (const G4DNAMeltonAttachmentModel &)=delete
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SetDissociationFlag (G4bool)
 
G4bool GetDissociationFlag ()
 
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 43 of file G4DNAMeltonAttachmentModel.hh.

Constructor & Destructor Documentation

◆ G4DNAMeltonAttachmentModel() [1/2]

G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNAMeltonAttachmentModel" )

Definition at line 43 of file G4DNAMeltonAttachmentModel.cc.

44 :
45G4VEmModel(nam)
46{
47 fpWaterDensity = nullptr;
48
49 SetLowEnergyLimit(4.*eV);
50 SetHighEnergyLimit(13.*eV);
51
52 verboseLevel = 0;
53 // Verbosity scale:
54 // 0 = nothing
55 // 1 = warning for energy non-conservation
56 // 2 = details of energy budget
57 // 3 = calculation of cross sections, file openings, sampling of atoms
58 // 4 = entering in methods
59
60#ifdef MELTON_VERBOSE
61 if (verboseLevel > 0)
62 {
63 G4cout << "Melton Attachment model is constructed "
64 << G4endl
65 << "Energy range: "
66 << LowEnergyLimit() / eV << " eV - "
67 << HighEnergyLimit() / eV << " eV"
68 << G4endl;
69 }
70#endif
71
73 fDissociationFlag = true;
74 fData = nullptr;
75
76 // Selection of stationary mode
77
78 statCode = false;
79}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetHighEnergyLimit(G4double)
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4DNAMeltonAttachmentModel()

G4DNAMeltonAttachmentModel::~G4DNAMeltonAttachmentModel ( )
override

Definition at line 83 of file G4DNAMeltonAttachmentModel.cc.

84{
85 delete fData;
86}

◆ G4DNAMeltonAttachmentModel() [2/2]

G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel ( const G4DNAMeltonAttachmentModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 180 of file G4DNAMeltonAttachmentModel.cc.

185{
186#ifdef MELTON_VERBOSE
187 if (verboseLevel > 3)
188 G4cout
189 << "Calling CrossSectionPerVolume() of G4DNAMeltonAttachmentModel"
190 << G4endl;
191#endif
192
193 // Calculate total cross section for model
194
195 G4double sigma = 0.;
196
197 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
198
199 if (ekin >= LowEnergyLimit() && ekin <= HighEnergyLimit())
200 sigma = fData->FindValue(ekin);
201
202#ifdef MELTON_VERBOSE
203 if (verboseLevel > 2)
204 {
205 G4cout << "__________________________________" << G4endl;
206 G4cout << "=== G4DNAMeltonAttachmentModel - XS INFO START" << G4endl;
207 G4cout << "--- Kinetic energy(eV)=" << ekin/eV
208 << " particle : " << particleDefinition->GetParticleName()
209 << G4endl;
210 G4cout << "--- Cross section per water molecule (cm^2)="
211 << sigma/cm/cm << G4endl;
212 G4cout << "--- Cross section per water molecule (cm^-1)="
213 << sigma*waterDensity/(1./cm) << G4endl;
214 G4cout << "--- G4DNAMeltonAttachmentModel - XS INFO END" << G4endl;
215 }
216#endif
217
218 return sigma*waterDensity;
219}
double G4double
Definition G4Types.hh:83
G4double FindValue(G4double e, G4int componentId=0) const override
std::size_t GetIndex() const

◆ GetDissociationFlag()

G4bool G4DNAMeltonAttachmentModel::GetDissociationFlag ( )
inline

Definition at line 96 of file G4DNAMeltonAttachmentModel.hh.

97{
98 return fDissociationFlag;
99}

◆ Initialise()

void G4DNAMeltonAttachmentModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 90 of file G4DNAMeltonAttachmentModel.cc.

92{
93#ifdef MELTON_VERBOSE
94 if (verboseLevel > 3)
95 G4cout
96 << "Calling G4DNAMeltonAttachmentModel::Initialise()" << G4endl;
97#endif
98
99 // Only electron
100
101 if(particle->GetParticleName() != "e-")
102 {
103 G4Exception("G4DNAMeltonAttachmentModel::Initialise",
104 "em0002",
106 "Model not applicable to particle type.");
107 }
108
109 // Energy limits
110
111 if (LowEnergyLimit() < 4.*eV)
112 {
114 errMsg << "G4DNAMeltonAttachmentModel: low energy limit increased from " <<
115 LowEnergyLimit()/eV << " eV to " << 4. << " eV" << G4endl;
116
117 G4Exception("G4DNAMeltonAttachmentModel::Initialise",
118 "Melton_LowerEBoundary",
120 errMsg);
121
122 SetLowEnergyLimit(4*eV);
123 }
124
125 if (HighEnergyLimit() > 13.*eV)
126 {
128 errMsg << "G4DNAMeltonAttachmentModel: high energy limit decreased from " <<
129 HighEnergyLimit()/eV << " eV to " << 13. << " eV" << G4endl;
130
131 G4Exception("G4DNAMeltonAttachmentModel::Initialise",
132 "Melton_HigherEBoundary",
134 errMsg);
135
136 SetHighEnergyLimit(13.*eV);
137 }
138
139 // Reading of data files
140
141 G4double scaleFactor = 1e-18*cm2;
142
143 // For total cross section
144 G4String fileElectron("dna/sigma_attachment_e_melton");
145
147 eV, scaleFactor);
148 fData->LoadData(fileElectron);
149
150
151#ifdef MELTON_VERBOSE
152 if( verboseLevel >0)
153 {
154 if (verboseLevel > 2)
155 {
156 G4cout << "Loaded cross section data for Melton Attachment model" << G4endl;
157 }
158
159 G4cout << "Melton Attachment model is initialized " << G4endl
160 << "Energy range: "
161 << LowEnergyLimit() / eV << " eV - "
162 << HighEnergyLimit() / eV << " eV"
163 << G4endl;
164 }
165#endif
166
167 // Initialize water density pointer
168 fpWaterDensity = G4DNAMolecularMaterial::Instance()->
169 GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER"));
170
171 if (isInitialised) return;
172
174 isInitialised = true;
175}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4bool LoadData(const G4String &argFileName) override
static G4DNAMolecularMaterial * Instance()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
const G4String & GetParticleName() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()

◆ operator=()

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

◆ SampleSecondaries()

void G4DNAMeltonAttachmentModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * ,
const G4MaterialCutsCouple * ,
const G4DynamicParticle * aDynamicElectron,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 224 of file G4DNAMeltonAttachmentModel.cc.

230{
231
232#ifdef MELTON_VERBOSE
233 if (verboseLevel > 3)
234 G4cout
235 << "Calling SampleSecondaries() of G4DNAMeltonAttachmentModel" << G4endl;
236#endif
237
238 // Electron is killed
239
240 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
241
242 if (!statCode)
243 {
247 }
248
249 else
250 {
253 }
254
255 if(fDissociationFlag)
256 {
258 CreateWaterMolecule(eDissociativeAttachment,
259 -1,
261 }
262 return;
263}
@ eDissociativeAttachment
@ fStopAndKill
static G4DNAChemistryManager * Instance()
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeTrackStatus(G4TrackStatus status)
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNAMeltonAttachmentModel::SelectStationary ( G4bool input)
inline

Definition at line 103 of file G4DNAMeltonAttachmentModel.hh.

104{
105 statCode = input;
106}

◆ SetDissociationFlag()

void G4DNAMeltonAttachmentModel::SetDissociationFlag ( G4bool flag)
inline

Definition at line 91 of file G4DNAMeltonAttachmentModel.hh.

92{
93 fDissociationFlag = flag;
94}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAMeltonAttachmentModel::fParticleChangeForGamma
protected

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