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

#include <G4VXTRenergyLoss.hh>

+ Inheritance diagram for G4VXTRenergyLoss:

Public Member Functions

 G4VXTRenergyLoss (G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="XTRenergyLoss", G4ProcessType type=fElectromagnetic)
 
virtual ~G4VXTRenergyLoss ()
 
virtual G4double GetStackFactor (G4double energy, G4double gamma, G4double varAngle)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void BuildEnergyTable ()
 
void BuildAngleForEnergyBank ()
 
void BuildTable ()
 
void BuildAngleTable ()
 
void BuildGlobalAngleTable ()
 
G4complex OneInterfaceXTRdEdx (G4double energy, G4double gamma, G4double varAngle)
 
G4double SpectralAngleXTRdEdx (G4double varAngle)
 
virtual G4double SpectralXTRdEdx (G4double energy)
 
G4double AngleSpectralXTRdEdx (G4double energy)
 
G4double AngleXTRdEdx (G4double varAngle)
 
G4double OneBoundaryXTRNdensity (G4double energy, G4double gamma, G4double varAngle) const
 
G4double XTRNSpectralAngleDensity (G4double varAngle)
 
G4double XTRNSpectralDensity (G4double energy)
 
G4double XTRNAngleSpectralDensity (G4double energy)
 
G4double XTRNAngleDensity (G4double varAngle)
 
void GetNumberOfPhotons ()
 
G4double GetPlateFormationZone (G4double, G4double, G4double)
 
G4complex GetPlateComplexFZ (G4double, G4double, G4double)
 
void ComputePlatePhotoAbsCof ()
 
G4double GetPlateLinearPhotoAbs (G4double)
 
void GetPlateZmuProduct ()
 
G4double GetPlateZmuProduct (G4double, G4double, G4double)
 
G4double GetGasFormationZone (G4double, G4double, G4double)
 
G4complex GetGasComplexFZ (G4double, G4double, G4double)
 
void ComputeGasPhotoAbsCof ()
 
G4double GetGasLinearPhotoAbs (G4double)
 
void GetGasZmuProduct ()
 
G4double GetGasZmuProduct (G4double, G4double, G4double)
 
G4double GetPlateCompton (G4double)
 
G4double GetGasCompton (G4double)
 
G4double GetComptonPerAtom (G4double, G4double)
 
G4double GetXTRrandomEnergy (G4double scaledTkin, G4int iTkin)
 
G4double GetXTRenergy (G4int iPlace, G4double position, G4int iTransfer)
 
G4double GetRandomAngle (G4double energyXTR, G4int iTkin)
 
G4double GetAngleXTR (G4int iTR, G4double position, G4int iAngle)
 
G4double GetGamma ()
 
G4double GetEnergy ()
 
G4double GetVarAngle ()
 
void SetGamma (G4double gamma)
 
void SetEnergy (G4double energy)
 
void SetVarAngle (G4double varAngle)
 
void SetAngleRadDistr (G4bool pAngleRadDistr)
 
void SetCompton (G4bool pC)
 
G4PhysicsLogVectorGetProtonVector ()
 
G4int GetTotBin ()
 
G4PhysicsFreeVectorGetAngleVector (G4double energy, G4int n)
 
- 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 &)
 

Protected Attributes

G4ParticleDefinitionfPtrGamma
 
G4doublefGammaCutInKineticEnergy
 
G4double fGammaTkinCut
 
G4LogicalVolumefEnvelope
 
G4PhysicsTablefAngleDistrTable
 
G4PhysicsTablefEnergyDistrTable
 
G4PhysicsLogVectorfProtonEnergyVector
 
G4PhysicsLogVectorfXTREnergyVector
 
G4double fTheMinEnergyTR
 
G4double fTheMaxEnergyTR
 
G4double fMinEnergyTR
 
G4double fMaxEnergyTR
 
G4double fTheMaxAngle
 
G4double fTheMinAngle
 
G4double fMaxThetaTR
 
G4int fBinTR
 
G4double fMinProtonTkin
 
G4double fMaxProtonTkin
 
G4int fTotBin
 
G4double fGamma
 
G4double fEnergy
 
G4double fVarAngle
 
G4double fLambda
 
G4double fPlasmaCof
 
G4double fCofTR
 
G4bool fExitFlux
 
G4bool fAngleRadDistr
 
G4bool fCompton
 
G4double fSigma1
 
G4double fSigma2
 
G4int fMatIndex1
 
G4int fMatIndex2
 
G4int fPlateNumber
 
G4double fTotalDist
 
G4double fPlateThick
 
G4double fGasThick
 
G4double fAlphaPlate
 
G4double fAlphaGas
 
G4SandiaTablefPlatePhotoAbsCof
 
G4SandiaTablefGasPhotoAbsCof
 
G4ParticleChange fParticleChange
 
G4PhysicsTablefAngleForEnergyTable
 
std::vector< G4PhysicsTable * > fAngleBank
 
- 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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
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 73 of file G4VXTRenergyLoss.hh.

Constructor & Destructor Documentation

◆ G4VXTRenergyLoss()

G4VXTRenergyLoss::G4VXTRenergyLoss ( G4LogicalVolume anEnvelope,
G4Material foilMat,
G4Material gasMat,
G4double  a,
G4double  b,
G4int  n,
const G4String processName = "XTRenergyLoss",
G4ProcessType  type = fElectromagnetic 
)
explicit

Definition at line 63 of file G4VXTRenergyLoss.cc.

67 :
68 G4VDiscreteProcess(processName, type),
70 fGammaTkinCut(0.0),
71 fAngleDistrTable(nullptr),
72 fEnergyDistrTable(nullptr),
73 fPlatePhotoAbsCof(nullptr),
74 fGasPhotoAbsCof(nullptr),
76{
77 verboseLevel = 1;
79
80 fPtrGamma = nullptr;
83
84 // Initialization of local constants
85 fTheMinEnergyTR = 1.0*keV;
86 fTheMaxEnergyTR = 100.0*keV;
87
88 fTheMaxAngle = 1.0e-2; // 100 mrad
89
90 fTheMinAngle = 2.5e-5;
91
92 fBinTR = 200; // 100; // 50;
93
94 fMinProtonTkin = 100.0*GeV;
95 fMaxProtonTkin = 100.0*TeV;
96 fTotBin = 50; // 100; //
97
98 // Proton energy vector initialization
101 fTotBin );
102
105 fBinTR );
106
107 fPlasmaCof = 4.0*pi*fine_structure_const*hbarc*hbarc*hbarc/electron_mass_c2;
108
109 fCofTR = fine_structure_const/pi;
110
111 fEnvelope = anEnvelope;
112
113 fPlateNumber = n;
114 if(verboseLevel > 0)
115 G4cout<<"### G4VXTRenergyLoss: the number of TR radiator plates = "
117 if(fPlateNumber == 0)
118 {
119 G4Exception("G4VXTRenergyLoss::G4VXTRenergyLoss()","VXTRELoss01",
120 FatalException,"No plates in X-ray TR radiator");
121 }
122 // default is XTR dEdx, not flux after radiator
123
124 fExitFlux = false;
125 fAngleRadDistr = true; // false;
126 fCompton = false;
127
129 // Mean thicknesses of plates and gas gaps
130
131 fPlateThick = a;
132 fGasThick = b;
134 if(verboseLevel > 0)
135 G4cout<<"total radiator thickness = "<<fTotalDist/cm<<" cm"<<G4endl;
136
137 // index of plate material
138 fMatIndex1 = foilMat->GetIndex();
139 if(verboseLevel > 0)
140 G4cout<<"plate material = "<<foilMat->GetName()<<G4endl;
141
142 // index of gas material
143 fMatIndex2 = gasMat->GetIndex();
144 if(verboseLevel > 0)
145 G4cout<<"gas material = "<<gasMat->GetName()<<G4endl;
146
147 // plasma energy squared for plate material
148
150 // fSigma1 = (20.9*eV)*(20.9*eV);
151 if(verboseLevel > 0)
152 G4cout<<"plate plasma energy = "<<std::sqrt(fSigma1)/eV<<" eV"<<G4endl;
153
154 // plasma energy squared for gas material
155
157 if(verboseLevel > 0)
158 G4cout<<"gas plasma energy = "<<std::sqrt(fSigma2)/eV<<" eV"<<G4endl;
159
160 // Compute cofs for preparation of linear photo absorption
161
164
166}
@ fTransitionRadiation
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetElectronDensity() const
Definition: G4Material.hh:215
const G4String & GetName() const
Definition: G4Material.hh:175
size_t GetIndex() const
Definition: G4Material.hh:258
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
G4PhysicsTable * fEnergyDistrTable
G4PhysicsTable * fAngleForEnergyTable
G4PhysicsLogVector * fProtonEnergyVector
G4PhysicsLogVector * fXTREnergyVector
G4SandiaTable * fPlatePhotoAbsCof
G4PhysicsTable * fAngleDistrTable
G4SandiaTable * fGasPhotoAbsCof
G4LogicalVolume * fEnvelope
G4ParticleDefinition * fPtrGamma
G4ParticleChange fParticleChange
G4double * fGammaCutInKineticEnergy
const G4double pi
#define DBL_MAX
Definition: templates.hh:62

