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

#include <G4DNACPA100ExcitationModel.hh>

+ Inheritance diagram for G4DNACPA100ExcitationModel:

Public Member Functions

 G4DNACPA100ExcitationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNACPA100ExcitationModel")
 
virtual ~G4DNACPA100ExcitationModel ()
 
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 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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

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 56 of file G4DNACPA100ExcitationModel.hh.

Constructor & Destructor Documentation

◆ G4DNACPA100ExcitationModel()

G4DNACPA100ExcitationModel::G4DNACPA100ExcitationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNACPA100ExcitationModel" 
)

Definition at line 53 of file G4DNACPA100ExcitationModel.cc.

55:G4VEmModel(nam),isInitialised(false)
56{
57 fpMolWaterDensity = 0;
58
59 SetLowEnergyLimit(11*eV);
60 SetHighEnergyLimit(255955*eV);
61
62 verboseLevel= 0;
63 // Verbosity scale:
64 // 0 = nothing
65 // 1 = warning for energy non-conservation
66 // 2 = details of energy budget
67 // 3 = calculation of cross sections, file openings, sampling of atoms
68 // 4 = entering in methods
69
70 if( verboseLevel>0 )
71 {
72 G4cout << "CPA100 excitation model is constructed " << G4endl;
73 }
75
76 // Selection of stationary mode
77
78 statCode = false;
79}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764

◆ ~G4DNACPA100ExcitationModel()

G4DNACPA100ExcitationModel::~G4DNACPA100ExcitationModel ( )
virtual

Definition at line 83 of file G4DNACPA100ExcitationModel.cc.

84{
85 // Cross section
86
87 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
88 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
89 {
90 G4DNACrossSectionDataSet* table = pos->second;
91 delete table;
92 }
93
94}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 154 of file G4DNACPA100ExcitationModel.cc.

159{
160
161 if (verboseLevel > 3)
162 G4cout << "Calling CrossSectionPerVolume() of G4DNACPA100ExcitationModel" << G4endl;
163
164 if (particleDefinition != G4Electron::ElectronDefinition()) return 0;
165
166 // Calculate total cross section for model
167
168 G4double sigma=0;
169
170 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
171
172 const G4String& particleName = particleDefinition->GetParticleName();
173
174 if (ekin >= LowEnergyLimit() && ekin <= HighEnergyLimit())
175 {
176 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
177 pos = tableData.find(particleName);
178
179 if (pos != tableData.end())
180 {
181 G4DNACrossSectionDataSet* table = pos->second;
182 if (table != 0)
183 {
184 sigma = table->FindValue(ekin);
185 }
186 }
187 else
188 {
189 G4Exception("G4DNACPA100ExcitationModel::CrossSectionPerVolume","em0002",
190 FatalException,"Model not applicable to particle type.");
191 }
192 }
193
194 if (verboseLevel > 2)
195 {
196 G4cout << "__________________________________" << G4endl;
197 G4cout << "G4DNACPA100ExcitationModel - XS INFO START" << G4endl;
198 G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
199 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
200 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
201 // G4cout << " - Cross section per water molecule (cm^-1)="
202 // << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
203 G4cout << "G4DNACPA100ExcitationModel - XS INFO END" << G4endl;
204 }
205
206 return sigma*waterDensity;
207
208}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
virtual G4double FindValue(G4double e, G4int componentId=0) const
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
size_t GetIndex() const
Definition: G4Material.hh:258
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 98 of file G4DNACPA100ExcitationModel.cc.

100{
101
102 if (verboseLevel > 3)
103 G4cout << "Calling G4DNACPA100ExcitationModel::Initialise()" << G4endl;
104
105 G4String fileElectron("dna/sigma_excitation_e_cpa100");
106
107 G4double scaleFactor = 1.e-20 *m*m;
108
109 // *** ELECTRON
110
113 electron = electronDef->GetParticleName();
114
115 tableFile[electron] = fileElectron;
116
117 // Cross section
118
120 = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV, scaleFactor );
121
122 /*
123 G4DNACrossSectionDataSet* tableE =
124 new G4DNACrossSectionDataSet(new G4DNACPA100LogLogInterpolation, eV, scaleFactor );
125 */
126
127 tableE->LoadData(fileElectron);
128
129 tableData[electron] = tableE;
130
131 //
132
133 if( verboseLevel>0 )
134 {
135 G4cout << "CPA100 excitation model is initialized " << G4endl
136 << "Energy range: "
137 << LowEnergyLimit() / eV << " eV - "
138 << HighEnergyLimit() / keV << " keV for "
139 << particle->GetParticleName()
140 << G4endl;
141 }
142
143 // Initialize water density pointer
144 fpMolWaterDensity =
146
147 if (isInitialised) return;
149 isInitialised = true;
150}
virtual G4bool LoadData(const G4String &argFileName)
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()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:651
const G4String & GetParticleName() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 212 of file G4DNACPA100ExcitationModel.cc.

