Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Attributes
G4MuElecElasticModel Class Reference

#include <G4MuElecElasticModel.hh>

Inheritance diagram for G4MuElecElasticModel:
G4VEmModel

Public Member Functions

 G4MuElecElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="MuElecElasticModel")
 
virtual ~G4MuElecElasticModel ()
 
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 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 *, 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 *)
 
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)
 
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=0)
 
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
 

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 *)
 

Detailed Description

Definition at line 51 of file G4MuElecElasticModel.hh.

Constructor & Destructor Documentation

G4MuElecElasticModel::G4MuElecElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "MuElecElasticModel" 
)

Definition at line 48 of file G4MuElecElasticModel.cc.

References python.hepunit::eV, G4NistManager::FindOrBuildMaterial(), fParticleChangeForGamma, G4cout, G4endl, G4NistManager::Instance(), python.hepunit::keV, python.hepunit::MeV, G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().

50 :G4VEmModel(nam),isInitialised(false)
51 {
52 
53  G4cout << G4endl;
54  G4cout << "*******************************************************************************" << G4endl;
55  G4cout << "*******************************************************************************" << G4endl;
56  G4cout << " The name of the class G4MuElecElasticModel is changed to G4MicroElecElasticModel. " << G4endl;
57  G4cout << " The obsolete class will be REMOVED with the next release of Geant4. " << G4endl;
58  G4cout << "*******************************************************************************" << G4endl;
59  G4cout << "*******************************************************************************" << G4endl;
60  G4cout << G4endl;
61 
62  nistSi = G4NistManager::Instance()->FindOrBuildMaterial("G4_Si");
63 
64  killBelowEnergy = 16.7 * eV; // Minimum e- energy for energy loss by excitation
65  lowEnergyLimit = 0 * eV;
66  lowEnergyLimitOfModel = 5 * eV; // The model lower energy is 5 eV
67  highEnergyLimit = 100. * MeV;
68  SetLowEnergyLimit(lowEnergyLimit);
69  SetHighEnergyLimit(highEnergyLimit);
70 
71  verboseLevel= 0;
72  // Verbosity scale:
73  // 0 = nothing
74  // 1 = warning for energy non-conservation
75  // 2 = details of energy budget
76  // 3 = calculation of cross sections, file openings, sampling of atoms
77  // 4 = entering in methods
78 
79  if( verboseLevel>0 )
80  {
81  G4cout << "MuElec Elastic model is constructed " << G4endl
82  << "Energy range: "
83  << lowEnergyLimit / eV << " eV - "
84  << highEnergyLimit / keV << " keV"
85  << G4endl;
86  }
88 }
G4Material * FindOrBuildMaterial(const G4String &name, G4bool isotopes=true, G4bool warning=false)
static G4NistManager * Instance()
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4GLOB_DLL std::ostream G4cout
G4ParticleChangeForGamma * fParticleChangeForGamma
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4MuElecElasticModel::~G4MuElecElasticModel ( )
virtual

Definition at line 92 of file G4MuElecElasticModel.cc.