◆ ~G4VXTRenergyLoss()

G4VXTRenergyLoss::~G4VXTRenergyLoss ( )
virtual

Definition at line 170 of file G4VXTRenergyLoss.cc.

171{
172 delete fProtonEnergyVector;
173 delete fXTREnergyVector;
174 if(fEnergyDistrTable) {
176 delete fEnergyDistrTable;
177 }
178 if(fAngleRadDistr) {
180 delete fAngleDistrTable;
181 }
185 }
186}
void clearAndDestroy()

Member Function Documentation

◆ AngleSpectralXTRdEdx()

G4double G4VXTRenergyLoss::AngleSpectralXTRdEdx ( G4double  energy)

Definition at line 948 of file G4VXTRenergyLoss.cc.

949{
950 G4double result = GetStackFactor(energy,fGamma,fVarAngle);
951 if(result < 0) result = 0.0;
952 return result;
953}
double G4double
Definition: G4Types.hh:83
virtual G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle)

◆ AngleXTRdEdx()

G4double G4VXTRenergyLoss::AngleXTRdEdx ( G4double  varAngle)

Definition at line 959 of file G4VXTRenergyLoss.cc.

960{
961 // G4cout<<"angle2 = "<<varAngle<<"; fGamma = "<<fGamma<<G4endl;
962
963 G4double result;
964 G4double sum = 0., tmp1, tmp2, tmp=0., cof1, cof2, cofMin, cofPHC, energy1, energy2;
965 G4int k, kMax, kMin, i;
966
967 cofPHC = twopi*hbarc;
968
969 cof1 = (fPlateThick + fGasThick)*(1./fGamma/fGamma + varAngle);
971
972 // G4cout<<"cof1 = "<<cof1<<"; cof2 = "<<cof2<<"; cofPHC = "<<cofPHC<<G4endl;
973
974 cofMin = std::sqrt(cof1*cof2);
975 cofMin /= cofPHC;
976
977 kMin = G4int(cofMin);
978 if (cofMin > kMin) kMin++;
979
980 kMax = kMin + 9; // 9; // kMin + G4int(tmp);
981
982 // G4cout<<"cofMin = "<<cofMin<<"; kMin = "<<kMin<<"; kMax = "<<kMax<<G4endl;
983
984 for( k = kMin; k <= kMax; k++ )
985 {
986 tmp1 = cofPHC*k;
987 tmp2 = std::sqrt(tmp1*tmp1-cof1*cof2);
988 energy1 = (tmp1+tmp2)/cof1;
989 energy2 = (tmp1-tmp2)/cof1;
990
991 for(i = 0; i < 2; i++)
992 {
993 if( i == 0 )
994 {
995 if (energy1 > fTheMaxEnergyTR || energy1 < fTheMinEnergyTR) continue;
996
997 tmp1 = ( energy1*energy1*(1./fGamma/fGamma + varAngle) + fSigma1 )
998 * fPlateThick/(4*hbarc*energy1);
999 tmp2 = std::sin(tmp1);
1000 tmp = energy1*tmp2*tmp2;
1001 tmp2 = fPlateThick/(4.*tmp1);
1002 tmp1 = hbarc*energy1/( energy1*energy1*(1./fGamma/fGamma + varAngle) + fSigma2 );
1003 tmp *= (tmp1-tmp2)*(tmp1-tmp2);
1004 tmp1 = cof1/(4.*hbarc) - cof2/(4.*hbarc*energy1*energy1);
1005 tmp2 = std::abs(tmp1);
1006
1007 if(tmp2 > 0.) tmp /= tmp2;
1008 else continue;
1009 }
1010 else
1011 {
1012 if (energy2 > fTheMaxEnergyTR || energy2 < fTheMinEnergyTR) continue;
1013
1014 tmp1 = ( energy2*energy2*(1./fGamma/fGamma + varAngle) + fSigma1 )
1015 * fPlateThick/(4.*hbarc*energy2);
1016 tmp2 = std::sin(tmp1);
1017 tmp = energy2*tmp2*tmp2;
1018 tmp2 = fPlateThick/(4.*tmp1);
1019 tmp1 = hbarc*energy2/( energy2*energy2*(1./fGamma/fGamma + varAngle) + fSigma2 );
1020 tmp *= (tmp1-tmp2)*(tmp1-tmp2);
1021 tmp1 = cof1/(4.*hbarc) - cof2/(4.*hbarc*energy2*energy2);
1022 tmp2 = std::abs(tmp1);
1023
1024 if(tmp2 > 0.) tmp /= tmp2;
1025 else continue;
1026 }
1027 sum += tmp;
1028 }
1029 // G4cout<<"k = "<<k<<"; energy1 = "<<energy1/keV<<" keV; energy2 = "<<energy2/keV
1030 // <<" keV; tmp = "<<tmp<<"; sum = "<<sum<<G4endl;
1031 }
1032 result = 4.*pi*fPlateNumber*sum*varAngle;
1033 result /= hbarc*hbarc;
1034
1035 // old code based on general numeric integration
1036 // fVarAngle = varAngle;
1037 // G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
1038 // result = integral.Legendre10(this,&G4VXTRenergyLoss::AngleSpectralXTRdEdx,
1039 // fMinEnergyTR,fMaxEnergyTR);
1040 return result;
1041}
int G4int
Definition: G4Types.hh:85

Referenced by BuildGlobalAngleTable().

◆ BuildAngleForEnergyBank()

void G4VXTRenergyLoss::BuildAngleForEnergyBank ( )

Definition at line 392 of file G4VXTRenergyLoss.cc.

393{
394 if( this->GetProcessName() == "TranspRegXTRadiator" ||
395 this->GetProcessName() == "TranspRegXTRmodel" ||
396 this->GetProcessName() == "RegularXTRadiator" ||
397 this->GetProcessName() == "RegularXTRmodel" )
398 {
400 return;
401 }
402 G4int i, iTkin, iTR;
403 G4double angleSum = 0.0;
404
405
406 fGammaTkinCut = 0.0;
407
408 // setting of min/max TR energies
409
412
415
418 fBinTR );
419
421
422 G4cout.precision(4);
423 G4Timer timer;
424 timer.Start();
425
426 for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
427 {
428
429 fGamma = 1.0 + (fProtonEnergyVector->
430 GetLowEdgeEnergy(iTkin)/proton_mass_c2);
431
432 fMaxThetaTR = 25*2500.0/(fGamma*fGamma) ; // theta^2
433
434 fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
435
438
440
441 for( iTR = 0; iTR < fBinTR; iTR++ )
442 {
443 angleSum = 0.0;
444 fEnergy = energyVector->GetLowEdgeEnergy(iTR);
445 G4PhysicsLinearVector* angleVector = new G4PhysicsLinearVector(0.0,
447 fBinTR );
448
449 angleVector ->PutValue(fBinTR - 1, angleSum);
450
451 for( i = fBinTR - 2; i >= 0; i-- )
452 {
453 // Legendre96 or Legendre10
454
455 angleSum += integral.Legendre10(
457 angleVector->GetLowEdgeEnergy(i),
458 angleVector->GetLowEdgeEnergy(i+1) );
459
460 angleVector ->PutValue(i, angleSum);
461 }
462 fAngleForEnergyTable->insertAt(iTR, angleVector);
463 }
465 }
466 timer.Stop();
467 G4cout.precision(6);
468 if(verboseLevel > 0)
469 {
470 G4cout<<G4endl;
471 G4cout<<"total time for build X-ray TR angle for energy loss tables = "
472 <<timer.GetUserElapsed()<<" s"<<G4endl;
473 }
474 fGamma = 0.;
475 delete energyVector;
476}
void insertAt(std::size_t, G4PhysicsVector *)
G4double GetLowEdgeEnergy(std::size_t binNumber) const
void PutValue(std::size_t index, G4double theValue)
void Stop()
G4double GetUserElapsed() const
Definition: G4Timer.cc:143
void Start()
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
G4double SpectralAngleXTRdEdx(G4double varAngle)
std::vector< G4PhysicsTable * > fAngleBank

Referenced by BuildPhysicsTable().

◆ BuildAngleTable()

void G4VXTRenergyLoss::BuildAngleTable ( )

Definition at line 483 of file G4VXTRenergyLoss.cc.

