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

#include <G4ForwardXrayTR.hh>

Inheritance diagram for G4ForwardXrayTR:
G4TransitionRadiation G4VDiscreteProcess G4VProcess

Public Member Functions

 G4ForwardXrayTR (const G4String &matName1, const G4String &matName2, const G4String &processName="XrayTR")
 
 G4ForwardXrayTR (const G4String &processName="XrayTR")
 
virtual ~G4ForwardXrayTR ()
 
void BuildXrayTRtables ()
 
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
G4double GetEnergyTR (G4int iMat, G4int jMat, G4int iTkin) const
 
G4double GetThetaTR (G4int iMat, G4int jMat, G4int iTkin) const
 
G4double SpectralAngleTRdensity (G4double energy, G4double varAngle) const
 
G4double AngleDensity (G4double energy, G4double varAngle) const
 
G4double EnergyInterval (G4double energy1, G4double energy2, G4double varAngle) const
 
G4double AngleSum (G4double varAngle1, G4double varAngle2) const
 
G4double SpectralDensity (G4double energy, G4double x) const
 
G4double AngleInterval (G4double energy, G4double varAngle1, G4double varAngle2) const
 
G4double EnergySum (G4double energy1, G4double energy2) const
 
G4PhysicsTableGetAngleDistrTable ()
 
G4PhysicsTableGetEnergyDistrTable ()
 
- Public Member Functions inherited from G4TransitionRadiation
 G4TransitionRadiation (const G4String &processName="TR", G4ProcessType type=fElectromagnetic)
 
virtual ~G4TransitionRadiation ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
G4double IntegralOverEnergy (G4double energy1, G4double energy2, G4double varAngle) const
 
G4double IntegralOverAngle (G4double energy, G4double varAngle1, G4double varAngle2) const
 
G4double AngleIntegralDistribution (G4double varAngle1, G4double varAngle2) const
 
G4double EnergyIntegralDistribution (G4double energy1, G4double energy2) const
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Static Public Member Functions

static G4int GetSympsonNumber ()
 
static G4int GetBinTR ()
 
static G4double GetMinProtonTkin ()
 
static G4double GetMaxProtonTkin ()
 
static G4int GetTotBin ()
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Protected Attributes

G4ParticleDefinitionfPtrGamma
 
const std::vector< G4double > * fGammaCutInKineticEnergy
 
G4double fGammaTkinCut
 
G4PhysicsTablefAngleDistrTable
 
G4PhysicsTablefEnergyDistrTable
 
G4PhysicsLogVectorfProtonEnergyVector
 
G4double fMinEnergyTR
 
G4double fMaxEnergyTR
 
G4double fMaxThetaTR
 
G4double fGamma
 
G4double fSigma1
 
G4double fSigma2
 
- Protected Attributes inherited from G4TransitionRadiation
G4int fMatIndex1
 
G4int fMatIndex2
 
G4double fGamma
 
G4double fEnergy
 
G4double fVarAngle
 
G4double fMinEnergy
 
G4double fMaxEnergy
 
G4double fMaxTheta
 
G4double fSigma1
 
G4double fSigma2
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Static Protected Attributes

static G4int fSympsonNumber = 100
 
static G4double fTheMinEnergyTR = 1.0*keV
 
static G4double fTheMaxEnergyTR = 100.0*keV
 
static G4double fTheMaxAngle = 1.0e-3
 
static G4double fTheMinAngle = 5.0e-6
 
static G4int fBinTR = 50
 
static G4double fMinProtonTkin = 100.0*GeV
 
static G4double fMaxProtonTkin = 100.0*TeV
 
static G4int fTotBin = 50
 
static G4double fPlasmaCof
 
static G4double fCofTR = fine_structure_const/pi
 
- Static Protected Attributes inherited from G4TransitionRadiation
static const G4int fSympsonNumber = 100
 
static const G4int fGammaNumber = 15
 
static const G4int fPointNumber = 100
 

Additional Inherited Members

- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 66 of file G4ForwardXrayTR.hh.

Constructor & Destructor Documentation

G4ForwardXrayTR::G4ForwardXrayTR ( const G4String matName1,
const G4String matName2,
const G4String processName = "XrayTR" 
)

Definition at line 83 of file G4ForwardXrayTR.cc.