217{
218
219 if (verboseLevel > 3)
220 G4cout << "Calling SampleSecondaries() of G4DNACPA100ExcitationModel" << G4endl;
221
222 G4double k = aDynamicParticle->GetKineticEnergy();
223
224 const G4String& particleName = aDynamicParticle->GetDefinition()->GetParticleName();
225
226 G4int level = RandomSelect(k,particleName);
227 G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
228 G4double newEnergy = k - excitationEnergy;
229
230 if (newEnergy > 0)
231 {
232 // fParticleChangeForGamma->ProposeMomentumDirection(aDynamicParticle->GetMomentumDirection());
233
234 // We take into account direction change as described page 87 (II.92) in thesis by S. Edel
235
236 G4double cosTheta =
237
238 (excitationEnergy/k) / (1. + (k/(2*electron_mass_c2))*(1.-excitationEnergy/k) );
239
240 cosTheta = std::sqrt(1.-cosTheta);
241
242 G4double phi = 2. * pi * G4UniformRand();
243
244 G4ThreeVector zVers = aDynamicParticle->GetMomentumDirection();
245
246 //G4ThreeVector xVers = zVers.orthogonal();
247 //G4ThreeVector yVers = zVers.cross(xVers);
248 //G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
249 //G4double yDir = xDir;
250 //xDir *= std::cos(phi);
251 //yDir *= std::sin(phi);
252 // G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
253
254 // Computation of scattering angles (from Subroutine DIRAN in CPA100)
255
256 G4double CT1, ST1, CF1, SF1, CT2, ST2, CF2, SF2;
257 G4double sinTheta = std::sqrt (1-cosTheta*cosTheta);
258
259 CT1=0;
260 ST1=0;
261 CF1=0;
262 SF1=0;
263 CT2=0;
264 ST2=0;
265 CF2=0;
266 SF2=0;
267
268 CT1 = zVers.z();
269 ST1=std::sqrt(1.-CT1*CT1);
270
271 if (ST1!=0) CF1 = zVers.x()/ST1; else CF1 = std::cos(2. * pi * G4UniformRand());
272 if (ST1!=0) SF1 = zVers.y()/ST1; else SF1 = std::sqrt(1.-CF1*CF1);
273
274 G4double A3, A4, A5, A2, A1;
275 A3=0;
276 A4=0;
277 A5=0;
278 A2=0;
279 A1=0;
280
281 A3 = sinTheta*std::cos(phi);
282 A4 = A3*CT1 + ST1*cosTheta;
283 A5 = sinTheta * std::sin(phi);
284 A2 = A4 * SF1 + A5 * CF1;
285 A1 = A4 * CF1 - A5 * SF1;
286
287 CT2 = CT1*cosTheta - ST1*A3;
288 ST2 = std::sqrt(1.-CT2*CT2);
289
290 if (ST2==0) ST2=1E-6;
291 CF2 = A1/ST2;
292 SF2 = A2/ST2;
293
294 /*
295 G4cout << "CT1=" << CT1 << G4endl;
296 G4cout << "ST1=" << ST1 << G4endl;
297 G4cout << "CF1=" << CF1 << G4endl;
298 G4cout << "SF1=" << SF1 << G4endl;
299 G4cout << "cosTheta=" << cosTheta << G4endl;
300 G4cout << "sinTheta=" << sinTheta << G4endl;
301 G4cout << "cosPhi=" << std::cos(phi) << G4endl;
302 G4cout << "sinPhi=" << std::sin(phi) << G4endl;
303 G4cout << "CT2=" << CT2 << G4endl;
304 G4cout << "ST2=" << ST2 << G4endl;
305 G4cout << "CF2=" << CF2 << G4endl;
306 G4cout << "SF2=" << SF2 << G4endl;
307 */
308
309 G4ThreeVector zPrimeVers(ST2*CF2,ST2*SF2,CT2);
310
311 //
312
314
315 //
316
317 if (!statCode) fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
319
321 }
322
323 // Chemistry
324
325 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
327 level,
328 theIncomingTrack);
329}
@ eExcitedMolecule
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
double z() const
double x() const
double y() const
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
const G4double pi

◆ SelectStationary()

void G4DNACPA100ExcitationModel::SelectStationary ( G4bool  input)
inline

Definition at line 124 of file G4DNACPA100ExcitationModel.hh.

125{
126 statCode = input;
127}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNACPA100ExcitationModel::fParticleChangeForGamma
protected

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