484{
485 G4int iTkin, iTR;
487
488 fGammaTkinCut = 0.0;
489
490 // setting of min/max TR energies
491
494
497
498 G4cout.precision(4);
499 G4Timer timer;
500 timer.Start();
501 if(verboseLevel > 0)
502 {
503 G4cout<<G4endl;
504 G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
505 G4cout<<G4endl;
506 }
507 for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
508 {
509
510 fGamma = 1.0 + (fProtonEnergyVector->
511 GetLowEdgeEnergy(iTkin)/proton_mass_c2);
512
513 fMaxThetaTR = 25*2500.0/(fGamma*fGamma); // theta^2
514
515 fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
516
518 else
519 {
521 }
522
524
525 for( iTR = 0; iTR < fBinTR; iTR++ )
526 {
527 // energy = fMinEnergyTR*(iTR+1);
528
530
531 G4PhysicsFreeVector* angleVector = GetAngleVector(energy,fBinTR);
532 // G4cout<<G4endl;
533
534 fAngleForEnergyTable->insertAt(iTR,angleVector);
535 }
536
538 }
539 timer.Stop();
540 G4cout.precision(6);
541 if(verboseLevel > 0)
542 {
543 G4cout<<G4endl;
544 G4cout<<"total time for build XTR angle for given energy tables = "
545 <<timer.GetUserElapsed()<<" s"<<G4endl;
546 }
547 fGamma = 0.;
548
549 return;
550}
G4PhysicsFreeVector * GetAngleVector(G4double energy, G4int n)
G4double energy(const ThreeVector &p, const G4double m)

Referenced by BuildAngleForEnergyBank().

◆ BuildEnergyTable()

void G4VXTRenergyLoss::BuildEnergyTable ( )

Definition at line 300 of file G4VXTRenergyLoss.cc.

301{
302 G4int iTkin, iTR, iPlace;
303 G4double radiatorCof = 1.0; // for tuning of XTR yield
304 G4double energySum = 0.0;
305
308
309 fGammaTkinCut = 0.0;
310
311 // setting of min/max TR energies
312
315
318
320
321 G4cout.precision(4);
322 G4Timer timer;
323 timer.Start();
324
325 if(verboseLevel > 0)
326 {
327 G4cout<<G4endl;
328 G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
329 G4cout<<G4endl;
330 }
331 for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
332 {
335 fBinTR );
336
337 fGamma = 1.0 + (fProtonEnergyVector->
338 GetLowEdgeEnergy(iTkin)/proton_mass_c2);
339
340 fMaxThetaTR = 25.*2500.0/(fGamma*fGamma) ; // theta^2
341
342 fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
343
346
347 energySum = 0.0;
348
349 energyVector->PutValue(fBinTR-1,energySum);
350
351 for( iTR = fBinTR - 2; iTR >= 0; iTR-- )
352 {
353 // Legendre96 or Legendre10
354
355 energySum += radiatorCof*fCofTR*integral.Legendre10(
357 energyVector->GetLowEdgeEnergy(iTR),
358 energyVector->GetLowEdgeEnergy(iTR+1) );
359
360 energyVector->PutValue(iTR,energySum/fTotalDist);
361 }
362 iPlace = iTkin;
363 fEnergyDistrTable->insertAt(iPlace,energyVector);
364
365 if(verboseLevel > 0)
366 {
367 G4cout
368 // <<iTkin<<"\t"
369 // <<"fGamma = "
370 <<fGamma<<"\t" // <<" fMaxThetaTR = "<<fMaxThetaTR
371 // <<"sumN = "
372 <<energySum // <<"; sumA = "<<angleSum
373 <<G4endl;
374 }
375 }
376 timer.Stop();
377 G4cout.precision(6);
378 if(verboseLevel > 0)
379 {
380 G4cout<<G4endl;
381 G4cout<<"total time for build X-ray TR energy loss tables = "
382 <<timer.GetUserElapsed()<<" s"<<G4endl;
383 }
384 fGamma = 0.;
385 return;
386}
virtual G4double SpectralXTRdEdx(G4double energy)

Referenced by BuildPhysicsTable().

◆ BuildGlobalAngleTable()

void G4VXTRenergyLoss::BuildGlobalAngleTable ( )

Definition at line 638 of file G4VXTRenergyLoss.cc.

639{
640 G4int iTkin, iTR, iPlace;
641 G4double radiatorCof = 1.0; // for tuning of XTR yield
642 G4double angleSum;
644
645 fGammaTkinCut = 0.0;
646
647 // setting of min/max TR energies
648
651
654
655 G4cout.precision(4);
656 G4Timer timer;
657 timer.Start();
658 if(verboseLevel > 0) {
659 G4cout<<G4endl;
660 G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
661 G4cout<<G4endl;
662 }
663 for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
664 {
665
666 fGamma = 1.0 + (fProtonEnergyVector->
667 GetLowEdgeEnergy(iTkin)/proton_mass_c2);
668
669 fMaxThetaTR = 25.0/(fGamma*fGamma); // theta^2
670
671 fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
672
674 else
675 {
677 }
678 G4PhysicsLinearVector* angleVector = new G4PhysicsLinearVector(0.0,
680 fBinTR );
681
682 angleSum = 0.0;
683
685
686
687 angleVector->PutValue(fBinTR-1,angleSum);
688
689 for( iTR = fBinTR - 2; iTR >= 0; iTR-- )
690 {
691
692 angleSum += radiatorCof*fCofTR*integral.Legendre96(
694 angleVector->GetLowEdgeEnergy(iTR),
695 angleVector->GetLowEdgeEnergy(iTR+1) );
696
697 angleVector ->PutValue(iTR,angleSum);
698 }
699 if(verboseLevel > 1) {
700 G4cout
701 // <<iTkin<<"\t"
702 // <<"fGamma = "
703 <<fGamma<<"\t" // <<" fMaxThetaTR = "<<fMaxThetaTR
704 // <<"sumN = "<<energySum // <<"; sumA = "
705 <<angleSum
706 <<G4endl;
707 }
708 iPlace = iTkin;
709 fAngleDistrTable->insertAt(iPlace,angleVector);
710 }
711 timer.Stop();
712 G4cout.precision(6);
713 if(verboseLevel > 0) {
714 G4cout<<G4endl;
715 G4cout<<"total time for build X-ray TR angle tables = "
716 <<timer.GetUserElapsed()<<" s"<<G4endl;
717 }
718 fGamma = 0.;
719
720 return;
721}
G4double AngleXTRdEdx(G4double varAngle)

◆ BuildPhysicsTable()

void G4VXTRenergyLoss::BuildPhysicsTable ( const G4ParticleDefinition pd)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 275 of file G4VXTRenergyLoss.cc.

276{
277 if(pd.GetPDGCharge() == 0.)
278 {
279 G4Exception("G4VXTRenergyLoss::BuildPhysicsTable", "Notification", JustWarning,
280 "XTR initialisation for neutral particle ?!" );
281 }
283
284 if ( fAngleRadDistr )
285 {
286 if(verboseLevel > 0)
287 {
288 G4cout<<"Build angle for energy distribution according the current radiator"
289 <<G4endl;
290 }
292 }
293}
@ JustWarning
G4double GetPDGCharge() const

◆ BuildTable()

void G4VXTRenergyLoss::BuildTable ( )
inline

Definition at line 101 of file G4VXTRenergyLoss.hh.

101{} ;

◆ ComputeGasPhotoAbsCof()

void G4VXTRenergyLoss::ComputeGasPhotoAbsCof ( )

Definition at line 1161 of file G4VXTRenergyLoss.cc.

1162{
1163 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1164 const G4Material* mat = (*theMaterialTable)[fMatIndex2];
1166 return;
1167}
std::vector< G4Material * > G4MaterialTable
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:227
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:637

Referenced by G4VXTRenergyLoss().

◆ ComputePlatePhotoAbsCof()

void G4VXTRenergyLoss::ComputePlatePhotoAbsCof ( )

Definition at line 1086 of file G4VXTRenergyLoss.cc.

1087{
1088 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1089 const G4Material* mat = (*theMaterialTable)[fMatIndex1];
1091
1092 return;
1093}

Referenced by G4VXTRenergyLoss().

◆ GetAngleVector()

G4PhysicsFreeVector * G4VXTRenergyLoss::GetAngleVector ( G4double  energy,
G4int  n 
)

Definition at line 556 of file G4VXTRenergyLoss.cc.