References BuildXrayTRtables(), fAngleDistrTable, fEnergyDistrTable, fGamma, fGammaCutInKineticEnergy, fGammaTkinCut, G4TransitionRadiation::fMatIndex1, G4TransitionRadiation::fMatIndex2, fMaxEnergyTR, fMaxProtonTkin, fMaxThetaTR, fMinEnergyTR, fMinProtonTkin, fProtonEnergyVector, fPtrGamma, fSigma1, fSigma2, fTotBin, G4Exception(), G4MaterialCutsCouple::GetIndex(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4Material::GetName(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), and JustWarning.

86  : G4TransitionRadiation(processName)
87 {
88  fPtrGamma = 0;
91  fSigma1 = fSigma2 = 0.0;
92  fAngleDistrTable = 0;
94  fMatIndex1 = fMatIndex2 = 0;
95 
96  // Proton energy vector initialization
97  //
100  G4int iMat;
101  const G4ProductionCutsTable* theCoupleTable=
103  G4int numOfCouples = theCoupleTable->GetTableSize();
104 
105  G4bool build = true;
106 
107  for(iMat=0;iMat<numOfCouples;iMat++) // check first material name
108  {
109  const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(iMat);
110  if( matName1 == couple->GetMaterial()->GetName() )
111  {
112  fMatIndex1 = couple->GetIndex();
113  break;
114  }
115  }
116  if(iMat == numOfCouples)
117  {
118  G4Exception("G4ForwardXrayTR::G4ForwardXrayTR", "ForwardXrayTR01", JustWarning, "Invalid first material name in G4ForwardXrayTR constructor!");
119  build = false;
120  }
121 
122  if(build) {
123  for(iMat=0;iMat<numOfCouples;iMat++) // check second material name
124  {
125  const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(iMat);
126  if( matName2 == couple->GetMaterial()->GetName() )
127  {
128  fMatIndex2 = couple->GetIndex();
129  break;
130  }
131  }
132  if(iMat == numOfCouples)
133  {
134  G4Exception("G4ForwardXrayTR::G4ForwardXrayTR", "ForwardXrayTR02", JustWarning, "Invalid second material name in G4ForwardXrayTR constructor!");
135  build = false;
136  }
137  }
138  // G4cout<<"G4ForwardXray constructor is called"<<G4endl;
139  if(build) { BuildXrayTRtables(); }
140 }
static G4double fMaxProtonTkin
G4TransitionRadiation(const G4String &processName="TR", G4ProcessType type=fElectromagnetic)
static G4double fMinProtonTkin
const G4String & GetName() const
Definition: G4Material.hh:176
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * fAngleDistrTable
bool G4bool
Definition: G4Types.hh:79
static G4int fTotBin
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4PhysicsTable * fEnergyDistrTable
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
G4ParticleDefinition * fPtrGamma
G4PhysicsLogVector * fProtonEnergyVector
const std::vector< G4double > * fGammaCutInKineticEnergy
const G4Material * GetMaterial() const
G4ForwardXrayTR::G4ForwardXrayTR ( const G4String processName = "XrayTR")

Definition at line 147 of file G4ForwardXrayTR.cc.

References fAngleDistrTable, fEnergyDistrTable, fGamma, fGammaCutInKineticEnergy, fGammaTkinCut, G4TransitionRadiation::fMatIndex1, G4TransitionRadiation::fMatIndex2, fMaxEnergyTR, fMaxProtonTkin, fMaxThetaTR, fMinEnergyTR, fMinProtonTkin, fProtonEnergyVector, fPtrGamma, fSigma1, fSigma2, and fTotBin.

148  : G4TransitionRadiation(processName)
149 {
150  fPtrGamma = 0;
153  fSigma1 = fSigma2 = 0.0;
154  fAngleDistrTable = 0;
155  fEnergyDistrTable = 0;
156  fMatIndex1 = fMatIndex2 = 0;
157 
158  // Proton energy vector initialization
159  //
162 }
static G4double fMaxProtonTkin
G4TransitionRadiation(const G4String &processName="TR", G4ProcessType type=fElectromagnetic)
static G4double fMinProtonTkin
G4PhysicsTable * fAngleDistrTable
static G4int fTotBin
G4PhysicsTable * fEnergyDistrTable
G4ParticleDefinition * fPtrGamma
G4PhysicsLogVector * fProtonEnergyVector
const std::vector< G4double > * fGammaCutInKineticEnergy
G4ForwardXrayTR::~G4ForwardXrayTR ( )
virtual

Definition at line 170 of file G4ForwardXrayTR.cc.

References fAngleDistrTable, fEnergyDistrTable, and fProtonEnergyVector.

171 {
172  delete fAngleDistrTable;
173  delete fEnergyDistrTable;
174  delete fProtonEnergyVector;
175 }
G4PhysicsTable * fAngleDistrTable
G4PhysicsTable * fEnergyDistrTable
G4PhysicsLogVector * fProtonEnergyVector

Member Function Documentation

G4double G4ForwardXrayTR::AngleDensity ( G4double  energy,
G4double  varAngle 
) const

Definition at line 340 of file G4ForwardXrayTR.cc.

References test::c, energy(), fGamma, fSigma1, fSigma2, and test::x.

Referenced by EnergyInterval().

342 {
343  G4double x, x2, /*a, b,*/ c, d, f, a2, b2, a4, b4;
344  G4double cof1, cof2, cof3;
345  x = 1.0/energy;
346  x2 = x*x;
347  c = 1.0/fSigma1;
348  d = 1.0/fSigma2;
349  f = (varAngle + 1.0/(fGamma*fGamma));
350  a2 = c*f;
351  b2 = d*f;
352  a4 = a2*a2;
353  b4 = b2*b2;
354  //a = std::sqrt(a2);
355  // b = std::sqrt(b2);
356  cof1 = c*c*(0.5/(a2*(x2 +a2)) +0.5*std::log(x2/(x2 +a2))/a4);
357  cof3 = d*d*(0.5/(b2*(x2 +b2)) +0.5*std::log(x2/(x2 +b2))/b4);
358  cof2 = -c*d*(std::log(x2/(x2 +b2))/b2 - std::log(x2/(x2 +a2))/a2)/(a2 - b2) ;
359  return -varAngle*(cof1 + cof2 + cof3);
360 }
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
double G4double
Definition: G4Types.hh:76
G4double G4ForwardXrayTR::AngleInterval ( G4double  energy,
G4double  varAngle1,
G4double  varAngle2 
) const

Definition at line 427 of file G4ForwardXrayTR.cc.

References SpectralDensity().

Referenced by EnergySum().

430 {
431  return SpectralDensity(energy,varAngle2)
432  - SpectralDensity(energy,varAngle1);
433 }
G4double SpectralDensity(G4double energy, G4double x) const
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4double G4ForwardXrayTR::AngleSum ( G4double  varAngle1,
G4double  varAngle2 
) const

Definition at line 382 of file G4ForwardXrayTR.cc.

References EnergyInterval(), fMaxEnergyTR, fMinEnergyTR, and fSympsonNumber.

Referenced by BuildXrayTRtables().

384 {
385  G4int i;
386  G4double h , sumEven = 0.0 , sumOdd = 0.0;
387  h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber;
388  for(i=1;i<fSympsonNumber;i++)
389  {
390  sumEven += EnergyInterval(fMinEnergyTR,fMaxEnergyTR,varAngle1 + 2*i*h );
392  varAngle1 + (2*i - 1)*h );
393  }
395  varAngle1 + (2*fSympsonNumber - 1)*h );
396 
397  return h*(EnergyInterval(fMinEnergyTR,fMaxEnergyTR,varAngle1)
399  + 4.0*sumOdd + 2.0*sumEven )/3.0;
400 }
int G4int
Definition: G4Types.hh:78
static G4int fSympsonNumber
double G4double
Definition: G4Types.hh:76
G4double EnergyInterval(G4double energy1, G4double energy2, G4double varAngle) const
void G4ForwardXrayTR::BuildXrayTRtables ( )

