Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ForwardXrayTR Class Reference

#include <G4ForwardXrayTR.hh>

+ Inheritance diagram for G4ForwardXrayTR:

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) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
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 override
 
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) override
 
virtual G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
virtual G4double SpectralAngleTRdensity (G4double energy, G4double varAngle) const =0
 
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 &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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 ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 G4bool IsApplicable (const G4ParticleDefinition &)
 
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
 
virtual void ProcessDescription (std::ostream &outfile) 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 = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

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

virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 65 of file G4ForwardXrayTR.hh.

Constructor & Destructor Documentation

◆ G4ForwardXrayTR() [1/2]

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

Definition at line 81 of file G4ForwardXrayTR.cc.

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

◆ G4ForwardXrayTR() [2/2]

G4ForwardXrayTR::G4ForwardXrayTR ( const G4String processName = "XrayTR")
explicit

Definition at line 145 of file G4ForwardXrayTR.cc.

147 : G4TransitionRadiation(processName)
148{
149 fPtrGamma = nullptr;
150 fGammaCutInKineticEnergy = nullptr;
152 fSigma1 = fSigma2 = 0.0;
153 fAngleDistrTable = nullptr;
154 fEnergyDistrTable = nullptr;
156
157 // Proton energy vector initialization
158 //
161}

◆ ~G4ForwardXrayTR()

G4ForwardXrayTR::~G4ForwardXrayTR ( )
virtual

Definition at line 169 of file G4ForwardXrayTR.cc.

170{
171 delete fAngleDistrTable;
172 delete fEnergyDistrTable;
173 delete fProtonEnergyVector;
174}

Member Function Documentation

◆ AngleDensity()

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

Definition at line 339 of file G4ForwardXrayTR.cc.

341{
342 G4double x, x2, /*a, b,*/ c, d, f, a2, b2, a4, b4;
343 G4double cof1, cof2, cof3;
344 x = 1.0/energy;
345 x2 = x*x;
346 c = 1.0/fSigma1;
347 d = 1.0/fSigma2;
348 f = (varAngle + 1.0/(fGamma*fGamma));
349 a2 = c*f;
350 b2 = d*f;
351 a4 = a2*a2;
352 b4 = b2*b2;
353 //a = std::sqrt(a2);
354 // b = std::sqrt(b2);
355 cof1 = c*c*(0.5/(a2*(x2 +a2)) +0.5*std::log(x2/(x2 +a2))/a4);
356 cof3 = d*d*(0.5/(b2*(x2 +b2)) +0.5*std::log(x2/(x2 +b2))/b4);
357 cof2 = -c*d*(std::log(x2/(x2 +b2))/b2 - std::log(x2/(x2 +a2))/a2)/(a2 - b2) ;
358 return -varAngle*(cof1 + cof2 + cof3);
359}
double G4double
Definition: G4Types.hh:83
G4double energy(const ThreeVector &p, const G4double m)

Referenced by EnergyInterval().

◆ AngleInterval()

G4double G4ForwardXrayTR::AngleInterval ( G4double  energy,
G4double  varAngle1,
G4double  varAngle2 
) const

Definition at line 426 of file G4ForwardXrayTR.cc.

429{
430 return SpectralDensity(energy,varAngle2)
431 - SpectralDensity(energy,varAngle1);
432}
G4double SpectralDensity(G4double energy, G4double x) const

Referenced by EnergySum().

◆ AngleSum()

G4double G4ForwardXrayTR::AngleSum ( G4double  varAngle1,
G4double  varAngle2 
) const

Definition at line 381 of file G4ForwardXrayTR.cc.

383{
384 G4int i;
385 G4double h , sumEven = 0.0 , sumOdd = 0.0;
386 h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber;
387 for(i=1;i<fSympsonNumber;i++)
388 {
389 sumEven += EnergyInterval(fMinEnergyTR,fMaxEnergyTR,varAngle1 + 2*i*h );
391 varAngle1 + (2*i - 1)*h );
392 }
394 varAngle1 + (2*fSympsonNumber - 1)*h );
395
396 return h*(EnergyInterval(fMinEnergyTR,fMaxEnergyTR,varAngle1)
398 + 4.0*sumOdd + 2.0*sumEven )/3.0;
399}
static G4int fSympsonNumber
G4double EnergyInterval(G4double energy1, G4double energy2, G4double varAngle) const

Referenced by BuildXrayTRtables().

◆ BuildXrayTRtables()

void G4ForwardXrayTR::BuildXrayTRtables ( )

Definition at line 188 of file G4ForwardXrayTR.cc.

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

Referenced by G4ForwardXrayTR().

◆ EnergyInterval()