557{
558 G4double theta=0., result, tmp=0., cof1, cof2, cofMin, cofPHC, angleSum = 0.;
559 G4int iTheta, k, /*kMax,*/ kMin;
560
561 G4PhysicsFreeVector* angleVector = new G4PhysicsFreeVector(n);
562
563 cofPHC = 4.*pi*hbarc;
564 tmp = (fSigma1 - fSigma2)/cofPHC/energy;
565 cof1 = fPlateThick*tmp;
566 cof2 = fGasThick*tmp;
567
569 cofMin += (fPlateThick*fSigma1 + fGasThick*fSigma2)/energy;
570 cofMin /= cofPHC;
571
572 kMin = G4int(cofMin);
573 if (cofMin > kMin) kMin++;
574
575 //kMax = kMin + fBinTR -1;
576
577 if(verboseLevel > 2)
578 {
579 G4cout<<"n-1 = "<<n-1<<"; theta = "
580 <<std::sqrt(fMaxThetaTR)*fGamma<<"; tmp = "
581 <<0.
582 <<"; angleSum = "<<angleSum<<G4endl;
583 }
584 // angleVector->PutValue(n-1,fMaxThetaTR, angleSum);
585
586 for( iTheta = n - 1; iTheta >= 1; iTheta-- )
587 {
588
589 k = iTheta - 1 + kMin;
590
591 tmp = pi*fPlateThick*(k + cof2)/(fPlateThick + fGasThick);
592
593 result = (k - cof1)*(k - cof1)*(k + cof2)*(k + cof2);
594
595 tmp = std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
596
597 if( k == kMin && kMin == G4int(cofMin) )
598 {
599 angleSum += 0.5*tmp; // 0.5*std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
600 }
601 else if(iTheta == n-1);
602 else
603 {
604 angleSum += tmp; // std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
605 }
606 theta = std::abs(k-cofMin)*cofPHC/energy/(fPlateThick + fGasThick);
607
608 if(verboseLevel > 2)
609 {
610 G4cout<<"iTheta = "<<iTheta<<"; k = "<<k<<"; theta = "
611 <<std::sqrt(theta)*fGamma<<"; tmp = "
612 <<tmp // std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result
613 <<"; angleSum = "<<angleSum<<G4endl;
614 }
615 angleVector->PutValue( iTheta, theta, angleSum );
616 }
617 if (theta > 0.)
618 {
619 angleSum += 0.5*tmp;
620 theta = 0.;
621 }
622 if(verboseLevel > 2)
623 {
624 G4cout<<"iTheta = "<<iTheta<<"; theta = "
625 <<std::sqrt(theta)*fGamma<<"; tmp = "
626 <<tmp
627 <<"; angleSum = "<<angleSum<<G4endl;
628 }
629 angleVector->PutValue( iTheta, theta, angleSum );
630
631 return angleVector;
632}
void PutValue(std::size_t index, G4double energy, G4double dValue)

Referenced by BuildAngleTable().

◆ GetAngleXTR()

G4double G4VXTRenergyLoss::GetAngleXTR ( G4int  iTR,
G4double  position,
G4int  iAngle 
)

Definition at line 1597 of file G4VXTRenergyLoss.cc.

1600{
1601 G4double x1, x2, y1, y2, result;
1602
1603 if( iTransfer == 0 )
1604 {
1605 result = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1606 }
1607 else
1608 {
1609 y1 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer-1);
1610 y2 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer);
1611
1612 x1 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
1613 x2 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1614
1615 if ( x1 == x2 ) result = x2;
1616 else
1617 {
1618 if ( y1 == y2 ) result = x1 + (x2 - x1)*G4UniformRand();
1619 else
1620 {
1621 result = x1 + (position - y1)*(x2 - x1)/(y2 - y1);
1622 }
1623 }
1624 }
1625 return result;
1626}
#define G4UniformRand()
Definition: Randomize.hh:52

Referenced by GetRandomAngle().

◆ GetComptonPerAtom()

G4double G4VXTRenergyLoss::GetComptonPerAtom ( G4double  GammaEnergy,
G4double  Z 
)

Definition at line 1317 of file G4VXTRenergyLoss.cc.

1318{
1319 G4double CrossSection = 0.0;
1320 if ( Z < 0.9999 ) return CrossSection;
1321 if ( GammaEnergy < 0.1*keV ) return CrossSection;
1322 if ( GammaEnergy > (100.*GeV/Z) ) return CrossSection;
1323
1324 static const G4double a = 20.0 , b = 230.0 , c = 440.0;
1325
1326 static const G4double
1327 d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527 *barn, d4=-1.9798e+1*barn,
1328 e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
1329 f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
1330
1331 G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
1332 p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
1333
1334 G4double T0 = 15.0*keV;
1335 if (Z < 1.5) T0 = 40.0*keV;
1336
1337 G4double X = std::max(GammaEnergy, T0) / electron_mass_c2;
1338 CrossSection = p1Z*std::log(1.+2.*X)/X
1339 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
1340
1341 // modification for low energy. (special case for Hydrogen)
1342
1343 if (GammaEnergy < T0)
1344 {
1345 G4double dT0 = 1.*keV;
1346 X = (T0+dT0) / electron_mass_c2;
1347 G4double sigma = p1Z*std::log(1.+2.*X)/X
1348 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
1349 G4double c1 = -T0*(sigma-CrossSection)/(CrossSection*dT0);
1350 G4double c2 = 0.150;
1351 if (Z > 1.5) c2 = 0.375-0.0556*std::log(Z);
1352 G4double y = std::log(GammaEnergy/T0);
1353 CrossSection *= std::exp(-y*(c1+c2*y));
1354 }
1355 // G4cout << "e= " << GammaEnergy << " Z= " << Z << " cross= " << CrossSection << G4endl;
1356 return CrossSection;
1357}

Referenced by GetGasCompton(), and GetPlateCompton().

◆ GetEnergy()

G4double G4VXTRenergyLoss::GetEnergy ( )
inline

Definition at line 164 of file G4VXTRenergyLoss.hh.

164{return fEnergy;};

◆ GetGamma()

G4double G4VXTRenergyLoss::GetGamma ( )
inline

Definition at line 163 of file G4VXTRenergyLoss.hh.

163{return fGamma;};

◆ GetGasComplexFZ()

G4complex G4VXTRenergyLoss::GetGasComplexFZ ( G4double  omega,
G4double  gamma,
G4double  varAngle 
)

Definition at line 1137 of file G4VXTRenergyLoss.cc.

1140{
1141 G4double cof, length,delta, real_v, image_v;
1142
1143 length = 0.5*GetGasFormationZone(omega,gamma,varAngle);
1144 delta = length*GetGasLinearPhotoAbs(omega);
1145 cof = 1.0/(1.0 + delta*delta);
1146
1147 real_v = length*cof;
1148 image_v = real_v*delta;
1149
1150 G4complex zone(real_v,image_v);
1151 return zone;
1152}
std::complex< G4double > G4complex
Definition: G4Types.hh:88
G4double GetGasFormationZone(G4double, G4double, G4double)
G4double GetGasLinearPhotoAbs(G4double)

Referenced by G4StrawTubeXTRadiator::GetStackFactor(), and OneInterfaceXTRdEdx().

◆ GetGasCompton()

G4double G4VXTRenergyLoss::GetGasCompton ( G4double  omega)

Definition at line 1293 of file G4VXTRenergyLoss.cc.

1294{
1295 G4int i, numberOfElements;
1296 G4double xSection = 0., nowZ, sumZ = 0.;
1297
1298 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1299 numberOfElements = (*theMaterialTable)[fMatIndex2]->GetNumberOfElements();
1300
1301 for( i = 0; i < numberOfElements; i++ )
1302 {
1303 nowZ = (*theMaterialTable)[fMatIndex2]->GetElement(i)->GetZ();
1304 sumZ += nowZ;
1305 xSection += GetComptonPerAtom(omega,nowZ); // *nowZ;
1306 }
1307 xSection /= sumZ;
1308 xSection *= (*theMaterialTable)[fMatIndex2]->GetElectronDensity();
1309 return xSection;
1310}
G4double GetComptonPerAtom(G4double, G4double)

Referenced by G4XTRTransparentRegRadModel::SpectralXTRdEdx().

◆ GetGasFormationZone()

◆ GetGasLinearPhotoAbs()

G4double G4VXTRenergyLoss::GetGasLinearPhotoAbs ( G4double  omega)

Definition at line 1174 of file G4VXTRenergyLoss.cc.

1175{
1176 G4double omega2, omega3, omega4;
1177
1178 omega2 = omega*omega;
1179 omega3 = omega2*omega;
1180 omega4 = omega2*omega2;
1181
1182 const G4double* SandiaCof = fGasPhotoAbsCof->GetSandiaCofForMaterial(omega);
1183 G4double cross = SandiaCof[0]/omega + SandiaCof[1]/omega2 +
1184 SandiaCof[2]/omega3 + SandiaCof[3]/omega4;
1185 return cross;
1186
1187}
G4double GetSandiaCofForMaterial(G4int, G4int) const

Referenced by GetGasComplexFZ(), GetGasZmuProduct(), G4GammaXTRadiator::GetStackFactor(), G4RegularXTRadiator::GetStackFactor(), G4StrawTubeXTRadiator::GetStackFactor(), G4TransparentRegXTRadiator::GetStackFactor(), G4XTRGammaRadModel::GetStackFactor(), G4XTRRegularRadModel::GetStackFactor(), G4XTRTransparentRegRadModel::GetStackFactor(), G4RegularXTRadiator::SpectralXTRdEdx(), G4XTRRegularRadModel::SpectralXTRdEdx(), and G4XTRTransparentRegRadModel::SpectralXTRdEdx().