Definition at line 189 of file G4ForwardXrayTR.cc.

References AngleSum(), EnergySum(), fAngleDistrTable, fBinTR, fCofTR, fEnergyDistrTable, fGamma, fGammaCutInKineticEnergy, fGammaTkinCut, G4TransitionRadiation::fMatIndex1, G4TransitionRadiation::fMatIndex2, fMaxEnergyTR, fMaxThetaTR, fMinEnergyTR, fPlasmaCof, fProtonEnergyVector, fSigma1, fSigma2, fTheMaxAngle, fTheMaxEnergyTR, fTheMinAngle, fTheMinEnergyTR, fTotBin, G4Material::GetElectronDensity(), G4ProductionCutsTable::GetEnergyCutsVector(), G4PhysicsVector::GetLowEdgeEnergy(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), idxG4GammaCut, G4PhysicsTable::insertAt(), python.hepunit::proton_mass_c2, and G4PhysicsVector::PutValue().

Referenced by G4ForwardXrayTR().

190 {
191  G4int iMat, jMat, iTkin, iTR, iPlace;
192  const G4ProductionCutsTable* theCoupleTable=
194  G4int numOfCouples = theCoupleTable->GetTableSize();
195 
197 
200 
201 
202  for(iMat=0;iMat<numOfCouples;iMat++) // loop over pairs of different materials
203  {
204  if( iMat != fMatIndex1 && iMat != fMatIndex2 ) continue;
205 
206  for(jMat=0;jMat<numOfCouples;jMat++) // transition iMat -> jMat !!!
207  {
208  if( iMat == jMat || ( jMat != fMatIndex1 && jMat != fMatIndex2 ) )
209  {
210  continue;
211  }
212  else
213  {
214  const G4MaterialCutsCouple* iCouple = theCoupleTable->GetMaterialCutsCouple(iMat);
215  const G4MaterialCutsCouple* jCouple = theCoupleTable->GetMaterialCutsCouple(jMat);
216  const G4Material* mat1 = iCouple->GetMaterial();
217  const G4Material* mat2 = jCouple->GetMaterial();
218 
221 
222  // fGammaTkinCut = fGammaCutInKineticEnergy[jMat]; // TR photon in jMat !
223 
224  fGammaTkinCut = 0.0;
225 
226  if(fGammaTkinCut > fTheMinEnergyTR) // setting of min/max TR energies
227  {
229  }
230  else
231  {
233  }
234  if(fGammaTkinCut > fTheMaxEnergyTR)
235  {
236  fMaxEnergyTR = 2.0*fGammaTkinCut; // usually very low TR rate
237  }
238  else
239  {
241  }
242  for(iTkin=0;iTkin<fTotBin;iTkin++) // Lorentz factor loop
243  {
245  energyVector = new G4PhysicsLogVector( fMinEnergyTR,
246  fMaxEnergyTR,
247  fBinTR );
248 
249  fGamma = 1.0 + (fProtonEnergyVector->
250  GetLowEdgeEnergy(iTkin)/proton_mass_c2);
251 
252  fMaxThetaTR = 10000.0/(fGamma*fGamma);
253 
255  {
257  }
258  else
259  {
261  {
263  }
264  }
265  // G4cout<<G4endl<<"fGamma = "<<fGamma<<" fMaxThetaTR = "<<fMaxThetaTR<<G4endl;
267  angleVector = new G4PhysicsLinearVector( 0.0,
268  fMaxThetaTR,
269  fBinTR );
270  G4double energySum = 0.0;
271  G4double angleSum = 0.0;
272 
273  energyVector->PutValue(fBinTR-1,energySum);
274  angleVector->PutValue(fBinTR-1,angleSum);
275 
276  for(iTR=fBinTR-2;iTR>=0;iTR--)
277  {
278  energySum += fCofTR*EnergySum(energyVector->GetLowEdgeEnergy(iTR),
279  energyVector->GetLowEdgeEnergy(iTR+1));
280 
281  angleSum += fCofTR*AngleSum(angleVector->GetLowEdgeEnergy(iTR),
282  angleVector->GetLowEdgeEnergy(iTR+1));
283 
284  energyVector->PutValue(iTR,energySum);
285  angleVector ->PutValue(iTR,angleSum);
286  }
287  // G4cout<<"sumE = "<<energySum<<"; sumA = "<<angleSum<<G4endl;
288 
289  if(jMat < iMat)
290  {
291  iPlace = fTotBin+iTkin; // (iMat*(numOfMat-1)+jMat)*
292  }
293  else // jMat > iMat right part of matrices (jMat-1) !
294  {
295  iPlace = iTkin; // (iMat*(numOfMat-1)+jMat-1)*fTotBin+
296  }
297  fEnergyDistrTable->insertAt(iPlace,energyVector);
298  fAngleDistrTable->insertAt(iPlace,angleVector);
299  } // iTkin
300  } // jMat != iMat
301  } // jMat
302  } // iMat
303  // G4cout<<"G4ForwardXrayTR::BuildXrayTRtables have been called"<<G4endl;
304 }
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4int fBinTR
static G4double fTheMaxEnergyTR
G4double AngleSum(G4double varAngle1, G4double varAngle2) const
G4double GetLowEdgeEnergy(size_t binNumber) const
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * fAngleDistrTable
static G4double fPlasmaCof
static G4double fTheMinAngle
G4double GetElectronDensity() const
Definition: G4Material.hh:215
void PutValue(size_t index, G4double theValue)
float proton_mass_c2
Definition: hepunit.py:275
static G4int fTotBin
G4PhysicsTable * fEnergyDistrTable
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
G4double EnergySum(G4double energy1, G4double energy2) const
void insertAt(size_t, G4PhysicsVector *)
static G4double fCofTR
static G4double fTheMinEnergyTR
static G4double fTheMaxAngle
double G4double
Definition: G4Types.hh:76
G4PhysicsLogVector * fProtonEnergyVector
const std::vector< G4double > * fGammaCutInKineticEnergy
const G4Material * GetMaterial() const
G4double G4ForwardXrayTR::EnergyInterval ( G4double  energy1,
G4double  energy2,
G4double  varAngle 
) const