93 {
94  // For total cross section
95 
96  std::map< G4String,G4MuElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
97  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
98  {
99  G4MuElecCrossSectionDataSet* table = pos->second;
100  delete table;
101  }
102 
103  // For final state
104 
105  eVecm.clear();
106 
107 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 217 of file G4MuElecElasticModel.cc.

References python.hepunit::cm, density, python.hepunit::eV, FatalException, G4MuElecCrossSectionDataSet::FindValue(), G4cout, G4endl, G4Exception(), G4Material::GetBaseMaterial(), G4ParticleDefinition::GetParticleName(), and G4Material::GetTotNbOfAtomsPerVolume().

222 {
223  if (verboseLevel > 3)
224  G4cout << "Calling CrossSectionPerVolume() of G4MuElecElasticModel" << G4endl;
225 
226  // Calculate total cross section for model
227 
228  G4double sigma=0;
229 
231 
232  if (material == nistSi || material->GetBaseMaterial() == nistSi)
233  {
234  const G4String& particleName = p->GetParticleName();
235 
236  if (ekin < highEnergyLimit)
237  {
238  //SI : XS must not be zero otherwise sampling of secondaries method ignored
239  if (ekin < lowEnergyLimitOfModel) ekin = lowEnergyLimitOfModel;
240  //
241 
242  std::map< G4String,G4MuElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
243  pos = tableData.find(particleName);
244 
245  if (pos != tableData.end())
246  {
247  G4MuElecCrossSectionDataSet* table = pos->second;
248  if (table != 0)
249  {
250  sigma = table->FindValue(ekin);
251  }
252  }
253  else
254  {
255  G4Exception("G4MuElecElasticModel::ComputeCrossSectionPerVolume","em0002",FatalException,"Model not applicable to particle type.");
256  }
257  }
258 
259  if (verboseLevel > 3)
260  {
261  G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl;
262  G4cout << " - Cross section per Si atom (cm^2)=" << sigma/cm/cm << G4endl;
263  G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density/(1./cm) << G4endl;
264  }
265 
266  }
267 
268  return sigma*density;
269 }
virtual G4double FindValue(G4double e, G4int componentId=0) const
const G4String & GetParticleName() const
G4double density
Definition: TRTMaterials.hh:39
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:207
const G4Material * GetBaseMaterial() const
Definition: G4Material.hh:231
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4MuElecElasticModel::GetKillBelowThreshold ( )
inline

Definition at line 76 of file G4MuElecElasticModel.hh.

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

Implements G4VEmModel.

Definition at line 111 of file G4MuElecElasticModel.cc.

References python.hepunit::cm, G4InuclParticleNames::electron, G4Electron::ElectronDefinition(), python.hepunit::eV, FatalException, fParticleChangeForGamma, G4cout, G4endl, G4Exception(), G4VEmModel::GetParticleChangeForGamma(), G4ParticleDefinition::GetParticleName(), G4VEmModel::HighEnergyLimit(), python.hepunit::keV, G4MuElecCrossSectionDataSet::LoadData(), G4VEmModel::LowEnergyLimit(), python.hepunit::MeV, G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().

113 {
114 
115  if (verboseLevel > 3)
116  G4cout << "Calling G4MuElecElasticModel::Initialise()" << G4endl;
117 
118  // Energy limits
119 
120  if (LowEnergyLimit() < lowEnergyLimit)
121  {
122  G4cout << "G4MuElecElasticModel: low energy limit increased from " <<
123  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
124  SetLowEnergyLimit(lowEnergyLimit);
125  }
126 
127  if (HighEnergyLimit() > highEnergyLimit)
128  {
129  G4cout << "G4MuElecElasticModel: high energy limit decreased from " <<
130  HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
131  SetHighEnergyLimit(highEnergyLimit);
132  }
133 
134  // Reading of data files
135 
136  G4double scaleFactor = 1e-18 * cm * cm;
137 
138  G4String fileElectron("microelec/sigma_elastic_e_Si");
139 
142 
143  // For total cross section
144 
145  electron = electronDef->GetParticleName();
146 
147  tableFile[electron] = fileElectron;
148 
150  tableE->LoadData(fileElectron);
151  tableData[electron] = tableE;
152 
153  // For final state
154 
155  char *path = getenv("G4LEDATA");
156 
157  if (!path)
158  {
159  G4Exception("G4MuElecElasticModel::Initialise","em0006",FatalException,"G4LEDATA environment variable not set.");
160  return;
161  }
162 
163  std::ostringstream eFullFileName;
164  eFullFileName << path << "/microelec/sigmadiff_elastic_e_Si.dat";
165  std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
166 
167  if (!eDiffCrossSection)
168  G4Exception("G4MuElecElasticModel::Initialise","em0003",FatalException,"Missing data file: /microelec/sigmadiff_elastic_e_Si.dat");
169 
170  eTdummyVec.push_back(0.);
171 
172  while(!eDiffCrossSection.eof())
173  {
174  double tDummy;
175  double eDummy;
176  eDiffCrossSection>>tDummy>>eDummy;
177 
178  // SI : mandatory eVecm initialization
179  if (tDummy != eTdummyVec.back())
180  {
181  eTdummyVec.push_back(tDummy);
182  eVecm[tDummy].push_back(0.);
183  }
184 
185  eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy];
186 
187  // SI : only if not end of file reached !
188  if (!eDiffCrossSection.eof()) eDiffCrossSectionData[tDummy][eDummy]*=scaleFactor;
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 MuElec Elastic model" << G4endl;
198 
199  if( verboseLevel>0 )
200  {
201  G4cout << "MuElec Elastic model is initialized " << G4endl
202  << "Energy range: "
203  << LowEnergyLimit() / eV << " eV - "
204  << HighEnergyLimit() / keV << " keV"
205  << G4endl;
206  }
207 
208  if (isInitialised) { return; }
210  isInitialised = true;
211 
212  // InitialiseElementSelectors(particle,cuts);
213 }
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
virtual G4bool LoadData(const G4String &argFileName)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4GLOB_DLL std::ostream G4cout
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
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4MuElecElasticModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicElectron,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 273 of file G4MuElecElasticModel.cc.

References CLHEP::Hep3Vector::cross(), fParticleChangeForGamma, fStopAndKill, G4cout, G4endl, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), CLHEP::Hep3Vector::orthogonal(), python.hepunit::pi, G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), and CLHEP::Hep3Vector::unit().

278 {
279 
280  if (verboseLevel > 3)
281  G4cout << "Calling SampleSecondaries() of G4MuElecElasticModel" << G4endl;
282 
283  G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
284 
285  if (electronEnergy0 < killBelowEnergy)
286  {
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:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector orthogonal() const
G4ParticleChangeForGamma * fParticleChangeForGamma
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void G4MuElecElasticModel::SetKillBelowThreshold ( G4double  threshold)
inline

Definition at line 138 of file G4MuElecElasticModel.hh.

References G4Exception(), and JustWarning.

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

Field Documentation

G4ParticleChangeForGamma* G4MuElecElasticModel::fParticleChangeForGamma
protected

Definition at line 80 of file G4MuElecElasticModel.hh.

Referenced by G4MuElecElasticModel(), Initialise(), and SampleSecondaries().


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