◆ GetGasZmuProduct() [1/2]

void G4VXTRenergyLoss::GetGasZmuProduct ( )

Definition at line 1244 of file G4VXTRenergyLoss.cc.

1245{
1246 std::ofstream outGas("gasZmu.dat", std::ios::out );
1247 outGas.setf( std::ios::scientific, std::ios::floatfield );
1248 G4int i;
1249 G4double omega, varAngle, gamma;
1250 gamma = 10000.;
1251 varAngle = 1/gamma/gamma;
1252 if(verboseLevel > 0)
1253 G4cout<<"energy, keV"<<"\t"<<"Zmu for gas"<<G4endl;
1254 for(i=0;i<100;i++)
1255 {
1256 omega = (1.0 + i)*keV;
1257 if(verboseLevel > 1)
1258 G4cout<<omega/keV<<"\t"<<GetGasZmuProduct(omega,gamma,varAngle)<<"\t";
1259 if(verboseLevel > 0)
1260 outGas<<omega/keV<<"\t\t"<<GetGasZmuProduct(omega,gamma,varAngle)<<G4endl;
1261 }
1262 return;
1263}

Referenced by GetGasZmuProduct().

◆ GetGasZmuProduct() [2/2]

G4double G4VXTRenergyLoss::GetGasZmuProduct ( G4double  omega,
G4double  gamma,
G4double  varAngle 
)

Definition at line 1233 of file G4VXTRenergyLoss.cc.

1236{
1237 return GetGasFormationZone(omega,gamma,varAngle)*GetGasLinearPhotoAbs(omega);
1238}

◆ GetMeanFreePath()

G4double G4VXTRenergyLoss::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
overridevirtual

Implements G4VDiscreteProcess.

Definition at line 202 of file G4VXTRenergyLoss.cc.

205{
206 G4int iTkin, iPlace;
207 G4double lambda, sigma, kinEnergy, mass, gamma;
208 G4double charge, chargeSq, massRatio, TkinScaled;
209 G4double E1,E2,W,W1,W2;
210
212
213 if( aTrack.GetVolume()->GetLogicalVolume() != fEnvelope ) lambda = DBL_MAX;
214 else
215 {
216 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
217 kinEnergy = aParticle->GetKineticEnergy();
218 mass = aParticle->GetDefinition()->GetPDGMass();
219 gamma = 1.0 + kinEnergy/mass;
220 if(verboseLevel > 1)
221 {
222 G4cout<<" gamma = "<<gamma<<"; fGamma = "<<fGamma<<G4endl;
223 }
224
225 if ( std::fabs( gamma - fGamma ) < 0.05*gamma ) lambda = fLambda;
226 else
227 {
228 charge = aParticle->GetDefinition()->GetPDGCharge();
229 chargeSq = charge*charge;
230 massRatio = proton_mass_c2/mass;
231 TkinScaled = kinEnergy*massRatio;
232
233 for(iTkin = 0; iTkin < fTotBin; iTkin++)
234 {
235 if( TkinScaled < fProtonEnergyVector->GetLowEdgeEnergy(iTkin)) break;
236 }
237 iPlace = iTkin - 1;
238
239 if(iTkin == 0) lambda = DBL_MAX; // Tkin is too small, neglect of TR photon generation
240 else // general case: Tkin between two vectors of the material
241 {
242 if(iTkin == fTotBin)
243 {
244 sigma = (*(*fEnergyDistrTable)(iPlace))(0)*chargeSq;
245 }
246 else
247 {
248 E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1);
250 W = 1.0/(E2 - E1);
251 W1 = (E2 - TkinScaled)*W;
252 W2 = (TkinScaled - E1)*W;
253 sigma = ( (*(*fEnergyDistrTable)(iPlace ))(0)*W1 +
254 (*(*fEnergyDistrTable)(iPlace+1))(0)*W2 )*chargeSq;
255
256 }
257 if (sigma < DBL_MIN) lambda = DBL_MAX;
258 else lambda = 1./sigma;
259 fLambda = lambda;
260 fGamma = gamma;
261 if(verboseLevel > 1)
262 {
263 G4cout<<" lambda = "<<lambda/mm<<" mm"<<G4endl;
264 }
265 }
266 }
267 }
268 return lambda;
269}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4VPhysicalVolume * GetVolume() const
const G4DynamicParticle * GetDynamicParticle() const
G4LogicalVolume * GetLogicalVolume() const
#define DBL_MIN
Definition: templates.hh:54

◆ GetNumberOfPhotons()

void G4VXTRenergyLoss::GetNumberOfPhotons ( )

Definition at line 1449 of file G4VXTRenergyLoss.cc.

1450{
1451 G4int iTkin;
1452 G4double gamma, numberE;
1453
1454 std::ofstream outEn("numberE.dat", std::ios::out );
1455 outEn.setf( std::ios::scientific, std::ios::floatfield );
1456
1457 std::ofstream outAng("numberAng.dat", std::ios::out );
1458 outAng.setf( std::ios::scientific, std::ios::floatfield );
1459
1460 for(iTkin=0;iTkin<fTotBin;iTkin++) // Lorentz factor loop
1461 {
1462 gamma = 1.0 + (fProtonEnergyVector->
1463 GetLowEdgeEnergy(iTkin)/proton_mass_c2);
1464 numberE = (*(*fEnergyDistrTable)(iTkin))(0);
1465 // numberA = (*(*fAngleDistrTable)(iTkin))(0);
1466 if(verboseLevel > 1)
1467 G4cout<<gamma<<"\t\t"<<numberE<<"\t" // <<numberA
1468 <<G4endl;
1469 if(verboseLevel > 0)
1470 outEn<<gamma<<"\t\t"<<numberE<<G4endl;
1471 }
1472 return;
1473}

◆ GetPlateComplexFZ()

G4complex G4VXTRenergyLoss::GetPlateComplexFZ ( G4double  omega,
G4double  gamma,
G4double  varAngle 
)

Definition at line 1064 of file G4VXTRenergyLoss.cc.

1067{
1068 G4double cof, length,delta, real_v, image_v;
1069
1070 length = 0.5*GetPlateFormationZone(omega,gamma,varAngle);
1071 delta = length*GetPlateLinearPhotoAbs(omega);
1072 cof = 1.0/(1.0 + delta*delta);
1073
1074 real_v = length*cof;
1075 image_v = real_v*delta;
1076
1077 G4complex zone(real_v,image_v);
1078 return zone;
1079}
G4double GetPlateLinearPhotoAbs(G4double)
G4double GetPlateFormationZone(G4double, G4double, G4double)

Referenced by G4StrawTubeXTRadiator::GetStackFactor(), and OneInterfaceXTRdEdx().

◆ GetPlateCompton()

G4double G4VXTRenergyLoss::GetPlateCompton ( G4double  omega)

Definition at line 1269 of file G4VXTRenergyLoss.cc.

1270{
1271 G4int i, numberOfElements;
1272 G4double xSection = 0., nowZ, sumZ = 0.;
1273
1274 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1275 numberOfElements = (*theMaterialTable)[fMatIndex1]->GetNumberOfElements();
1276
1277 for( i = 0; i < numberOfElements; i++ )
1278 {
1279 nowZ = (*theMaterialTable)[fMatIndex1]->GetElement(i)->GetZ();
1280 sumZ += nowZ;
1281 xSection += GetComptonPerAtom(omega,nowZ); // *nowZ;
1282 }
1283 xSection /= sumZ;
1284 xSection *= (*theMaterialTable)[fMatIndex1]->GetElectronDensity();
1285 return xSection;
1286}

Referenced by G4XTRTransparentRegRadModel::SpectralXTRdEdx().

◆ GetPlateFormationZone()

G4double G4VXTRenergyLoss::GetPlateFormationZone ( G4double  omega,
G4double  gamma,
G4double  varAngle 
)

◆ GetPlateLinearPhotoAbs()

G4double G4VXTRenergyLoss::GetPlateLinearPhotoAbs ( G4double  omega)

Definition at line 1102 of file G4VXTRenergyLoss.cc.

1103{
1104 // G4int i;
1105 G4double omega2, omega3, omega4;
1106
1107 omega2 = omega*omega;
1108 omega3 = omega2*omega;
1109 omega4 = omega2*omega2;
1110
1111 const G4double* SandiaCof = fPlatePhotoAbsCof->GetSandiaCofForMaterial(omega);
1112 G4double cross = SandiaCof[0]/omega + SandiaCof[1]/omega2 +
1113 SandiaCof[2]/omega3 + SandiaCof[3]/omega4;
1114 return cross;
1115}