Definition at line 368 of file G4ForwardXrayTR.cc.

References AngleDensity().

Referenced by AngleSum().

371 {
372  return AngleDensity(energy2,varAngle)
373  - AngleDensity(energy1,varAngle);
374 }
G4double AngleDensity(G4double energy, G4double varAngle) const
G4double G4ForwardXrayTR::EnergySum ( G4double  energy1,
G4double  energy2 
) const

Definition at line 441 of file G4ForwardXrayTR.cc.

References AngleInterval(), fMaxThetaTR, and fSympsonNumber.

Referenced by BuildXrayTRtables().

443 {
444  G4int i;
445  G4double h , sumEven = 0.0 , sumOdd = 0.0;
446  h = 0.5*(energy2 - energy1)/fSympsonNumber;
447  for(i=1;i<fSympsonNumber;i++)
448  {
449  sumEven += AngleInterval(energy1 + 2*i*h,0.0,fMaxThetaTR);
450  sumOdd += AngleInterval(energy1 + (2*i - 1)*h,0.0,fMaxThetaTR);
451  }
452  sumOdd += AngleInterval(energy1 + (2*fSympsonNumber - 1)*h,
453  0.0,fMaxThetaTR);
454 
455  return h*( AngleInterval(energy1,0.0,fMaxThetaTR)
456  + AngleInterval(energy2,0.0,fMaxThetaTR)
457  + 4.0*sumOdd + 2.0*sumEven )/3.0;
458 }
G4double AngleInterval(G4double energy, G4double varAngle1, G4double varAngle2) const
int G4int
Definition: G4Types.hh:78
static G4int fSympsonNumber
double G4double
Definition: G4Types.hh:76
G4PhysicsTable* G4ForwardXrayTR::GetAngleDistrTable ( )
inline

Definition at line 128 of file G4ForwardXrayTR.hh.

References fAngleDistrTable.

128 { return fAngleDistrTable; };
G4PhysicsTable * fAngleDistrTable
static G4int G4ForwardXrayTR::GetBinTR ( )
inlinestatic

Definition at line 132 of file G4ForwardXrayTR.hh.

References fBinTR.

132 { return fBinTR; };
static G4int fBinTR
G4PhysicsTable* G4ForwardXrayTR::GetEnergyDistrTable ( )
inline

Definition at line 129 of file G4ForwardXrayTR.hh.

References fEnergyDistrTable.

129 { return fEnergyDistrTable; };
G4PhysicsTable * fEnergyDistrTable
G4double G4ForwardXrayTR::GetEnergyTR ( G4int  iMat,
G4int  jMat,
G4int  iTkin 
) const

Definition at line 696 of file G4ForwardXrayTR.cc.

