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

#include <G4LivermoreRayleighModel.hh>

Inheritance diagram for G4LivermoreRayleighModel:
G4VEmModel

Public Member Functions

 G4LivermoreRayleighModel ()
 
virtual ~G4LivermoreRayleighModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetLowEnergyThreshold (G4double)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
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 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
 

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 *)
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Detailed Description

Definition at line 39 of file G4LivermoreRayleighModel.hh.

Constructor & Destructor Documentation

G4LivermoreRayleighModel::G4LivermoreRayleighModel ( )

Definition at line 43 of file G4LivermoreRayleighModel.cc.

References python.hepunit::eV, G4cout, G4endl, and G4VEmModel::SetAngularDistribution().

44  :G4VEmModel("LivermoreRayleigh"),isInitialised(false)
45 {
46  fParticleChange = 0;
47  lowEnergyLimit = 10 * eV;
48 
50 
51  verboseLevel= 0;
52  // Verbosity scale for debugging purposes:
53  // 0 = nothing
54  // 1 = calculation of cross sections, file openings...
55  // 2 = entering in methods
56 
57  if(verboseLevel > 0)
58  {
59  G4cout << "G4LivermoreRayleighModel is constructed " << G4endl;
60  }
61 
62 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4GLOB_DLL std::ostream G4cout
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:585
#define G4endl
Definition: G4ios.hh:61
G4LivermoreRayleighModel::~G4LivermoreRayleighModel ( )
virtual

Definition at line 66 of file G4LivermoreRayleighModel.cc.

67 {}

Member Function Documentation

G4double G4LivermoreRayleighModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 184 of file G4LivermoreRayleighModel.cc.

References G4PhysicsVector::Energy(), G4cout, G4endl, G4lrint(), G4PhysicsVector::GetVectorLength(), InitialiseForElement(), python.hepunit::MeV, n, and G4PhysicsVector::Value().

189 {
190  if (verboseLevel > 1)
191  {
192  G4cout << "G4LivermoreRayleighModel::ComputeCrossSectionPerAtom()"
193  << G4endl;
194  }
195 
196  if(GammaEnergy < lowEnergyLimit) { return 0.0; }
197 
198  G4double xs = 0.0;
199 
200  G4int intZ = G4lrint(Z);
201 
202  if(intZ < 1 || intZ > maxZ) { return xs; }
203 
204  G4LPhysicsFreeVector* pv = dataCS[intZ];
205 
206  // if element was not initialised
207  // do initialisation safely for MT mode
208  if(!pv) {
209  InitialiseForElement(0, intZ);
210  pv = dataCS[intZ];
211  if(!pv) { return xs; }
212  }
213 
214  G4int n = pv->GetVectorLength() - 1;
215  G4double e = GammaEnergy/MeV;
216  if(e >= pv->Energy(n)) {
217  xs = (*pv)[n]/(e*e);
218  } else if(e >= pv->Energy(0)) {
219  xs = pv->Value(e)/(e*e);
220  }
221 
222  if(verboseLevel > 0)
223  {
224  G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
225  << e << G4endl;
226  G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
227  G4cout << " -> first E*E*cs value in CS data file (iu) =" << (*pv)[0]
228  << G4endl;
229  G4cout << " -> last E*E*cs value in CS data file (iu) =" << (*pv)[n]
230  << G4endl;
231  G4cout << "*********************************************************"
232  << G4endl;
233  }
234  return xs;
235 }
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
const G4int n
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
int G4lrint(double ad)
Definition: templates.hh:163
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void G4LivermoreRayleighModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 71 of file G4LivermoreRayleighModel.cc.

References python.hepunit::eV, G4cout, G4endl, G4lrint(), G4Material::GetElementVector(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4Material::GetNumberOfElements(), G4VEmModel::GetParticleChangeForGamma(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), python.hepunit::GeV, G4VEmModel::HighEnergyLimit(), G4VEmModel::InitialiseElementSelectors(), G4VEmModel::IsMaster(), G4VEmModel::LowEnergyLimit(), and eplot::material.

73 {
74  if (verboseLevel > 1)
75  {
76  G4cout << "Calling Initialise() of G4LivermoreRayleighModel." << G4endl
77  << "Energy range: "
78  << LowEnergyLimit() / eV << " eV - "
79  << HighEnergyLimit() / GeV << " GeV"
80  << G4endl;
81  }
82 
83  if(IsMaster()) {
84 
85  // Initialise element selector
86  InitialiseElementSelectors(particle, cuts);
87 
88  // Access to elements
89  char* path = getenv("G4LEDATA");
90  G4ProductionCutsTable* theCoupleTable =
92  G4int numOfCouples = theCoupleTable->GetTableSize();
93 
94  for(G4int i=0; i<numOfCouples; ++i)
95  {
96  const G4MaterialCutsCouple* couple =
97  theCoupleTable->GetMaterialCutsCouple(i);
98  const G4Material* material = couple->GetMaterial();
99  const G4ElementVector* theElementVector = material->GetElementVector();
100  G4int nelm = material->GetNumberOfElements();
101 
102  for (G4int j=0; j<nelm; ++j)
103  {
104  G4int Z = G4lrint((*theElementVector)[j]->GetZ());
105  if(Z < 1) { Z = 1; }
106  else if(Z > maxZ) { Z = maxZ; }
107  if( (!dataCS[Z]) ) { ReadData(Z, path); }
108  }
109  }
110  }
111 
112  if(isInitialised) { return; }
113  fParticleChange = GetParticleChangeForGamma();
114  isInitialised = true;
115 
116 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
std::vector< G4Element * > G4ElementVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
G4GLOB_DLL std::ostream G4cout
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
int G4lrint(double ad)
Definition: templates.hh:163
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
const G4Material * GetMaterial() const
void G4LivermoreRayleighModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 280 of file G4LivermoreRayleighModel.cc.

References G4TemplateAutoLock< M, L, U >::unlock().

Referenced by ComputeCrossSectionPerAtom().

282 {
283  G4AutoLock l(&LivermoreRayleighModelMutex);
284  // G4cout << "G4LivermoreRayleighModel::InitialiseForElement Z= "
285  // << Z << G4endl;
286  if(!dataCS[Z]) { ReadData(Z); }
287  l.unlock();
288 }
void G4LivermoreRayleighModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 120 of file G4LivermoreRayleighModel.cc.

References G4VEmModel::GetElementSelectors(), and G4VEmModel::SetElementSelectors().

122 {
124 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
void G4LivermoreRayleighModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 239 of file G4LivermoreRayleighModel.cc.

References fStopAndKill, G4cout, G4endl, G4lrint(), G4VEmModel::GetAngularDistribution(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4Element::GetZ(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), G4VEmAngularDistribution::SampleDirection(), G4VEmModel::SelectRandomAtom(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().

244 {
245  if (verboseLevel > 1) {
246  G4cout << "Calling SampleSecondaries() of G4LivermoreRayleighModel"
247  << G4endl;
248  }
249  G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
250 
251  // absorption of low-energy gamma
252  if (photonEnergy0 <= lowEnergyLimit)
253  {
254  fParticleChange->ProposeTrackStatus(fStopAndKill);
255  fParticleChange->SetProposedKineticEnergy(0.);
256  fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
257  return ;
258  }
259 
260  // Select randomly one element in the current material
261  const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
262  const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
263  G4int Z = G4lrint(elm->GetZ());
264 
265  // Sample the angle of the scattered photon
266 
267  G4ThreeVector photonDirection =
268  GetAngularDistribution()->SampleDirection(aDynamicGamma,
269  photonEnergy0,
270  Z, couple->GetMaterial());
271  fParticleChange->ProposeMomentumDirection(photonDirection);
272 }
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:578
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4GLOB_DLL std::ostream G4cout
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, const G4Material *)=0
int G4lrint(double ad)
Definition: templates.hh:163
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510
const G4Material * GetMaterial() const
void G4LivermoreRayleighModel::SetLowEnergyThreshold ( G4double  val)
inline

Definition at line 90 of file G4LivermoreRayleighModel.hh.

91 {
92  lowEnergyLimit = val;
93 }

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