Referenced by GetPlateComplexFZ(), GetPlateZmuProduct(), G4GammaXTRadiator::GetStackFactor(), G4RegularXTRadiator::GetStackFactor(), G4StrawTubeXTRadiator::GetStackFactor(), G4TransparentRegXTRadiator::GetStackFactor(), G4XTRGammaRadModel::GetStackFactor(), G4XTRRegularRadModel::GetStackFactor(), G4XTRTransparentRegRadModel::GetStackFactor(), G4RegularXTRadiator::SpectralXTRdEdx(), G4XTRRegularRadModel::SpectralXTRdEdx(), and G4XTRTransparentRegRadModel::SpectralXTRdEdx().

◆ GetPlateZmuProduct() [1/2]

void G4VXTRenergyLoss::GetPlateZmuProduct ( )

Definition at line 1206 of file G4VXTRenergyLoss.cc.

1207{
1208 std::ofstream outPlate("plateZmu.dat", std::ios::out );
1209 outPlate.setf( std::ios::scientific, std::ios::floatfield );
1210
1211 G4int i;
1212 G4double omega, varAngle, gamma;
1213 gamma = 10000.;
1214 varAngle = 1/gamma/gamma;
1215 if(verboseLevel > 0)
1216 G4cout<<"energy, keV"<<"\t"<<"Zmu for plate"<<G4endl;
1217 for(i=0;i<100;i++)
1218 {
1219 omega = (1.0 + i)*keV;
1220 if(verboseLevel > 1)
1221 G4cout<<omega/keV<<"\t"<<GetPlateZmuProduct(omega,gamma,varAngle)<<"\t";
1222 if(verboseLevel > 0)
1223 outPlate<<omega/keV<<"\t\t"<<GetPlateZmuProduct(omega,gamma,varAngle)<<G4endl;
1224 }
1225 return;
1226}

Referenced by GetPlateZmuProduct().

◆ GetPlateZmuProduct() [2/2]

G4double G4VXTRenergyLoss::GetPlateZmuProduct ( G4double  omega,
G4double  gamma,
G4double  varAngle 
)

Definition at line 1194 of file G4VXTRenergyLoss.cc.

1197{
1198 return GetPlateFormationZone(omega,gamma,varAngle)
1199 * GetPlateLinearPhotoAbs(omega);
1200}

◆ GetProtonVector()

G4PhysicsLogVector * G4VXTRenergyLoss::GetProtonVector ( )
inline

Definition at line 173 of file G4VXTRenergyLoss.hh.

173{ return fProtonEnergyVector;};

◆ GetRandomAngle()

G4double G4VXTRenergyLoss::GetRandomAngle ( G4double  energyXTR,
G4int  iTkin 
)

Definition at line 1567 of file G4VXTRenergyLoss.cc.

1568{
1569 G4int iTR, iAngle;
1570 G4double position, angle;
1571
1572 if (iTkin == fTotBin) iTkin--;
1573
1575
1576 for( iTR = 0; iTR < fBinTR; iTR++ )
1577 {
1578 if( energyXTR < fXTREnergyVector->GetLowEdgeEnergy(iTR) ) break;
1579 }
1580 if (iTR == fBinTR) iTR--;
1581
1582 position = (*(*fAngleForEnergyTable)(iTR))(0)*G4UniformRand();
1583
1584 for( iAngle = 0;; iAngle++)
1585 {
1586 if( position >= (*(*fAngleForEnergyTable)(iTR))(iAngle) ) break;
1587 }
1588 angle = GetAngleXTR(iTR,position,iAngle);
1589 return angle;
1590}
G4double GetAngleXTR(G4int iTR, G4double position, G4int iAngle)
#define position
Definition: xmlparse.cc:622

Referenced by PostStepDoIt().

◆ GetStackFactor()

G4double G4VXTRenergyLoss::GetStackFactor ( G4double  energy,
G4double  gamma,
G4double  varAngle 
)
virtual

Reimplemented in G4GammaXTRadiator, G4RegularXTRadiator, G4StrawTubeXTRadiator, G4TransparentRegXTRadiator, G4XTRGammaRadModel, G4XTRRegularRadModel, and G4XTRTransparentRegRadModel.

Definition at line 1388 of file G4VXTRenergyLoss.cc.

1390{
1391 // return stack factor corresponding to one interface
1392
1393 return std::real( OneInterfaceXTRdEdx(energy,gamma,varAngle) );
1394}
G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma, G4double varAngle)

Referenced by AngleSpectralXTRdEdx(), SpectralAngleXTRdEdx(), XTRNAngleSpectralDensity(), and XTRNSpectralAngleDensity().

◆ GetTotBin()

G4int G4VXTRenergyLoss::GetTotBin ( )
inline

Definition at line 174 of file G4VXTRenergyLoss.hh.

174{return fTotBin;};

◆ GetVarAngle()

G4double G4VXTRenergyLoss::GetVarAngle ( )
inline

Definition at line 165 of file G4VXTRenergyLoss.hh.

165{return fVarAngle;};

◆ GetXTRenergy()

G4double G4VXTRenergyLoss::GetXTRenergy ( G4int  iPlace,
G4double  position,
G4int  iTransfer 
)

Definition at line 1531 of file G4VXTRenergyLoss.cc.

1534{
1535 G4double x1, x2, y1, y2, result;
1536
1537 if(iTransfer == 0)
1538 {
1539 result = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1540 }
1541 else
1542 {
1543 y1 = (*(*fEnergyDistrTable)(iPlace))(iTransfer-1);
1544 y2 = (*(*fEnergyDistrTable)(iPlace))(iTransfer);
1545
1546 x1 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
1547 x2 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1548
1549 if ( x1 == x2 ) result = x2;
1550 else
1551 {
1552 if ( y1 == y2 ) result = x1 + (x2 - x1)*G4UniformRand();
1553 else
1554 {
1555 // result = x1 + (position - y1)*(x2 - x1)/(y2 - y1);
1556 result = x1 + (x2 - x1)*G4UniformRand();
1557 }
1558 }
1559 }
1560 return result;
1561}

Referenced by GetXTRrandomEnergy().

◆ GetXTRrandomEnergy()

G4double G4VXTRenergyLoss::GetXTRrandomEnergy ( G4double  scaledTkin,
G4int  iTkin 
)

Definition at line 1480 of file G4VXTRenergyLoss.cc.

1481{
1482 G4int iTransfer, iPlace;
1483 G4double transfer = 0.0, position, E1, E2, W1, W2, W;
1484
1485 iPlace = iTkin - 1;
1486
1487 // G4cout<<"iPlace = "<<iPlace<<endl;
1488
1489 if(iTkin == fTotBin) // relativistic plato, try from left
1490 {
1491 position = (*(*fEnergyDistrTable)(iPlace))(0)*G4UniformRand();
1492
1493 for(iTransfer=0;;iTransfer++)
1494 {
1495 if(position >= (*(*fEnergyDistrTable)(iPlace))(iTransfer)) break;
1496 }
1497 transfer = GetXTRenergy(iPlace,position,iTransfer);
1498 }
1499 else
1500 {
1501 E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1);
1503 W = 1.0/(E2 - E1);
1504 W1 = (E2 - scaledTkin)*W;
1505 W2 = (scaledTkin - E1)*W;
1506
1507 position =( (*(*fEnergyDistrTable)(iPlace))(0)*W1 +
1508 (*(*fEnergyDistrTable)(iPlace+1))(0)*W2 )*G4UniformRand();
1509
1510 // G4cout<<position<<"\t";
1511
1512 for(iTransfer=0;;iTransfer++)
1513 {
1514 if( position >=
1515 ( (*(*fEnergyDistrTable)(iPlace))(iTransfer)*W1 +
1516 (*(*fEnergyDistrTable)(iPlace+1))(iTransfer)*W2) ) break;
1517 }
1518 transfer = GetXTRenergy(iPlace,position,iTransfer);
1519
1520 }
1521 // G4cout<<"XTR transfer = "<<transfer/keV<<" keV"<<endl;
1522 if(transfer < 0.0 ) transfer = 0.0;
1523 return transfer;
1524}
G4double GetXTRenergy(G4int iPlace, G4double position, G4int iTransfer)

Referenced by PostStepDoIt().

◆ IsApplicable()

G4bool G4VXTRenergyLoss::IsApplicable ( const G4ParticleDefinition particle)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 193 of file G4VXTRenergyLoss.cc.

194{
195 return ( particle.GetPDGCharge() != 0.0 );
196}

◆ OneBoundaryXTRNdensity()

G4double G4VXTRenergyLoss::OneBoundaryXTRNdensity ( G4double  energy,
G4double  gamma,
G4double  varAngle 
) const

Definition at line 1371 of file G4VXTRenergyLoss.cc.

