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

#include <G4MicroElecElasticModel.hh>

Inheritance diagram for G4MicroElecElasticModel:
Collaboration diagram for G4MicroElecElasticModel:

Public Member Functions

 G4MicroElecElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="MicroElecElasticModel")
 
virtual ~G4MicroElecElasticModel ()
 
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 SetKillBelowThreshold (G4double threshold)
 
G4double GetKillBelowThreshold ()
 
- 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 51 of file G4MicroElecElasticModel.hh.

Constructor & Destructor Documentation

G4MicroElecElasticModel::G4MicroElecElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "MicroElecElasticModel" 
)

Definition at line 49 of file G4MicroElecElasticModel.cc.

51 :G4VEmModel(nam),isInitialised(false)
52 {
53  nistSi = G4NistManager::Instance()->FindOrBuildMaterial("G4_Si");
54 
55  killBelowEnergy = 16.7 * eV; // Minimum e- energy for energy loss by excitation
56  lowEnergyLimit = 0 * eV;
57  lowEnergyLimitOfModel = 5 * eV; // The model lower energy is 5 eV
58  highEnergyLimit = 100. * MeV;
59  SetLowEnergyLimit(lowEnergyLimit);
60  SetHighEnergyLimit(highEnergyLimit);
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 << "MicroElec Elastic model is constructed " << G4endl
73  << "Energy range: "
74  << lowEnergyLimit / eV << " eV - "
75  << highEnergyLimit / MeV << " MeV"
76  << G4endl;
77  }
79 }
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
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
G4ParticleChangeForGamma * fParticleChangeForGamma
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731

Here is the call graph for this function:

G4MicroElecElasticModel::~G4MicroElecElasticModel ( )
virtual

Definition at line 83 of file G4MicroElecElasticModel.cc.

84 {
85  // For total cross section
86 
87  std::map< G4String,G4MicroElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
88  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
89  {
90  G4MicroElecCrossSectionDataSet* table = pos->second;
91  delete table;
92  }
93 
94  // For final state
95 
96  eVecm.clear();
97 
98 }
static const G4double pos

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 216 of file G4MicroElecElasticModel.cc.

221 {
222  if (verboseLevel > 3)
223  G4cout << "Calling CrossSectionPerVolume() of G4MicroElecElasticModel" << G4endl;
224 
225  // Calculate total cross section for model
226 
227  G4double sigma=0;
228 
229  G4double density = material->GetTotNbOfAtomsPerVolume();
230 
231  if (material == nistSi || material->GetBaseMaterial() == nistSi)
232  {
233  const G4String& particleName = p->GetParticleName();
234 
235  if (ekin < highEnergyLimit)
236  {
237  //SI : XS must not be zero otherwise sampling of secondaries method ignored
238  if (ekin < killBelowEnergy) return DBL_MAX;
239  //
240 
241  std::map< G4String,G4MicroElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
242  pos = tableData.find(particleName);
243 
244  if (pos != tableData.end())
245  {
246  G4MicroElecCrossSectionDataSet* table = pos->second;
247  if (table != 0)
248  {
249  sigma = table->FindValue(ekin);
250  }
251  }
252  else
253  {
254  G4Exception("G4MicroElecElasticModel::ComputeCrossSectionPerVolume","em0002",FatalException,"Model not applicable to particle type.");
255  }
256  }
257 
258  if (verboseLevel > 3)
259  {
260  G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl;
261  G4cout << " - Cross section per Si atom (cm^2)=" << sigma/cm/cm << G4endl;
262  G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density/(1./cm) << G4endl;
263  }
264 
265  }
266 
267  return sigma*density;
268 }
const G4String & GetParticleName() const
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
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:209
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:233
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
static const G4double pos

Here is the call graph for this function:

G4double G4MicroElecElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 76 of file G4MicroElecElasticModel.hh.

76 { return killBelowEnergy; }
void G4MicroElecElasticModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 102 of file G4MicroElecElasticModel.cc.

