Geant4 9.6.0
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=0, const G4String &nam="DNAMeltonAttachmentModel")
 
virtual ~G4DNAMeltonAttachmentModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetDissociationFlag (G4bool)
 
G4bool GetDissociationFlag ()
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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 44 of file G4DNAMeltonAttachmentModel.hh.

Constructor & Destructor Documentation

◆ G4DNAMeltonAttachmentModel()

G4DNAMeltonAttachmentModel::G4DNAMeltonAttachmentModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAMeltonAttachmentModel" 
)

Definition at line 42 of file G4DNAMeltonAttachmentModel.cc.

44 :G4VEmModel(nam),isInitialised(false)
45{
46// nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
47 fpWaterDensity = 0;
48
49 lowEnergyLimit = 4 * eV;
50 lowEnergyLimitOfModel = 4 * eV;
51 highEnergyLimit = 13 * eV;
52 SetLowEnergyLimit(lowEnergyLimit);
53 SetHighEnergyLimit(highEnergyLimit);
54
55 verboseLevel= 0;
56 // Verbosity scale:
57 // 0 = nothing
58 // 1 = warning for energy non-conservation
59 // 2 = details of energy budget
60 // 3 = calculation of cross sections, file openings, sampling of atoms
61 // 4 = entering in methods
62
63 if( verboseLevel>0 )
64 {
65 G4cout << "Melton Attachment model is constructed " << G4endl
66 << "Energy range: "
67 << lowEnergyLimit / eV << " eV - "
68 << highEnergyLimit / eV << " eV"
69 << G4endl;
70 }
72 fDissociationFlag = true;
73}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592

◆ ~G4DNAMeltonAttachmentModel()

G4DNAMeltonAttachmentModel::~G4DNAMeltonAttachmentModel ( )
virtual

Definition at line 77 of file G4DNAMeltonAttachmentModel.cc.

78{
79 // For total cross section
80
81 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
82
83 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
84 {
85 G4DNACrossSectionDataSet* table = pos->second;
86 delete table;
87 }
88
89 // For final state
90
91}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 163 of file G4DNAMeltonAttachmentModel.cc.

168{
169 if (verboseLevel > 3)
170 G4cout << "Calling CrossSectionPerVolume() of G4DNAMeltonAttachmentModel" << G4endl;
171
172 // Calculate total cross section for model
173
174 G4double sigma=0;
175
176 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
177
178 if(waterDensity!= 0.0)
179 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
180 {
181 const G4String& particleName = particleDefinition->GetParticleName();
182
183 if (ekin >= lowEnergyLimit && ekin < highEnergyLimit)
184 {
185
186 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
187 pos = tableData.find(particleName);
188
189 if (pos != tableData.end())
190 {
191 G4DNACrossSectionDataSet* table = pos->second;
192 if (table != 0)
193 {
194 sigma = table->FindValue(ekin);
195 }
196 }
197 else
198 {
199 G4Exception("G4DNAMeltonAttachmentModel::ComputeCrossSectionPerVolume","em0002",
200 FatalException,"Model not applicable to particle type.");
201 }
202 }
203
204 if (verboseLevel > 2)
205 {
206 G4cout << "__________________________________" << G4endl;
207 G4cout << "°°° G4DNAMeltonAttachmentModel - XS INFO START" << G4endl;
208 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
209 G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
210 G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
211 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
212 G4cout << "°°° G4DNAMeltonAttachmentModel - XS INFO END" << G4endl;
213 }
214
215 } // if water
216
217 return sigma*waterDensity;
218// return sigma*material->GetAtomicNumDensityVector()[1];
219}
@ FatalException
double G4double
Definition: G4Types.hh:64
virtual G4double FindValue(G4double e, G4int componentId=0) const
size_t GetIndex() const
Definition: G4Material.hh:261
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ GetDissociationFlag()

G4bool G4DNAMeltonAttachmentModel::GetDissociationFlag ( )
inline

Definition at line 106 of file G4DNAMeltonAttachmentModel.hh.

107{
108 return fDissociationFlag;
109}

◆ Initialise()

void G4DNAMeltonAttachmentModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 95 of file G4DNAMeltonAttachmentModel.cc.

97{
98
99 if (verboseLevel > 3)
100 G4cout << "Calling G4DNAMeltonAttachmentModel::Initialise()" << G4endl;
101
102 // Energy limits
103
104 if (LowEnergyLimit() < lowEnergyLimit)
105 {
106 G4cout << "G4DNAMeltonAttachmentModel: low energy limit increased from " <<
107 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
108 SetLowEnergyLimit(lowEnergyLimit);
109 }
110
111 if (HighEnergyLimit() > highEnergyLimit)
112 {
113 G4cout << "G4DNAMeltonAttachmentModel: high energy limit decreased from " <<
114 HighEnergyLimit()/eV << " eV to " << highEnergyLimit/eV << " eV" << G4endl;
115 SetHighEnergyLimit(highEnergyLimit);
116 }
117
118 // Reading of data files
119
120 G4double scaleFactor = 1e-18*cm*cm;
121
122 G4String fileElectron("dna/sigma_attachment_e_melton");
123
125 G4String electron;
126
127 // ELECTRON
128
129 // For total cross section
130
131 electron = electronDef->GetParticleName();
132
133 tableFile[electron] = fileElectron;
134
136 tableE->LoadData(fileElectron);
137 tableData[electron] = tableE;
138
139 //
140
141 if (verboseLevel > 2)
142 G4cout << "Loaded cross section data for Melton Attachment model" << G4endl;
143
144 if( verboseLevel>0 )
145 {
146 G4cout << "Melton Attachment model is initialized " << G4endl
147 << "Energy range: "
148 << LowEnergyLimit() / eV << " eV - "
149 << HighEnergyLimit() / eV << " eV"
150 << G4endl;
151 }
152 // Initialize water density pointer
154
155 if (isInitialised) { return; }
157 isInitialised = true;
158
159}
virtual G4bool LoadData(const G4String &argFileName)
static G4DNAMolecularMaterial * Instance()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
const G4String & GetParticleName() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 223 of file G4DNAMeltonAttachmentModel.cc.

228{
229
230 if (verboseLevel > 3)
231 G4cout << "Calling SampleSecondaries() of G4DNAMeltonAttachmentModel" << G4endl;
232
233 // Electron is killed
234
235 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
239
240 if(fDissociationFlag)
241 {
244 }
245 return ;
246}
@ eDissociativeAttachment
@ fStopAndKill
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
G4double GetKineticEnergy() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SetDissociationFlag()

void G4DNAMeltonAttachmentModel::SetDissociationFlag ( G4bool  flag)
inline

Definition at line 101 of file G4DNAMeltonAttachmentModel.hh.

102{
103 fDissociationFlag = flag;
104}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAMeltonAttachmentModel::fParticleChangeForGamma
protected

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