1373{
1374 G4double formationLength1, formationLength2;
1375 formationLength1 = 1.0/
1376 (1.0/(gamma*gamma)
1377 + fSigma1/(energy*energy)
1378 + varAngle);
1379 formationLength2 = 1.0/
1380 (1.0/(gamma*gamma)
1381 + fSigma2/(energy*energy)
1382 + varAngle);
1383 return (varAngle/energy)*(formationLength1 - formationLength2)
1384 *(formationLength1 - formationLength2);
1385
1386}

Referenced by XTRNAngleSpectralDensity(), and XTRNSpectralAngleDensity().

◆ OneInterfaceXTRdEdx()

G4complex G4VXTRenergyLoss::OneInterfaceXTRdEdx ( G4double  energy,
G4double  gamma,
G4double  varAngle 
)

Definition at line 871 of file G4VXTRenergyLoss.cc.

874{
875 G4complex Z1 = GetPlateComplexFZ(energy,gamma,varAngle);
876 G4complex Z2 = GetGasComplexFZ(energy,gamma,varAngle);
877
878 G4complex zOut = (Z1 - Z2)*(Z1 - Z2)
879 * (varAngle*energy/hbarc/hbarc);
880 return zOut;
881
882}
G4complex GetPlateComplexFZ(G4double, G4double, G4double)
G4complex GetGasComplexFZ(G4double, G4double, G4double)

Referenced by GetStackFactor(), G4GammaXTRadiator::GetStackFactor(), G4RegularXTRadiator::GetStackFactor(), G4TransparentRegXTRadiator::GetStackFactor(), G4XTRGammaRadModel::GetStackFactor(), G4XTRRegularRadModel::GetStackFactor(), and G4XTRTransparentRegRadModel::GetStackFactor().

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 729 of file G4VXTRenergyLoss.cc.

731{
732 G4int iTkin /*, iPlace*/;
733 G4double energyTR, theta,theta2, phi, dirX, dirY, dirZ;
734
735
737
738 if(verboseLevel > 1)
739 {
740 G4cout<<"Start of G4VXTRenergyLoss::PostStepDoIt "<<G4endl;
741 G4cout<<"name of current material = "
743 }
744 if( aTrack.GetVolume()->GetLogicalVolume() != fEnvelope )
745 {
746 if(verboseLevel > 0)
747 {
748 G4cout<<"Go out from G4VXTRenergyLoss::PostStepDoIt: wrong volume "<<G4endl;
749 }
750 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
751 }
752 else
753 {
754 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
755 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
756
757 // Now we are ready to Generate one TR photon
758
759 G4double kinEnergy = aParticle->GetKineticEnergy();
760 G4double mass = aParticle->GetDefinition()->GetPDGMass();
761 G4double gamma = 1.0 + kinEnergy/mass;
762
763 if(verboseLevel > 1 )
764 {
765 G4cout<<"gamma = "<<gamma<<G4endl;
766 }
767 G4double massRatio = proton_mass_c2/mass;
768 G4double TkinScaled = kinEnergy*massRatio;
769 G4ThreeVector position = pPostStepPoint->GetPosition();
770 G4ParticleMomentum direction = aParticle->GetMomentumDirection();
771 G4double startTime = pPostStepPoint->GetGlobalTime();
772
773 for( iTkin = 0; iTkin < fTotBin; iTkin++ )
774 {
775 if(TkinScaled < fProtonEnergyVector->GetLowEdgeEnergy(iTkin)) break;
776 }
777 //iPlace = iTkin - 1;
778
779 if(iTkin == 0) // Tkin is too small, neglect of TR photon generation
780 {
781 if( verboseLevel > 0)
782 {
783 G4cout<<"Go out from G4VXTRenergyLoss::PostStepDoIt:iTkin = "<<iTkin<<G4endl;
784 }
785 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
786 }
787 else // general case: Tkin between two vectors of the material
788 {
790
791 energyTR = GetXTRrandomEnergy(TkinScaled,iTkin);
792
793 if( verboseLevel > 1)
794 {
795 G4cout<<"energyTR = "<<energyTR/keV<<" keV"<<G4endl;
796 }
797 if ( fAngleRadDistr )
798 {
799 // theta = std::fabs(G4RandGauss::shoot(0.0,pi/gamma));
800
801 theta2 = GetRandomAngle(energyTR,iTkin);
802
803 if( theta2 > 0.) theta = std::sqrt(theta2);
804 else theta = 0.; // theta2;
805 }
806 else theta = std::fabs(G4RandGauss::shoot(0.0,pi/gamma));
807
808 // theta = 0.; // check no spread
809
810 if( theta >= 0.1 ) theta = 0.1;
811
812 // G4cout<<" : theta = "<<theta<<endl;
813
814 phi = twopi*G4UniformRand();
815
816 dirX = std::sin(theta)*std::cos(phi);
817 dirY = std::sin(theta)*std::sin(phi);
818 dirZ = std::cos(theta);
819
820 G4ThreeVector directionTR(dirX,dirY,dirZ);
821 directionTR.rotateUz(direction);
822 directionTR.unit();
823
825 directionTR, energyTR);
826
827 // A XTR photon is set on the particle track inside the radiator
828 // and is moved to the G4Envelope surface for standard X-ray TR models
829 // only. The case of fExitFlux=true
830
831 if( fExitFlux )
832 {
833 const G4RotationMatrix* rotM = pPostStepPoint->GetTouchable()->GetRotation();
834 G4ThreeVector transl = pPostStepPoint->GetTouchable()->GetTranslation();
835 G4AffineTransform transform = G4AffineTransform(rotM,transl);
836 transform.Invert();
837 G4ThreeVector localP = transform.TransformPoint(position);
838 G4ThreeVector localV = transform.TransformAxis(directionTR);
839
840 G4double distance = fEnvelope->GetSolid()->DistanceToOut(localP, localV);
841 if(verboseLevel > 1)
842 {
843 G4cout<<"distance to exit = "<<distance/mm<<" mm"<<G4endl;
844 }
845 position += distance*directionTR;
846 startTime += distance/c_light;
847 }
848 G4Track* aSecondaryTrack = new G4Track( aPhotonTR,
849 startTime, position );
850 aSecondaryTrack->SetTouchableHandle(
852 aSecondaryTrack->SetParentID( aTrack.GetTrackID() );
853
854 fParticleChange.AddSecondary(aSecondaryTrack);
855 fParticleChange.ProposeEnergy(kinEnergy);
856 }
857 }
858 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
859}
G4AffineTransform & Invert()
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
const G4ThreeVector & GetMomentumDirection() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
G4VSolid * GetSolid() const
G4Material * GetMaterial() const
void AddSecondary(G4Track *aSecondary)
void ProposeEnergy(G4double finalEnergy)
virtual void Initialize(const G4Track &)
const G4VTouchable * GetTouchable() const
G4double GetGlobalTime() const
const G4ThreeVector & GetPosition() const
const G4TouchableHandle & GetTouchableHandle() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
void SetParentID(const G4int aValue)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void SetNumberOfSecondaries(G4int totSecondaries)
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
virtual const G4ThreeVector & GetTranslation(G4int depth=0) const =0
virtual const G4RotationMatrix * GetRotation(G4int depth=0) const =0
G4double GetRandomAngle(G4double energyXTR, G4int iTkin)
G4double GetXTRrandomEnergy(G4double scaledTkin, G4int iTkin)

◆ SetAngleRadDistr()

void G4VXTRenergyLoss::SetAngleRadDistr ( G4bool  pAngleRadDistr)
inline

Definition at line 170 of file G4VXTRenergyLoss.hh.

170{fAngleRadDistr=pAngleRadDistr;};

◆ SetCompton()

void G4VXTRenergyLoss::SetCompton ( G4bool  pC)
inline

Definition at line 171 of file G4VXTRenergyLoss.hh.

171{fCompton=pC;};

◆ SetEnergy()

void G4VXTRenergyLoss::SetEnergy ( G4double  energy)
inline

Definition at line 168 of file G4VXTRenergyLoss.hh.

168{fEnergy = energy;};

◆ SetGamma()

void G4VXTRenergyLoss::SetGamma ( G4double  gamma)
inline

Definition at line 167 of file G4VXTRenergyLoss.hh.

167{fGamma = gamma;};

◆ SetVarAngle()

void G4VXTRenergyLoss::SetVarAngle ( G4double  varAngle)
inline

Definition at line 169 of file G4VXTRenergyLoss.hh.

169{fVarAngle = varAngle;};

◆ SpectralAngleXTRdEdx()

G4double G4VXTRenergyLoss::SpectralAngleXTRdEdx ( G4double  varAngle)

Definition at line 890 of file G4VXTRenergyLoss.cc.

891{
892 G4double result = GetStackFactor(fEnergy,fGamma,varAngle);
893 if(result < 0.0) result = 0.0;
894 return result;
895}

Referenced by BuildAngleForEnergyBank(), and SpectralXTRdEdx().

◆ SpectralXTRdEdx()

G4double G4VXTRenergyLoss::SpectralXTRdEdx ( G4double  energy)
virtual

