Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNAEmfietzoglouExcitationModel Class Reference

#include <G4DNAEmfietzoglouExcitationModel.hh>

Inheritance diagram for G4DNAEmfietzoglouExcitationModel:
Collaboration diagram for G4DNAEmfietzoglouExcitationModel:

Public Member Functions

 G4DNAEmfietzoglouExcitationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAEmfietzoglouExcitationModel")
 
virtual ~G4DNAEmfietzoglouExcitationModel ()
 
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 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 *, 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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

Protected Attributes

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

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 *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 50 of file G4DNAEmfietzoglouExcitationModel.hh.

Constructor & Destructor Documentation

G4DNAEmfietzoglouExcitationModel::G4DNAEmfietzoglouExcitationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAEmfietzoglouExcitationModel" 
)

Definition at line 47 of file G4DNAEmfietzoglouExcitationModel.cc.

49  :G4VEmModel(nam),isInitialised(false)
50 {
51 
52  // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
53  fpMolWaterDensity = 0;
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 << "Emfietzoglou excitation model is constructed " << G4endl;
66  }
68  SetLowEnergyLimit(8. * eV);
69  SetHighEnergyLimit(10. * keV);
70 
71  // Selection of stationary mode
72 
73  statCode = false;
74 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4DNAEmfietzoglouExcitationModel::~G4DNAEmfietzoglouExcitationModel ( )
virtual

Definition at line 78 of file G4DNAEmfietzoglouExcitationModel.cc.

79 {
80  // Cross section
81 
82  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
83  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
84  {
85  G4DNACrossSectionDataSet* table = pos->second;
86  delete table;
87  }
88 
89 }
static const G4double pos

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 152 of file G4DNAEmfietzoglouExcitationModel.cc.

157 {
158  if (verboseLevel > 3)
159  G4cout << "Calling CrossSectionPerVolume() of G4DNAEmfietzoglouExcitationModel" << G4endl;
160 
161  if (particleDefinition != G4Electron::ElectronDefinition()) return 0;
162 
163  // Calculate total cross section for model
164 
165  G4double lowLim = 0;
166  G4double highLim = 0;
167  G4double sigma=0;
168 
169  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
170 
171  if(waterDensity!= 0.0)
172  {
173  const G4String& particleName = particleDefinition->GetParticleName();
174 
175  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
176  pos1 = lowEnergyLimit.find(particleName);
177  if (pos1 != lowEnergyLimit.end())
178  {
179  lowLim = pos1->second;
180  }
181 
182  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
183  pos2 = highEnergyLimit.find(particleName);
184  if (pos2 != highEnergyLimit.end())
185  {
186  highLim = pos2->second;
187  }
188 
189  if (ekin >= lowLim && ekin < highLim)
190  {
191  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
192  pos = tableData.find(particleName);
193 
194  if (pos != tableData.end())
195  {
196  G4DNACrossSectionDataSet* table = pos->second;
197  if (table != 0)
198  {
199  sigma = table->FindValue(ekin);
200  }
201  }
202  else
203  {
204  G4Exception("G4DNAEmfietzoglouExcitationModel::CrossSectionPerVolume","em0002",
205  FatalException,"Model not applicable to particle type.");
206  }
207  }
208 
209  if (verboseLevel > 2)
210  {
211  G4cout << "__________________________________" << G4endl;
212  G4cout << "G4DNAEmfietzoglouExcitationModel - XS INFO START" << G4endl;
213  G4cout << "Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
214  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
215  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
216  // G4cout << " Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
217  G4cout << "G4DNAEmfietzoglouExcitationModel - XS INFO END" << G4endl;
218  }
219 
220  } // if (waterMaterial)
221 
222  return sigma*waterDensity;
223 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
size_t GetIndex() const
Definition: G4Material.hh:262
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
virtual G4double FindValue(G4double e, G4int componentId=0) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static const G4double pos

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 93 of file G4DNAEmfietzoglouExcitationModel.cc.

95 {
96 
97  if (verboseLevel > 3)
98  G4cout << "Calling G4DNAEmfietzoglouExcitationModel::Initialise()" << G4endl;
99 
100  G4String fileElectron("dna/sigma_excitation_e_emfietzoglou");
101 
103 
105 
106  G4double scaleFactor = (1.e-22 / 3.343) * m*m;
107 
108  // *** ELECTRON
109 
110  electron = electronDef->GetParticleName();
111 
112  tableFile[electron] = fileElectron;
113 
114  lowEnergyLimit[electron] = 8. * eV;
115  highEnergyLimit[electron] = 10. * keV;
116 
117  // Cross section
118 
120  tableE->LoadData(fileElectron);
121 
122  tableData[electron] = tableE;
123 
124  //
125 
126  if (particle==electronDef)
127  {
128  SetLowEnergyLimit(lowEnergyLimit[electron]);
129  SetHighEnergyLimit(highEnergyLimit[electron]);
130  }
131 
132  if( verboseLevel>0 )
133  {
134  G4cout << "Emfietzoglou excitation model is initialized " << G4endl
135  << "Energy range: "
136  << LowEnergyLimit() / eV << " eV - "
137  << HighEnergyLimit() / keV << " keV for "
138  << particle->GetParticleName()
139  << G4endl;
140  }
141 
142  // Initialize water density pointer
144 
145  if (isInitialised) { return; }
147  isInitialised = true;
148 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
virtual G4bool LoadData(const G4String &argFileName)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static constexpr double m
Definition: G4SIunits.hh:129
static constexpr double eV
Definition: G4SIunits.hh:215
static G4DNAMolecularMaterial * Instance()
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
static constexpr double keV
Definition: G4SIunits.hh:216
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 227 of file G4DNAEmfietzoglouExcitationModel.cc.

232 {
233 
234  if (verboseLevel > 3)
235  G4cout << "Calling SampleSecondaries() of G4DNAEmfietzoglouExcitationModel" << G4endl;
236 
237  G4double k = aDynamicParticle->GetKineticEnergy();
238 
239  const G4String& particleName = aDynamicParticle->GetDefinition()->GetParticleName();
240 
241  G4int level = RandomSelect(k,particleName);
242  G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
243  G4double newEnergy = k - excitationEnergy;
244 
245  if (newEnergy > 0)
246  {
248 
249  if (!statCode) fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy);
251 
253  }
254 
255  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
257  level,
258  theIncomingTrack);
259 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4DNAEmfietzoglouExcitationModel::SelectStationary ( G4bool  input)
inline

Definition at line 120 of file G4DNAEmfietzoglouExcitationModel.hh.

121 {
122  statCode = input;
123 }

Member Data Documentation

G4ParticleChangeForGamma* G4DNAEmfietzoglouExcitationModel::fParticleChangeForGamma
protected

Definition at line 80 of file G4DNAEmfietzoglouExcitationModel.hh.


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