Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LivermorePolarizedRayleighModel Class Reference

#include <G4LivermorePolarizedRayleighModel.hh>

Inheritance diagram for G4LivermorePolarizedRayleighModel:
Collaboration diagram for G4LivermorePolarizedRayleighModel:

Public Member Functions

 G4LivermorePolarizedRayleighModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermorePolarizedRayleigh")
 
virtual ~G4LivermorePolarizedRayleighModel ()
 
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)
 
- 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 GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 SetFluctuationFlag (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)
 

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

Detailed Description

Definition at line 50 of file G4LivermorePolarizedRayleighModel.hh.

Constructor & Destructor Documentation

G4LivermorePolarizedRayleighModel::G4LivermorePolarizedRayleighModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LivermorePolarizedRayleigh" 
)

Definition at line 59 of file G4LivermorePolarizedRayleighModel.cc.

61  :G4VEmModel(nam),fParticleChange(0),isInitialised(false)
62 {
63  fParticleChange =0;
64  lowEnergyLimit = 250 * eV;
65  //SetLowEnergyLimit(lowEnergyLimit);
66  //SetHighEnergyLimit(highEnergyLimit);
67  //
68  verboseLevel= 0;
69  // Verbosity scale:
70  // 0 = nothing
71  // 1 = warning for energy non-conservation
72  // 2 = details of energy budget
73  // 3 = calculation of cross sections, file openings, sampling of atoms
74  // 4 = entering in methods
75 
76  if(verboseLevel > 0) {
77  G4cout << "Livermore Polarized Rayleigh is constructed " << G4endl
78  << "Energy range: "
79  << LowEnergyLimit() / eV << " eV - "
80  << HighEnergyLimit() / GeV << " GeV"
81  << G4endl;
82  }
83 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

G4LivermorePolarizedRayleighModel::~G4LivermorePolarizedRayleighModel ( )
virtual

Definition at line 87 of file G4LivermorePolarizedRayleighModel.cc.

88 {
89  if(IsMaster()) {
90  for(G4int i=0; i<maxZ; ++i) {
91  if(dataCS[i]) {
92  delete dataCS[i];
93  dataCS[i] = 0;
94  }
95  }
96  delete formFactorData;
97  formFactorData = 0;
98 
99  }
100 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
int G4int
Definition: G4Types.hh:78

Here is the call graph for this function:

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 227 of file G4LivermorePolarizedRayleighModel.cc.

232  {
233  if (verboseLevel > 1)
234  {
235  G4cout << "G4LivermoreRayleighModel::ComputeCrossSectionPerAtom()"
236  << G4endl;
237  }
238 
239  if(GammaEnergy < lowEnergyLimit) { return 0.0; }
240 
241  G4double xs = 0.0;
242 
243  G4int intZ = G4lrint(Z);
244 
245  if(intZ < 1 || intZ > maxZ) { return xs; }
246 
247  G4LPhysicsFreeVector* pv = dataCS[intZ];
248 
249  // if element was not initialised
250  // do initialisation safely for MT mode
251  if(!pv) {
252  InitialiseForElement(0, intZ);
253  pv = dataCS[intZ];
254  if(!pv) { return xs; }
255  }
256 
257  G4int n = pv->GetVectorLength() - 1;
258  G4double e = GammaEnergy/MeV;
259  if(e >= pv->Energy(n)) {
260  xs = (*pv)[n]/(e*e);
261  } else if(e >= pv->Energy(0)) {
262  xs = pv->Value(e)/(e*e);
263  }
264 
265  /* if(verboseLevel > 0)
266  {
267  G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
268  << e << G4endl;
269  G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
270  G4cout << " -> first E*E*cs value in CS data file (iu) =" << (*pv)[0]
271  << G4endl;
272  G4cout << " -> last E*E*cs value in CS data file (iu) =" << (*pv)[n]
273  << G4endl;
274  G4cout << "*********************************************************"
275  << G4endl;
276  }
277  */
278 
279  return xs;
280  }
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
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

void G4LivermorePolarizedRayleighModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 104 of file G4LivermorePolarizedRayleighModel.cc.

106 {
107 // Rayleigh process: The Quantum Theory of Radiation
108 // W. Heitler, Oxford at the Clarendon Press, Oxford (1954)
109 // Scattering function: A simple model of photon transport
110 // D.E. Cullen, Nucl. Instr. Meth. in Phys. Res. B 101 (1995) 499-510
111 // Polarization of the outcoming photon: Beam test of a prototype detector array for the PoGO astronomical hard X-ray/soft gamma-ray polarimeter
112 // T. Mizuno et al., Nucl. Instr. Meth. in Phys. Res. A 540 (2005) 158-168
113 
114  if (verboseLevel > 3)
115  G4cout << "Calling G4LivermorePolarizedRayleighModel::Initialise()" << G4endl;
116 
117 
118  if(IsMaster()) {
119 
120  // Form Factor
121 
122  G4VDataSetAlgorithm* ffInterpolation = new G4LogLogInterpolation;
123  G4String formFactorFile = "rayl/re-ff-";
124  formFactorData = new G4CompositeEMDataSet(ffInterpolation,1.,1.);
125  formFactorData->LoadData(formFactorFile);
126 
127  // Initialise element selector
128  InitialiseElementSelectors(particle, cuts);
129 
130  // Access to elements
131  char* path = getenv("G4LEDATA");
132  G4ProductionCutsTable* theCoupleTable =
134  G4int numOfCouples = theCoupleTable->GetTableSize();
135 
136  for(G4int i=0; i<numOfCouples; ++i)
137  {
138  const G4MaterialCutsCouple* couple =
139  theCoupleTable->GetMaterialCutsCouple(i);
140  const G4Material* material = couple->GetMaterial();
141  const G4ElementVector* theElementVector = material->GetElementVector();
142  G4int nelm = material->GetNumberOfElements();
143 
144  for (G4int j=0; j<nelm; ++j)
145  {
146  G4int Z = G4lrint((*theElementVector)[j]->GetZ());
147  if(Z < 1) { Z = 1; }
148  else if(Z > maxZ) { Z = maxZ; }
149  if( (!dataCS[Z]) ) { ReadData(Z, path); }
150  }
151  }
152  }
153 
154  if(isInitialised) { return; }
155  fParticleChange = GetParticleChangeForGamma();
156  isInitialised = true;
157 
158 }
std::vector< G4Element * > G4ElementVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
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
virtual G4bool LoadData(const G4String &fileName)=0
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
const G4Material * GetMaterial() const

Here is the call graph for this function:

void G4LivermorePolarizedRayleighModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 494 of file G4LivermorePolarizedRayleighModel.cc.

496  {
497  G4AutoLock l(&LivermorePolarizedRayleighModelMutex);
498  // G4cout << "G4LivermoreRayleighModel::InitialiseForElement Z= "
499  // << Z << G4endl;
500  if(!dataCS[Z]) { ReadData(Z); }
501  l.unlock();
502  }

Here is the call graph for this function:

Here is the caller graph for this function:

void G4LivermorePolarizedRayleighModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 163 of file G4LivermorePolarizedRayleighModel.cc.

165 {
167 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:809
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:817

Here is the call graph for this function:

void G4LivermorePolarizedRayleighModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 284 of file G4LivermorePolarizedRayleighModel.cc.

289 {
290  if (verboseLevel > 3)
291  G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedRayleighModel" << G4endl;
292 
293  G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
294 
295  if (photonEnergy0 <= lowEnergyLimit)
296  {
297  fParticleChange->ProposeTrackStatus(fStopAndKill);
298  fParticleChange->SetProposedKineticEnergy(0.);
299  fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
300  return ;
301  }
302 
303  G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
304 
305  // Select randomly one element in the current material
306  // G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy0);
307  const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
308  const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
309  G4int Z = (G4int)elm->GetZ();
310 
311  G4double outcomingPhotonCosTheta = GenerateCosTheta(photonEnergy0, Z);
312  G4double outcomingPhotonPhi = GeneratePhi(outcomingPhotonCosTheta);
313  G4double beta=GeneratePolarizationAngle();
314 
315  // incomingPhoton reference frame:
316  // z = versor parallel to the incomingPhotonDirection
317  // x = versor parallel to the incomingPhotonPolarization
318  // y = defined as z^x
319 
320  // outgoingPhoton reference frame:
321  // z' = versor parallel to the outgoingPhotonDirection
322  // x' = defined as x-x*z'z' normalized
323  // y' = defined as z'^x'
324 
325  G4ThreeVector z(aDynamicGamma->GetMomentumDirection().unit());
326  G4ThreeVector x(GetPhotonPolarization(*aDynamicGamma));
327  G4ThreeVector y(z.cross(x));
328 
329  // z' = std::cos(phi)*std::sin(theta) x + std::sin(phi)*std::sin(theta) y + std::cos(theta) z
330  G4double xDir;
331  G4double yDir;
332  G4double zDir;
333  zDir=outcomingPhotonCosTheta;
334  xDir=std::sqrt(1-outcomingPhotonCosTheta*outcomingPhotonCosTheta);
335  yDir=xDir;
336  xDir*=std::cos(outcomingPhotonPhi);
337  yDir*=std::sin(outcomingPhotonPhi);
338 
339  G4ThreeVector zPrime((xDir*x + yDir*y + zDir*z).unit());
340  G4ThreeVector xPrime(x.perpPart(zPrime).unit());
341  G4ThreeVector yPrime(zPrime.cross(xPrime));
342 
343  // outgoingPhotonPolarization is directed as x' std::cos(beta) + y' std::sin(beta)
344  G4ThreeVector outcomingPhotonPolarization(xPrime*std::cos(beta) + yPrime*std::sin(beta));
345 
346  fParticleChange->ProposeMomentumDirection(zPrime);
347  fParticleChange->ProposePolarization(outcomingPhotonPolarization);
348  fParticleChange->SetProposedKineticEnergy(photonEnergy0);
349 
350 }
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
G4ParticleDefinition * GetDefinition() const
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
void ProposePolarization(const G4ThreeVector &dir)
Hep3Vector unit() const
tuple z
Definition: test.py:28
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:544

Here is the call graph for this function:


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