Reimplemented in G4RegularXTRadiator, G4TransparentRegXTRadiator, G4XTRRegularRadModel, and G4XTRTransparentRegRadModel.

Definition at line 901 of file G4VXTRenergyLoss.cc.

902{
903 G4int i, iMax = 8;
904 G4double angleSum = 0.0;
905
906 G4double lim[8] = { 0.0, 0.01, 0.02, 0.05, 0.1, 0.2, 0.5, 1.0 };
907
908 for( i = 0; i < iMax; i++ ) lim[i] *= fMaxThetaTR;
909
911
912 fEnergy = energy;
913 /*
914 if( fAngleRadDistr && ( fEnergy == fEnergyForAngle ) )
915 {
916 fAngleVector ->PutValue(fBinTR - 1, angleSum);
917
918 for( i = fBinTR - 2; i >= 0; i-- )
919 {
920
921 angleSum += integral.Legendre10(
922 this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
923 fAngleVector->GetLowEdgeEnergy(i),
924 fAngleVector->GetLowEdgeEnergy(i+1) );
925
926 fAngleVector ->PutValue(i, angleSum);
927 }
928 }
929 else
930 */
931 {
932 for( i = 0; i < iMax-1; i++ )
933 {
934 angleSum += integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
935 lim[i],lim[i+1]);
936 // result += integral.Legendre10(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
937 // lim[i],lim[i+1]);
938 }
939 }
940 return angleSum;
941}

Referenced by BuildEnergyTable().

◆ XTRNAngleDensity()

G4double G4VXTRenergyLoss::XTRNAngleDensity ( G4double  varAngle)

Definition at line 1436 of file G4VXTRenergyLoss.cc.

1437{
1438 fVarAngle = varAngle;
1440 return integral.Legendre96(this,&G4VXTRenergyLoss::XTRNAngleSpectralDensity,
1442}
G4double XTRNAngleSpectralDensity(G4double energy)

◆ XTRNAngleSpectralDensity()

G4double G4VXTRenergyLoss::XTRNAngleSpectralDensity ( G4double  energy)

Definition at line 1426 of file G4VXTRenergyLoss.cc.

1427{
1430}
G4double OneBoundaryXTRNdensity(G4double energy, G4double gamma, G4double varAngle) const

Referenced by XTRNAngleDensity().

◆ XTRNSpectralAngleDensity()

G4double G4VXTRenergyLoss::XTRNSpectralAngleDensity ( G4double  varAngle)

Definition at line 1401 of file G4VXTRenergyLoss.cc.

1402{
1403 return OneBoundaryXTRNdensity(fEnergy,fGamma,varAngle)*
1404 GetStackFactor(fEnergy,fGamma,varAngle);
1405}

Referenced by XTRNSpectralDensity().

◆ XTRNSpectralDensity()

G4double G4VXTRenergyLoss::XTRNSpectralDensity ( G4double  energy)

Definition at line 1411 of file G4VXTRenergyLoss.cc.

1412{
1413 fEnergy = energy;
1415 return integral.Legendre96(this,&G4VXTRenergyLoss::XTRNSpectralAngleDensity,
1416 0.0,0.2*fMaxThetaTR) +
1417 integral.Legendre10(this,&G4VXTRenergyLoss::XTRNSpectralAngleDensity,
1419}
G4double XTRNSpectralAngleDensity(G4double varAngle)

Member Data Documentation

◆ fAlphaGas

◆ fAlphaPlate

◆ fAngleBank

std::vector<G4PhysicsTable*> G4VXTRenergyLoss::fAngleBank
protected

Definition at line 234 of file G4VXTRenergyLoss.hh.

Referenced by BuildAngleForEnergyBank(), BuildAngleTable(), and GetRandomAngle().

◆ fAngleDistrTable

G4PhysicsTable* G4VXTRenergyLoss::fAngleDistrTable
protected

Definition at line 185 of file G4VXTRenergyLoss.hh.

Referenced by BuildEnergyTable(), BuildGlobalAngleTable(), and ~G4VXTRenergyLoss().

◆ fAngleForEnergyTable

G4PhysicsTable* G4VXTRenergyLoss::fAngleForEnergyTable
protected

◆ fAngleRadDistr

G4bool G4VXTRenergyLoss::fAngleRadDistr
protected

◆ fBinTR

G4int G4VXTRenergyLoss::fBinTR
protected

◆ fCofTR

G4double G4VXTRenergyLoss::fCofTR
protected

Definition at line 209 of file G4VXTRenergyLoss.hh.

Referenced by BuildEnergyTable(), BuildGlobalAngleTable(), and G4VXTRenergyLoss().

◆ fCompton

G4bool G4VXTRenergyLoss::fCompton
protected

◆ fEnergy

◆ fEnergyDistrTable

G4PhysicsTable* G4VXTRenergyLoss::fEnergyDistrTable
protected

◆ fEnvelope

G4LogicalVolume* G4VXTRenergyLoss::fEnvelope
protected

Definition at line 184 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss(), GetMeanFreePath(), and PostStepDoIt().

◆ fExitFlux

◆ fGamma

◆ fGammaCutInKineticEnergy

G4double* G4VXTRenergyLoss::fGammaCutInKineticEnergy
protected

Definition at line 181 of file G4VXTRenergyLoss.hh.

◆ fGammaTkinCut

G4double G4VXTRenergyLoss::fGammaTkinCut
protected

◆ fGasPhotoAbsCof

G4SandiaTable* G4VXTRenergyLoss::fGasPhotoAbsCof
protected

Definition at line 229 of file G4VXTRenergyLoss.hh.

Referenced by ComputeGasPhotoAbsCof(), and GetGasLinearPhotoAbs().

◆ fGasThick

◆ fLambda

G4double G4VXTRenergyLoss::fLambda
protected

Definition at line 206 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss(), and GetMeanFreePath().

◆ fMatIndex1

G4int G4VXTRenergyLoss::fMatIndex1
protected

Definition at line 217 of file G4VXTRenergyLoss.hh.

Referenced by ComputePlatePhotoAbsCof(), G4VXTRenergyLoss(), and GetPlateCompton().

◆ fMatIndex2

G4int G4VXTRenergyLoss::fMatIndex2
protected

Definition at line 218 of file G4VXTRenergyLoss.hh.

Referenced by ComputeGasPhotoAbsCof(), G4VXTRenergyLoss(), and GetGasCompton().

◆ fMaxEnergyTR

G4double G4VXTRenergyLoss::fMaxEnergyTR
protected

◆ fMaxProtonTkin

G4double G4VXTRenergyLoss::fMaxProtonTkin
protected

Definition at line 201 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss().

◆ fMaxThetaTR

◆ fMinEnergyTR

G4double G4VXTRenergyLoss::fMinEnergyTR
protected

◆ fMinProtonTkin

G4double G4VXTRenergyLoss::fMinProtonTkin
protected

Definition at line 200 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss().

◆ fParticleChange

G4ParticleChange G4VXTRenergyLoss::fParticleChange
protected

Definition at line 231 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss(), and PostStepDoIt().

◆ fPlasmaCof

G4double G4VXTRenergyLoss::fPlasmaCof
protected

◆ fPlateNumber

◆ fPlatePhotoAbsCof

G4SandiaTable* G4VXTRenergyLoss::fPlatePhotoAbsCof
protected

Definition at line 227 of file G4VXTRenergyLoss.hh.

Referenced by ComputePlatePhotoAbsCof(), and GetPlateLinearPhotoAbs().

◆ fPlateThick

◆ fProtonEnergyVector

◆ fPtrGamma

G4ParticleDefinition* G4VXTRenergyLoss::fPtrGamma
protected

Definition at line 179 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss().

◆ fSigma1

◆ fSigma2

◆ fTheMaxAngle

G4double G4VXTRenergyLoss::fTheMaxAngle
protected

◆ fTheMaxEnergyTR

G4double G4VXTRenergyLoss::fTheMaxEnergyTR
protected

◆ fTheMinAngle

G4double G4VXTRenergyLoss::fTheMinAngle
protected

◆ fTheMinEnergyTR

G4double G4VXTRenergyLoss::fTheMinEnergyTR
protected

◆ fTotalDist

G4double G4VXTRenergyLoss::fTotalDist
protected

Definition at line 221 of file G4VXTRenergyLoss.hh.

Referenced by BuildEnergyTable(), and G4VXTRenergyLoss().

◆ fTotBin

◆ fVarAngle

G4double G4VXTRenergyLoss::fVarAngle
protected

◆ fXTREnergyVector

G4PhysicsLogVector* G4VXTRenergyLoss::fXTREnergyVector
protected

Definition at line 189 of file G4VXTRenergyLoss.hh.

Referenced by BuildAngleTable(), G4VXTRenergyLoss(), and ~G4VXTRenergyLoss().


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