References fBinTR, fTotBin, G4cout, G4endl, G4Poisson(), G4UniformRand, G4PhysicsVector::GetLowEdgeEnergy(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4ProductionCutsTable::GetProductionCutsTable(), G4Material::GetState(), G4ProductionCutsTable::GetTableSize(), kStateLiquid, and kStateSolid.

697 {
698  G4int iPlace, numOfTR, iTR, iTransfer;
699  G4double energyTR = 0.0; // return this value for no TR photons
700  G4double energyPos ;
701  G4double W1, W2;
702 
703  const G4ProductionCutsTable* theCoupleTable=
705  G4int numOfCouples = theCoupleTable->GetTableSize();
706 
707  // The case of equal or approximate (in terms of plasma energy) materials
708  // No TR photons ?!
709 
710  const G4MaterialCutsCouple* iCouple = theCoupleTable->GetMaterialCutsCouple(iMat);
711  const G4MaterialCutsCouple* jCouple = theCoupleTable->GetMaterialCutsCouple(jMat);
712  const G4Material* iMaterial = iCouple->GetMaterial();
713  const G4Material* jMaterial = jCouple->GetMaterial();
714 
715  if ( iMat == jMat
716 
717  || iMaterial->GetState() == jMaterial->GetState()
718 
719  ||(iMaterial->GetState() == kStateSolid && jMaterial->GetState() == kStateLiquid )
720 
721  ||(iMaterial->GetState() == kStateLiquid && jMaterial->GetState() == kStateSolid ) )
722 
723  {
724  return energyTR;
725  }
726 
727  if(jMat < iMat)
728  {
729  iPlace = (iMat*(numOfCouples - 1) + jMat)*fTotBin + iTkin - 1;
730  }
731  else
732  {
733  iPlace = (iMat*(numOfCouples - 1) + jMat - 1)*fTotBin + iTkin - 1;
734  }
735  G4PhysicsVector* energyVector1 = (*fEnergyDistrTable)(iPlace) ;
736  G4PhysicsVector* energyVector2 = (*fEnergyDistrTable)(iPlace + 1);
737 
738  if(iTkin == fTotBin) // TR plato, try from left
739  {
740  numOfTR = G4Poisson( (*energyVector1)(0) );
741  if(numOfTR == 0)
742  {
743  return energyTR;
744  }
745  else
746  {
747  for(iTR=0;iTR<numOfTR;iTR++)
748  {
749  energyPos = (*energyVector1)(0)*G4UniformRand();
750  for(iTransfer=0;iTransfer<fBinTR-1;iTransfer++)
751  {
752  if(energyPos >= (*energyVector1)(iTransfer)) break;
753  }
754  energyTR += energyVector1->GetLowEdgeEnergy(iTransfer);
755  }
756  }
757  }
758  else
759  {
760  if(iTkin == 0) // Tkin is too small, neglect of TR photon generation
761  {
762  return energyTR;
763  }
764  else // general case: Tkin between two vectors of the material
765  { // use trivial mean half/half
766  W1 = 0.5;
767  W2 = 0.5;
768  numOfTR = G4Poisson( (*energyVector1)(0)*W1 +
769  (*energyVector2)(0)*W2 );
770  if(numOfTR == 0)
771  {
772  return energyTR;
773  }
774  else
775  {
776  G4cout<<"It is still OK in GetEnergyTR(int,int,int)"<<G4endl;
777  for(iTR=0;iTR<numOfTR;iTR++)
778  {
779  energyPos = ((*energyVector1)(0)*W1+
780  (*energyVector2)(0)*W2)*G4UniformRand();
781  for(iTransfer=0;iTransfer<fBinTR-1;iTransfer++)
782  {
783  if(energyPos >= ((*energyVector1)(iTransfer)*W1+
784  (*energyVector2)(iTransfer)*W2)) break;
785  }
786  energyTR += (energyVector1->GetLowEdgeEnergy(iTransfer))*W1+
787  (energyVector2->GetLowEdgeEnergy(iTransfer))*W2;
788 
789  }
790  }
791  }
792  }
793 
794  return energyTR;
795 }
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
static G4int fBinTR
G4double GetLowEdgeEnergy(size_t binNumber) const
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
static G4int fTotBin
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
#define G4endl
Definition: G4ios.hh:61
G4State GetState() const
Definition: G4Material.hh:179
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const
static G4double G4ForwardXrayTR::GetMaxProtonTkin ( )
inlinestatic

Definition at line 135 of file G4ForwardXrayTR.hh.

References fMaxProtonTkin.

135 { return fMaxProtonTkin; };
static G4double fMaxProtonTkin
G4double G4ForwardXrayTR::GetMeanFreePath ( const G4Track ,
G4double  ,
G4ForceCondition condition 
)
virtual

Implements G4VDiscreteProcess.

Definition at line 177 of file G4ForwardXrayTR.cc.

References DBL_MAX, and Forced.

180 {
181  *condition = Forced;
182  return DBL_MAX; // so TR doesn't limit mean free path
183 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
static G4double G4ForwardXrayTR::GetMinProtonTkin ( )
inlinestatic

Definition at line 134 of file G4ForwardXrayTR.hh.

References fMinProtonTkin.

134 { return fMinProtonTkin; };
static G4double fMinProtonTkin
static G4int G4ForwardXrayTR::GetSympsonNumber ( )
inlinestatic

Definition at line 131 of file G4ForwardXrayTR.hh.

References fSympsonNumber.

131 { return fSympsonNumber; };
static G4int fSympsonNumber
G4double G4ForwardXrayTR::GetThetaTR ( G4int  iMat,
G4int  jMat,
G4int  iTkin 
) const

Definition at line 805 of file G4ForwardXrayTR.cc.

806 {
807  G4double theta = 0.0;
808 
809  return theta;
810 }
double G4double
Definition: G4Types.hh:76
static G4int G4ForwardXrayTR::GetTotBin ( )
inlinestatic

Definition at line 136 of file G4ForwardXrayTR.hh.

References fTotBin.

136 { return fTotBin; };
static G4int fTotBin
G4VParticleChange * G4ForwardXrayTR::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 466 of file G4ForwardXrayTR.cc.

References G4ParticleChange::AddSecondary(), G4VProcess::aParticleChange, fAngleDistrTable, fBinTR, fEnergyDistrTable, fGeomBoundary, G4TransitionRadiation::fMatIndex1, G4TransitionRadiation::fMatIndex2, fProtonEnergyVector, fTotBin, G4Poisson(), G4UniformRand, G4Gamma::Gamma(), G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4MaterialCutsCouple::GetIndex(), G4GeometryTolerance::GetInstance(), G4DynamicParticle::GetKineticEnergy(), G4PhysicsVector::GetLowEdgeEnergy(), G4MaterialCutsCouple::GetMaterial(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4StepPoint::GetPhysicalVolume(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4Material::GetState(), G4Track::GetStepLength(), G4StepPoint::GetStepStatus(), G4GeometryTolerance::GetSurfaceTolerance(), G4ParticleChange::Initialize(), kStateLiquid, kStateSolid, G4VDiscreteProcess::PostStepDoIt(), G4ParticleChange::ProposeEnergy(), python.hepunit::proton_mass_c2, CLHEP::Hep3Vector::rotateUz(), G4VParticleChange::SetNumberOfSecondaries(), and python.hepunit::twopi.

468 {
469  aParticleChange.Initialize(aTrack);
470  // G4cout<<"call G4ForwardXrayTR::PostStepDoIt"<<G4endl;
471  G4int iMat, jMat, iTkin, iPlace, numOfTR, iTR, iTransfer;
472 
473  G4double energyPos, anglePos, energyTR, theta, phi, dirX, dirY, dirZ;
474  G4double W, W1, W2, E1, E2;
475 
476  G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint();
477  G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
479 
480  if (pPostStepPoint->GetStepStatus() != fGeomBoundary)
481  {
482  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
483  }
484  if (aTrack.GetStepLength() <= tol)
485  {
486  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
487  }
488  // Come on boundary, so begin to try TR
489 
490  const G4MaterialCutsCouple* iCouple = pPreStepPoint ->GetPhysicalVolume()->
491  GetLogicalVolume()->GetMaterialCutsCouple();
492  const G4MaterialCutsCouple* jCouple = pPostStepPoint ->GetPhysicalVolume()->
493  GetLogicalVolume()->GetMaterialCutsCouple();
494  const G4Material* iMaterial = iCouple->GetMaterial();
495  const G4Material* jMaterial = jCouple->GetMaterial();
496  iMat = iCouple->GetIndex();
497  jMat = jCouple->GetIndex();
498 
499  // The case of equal or approximate (in terms of plasma energy) materials
500  // No TR photons ?!
501 
502  if ( iMat == jMat
503  || ( (fMatIndex1 >= 0 && fMatIndex1 >= 0)
504  && ( iMat != fMatIndex1 && iMat != fMatIndex2 )
505  && ( jMat != fMatIndex1 && jMat != fMatIndex2 ) )
506 
507  || iMaterial->GetState() == jMaterial->GetState()
508 
509  ||(iMaterial->GetState() == kStateSolid && jMaterial->GetState() == kStateLiquid )
510 
511  ||(iMaterial->GetState() == kStateLiquid && jMaterial->GetState() == kStateSolid ) )
512  {
513  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
514  }
515 
516  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
517  G4double charge = aParticle->GetDefinition()->GetPDGCharge();
518 
519  if(charge == 0.0) // Uncharged particle doesn't Generate TR photons
520  {
521  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
522  }
523  // Now we are ready to Generate TR photons
524 
525  G4double chargeSq = charge*charge;
526  G4double kinEnergy = aParticle->GetKineticEnergy();
527  G4double massRatio = proton_mass_c2/aParticle->GetDefinition()->GetPDGMass();
528  G4double TkinScaled = kinEnergy*massRatio;
529  for(iTkin=0;iTkin<fTotBin;iTkin++)
530  {
531  if(TkinScaled < fProtonEnergyVector->GetLowEdgeEnergy(iTkin)) // <= ?
532  {
533  break;
534  }
535  }
536  if(jMat < iMat)
537  {
538  iPlace = fTotBin + iTkin - 1; // (iMat*(numOfMat - 1) + jMat)*
539  }
540  else
541  {
542  iPlace = iTkin - 1; // (iMat*(numOfMat - 1) + jMat - 1)*fTotBin +
543  }
544  // G4PhysicsVector* energyVector1 = (*fEnergyDistrTable)(iPlace) ;
545  // G4PhysicsVector* energyVector2 = (*fEnergyDistrTable)(iPlace + 1);
546 
547  // G4PhysicsVector* angleVector1 = (*fAngleDistrTable)(iPlace) ;
548  // G4PhysicsVector* angleVector2 = (*fAngleDistrTable)(iPlace + 1) ;
549 
550  G4ParticleMomentum particleDir = aParticle->GetMomentumDirection();
551 
552  if(iTkin == fTotBin) // TR plato, try from left
553  {
554  // G4cout<<iTkin<<" mean TR number = "<<( (*(*fEnergyDistrTable)(iPlace))(0) +
555  // (*(*fAngleDistrTable)(iPlace))(0) )
556  // *chargeSq*0.5<<G4endl;
557 
558  numOfTR = G4Poisson( ( (*(*fEnergyDistrTable)(iPlace))(0) +
559  (*(*fAngleDistrTable)(iPlace))(0) )
560  *chargeSq*0.5 );
561  if(numOfTR == 0)
562  {
563  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
564  }
565  else
566  {
567  // G4cout<<"Number of X-ray TR photons = "<<numOfTR<<G4endl;
568 
570 
571  for(iTR=0;iTR<numOfTR;iTR++)
572  {
573  energyPos = (*(*fEnergyDistrTable)(iPlace))(0)*G4UniformRand();
574  for(iTransfer=0;iTransfer<fBinTR-1;iTransfer++)
575  {
576  if(energyPos >= (*(*fEnergyDistrTable)(iPlace))(iTransfer)) break;
577  }
578  energyTR = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
579 
580  // G4cout<<"energyTR = "<<energyTR/keV<<"keV"<<G4endl;
581 
582  kinEnergy -= energyTR;
583  aParticleChange.ProposeEnergy(kinEnergy);
584 
585  anglePos = (*(*fAngleDistrTable)(iPlace))(0)*G4UniformRand();
586  for(iTransfer=0;iTransfer<fBinTR-1;iTransfer++)
587  {
588  if(anglePos > (*(*fAngleDistrTable)(iPlace))(iTransfer)) break;
589  }
590  theta = std::sqrt((*fAngleDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1));
591 
592  // G4cout<<iTransfer<<" : theta = "<<theta<<G4endl;
593 
594  phi = twopi*G4UniformRand();
595  dirX = std::sin(theta)*std::cos(phi) ;
596  dirY = std::sin(theta)*std::sin(phi) ;
597  dirZ = std::cos(theta) ;
598  G4ThreeVector directionTR(dirX,dirY,dirZ);
599  directionTR.rotateUz(particleDir);
601  directionTR,
602  energyTR );
603  aParticleChange.AddSecondary(aPhotonTR);
604  }
605  }
606  }
607  else
608  {
609  if(iTkin == 0) // Tkin is too small, neglect of TR photon generation
610  {
611  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
612  }
613  else // general case: Tkin between two vectors of the material
614  {
615  E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1);
617  W = 1.0/(E2 - E1);
618  W1 = (E2 - TkinScaled)*W;
619  W2 = (TkinScaled - E1)*W;
620 
621  // G4cout<<iTkin<<" mean TR number = "<<(((*(*fEnergyDistrTable)(iPlace))(0)+
622  // (*(*fAngleDistrTable)(iPlace))(0))*W1 +
623  // ((*(*fEnergyDistrTable)(iPlace + 1))(0)+
624  // (*(*fAngleDistrTable)(iPlace + 1))(0))*W2)
625  // *chargeSq*0.5<<G4endl;
626 
627  numOfTR = G4Poisson((((*(*fEnergyDistrTable)(iPlace))(0)+
628  (*(*fAngleDistrTable)(iPlace))(0))*W1 +
629  ((*(*fEnergyDistrTable)(iPlace + 1))(0)+
630  (*(*fAngleDistrTable)(iPlace + 1))(0))*W2)
631  *chargeSq*0.5 );
632  if(numOfTR == 0)
633  {
634  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
635  }
636  else
637  {
638  // G4cout<<"Number of X-ray TR photons = "<<numOfTR<<G4endl;
639 
641  for(iTR=0;iTR<numOfTR;iTR++)
642  {
643  energyPos = ((*(*fEnergyDistrTable)(iPlace))(0)*W1+
644  (*(*fEnergyDistrTable)(iPlace + 1))(0)*W2)*G4UniformRand();
645  for(iTransfer=0;iTransfer<fBinTR-1;iTransfer++)
646  {
647  if(energyPos >= ((*(*fEnergyDistrTable)(iPlace))(iTransfer)*W1+
648  (*(*fEnergyDistrTable)(iPlace + 1))(iTransfer)*W2)) break;
649  }
650  energyTR = ((*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer))*W1+
651  ((*fEnergyDistrTable)(iPlace + 1)->GetLowEdgeEnergy(iTransfer))*W2;
652 
653  // G4cout<<"energyTR = "<<energyTR/keV<<"keV"<<G4endl;
654 
655  kinEnergy -= energyTR;
656  aParticleChange.ProposeEnergy(kinEnergy);
657 
658  anglePos = ((*(*fAngleDistrTable)(iPlace))(0)*W1+
659  (*(*fAngleDistrTable)(iPlace + 1))(0)*W2)*G4UniformRand();
660  for(iTransfer=0;iTransfer<fBinTR-1;iTransfer++)
661  {
662  if(anglePos > ((*(*fAngleDistrTable)(iPlace))(iTransfer)*W1+
663  (*(*fAngleDistrTable)(iPlace + 1))(iTransfer)*W2)) break;
664  }
665  theta = std::sqrt(((*fAngleDistrTable)(iPlace)->
666  GetLowEdgeEnergy(iTransfer-1))*W1+
667  ((*fAngleDistrTable)(iPlace + 1)->
668  GetLowEdgeEnergy(iTransfer-1))*W2);
669 
670  // G4cout<<iTransfer<<" : theta = "<<theta<<G4endl;
671 
672  phi = twopi*G4UniformRand();
673  dirX = std::sin(theta)*std::cos(phi) ;
674  dirY = std::sin(theta)*std::sin(phi) ;
675  dirZ = std::cos(theta) ;
676  G4ThreeVector directionTR(dirX,dirY,dirZ);
677  directionTR.rotateUz(particleDir);
679  directionTR,
680  energyTR );
681  aParticleChange.AddSecondary(aPhotonTR);
682  }
683  }
684  }
685  }
686  return &aParticleChange;
687 }
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
static G4int fBinTR
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4StepStatus GetStepStatus() const
G4double GetSurfaceTolerance() const
G4ParticleDefinition * GetDefinition() const
G4double GetLowEdgeEnergy(size_t binNumber) const
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * fAngleDistrTable
G4StepPoint * GetPreStepPoint() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4VPhysicalVolume * GetPhysicalVolume() const
const G4ThreeVector & GetMomentumDirection() const
float proton_mass_c2
Definition: hepunit.py:275
static G4int fTotBin
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4PhysicsTable * fEnergyDistrTable
virtual void Initialize(const G4Track &)
G4double GetPDGMass() const
void SetNumberOfSecondaries(G4int totSecondaries)
G4StepPoint * GetPostStepPoint() const
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void AddSecondary(G4Track *aSecondary)
G4State GetState() const
Definition: G4Material.hh:179
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
G4PhysicsLogVector * fProtonEnergyVector
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4double GetStepLength() const
static G4GeometryTolerance * GetInstance()
const G4Material * GetMaterial() const
G4double G4ForwardXrayTR::SpectralAngleTRdensity ( G4double  energy,
G4double  varAngle 
) const
virtual

Implements G4TransitionRadiation.

Definition at line 317 of file G4ForwardXrayTR.cc.

References fGamma, fSigma1, and fSigma2.

319 {
320  G4double formationLength1, formationLength2;
321  formationLength1 = 1.0/
322  (1.0/(fGamma*fGamma)
323  + fSigma1/(energy*energy)
324  + varAngle);
325  formationLength2 = 1.0/
326  (1.0/(fGamma*fGamma)
327  + fSigma2/(energy*energy)
328  + varAngle);
329  return (varAngle/energy)*(formationLength1 - formationLength2)
330  *(formationLength1 - formationLength2);
331 
332 }
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
double G4double
Definition: G4Types.hh:76
G4double G4ForwardXrayTR::SpectralDensity ( G4double  energy,
G4double  x 
) const

Definition at line 408 of file G4ForwardXrayTR.cc.

References test::a, test::b, energy(), fGamma, fSigma1, and fSigma2.

Referenced by AngleInterval().

410 {
411  G4double a, b;
412  a = 1.0/(fGamma*fGamma)
413  + fSigma1/(energy*energy) ;
414  b = 1.0/(fGamma*fGamma)
415  + fSigma2/(energy*energy) ;
416  return ( (a + b)*std::log((x + b)/(x + a))/(a - b)
417  + a/(x + a) + b/(x + b) )/energy;
418 
419 }
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
double G4double
Definition: G4Types.hh:76

Field Documentation

G4PhysicsTable* G4ForwardXrayTR::fAngleDistrTable
protected
G4int G4ForwardXrayTR::fBinTR = 50
staticprotected

Definition at line 163 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables(), GetBinTR(), GetEnergyTR(), and PostStepDoIt().

G4double G4ForwardXrayTR::fCofTR = fine_structure_const/pi
staticprotected

Definition at line 171 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

G4PhysicsTable* G4ForwardXrayTR::fEnergyDistrTable
protected
G4double G4ForwardXrayTR::fGamma
protected
const std::vector<G4double>* G4ForwardXrayTR::fGammaCutInKineticEnergy
protected

Definition at line 145 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables(), and G4ForwardXrayTR().

G4double G4ForwardXrayTR::fGammaTkinCut
protected

Definition at line 147 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables(), and G4ForwardXrayTR().

G4double G4ForwardXrayTR::fMaxEnergyTR
protected

Definition at line 159 of file G4ForwardXrayTR.hh.

Referenced by AngleSum(), BuildXrayTRtables(), and G4ForwardXrayTR().

G4double G4ForwardXrayTR::fMaxProtonTkin = 100.0*TeV
staticprotected

Definition at line 166 of file G4ForwardXrayTR.hh.

Referenced by G4ForwardXrayTR(), and GetMaxProtonTkin().

G4double G4ForwardXrayTR::fMaxThetaTR
protected

Definition at line 162 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables(), EnergySum(), and G4ForwardXrayTR().

G4double G4ForwardXrayTR::fMinEnergyTR
protected

Definition at line 158 of file G4ForwardXrayTR.hh.

Referenced by AngleSum(), BuildXrayTRtables(), and G4ForwardXrayTR().

G4double G4ForwardXrayTR::fMinProtonTkin = 100.0*GeV
staticprotected

Definition at line 165 of file G4ForwardXrayTR.hh.

Referenced by G4ForwardXrayTR(), and GetMinProtonTkin().

G4double G4ForwardXrayTR::fPlasmaCof
staticprotected
Initial value:

Definition at line 170 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

G4PhysicsLogVector* G4ForwardXrayTR::fProtonEnergyVector
protected
G4ParticleDefinition* G4ForwardXrayTR::fPtrGamma
protected

Definition at line 136 of file G4ForwardXrayTR.hh.

Referenced by G4ForwardXrayTR().

G4double G4ForwardXrayTR::fSigma1
protected
G4double G4ForwardXrayTR::fSigma2
protected
G4int G4ForwardXrayTR::fSympsonNumber = 100
staticprotected

Definition at line 154 of file G4ForwardXrayTR.hh.

Referenced by AngleSum(), EnergySum(), and GetSympsonNumber().

G4double G4ForwardXrayTR::fTheMaxAngle = 1.0e-3
staticprotected

Definition at line 160 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

G4double G4ForwardXrayTR::fTheMaxEnergyTR = 100.0*keV
staticprotected

Definition at line 157 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

G4double G4ForwardXrayTR::fTheMinAngle = 5.0e-6
staticprotected

Definition at line 161 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

G4double G4ForwardXrayTR::fTheMinEnergyTR = 1.0*keV
staticprotected

Definition at line 156 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

G4int G4ForwardXrayTR::fTotBin = 50
staticprotected

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