G4double G4ForwardXrayTR::EnergyInterval ( G4double  energy1,
G4double  energy2,
G4double  varAngle 
) const

Definition at line 367 of file G4ForwardXrayTR.cc.

370{
371 return AngleDensity(energy2,varAngle)
372 - AngleDensity(energy1,varAngle);
373}
G4double AngleDensity(G4double energy, G4double varAngle) const

Referenced by AngleSum().

◆ EnergySum()

G4double G4ForwardXrayTR::EnergySum ( G4double  energy1,
G4double  energy2 
) const

Definition at line 440 of file G4ForwardXrayTR.cc.

442{
443 G4int i;
444 G4double h , sumEven = 0.0 , sumOdd = 0.0;
445 h = 0.5*(energy2 - energy1)/fSympsonNumber;
446 for(i=1;i<fSympsonNumber;i++)
447 {
448 sumEven += AngleInterval(energy1 + 2*i*h,0.0,fMaxThetaTR);
449 sumOdd += AngleInterval(energy1 + (2*i - 1)*h,0.0,fMaxThetaTR);
450 }
451 sumOdd += AngleInterval(energy1 + (2*fSympsonNumber - 1)*h,
452 0.0,fMaxThetaTR);
453
454 return h*( AngleInterval(energy1,0.0,fMaxThetaTR)
455 + AngleInterval(energy2,0.0,fMaxThetaTR)
456 + 4.0*sumOdd + 2.0*sumEven )/3.0;
457}
G4double AngleInterval(G4double energy, G4double varAngle1, G4double varAngle2) const

Referenced by BuildXrayTRtables().

◆ GetAngleDistrTable()

G4PhysicsTable * G4ForwardXrayTR::GetAngleDistrTable ( )

◆ GetBinTR()

G4int G4ForwardXrayTR::GetBinTR ( )
static

Definition at line 813 of file G4ForwardXrayTR.cc.

814{
815 return fBinTR;
816}

◆ GetEnergyDistrTable()

G4PhysicsTable * G4ForwardXrayTR::GetEnergyDistrTable ( )

◆ GetEnergyTR()

G4double G4ForwardXrayTR::GetEnergyTR ( G4int  iMat,
G4int  jMat,
G4int  iTkin 
) const

Definition at line 695 of file G4ForwardXrayTR.cc.

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

◆ GetMaxProtonTkin()

G4double G4ForwardXrayTR::GetMaxProtonTkin ( )
static

Definition at line 823 of file G4ForwardXrayTR.cc.

824{
825 return fMaxProtonTkin;
826}

◆ GetMeanFreePath()

G4double G4ForwardXrayTR::GetMeanFreePath ( const G4Track ,
G4double  ,
G4ForceCondition condition 
)
overridevirtual

Reimplemented from G4TransitionRadiation.

Definition at line 176 of file G4ForwardXrayTR.cc.

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

◆ GetMinProtonTkin()

G4double G4ForwardXrayTR::GetMinProtonTkin ( )
static

Definition at line 818 of file G4ForwardXrayTR.cc.

819{
820 return fMinProtonTkin;
821}

◆ GetSympsonNumber()

G4int G4ForwardXrayTR::GetSympsonNumber ( )
static

Definition at line 808 of file G4ForwardXrayTR.cc.

809{
810 return fSympsonNumber;
811}

◆ GetThetaTR()

G4double G4ForwardXrayTR::GetThetaTR ( G4int  iMat,
G4int  jMat,
G4int  iTkin 
) const

Definition at line 803 of file G4ForwardXrayTR.cc.

804{
805 return 0.0;
806}

◆ GetTotBin()

G4int G4ForwardXrayTR::GetTotBin ( )
static

Definition at line 828 of file G4ForwardXrayTR.cc.

829{
830 return fTotBin;
831}

◆ PostStepDoIt()

G4VParticleChange * G4ForwardXrayTR::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
overridevirtual

Reimplemented from G4TransitionRadiation.

Definition at line 465 of file G4ForwardXrayTR.cc.

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

◆ SpectralAngleTRdensity()

G4double G4ForwardXrayTR::SpectralAngleTRdensity ( G4double  energy,
G4double  varAngle 
) const
overridevirtual

Implements G4TransitionRadiation.

Definition at line 316 of file G4ForwardXrayTR.cc.

318{
319 G4double formationLength1, formationLength2;
320 formationLength1 = 1.0/
321 (1.0/(fGamma*fGamma)
322 + fSigma1/(energy*energy)
323 + varAngle);
324 formationLength2 = 1.0/
325 (1.0/(fGamma*fGamma)
326 + fSigma2/(energy*energy)
327 + varAngle);
328 return (varAngle/energy)*(formationLength1 - formationLength2)
329 *(formationLength1 - formationLength2);
330
331}