104 {
105 
106  if (verboseLevel > 3)
107  G4cout << "Calling G4MicroElecElasticModel::Initialise()" << G4endl;
108 
109  // Energy limits
110 
111  if (LowEnergyLimit() < lowEnergyLimit)
112  {
113  G4cout << "G4MicroElecElasticModel: low energy limit increased from " <<
114  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
115  SetLowEnergyLimit(lowEnergyLimit);
116  }
117 
118  if (HighEnergyLimit() > highEnergyLimit)
119  {
120  G4cout << "G4MicroElecElasticModel: high energy limit decreased from " <<
121  HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
122  SetHighEnergyLimit(highEnergyLimit);
123  }
124 
125  // Reading of data files
126 
127  G4double scaleFactor = 1e-18 * cm * cm;
128 
129  G4String fileElectron("microelec/sigma_elastic_e_Si");
130 
133 
134  // For total cross section
135 
136  electron = electronDef->GetParticleName();
137 
138  tableFile[electron] = fileElectron;
139 
141  tableE->LoadData(fileElectron);
142  tableData[electron] = tableE;
143 
144  // For final state
145 
146  char *path = getenv("G4LEDATA");
147 
148  if (!path)
149  {
150  G4Exception("G4MicroElecElasticModel::Initialise","em0006",FatalException,"G4LEDATA environment variable not set.");
151  return;
152  }
153 
154  std::ostringstream eFullFileName;
155  eFullFileName << path << "/microelec/sigmadiff_cumulated_elastic_e_Si.dat";
156  std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
157 
158  if (!eDiffCrossSection)
159  G4Exception("G4MicroElecElasticModel::Initialise","em0003",FatalException,"Missing data file: /microelec/sigmadiff_cumulated_elastic_e_Si.dat");
160 
161 
162  // October 21th, 2014 - Melanie Raine
163  // Added clear for MT
164 
165  eTdummyVec.clear();
166  eVecm.clear();
167  eDiffCrossSectionData.clear();
168 
169  //
170 
171 
172  eTdummyVec.push_back(0.);
173 
174  while(!eDiffCrossSection.eof())
175  {
176  double tDummy;
177  double eDummy;
178  eDiffCrossSection>>tDummy>>eDummy;
179 
180  // SI : mandatory eVecm initialization
181 
182  if (tDummy != eTdummyVec.back())
183  {
184  eTdummyVec.push_back(tDummy);
185  eVecm[tDummy].push_back(0.);
186  }
187 
188  eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy];
189 
190  if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
191 
192  }
193 
194  // End final state
195 
196  if (verboseLevel > 2)
197  G4cout << "Loaded cross section files for MicroElec Elastic model" << G4endl;
198 
199  if( verboseLevel>0 )
200  {
201  G4cout << "MicroElec Elastic model is initialized " << G4endl
202  << "Energy range: "
203  << LowEnergyLimit() / eV << " eV - "
204  << HighEnergyLimit() / MeV << " MeV"
205  << G4endl;
206  }
207 
208  if (isInitialised) { return; }
210  isInitialised = true;
211 
212 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
virtual G4bool LoadData(const G4String &argFileName)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ParticleChangeForGamma * fParticleChangeForGamma
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 272 of file G4MicroElecElasticModel.cc.

277 {
278 
279  if (verboseLevel > 3)
280  G4cout << "Calling SampleSecondaries() of G4MicroElecElasticModel" << G4endl;
281 
282  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
283 
284  if (electronEnergy0 < killBelowEnergy)
285  {
289  return ;
290  }
291 
292  if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
293  {
294  G4double cosTheta = RandomizeCosTheta(electronEnergy0);
295 
296  G4double phi = 2. * pi * G4UniformRand();
297 
298  G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
299  G4ThreeVector xVers = zVers.orthogonal();
300  G4ThreeVector yVers = zVers.cross(xVers);
301 
302  G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
303  G4double yDir = xDir;
304  xDir *= std::cos(phi);
305  yDir *= std::sin(phi);
306 
307  G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
308 
310 
312  }
313 
314 }
G4double GetKineticEnergy() const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
G4ParticleChangeForGamma * fParticleChangeForGamma
Hep3Vector orthogonal() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double pi
Definition: G4SIunits.hh:75
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)

Here is the call graph for this function:

void G4MicroElecElasticModel::SetKillBelowThreshold ( G4double  threshold)
inline

Definition at line 140 of file G4MicroElecElasticModel.hh.

141 {
142  killBelowEnergy = threshold;
143 
144  if (threshold < 5*CLHEP::eV)
145  {
146  G4Exception ("*** WARNING : the G4MicroElecElasticModel class is not validated below 5 eV !","",JustWarning,"") ;
147  threshold = 5*CLHEP::eV;
148  }
149 
150 }
static constexpr double eV
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Here is the call graph for this function:

Member Data Documentation

G4ParticleChangeForGamma* G4MicroElecElasticModel::fParticleChangeForGamma
protected

Definition at line 80 of file G4MicroElecElasticModel.hh.


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