◆ SpectralDensity()

G4double G4ForwardXrayTR::SpectralDensity ( G4double  energy,
G4double  x 
) const

Definition at line 407 of file G4ForwardXrayTR.cc.

409{
410 G4double a, b;
411 a = 1.0/(fGamma*fGamma)
412 + fSigma1/(energy*energy) ;
413 b = 1.0/(fGamma*fGamma)
414 + fSigma2/(energy*energy) ;
415 return ( (a + b)*std::log((x + b)/(x + a))/(a - b)
416 + a/(x + a) + b/(x + b) )/energy;
417
418}

Referenced by AngleInterval().

Member Data Documentation

◆ fAngleDistrTable

G4PhysicsTable* G4ForwardXrayTR::fAngleDistrTable
protected

◆ fBinTR

G4int G4ForwardXrayTR::fBinTR = 50
staticprotected

Definition at line 162 of file G4ForwardXrayTR.hh.

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

◆ fCofTR

G4double G4ForwardXrayTR::fCofTR = fine_structure_const/pi
staticprotected

Definition at line 170 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

◆ fEnergyDistrTable

G4PhysicsTable* G4ForwardXrayTR::fEnergyDistrTable
protected

◆ fGamma

G4double G4ForwardXrayTR::fGamma
protected

◆ fGammaCutInKineticEnergy

const std::vector<G4double>* G4ForwardXrayTR::fGammaCutInKineticEnergy
protected

Definition at line 144 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables(), and G4ForwardXrayTR().

◆ fGammaTkinCut

G4double G4ForwardXrayTR::fGammaTkinCut
protected

Definition at line 146 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables(), and G4ForwardXrayTR().

◆ fMaxEnergyTR

G4double G4ForwardXrayTR::fMaxEnergyTR
protected

Definition at line 158 of file G4ForwardXrayTR.hh.

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

◆ fMaxProtonTkin

G4double G4ForwardXrayTR::fMaxProtonTkin = 100.0*TeV
staticprotected

Definition at line 165 of file G4ForwardXrayTR.hh.

Referenced by G4ForwardXrayTR(), and GetMaxProtonTkin().

◆ fMaxThetaTR

G4double G4ForwardXrayTR::fMaxThetaTR
protected

Definition at line 161 of file G4ForwardXrayTR.hh.

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

◆ fMinEnergyTR

G4double G4ForwardXrayTR::fMinEnergyTR
protected

Definition at line 157 of file G4ForwardXrayTR.hh.

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

◆ fMinProtonTkin

G4double G4ForwardXrayTR::fMinProtonTkin = 100.0*GeV
staticprotected

Definition at line 164 of file G4ForwardXrayTR.hh.

Referenced by G4ForwardXrayTR(), and GetMinProtonTkin().

◆ fPlasmaCof

G4double G4ForwardXrayTR::fPlasmaCof
staticprotected
Initial value:
= 4.0*pi*fine_structure_const*
hbarc*hbarc*hbarc/electron_mass_c2

Definition at line 169 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

◆ fProtonEnergyVector

G4PhysicsLogVector* G4ForwardXrayTR::fProtonEnergyVector
protected

◆ fPtrGamma

G4ParticleDefinition* G4ForwardXrayTR::fPtrGamma
protected

Definition at line 142 of file G4ForwardXrayTR.hh.

Referenced by G4ForwardXrayTR().

◆ fSigma1

G4double G4ForwardXrayTR::fSigma1
protected

◆ fSigma2

G4double G4ForwardXrayTR::fSigma2
protected

◆ fSympsonNumber

G4int G4ForwardXrayTR::fSympsonNumber = 100
staticprotected

Definition at line 153 of file G4ForwardXrayTR.hh.

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

◆ fTheMaxAngle

G4double G4ForwardXrayTR::fTheMaxAngle = 1.0e-3
staticprotected

Definition at line 159 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

◆ fTheMaxEnergyTR

G4double G4ForwardXrayTR::fTheMaxEnergyTR = 100.0*keV
staticprotected

Definition at line 156 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

◆ fTheMinAngle

G4double G4ForwardXrayTR::fTheMinAngle = 5.0e-6
staticprotected

Definition at line 160 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

◆ fTheMinEnergyTR

G4double G4ForwardXrayTR::fTheMinEnergyTR = 1.0*keV
staticprotected

Definition at line 155 of file G4ForwardXrayTR.hh.

Referenced by BuildXrayTRtables().

◆ fTotBin

G4int G4ForwardXrayTR::fTotBin = 50
staticprotected

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