Geant4 11.2.2
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 void ProcessDescription (std::ostream &) const override
 
virtual void DumpInfo () const override
 
 G4VXTRenergyLoss (G4VXTRenergyLoss &)=delete
 
G4VXTRenergyLossoperator= (const G4VXTRenergyLoss &right)=delete
 
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)
 
void SetGamma (G4double gamma)
 
G4double GetGamma ()
 
void SetEnergy (G4double energy)
 
G4double GetEnergy ()
 
void SetVarAngle (G4double varAngle)
 
G4double GetVarAngle ()
 
void SetCompton (G4bool pC)
 
G4bool GetCompton ()
 
G4int GetKrange ()
 
void SetKrange (G4int kk)
 
void SetAlphaGas (G4double ag)
 
G4double GetAlphaGas ()
 
void SetAlphaPlate (G4double ap)
 
G4double GetAlphaPlate ()
 
void SetTheMinEnergyTR (G4double minetr)
 
G4double GetTheMinEnergyTR ()
 
void SetTheMaxEnergyTR (G4double maxetr)
 
G4double GetTheMaxEnergyTR ()
 
void SetMinEnergyTR (G4double minetr)
 
G4double GetMinEnergyTR ()
 
void SetMaxEnergyTR (G4double maxetr)
 
G4double GetMaxEnergyTR ()
 
void SetTheMinAngle (G4double minang)
 
G4double GetTheMinAngle ()
 
void SetTheMaxAngle (G4double maxang)
 
G4double GetTheMaxAngle ()
 
void SetMinThetaTR (G4double minatr)
 
G4double GetMinThetaTR ()
 
void SetMaxThetaTR (G4double maxatr)
 
G4double GetMaxThetaTR ()
 
void SetFastAngle (G4bool fatr)
 
G4bool GetFastAngle ()
 
void SetAngleRadDistr (G4bool fatr)
 
G4bool GetAngleRadDistr ()
 
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 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 &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- 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
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void 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 const G4VProcessGetCreatorProcess () const
 
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
 
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

G4double fTheMinEnergyTR
 
G4double fTheMaxEnergyTR
 
G4double fTheMinAngle
 
G4double fTheMaxAngle
 
G4int fTotBin
 
G4int fBinTR
 
G4int fKrange
 
G4ParticleDefinitionfPtrGamma
 
G4doublefGammaCutInKineticEnergy
 
G4LogicalVolumefEnvelope
 
G4PhysicsTablefAngleDistrTable
 
G4PhysicsTablefEnergyDistrTable
 
G4PhysicsTablefAngleForEnergyTable
 
G4PhysicsLogVectorfProtonEnergyVector
 
G4PhysicsLogVectorfXTREnergyVector
 
G4SandiaTablefPlatePhotoAbsCof
 
G4SandiaTablefGasPhotoAbsCof
 
G4ParticleChange fParticleChange
 
std::vector< G4PhysicsTable * > fAngleBank
 
G4double fGammaTkinCut
 
G4double fMinEnergyTR
 
G4double fMaxEnergyTR
 
G4double fMinThetaTR
 
G4double fMaxThetaTR
 
G4double fTotalDist
 
G4double fPlateThick
 
G4double fGasThick
 
G4double fAlphaPlate
 
G4double fAlphaGas
 
G4double fGamma
 
G4double fEnergy
 
G4double fVarAngle
 
G4double fLambda
 
G4double fSigma1
 
G4double fSigma2
 
G4int fMatIndex1
 
G4int fMatIndex2
 
G4int fPlateNumber
 
G4bool fExitFlux
 
G4bool fFastAngle
 
G4bool fAngleRadDistr
 
G4bool fCompton
 
G4int secID = -1
 
- 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 constexpr G4double fMinProtonTkin = 100. * CLHEP::GeV
 
static constexpr G4double fMaxProtonTkin = 100. * CLHEP::TeV
 
static constexpr G4double fPlasmaCof
 
static constexpr G4double fCofTR = CLHEP::fine_structure_const / CLHEP::pi
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VDiscreteProcess
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 58 of file G4VXTRenergyLoss.hh.

Constructor & Destructor Documentation

◆ G4VXTRenergyLoss() [1/2]

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 59 of file G4VXTRenergyLoss.cc.

64 : G4VDiscreteProcess(processName, type)
66 , fAngleDistrTable(nullptr)
67 , fEnergyDistrTable(nullptr)
68 , fAngleForEnergyTable(nullptr)
69 , fPlatePhotoAbsCof(nullptr)
70 , fGasPhotoAbsCof(nullptr)
71 , fGammaTkinCut(0.0)
72{
73 verboseLevel = 1;
74 secID = G4PhysicsModelCatalog::GetModelID("model_XTRenergyLoss");
76
77 fPtrGamma = nullptr;
80 fAlphaPlate = 100.;
81 fAlphaGas = 40.;
82
83 fTheMinEnergyTR = CLHEP::keV * 1.; // 1.; //
84 fTheMaxEnergyTR = CLHEP::keV * 100.; // 40.; //
85
86 fTheMinAngle = 1.e-8; //
87 fTheMaxAngle = 4.e-4;
88
89 fTotBin = 50; // number of bins in log scale
90 fBinTR = 100; // number of bins in TR vectors
91 fKrange = 229;
92 // min/max angle2 in log-vectors
93
94 fMinThetaTR = 3.0e-9;
95 fMaxThetaTR = 1.0e-4;
96
97
98 // Proton energy vector initialization
101
104
105 fEnvelope = anEnvelope;
106
107 fPlateNumber = n;
108 if(verboseLevel > 0)
109 G4cout << "### G4VXTRenergyLoss: the number of TR radiator plates = "
110 << fPlateNumber << G4endl;
111 if(fPlateNumber == 0)
112 {
113 G4Exception("G4VXTRenergyLoss::G4VXTRenergyLoss()", "VXTRELoss01",
114 FatalException, "No plates in X-ray TR radiator");
115 }
116 // default is XTR dEdx, not flux after radiator
117 fExitFlux = false;
118 // default angle distribution according numerical integration
119 fFastAngle = false; // no angle according sum of delta-functions by default
120 fAngleRadDistr = true;
121 fCompton = false;
122
124
125 // Mean thicknesses of plates and gas gaps
126 fPlateThick = a;
127 fGasThick = b;
129 if(verboseLevel > 0)
130 G4cout << "total radiator thickness = " << fTotalDist / cm << " cm"
131 << G4endl;
132
133 // index of plate material
134 fMatIndex1 = (G4int)foilMat->GetIndex();
135 if(verboseLevel > 0)
136 G4cout << "plate material = " << foilMat->GetName() << G4endl;
137
138 // index of gas material
139 fMatIndex2 = (G4int)gasMat->GetIndex();
140 if(verboseLevel > 0)
141 G4cout << "gas material = " << gasMat->GetName() << G4endl;
142
143 // plasma energy squared for plate material
145 if(verboseLevel > 0)
146 G4cout << "plate plasma energy = " << std::sqrt(fSigma1) / eV << " eV"
147 << G4endl;
148
149 // plasma energy squared for gas material
151 if(verboseLevel > 0)
152 G4cout << "gas plasma energy = " << std::sqrt(fSigma2) / eV << " eV"
153 << G4endl;
154
155 // Compute cofs for preparation of linear photo absorption
158
160}
@ fTransitionRadiation
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetElectronDensity() const
std::size_t GetIndex() const
const G4String & GetName() const
static G4int GetModelID(const G4int modelIndex)
G4int verboseLevel
void SetProcessSubType(G4int)
G4VParticleChange * pParticleChange
G4PhysicsTable * fEnergyDistrTable
G4PhysicsTable * fAngleForEnergyTable
static constexpr G4double fMaxProtonTkin
G4PhysicsLogVector * fProtonEnergyVector
G4PhysicsLogVector * fXTREnergyVector
G4SandiaTable * fPlatePhotoAbsCof
G4PhysicsTable * fAngleDistrTable
G4SandiaTable * fGasPhotoAbsCof
G4LogicalVolume * fEnvelope
static constexpr G4double fPlasmaCof
static constexpr G4double fMinProtonTkin
G4ParticleDefinition * fPtrGamma
G4ParticleChange fParticleChange
G4double * fGammaCutInKineticEnergy
#define DBL_MAX
Definition templates.hh:62

Referenced by BuildAngleForEnergyBank(), BuildEnergyTable(), BuildGlobalAngleTable(), SpectralXTRdEdx(), XTRNAngleDensity(), and XTRNSpectralDensity().

◆ ~G4VXTRenergyLoss()

G4VXTRenergyLoss::~G4VXTRenergyLoss ( )
virtual

Definition at line 163 of file G4VXTRenergyLoss.cc.

164{
165 delete fProtonEnergyVector;
166 delete fXTREnergyVector;
168 {
170 delete fEnergyDistrTable;
171 }
173 {
175 delete fAngleDistrTable;
176 }
178 {
181 }
182}
void clearAndDestroy()

◆ G4VXTRenergyLoss() [2/2]

G4VXTRenergyLoss::G4VXTRenergyLoss ( G4VXTRenergyLoss & )
delete

Member Function Documentation

◆ AngleSpectralXTRdEdx()

G4double G4VXTRenergyLoss::AngleSpectralXTRdEdx ( G4double energy)

Definition at line 901 of file G4VXTRenergyLoss.cc.

902{
903 G4double result = GetStackFactor(energy, fGamma, fVarAngle);
904 if(result < 0)
905 result = 0.0;
906 return result;
907}
double G4double
Definition G4Types.hh:83
virtual G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle)

◆ AngleXTRdEdx()

G4double G4VXTRenergyLoss::AngleXTRdEdx ( G4double varAngle)

Definition at line 911 of file G4VXTRenergyLoss.cc.

912{
913 G4double result;
914 G4double sum = 0., tmp1, tmp2, tmp = 0., cof1, cof2, cofMin, cofPHC, energy1,
915 energy2;
916 G4int k, kMax, kMin, i;
917
918 cofPHC = twopi * hbarc;
919
920 cof1 = (fPlateThick + fGasThick) * (1. / fGamma / fGamma + varAngle);
922
923 cofMin = std::sqrt(cof1 * cof2);
924 cofMin /= cofPHC;
925
926 kMin = G4int(cofMin);
927 if(cofMin > kMin)
928 kMin++;
929
930 kMax = kMin + 9;
931
932 for(k = kMin; k <= kMax; ++k)
933 {
934 tmp1 = cofPHC * k;
935 tmp2 = std::sqrt(tmp1 * tmp1 - cof1 * cof2);
936 energy1 = (tmp1 + tmp2) / cof1;
937 energy2 = (tmp1 - tmp2) / cof1;
938
939 for(i = 0; i < 2; ++i)
940 {
941 if(i == 0)
942 {
943 if(energy1 > fTheMaxEnergyTR || energy1 < fTheMinEnergyTR)
944 continue;
945
946 tmp1 =
947 (energy1 * energy1 * (1. / fGamma / fGamma + varAngle) + fSigma1) *
948 fPlateThick / (4 * hbarc * energy1);
949 tmp2 = std::sin(tmp1);
950 tmp = energy1 * tmp2 * tmp2;
951 tmp2 = fPlateThick / (4. * tmp1);
952 tmp1 =
953 hbarc * energy1 /
954 (energy1 * energy1 * (1. / fGamma / fGamma + varAngle) + fSigma2);
955 tmp *= (tmp1 - tmp2) * (tmp1 - tmp2);
956 tmp1 = cof1 / (4. * hbarc) - cof2 / (4. * hbarc * energy1 * energy1);
957 tmp2 = std::abs(tmp1);
958
959 if(tmp2 > 0.)
960 tmp /= tmp2;
961 else
962 continue;
963 }
964 else
965 {
966 if(energy2 > fTheMaxEnergyTR || energy2 < fTheMinEnergyTR)
967 continue;
968
969 tmp1 =
970 (energy2 * energy2 * (1. / fGamma / fGamma + varAngle) + fSigma1) *
971 fPlateThick / (4. * hbarc * energy2);
972 tmp2 = std::sin(tmp1);
973 tmp = energy2 * tmp2 * tmp2;
974 tmp2 = fPlateThick / (4. * tmp1);
975 tmp1 =
976 hbarc * energy2 /
977 (energy2 * energy2 * (1. / fGamma / fGamma + varAngle) + fSigma2);
978 tmp *= (tmp1 - tmp2) * (tmp1 - tmp2);
979 tmp1 = cof1 / (4. * hbarc) - cof2 / (4. * hbarc * energy2 * energy2);
980 tmp2 = std::abs(tmp1);
981
982 if(tmp2 > 0.)
983 tmp /= tmp2;
984 else
985 continue;
986 }
987 sum += tmp;
988 }
989 }
990 result = 4. * pi * fPlateNumber * sum * varAngle;
991 result /= hbarc * hbarc;
992
993 return result;
994}

Referenced by BuildGlobalAngleTable().

◆ BuildAngleForEnergyBank()

void G4VXTRenergyLoss::BuildAngleForEnergyBank ( )

Definition at line 390 of file G4VXTRenergyLoss.cc.

391{
392
393 if( ( this->GetProcessName() == "TranspRegXTRadiator" ||
394 this->GetProcessName() == "TranspRegXTRmodel" ||
395 this->GetProcessName() == "RegularXTRadiator" ||
396 this->GetProcessName() == "RegularXTRmodel" ) && fFastAngle ) // ffastAngle=true!
397 {
398 BuildAngleTable(); // by sum of delta-functions
399 return;
400 }
401 G4int i, iTkin, iTR;
402 G4double angleSum = 0.0;
403
404 fGammaTkinCut = 0.0;
405
406 // setting of min/max TR energies
409 else
411
414 else
416
417 auto energyVector =
419
421 integral;
422
423 G4cout.precision(4);
424 G4Timer timer;
425 timer.Start();
426
427 for(iTkin = 0; iTkin < fTotBin; ++iTkin) // Lorentz factor loop
428 {
429 fGamma =
430 1.0 + (fProtonEnergyVector->GetLowEdgeEnergy(iTkin) / proton_mass_c2);
431
434 else if(fMaxThetaTR < fTheMinAngle)
436
438
439 for(iTR = 0; iTR < fBinTR; ++iTR)
440 {
441 angleSum = 0.0;
442 fEnergy = energyVector->GetLowEdgeEnergy(iTR);
443
444 // log-vector to increase number of thin bins for small angles
445 auto angleVector = new G4PhysicsLogVector(fMinThetaTR, fMaxThetaTR, fBinTR);
446
447
448
449 angleVector->PutValue(fBinTR - 1, angleSum);
450
451 for(i = fBinTR - 2; i >= 0; --i)
452 {
453 // Legendre96 or Legendre10
454
455 angleSum +=
456 integral.Legendre10(this, &G4VXTRenergyLoss::SpectralAngleXTRdEdx,
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(const std::size_t index) const
void Stop()
G4double GetUserElapsed() const
Definition G4Timer.cc:135
void Start()
const G4String & GetProcessName() const
G4VXTRenergyLoss(G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="XTRenergyLoss", G4ProcessType type=fElectromagnetic)
G4double SpectralAngleXTRdEdx(G4double varAngle)
std::vector< G4PhysicsTable * > fAngleBank

Referenced by BuildPhysicsTable().

◆ BuildAngleTable()

void G4VXTRenergyLoss::BuildAngleTable ( )

Definition at line 481 of file G4VXTRenergyLoss.cc.

482{
483 G4int iTkin, iTR;
485
486 fGammaTkinCut = 0.0;
487
488 // setting of min/max TR energies
491 else
493
496 else
498
499 G4cout.precision(4);
500 G4Timer timer;
501 timer.Start();
502 if(verboseLevel > 0)
503 {
504 G4cout << G4endl << "Lorentz Factor" << "\t"
505 << "XTR photon number" << G4endl << G4endl;
506 }
507 for(iTkin = 0; iTkin < fTotBin; ++iTkin) // Lorentz factor loop
508 {
509 fGamma =
510 1.0 + (fProtonEnergyVector->GetLowEdgeEnergy(iTkin) / proton_mass_c2);
511
512 // fMaxThetaTR = 25. * 2500.0 / (fGamma * fGamma); // theta^2
513
516 else
517 {
520 }
521
523
524 for(iTR = 0; iTR < fBinTR; ++iTR)
525 {
527
528 G4PhysicsFreeVector* angleVector = GetAngleVector(energy, fBinTR);
529
530 fAngleForEnergyTable->insertAt(iTR, angleVector);
531 }
533 }
534 timer.Stop();
535 G4cout.precision(6);
536 if(verboseLevel > 0)
537 {
538 G4cout << G4endl;
539 G4cout << "total time for build XTR angle for given energy tables = "
540 << timer.GetUserElapsed() << " s" << G4endl;
541 }
542 fGamma = 0.;
543
544 return;
545}
G4PhysicsFreeVector * GetAngleVector(G4double energy, G4int n)
G4double energy(const ThreeVector &p, const G4double m)

Referenced by BuildAngleForEnergyBank().

◆ BuildEnergyTable()

void G4VXTRenergyLoss::BuildEnergyTable ( )

Definition at line 299 of file G4VXTRenergyLoss.cc.

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

Referenced by BuildPhysicsTable().

◆ BuildGlobalAngleTable()

void G4VXTRenergyLoss::BuildGlobalAngleTable ( )

Definition at line 624 of file G4VXTRenergyLoss.cc.

625{
626 G4int iTkin, iTR, iPlace;
627 G4double radiatorCof = 1.0; // for tuning of XTR yield
628 G4double angleSum;
630
631 fGammaTkinCut = 0.0;
632
633 // setting of min/max TR energies
636 else
638
641 else
643
644 G4cout.precision(4);
645 G4Timer timer;
646 timer.Start();
647 if(verboseLevel > 0)
648 {
649 G4cout << G4endl;
650 G4cout << "Lorentz Factor"
651 << "\t"
652 << "XTR photon number" << G4endl;
653 G4cout << G4endl;
654 }
655 for(iTkin = 0; iTkin < fTotBin; ++iTkin) // Lorentz factor loop
656 {
657 fGamma =
658 1.0 + (fProtonEnergyVector->GetLowEdgeEnergy(iTkin) / proton_mass_c2);
659
660 // fMaxThetaTR = 25.0 / (fGamma * fGamma); // theta^2
661 // fMaxThetaTR = 1.e-4; // theta^2
662
665 else
666 {
669 }
670 auto angleVector =
671 // G4PhysicsLogVector* angleVector =
673 // new G4PhysicsLogVector(1.e-8, fMaxThetaTR, fBinTR);
674
675 angleSum = 0.0;
676
678 integral;
679
680 angleVector->PutValue(fBinTR - 1, angleSum);
681
682 for(iTR = fBinTR - 2; iTR >= 0; --iTR)
683 {
684 angleSum += radiatorCof * fCofTR *
685 integral.Legendre96(this, &G4VXTRenergyLoss::AngleXTRdEdx,
686 angleVector->GetLowEdgeEnergy(iTR),
687 angleVector->GetLowEdgeEnergy(iTR + 1));
688
689 angleVector->PutValue(iTR, angleSum);
690 }
691 if(verboseLevel > 1)
692 {
693 G4cout << fGamma << "\t" << angleSum << G4endl;
694 }
695 iPlace = iTkin;
696 fAngleDistrTable->insertAt(iPlace, angleVector);
697 }
698 timer.Stop();
699 G4cout.precision(6);
700 if(verboseLevel > 0)
701 {
702 G4cout << G4endl;
703 G4cout << "total time for build X-ray TR angle tables = "
704 << timer.GetUserElapsed() << " s" << G4endl;
705 }
706 fGamma = 0.;
707
708 return;
709}
G4double AngleXTRdEdx(G4double varAngle)

◆ BuildPhysicsTable()

void G4VXTRenergyLoss::BuildPhysicsTable ( const G4ParticleDefinition & pd)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 276 of file G4VXTRenergyLoss.cc.

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

◆ BuildTable()

void G4VXTRenergyLoss::BuildTable ( )
inline

Definition at line 90 of file G4VXTRenergyLoss.hh.

90{};

◆ ComputeGasPhotoAbsCof()

void G4VXTRenergyLoss::ComputeGasPhotoAbsCof ( )

Definition at line 1086 of file G4VXTRenergyLoss.cc.

1087{
1088 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1089 const G4Material* mat = (*theMaterialTable)[fMatIndex2];
1091 return;
1092}
std::vector< G4Material * > G4MaterialTable
G4SandiaTable * GetSandiaTable() const
static G4MaterialTable * GetMaterialTable()

Referenced by G4VXTRenergyLoss().

◆ ComputePlatePhotoAbsCof()

void G4VXTRenergyLoss::ComputePlatePhotoAbsCof ( )

Definition at line 1028 of file G4VXTRenergyLoss.cc.

1029{
1030 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1031 const G4Material* mat = (*theMaterialTable)[fMatIndex1];
1033
1034 return;
1035}

Referenced by G4VXTRenergyLoss().

◆ DumpInfo()

virtual void G4VXTRenergyLoss::DumpInfo ( ) const
inlineoverridevirtual

◆ GetAlphaGas()

G4double G4VXTRenergyLoss::GetAlphaGas ( )
inline

Definition at line 159 of file G4VXTRenergyLoss.hh.

159{ return fAlphaGas; };

◆ GetAlphaPlate()

G4double G4VXTRenergyLoss::GetAlphaPlate ( )
inline

Definition at line 161 of file G4VXTRenergyLoss.hh.

161{ return fAlphaPlate; };

◆ GetAngleRadDistr()

G4bool G4VXTRenergyLoss::GetAngleRadDistr ( )
inline

Definition at line 188 of file G4VXTRenergyLoss.hh.

188{ return fAngleRadDistr; };

◆ GetAngleVector()

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

Definition at line 549 of file G4VXTRenergyLoss.cc.

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

Referenced by BuildAngleTable().

◆ GetAngleXTR()

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

Definition at line 1492 of file G4VXTRenergyLoss.cc.

1494{
1495 G4double x1, x2, y1, y2, result;
1496
1497 if( iTransfer == 0 )
1498 // if( iTransfer == 1 ) // ATLAS TB
1499 {
1500 result = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1501 }
1502 else
1503 {
1504 y1 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer - 1);
1505 y2 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer);
1506
1507 x1 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer - 1);
1508 x2 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1509
1510 if(x1 == x2) result = x2;
1511 else
1512 {
1513 if( y1 == y2 ) result = x1 + (x2 - x1) * G4UniformRand();
1514 else
1515 {
1516 result = x1 + (position - y1) * (x2 - x1) / (y2 - y1);
1517 // result = x1 + 0.1*(position - y1) * (x2 - x1) / (y2 - y1); // ATLAS TB
1518 // result = x1 + 0.05*(position - y1) * (x2 - x1) / (y2 - y1); // ATLAS TB
1519 }
1520 }
1521 }
1522 return result;
1523}
#define G4UniformRand()
Definition Randomize.hh:52

Referenced by GetRandomAngle().

◆ GetCompton()

G4bool G4VXTRenergyLoss::GetCompton ( )
inline

Definition at line 152 of file G4VXTRenergyLoss.hh.

152{ return fCompton; };

◆ GetComptonPerAtom()

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

Definition at line 1228 of file G4VXTRenergyLoss.cc.

1229{
1230 G4double CrossSection = 0.0;
1231 if(Z < 0.9999)
1232 return CrossSection;
1233 if(GammaEnergy < 0.1 * keV)
1234 return CrossSection;
1235 if(GammaEnergy > (100. * GeV / Z))
1236 return CrossSection;
1237
1238 static constexpr G4double a = 20.0;
1239 static constexpr G4double b = 230.0;
1240 static constexpr G4double c = 440.0;
1241
1242 static constexpr G4double d1 = 2.7965e-1 * barn, d2 = -1.8300e-1 * barn,
1243 d3 = 6.7527 * barn, d4 = -1.9798e+1 * barn,
1244 e1 = 1.9756e-5 * barn, e2 = -1.0205e-2 * barn,
1245 e3 = -7.3913e-2 * barn, e4 = 2.7079e-2 * barn,
1246 f1 = -3.9178e-7 * barn, f2 = 6.8241e-5 * barn,
1247 f3 = 6.0480e-5 * barn, f4 = 3.0274e-4 * barn;
1248
1249 G4double p1Z = Z * (d1 + e1 * Z + f1 * Z * Z);
1250 G4double p2Z = Z * (d2 + e2 * Z + f2 * Z * Z);
1251 G4double p3Z = Z * (d3 + e3 * Z + f3 * Z * Z);
1252 G4double p4Z = Z * (d4 + e4 * Z + f4 * Z * Z);
1253
1254 G4double T0 = 15.0 * keV;
1255 if(Z < 1.5)
1256 T0 = 40.0 * keV;
1257
1258 G4double X = std::max(GammaEnergy, T0) / electron_mass_c2;
1259 CrossSection =
1260 p1Z * std::log(1. + 2. * X) / X +
1261 (p2Z + p3Z * X + p4Z * X * X) / (1. + a * X + b * X * X + c * X * X * X);
1262
1263 // modification for low energy. (special case for Hydrogen)
1264 if(GammaEnergy < T0)
1265 {
1266 G4double dT0 = 1. * keV;
1267 X = (T0 + dT0) / electron_mass_c2;
1268 G4double sigma =
1269 p1Z * std::log(1. + 2. * X) / X +
1270 (p2Z + p3Z * X + p4Z * X * X) / (1. + a * X + b * X * X + c * X * X * X);
1271 G4double c1 = -T0 * (sigma - CrossSection) / (CrossSection * dT0);
1272 G4double c2 = 0.150;
1273 if(Z > 1.5)
1274 c2 = 0.375 - 0.0556 * std::log(Z);
1275 G4double y = std::log(GammaEnergy / T0);
1276 CrossSection *= std::exp(-y * (c1 + c2 * y));
1277 }
1278 return CrossSection;
1279}

Referenced by GetGasCompton(), and GetPlateCompton().

◆ GetEnergy()

G4double G4VXTRenergyLoss::GetEnergy ( )
inline

Definition at line 148 of file G4VXTRenergyLoss.hh.

148{ return fEnergy; };

◆ GetFastAngle()

G4bool G4VXTRenergyLoss::GetFastAngle ( )
inline

Definition at line 186 of file G4VXTRenergyLoss.hh.

186{ return fFastAngle; };

◆ GetGamma()

G4double G4VXTRenergyLoss::GetGamma ( )
inline

Definition at line 146 of file G4VXTRenergyLoss.hh.

146{ return fGamma; };

◆ GetGasComplexFZ()

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

Definition at line 1067 of file G4VXTRenergyLoss.cc.

1069{
1070 G4double cof, length, delta, real_v, image_v;
1071
1072 length = 0.5 * GetGasFormationZone(omega, gamma, varAngle);
1073 delta = length * GetGasLinearPhotoAbs(omega);
1074 cof = 1.0 / (1.0 + delta * delta);
1075
1076 real_v = length * cof;
1077 image_v = real_v * delta;
1078
1079 G4complex zone(real_v, image_v);
1080 return zone;
1081}
std::complex< G4double > G4complex
Definition G4Types.hh:88
G4double GetGasFormationZone(G4double, G4double, G4double)
G4double GetGasLinearPhotoAbs(G4double)

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

◆ GetGasCompton()

G4double G4VXTRenergyLoss::GetGasCompton ( G4double omega)

Definition at line 1206 of file G4VXTRenergyLoss.cc.

1207{
1208 G4int i, numberOfElements;
1209 G4double xSection = 0., nowZ, sumZ = 0.;
1210
1211 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1212 numberOfElements = (G4int)(*theMaterialTable)[fMatIndex2]->GetNumberOfElements();
1213
1214 for(i = 0; i < numberOfElements; ++i)
1215 {
1216 nowZ = (*theMaterialTable)[fMatIndex2]->GetElement(i)->GetZ();
1217 sumZ += nowZ;
1218 xSection += GetComptonPerAtom(omega, nowZ);
1219 }
1220 xSection /= sumZ;
1221 xSection *= (*theMaterialTable)[fMatIndex2]->GetElectronDensity();
1222 return xSection;
1223}
G4double GetComptonPerAtom(G4double, G4double)

Referenced by G4XTRTransparentRegRadModel::SpectralXTRdEdx().

◆ GetGasFormationZone()

◆ GetGasLinearPhotoAbs()

G4double G4VXTRenergyLoss::GetGasLinearPhotoAbs ( G4double omega)

Definition at line 1097 of file G4VXTRenergyLoss.cc.

1098{
1099 G4double omega2, omega3, omega4;
1100
1101 omega2 = omega * omega;
1102 omega3 = omega2 * omega;
1103 omega4 = omega2 * omega2;
1104
1105 const G4double* SandiaCof = fGasPhotoAbsCof->GetSandiaCofForMaterial(omega);
1106 G4double cross = SandiaCof[0] / omega + SandiaCof[1] / omega2 +
1107 SandiaCof[2] / omega3 + SandiaCof[3] / omega4;
1108 return cross;
1109}
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 1160 of file G4VXTRenergyLoss.cc.

1161{
1162 std::ofstream outGas("gasZmu.dat", std::ios::out);
1163 outGas.setf(std::ios::scientific, std::ios::floatfield);
1164 G4int i;
1165 G4double omega, varAngle, gamma;
1166 gamma = 10000.;
1167 varAngle = 1 / gamma / gamma;
1168 if(verboseLevel > 0)
1169 G4cout << "energy, keV" << "\t" << "Zmu for gas" << G4endl;
1170 for(i = 0; i < 100; ++i)
1171 {
1172 omega = (1.0 + i) * keV;
1173 if(verboseLevel > 1)
1174 G4cout << omega / keV << "\t" << GetGasZmuProduct(omega, gamma, varAngle)
1175 << "\t";
1176 if(verboseLevel > 0)
1177 outGas << omega / keV << "\t\t"
1178 << GetGasZmuProduct(omega, gamma, varAngle) << G4endl;
1179 }
1180 return;
1181}

Referenced by GetGasZmuProduct().

◆ GetGasZmuProduct() [2/2]

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

Definition at line 1150 of file G4VXTRenergyLoss.cc.

1152{
1153 return GetGasFormationZone(omega, gamma, varAngle) *
1154 GetGasLinearPhotoAbs(omega);
1155}

◆ GetKrange()

G4int G4VXTRenergyLoss::GetKrange ( )
inline

Definition at line 154 of file G4VXTRenergyLoss.hh.

154{ return fKrange;};

◆ GetMaxEnergyTR()

G4double G4VXTRenergyLoss::GetMaxEnergyTR ( )
inline

Definition at line 171 of file G4VXTRenergyLoss.hh.

171{ return fMaxEnergyTR; };

◆ GetMaxThetaTR()

G4double G4VXTRenergyLoss::GetMaxThetaTR ( )
inline

Definition at line 181 of file G4VXTRenergyLoss.hh.

181{ return fMaxThetaTR; };

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 200 of file G4VXTRenergyLoss.cc.

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

◆ GetMinEnergyTR()

G4double G4VXTRenergyLoss::GetMinEnergyTR ( )
inline

Definition at line 169 of file G4VXTRenergyLoss.hh.

169{ return fMinEnergyTR; };

◆ GetMinThetaTR()

G4double G4VXTRenergyLoss::GetMinThetaTR ( )
inline

Definition at line 179 of file G4VXTRenergyLoss.hh.

179{ return fMinThetaTR; };

◆ GetNumberOfPhotons()

void G4VXTRenergyLoss::GetNumberOfPhotons ( )

Definition at line 1350 of file G4VXTRenergyLoss.cc.

1351{
1352 G4int iTkin;
1353 G4double gamma, numberE;
1354
1355 std::ofstream outEn("numberE.dat", std::ios::out);
1356 outEn.setf(std::ios::scientific, std::ios::floatfield);
1357
1358 std::ofstream outAng("numberAng.dat", std::ios::out);
1359 outAng.setf(std::ios::scientific, std::ios::floatfield);
1360
1361 for(iTkin = 0; iTkin < fTotBin; ++iTkin) // Lorentz factor loop
1362 {
1363 gamma =
1364 1.0 + (fProtonEnergyVector->GetLowEdgeEnergy(iTkin) / proton_mass_c2);
1365 numberE = (*(*fEnergyDistrTable)(iTkin))(0);
1366 if(verboseLevel > 1)
1367 G4cout << gamma << "\t\t" << numberE << "\t" << G4endl;
1368 if(verboseLevel > 0)
1369 outEn << gamma << "\t\t" << numberE << G4endl;
1370 }
1371 return;
1372}

◆ GetPlateComplexFZ()

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

Definition at line 1009 of file G4VXTRenergyLoss.cc.

1011{
1012 G4double cof, length, delta, real_v, image_v;
1013
1014 length = 0.5 * GetPlateFormationZone(omega, gamma, varAngle);
1015 delta = length * GetPlateLinearPhotoAbs(omega);
1016 cof = 1.0 / (1.0 + delta * delta);
1017
1018 real_v = length * cof;
1019 image_v = real_v * delta;
1020
1021 G4complex zone(real_v, image_v);
1022 return zone;
1023}
G4double GetPlateLinearPhotoAbs(G4double)
G4double GetPlateFormationZone(G4double, G4double, G4double)

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

◆ GetPlateCompton()

G4double G4VXTRenergyLoss::GetPlateCompton ( G4double omega)

Definition at line 1185 of file G4VXTRenergyLoss.cc.

1186{
1187 G4int i, numberOfElements;
1188 G4double xSection = 0., nowZ, sumZ = 0.;
1189
1190 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1191 numberOfElements = (G4int)(*theMaterialTable)[fMatIndex1]->GetNumberOfElements();
1192
1193 for(i = 0; i < numberOfElements; ++i)
1194 {
1195 nowZ = (*theMaterialTable)[fMatIndex1]->GetElement(i)->GetZ();
1196 sumZ += nowZ;
1197 xSection += GetComptonPerAtom(omega, nowZ);
1198 }
1199 xSection /= sumZ;
1200 xSection *= (*theMaterialTable)[fMatIndex1]->GetElectronDensity();
1201 return xSection;
1202}

Referenced by G4XTRTransparentRegRadModel::SpectralXTRdEdx().

◆ GetPlateFormationZone()

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

◆ GetPlateLinearPhotoAbs()

G4double G4VXTRenergyLoss::GetPlateLinearPhotoAbs ( G4double omega)

Definition at line 1040 of file G4VXTRenergyLoss.cc.

1041{
1042 G4double omega2, omega3, omega4;
1043
1044 omega2 = omega * omega;
1045 omega3 = omega2 * omega;
1046 omega4 = omega2 * omega2;
1047
1048 const G4double* SandiaCof = fPlatePhotoAbsCof->GetSandiaCofForMaterial(omega);
1049 G4double cross = SandiaCof[0] / omega + SandiaCof[1] / omega2 +
1050 SandiaCof[2] / omega3 + SandiaCof[3] / omega4;
1051 return cross;
1052}

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 1123 of file G4VXTRenergyLoss.cc.

1124{
1125 std::ofstream outPlate("plateZmu.dat", std::ios::out);
1126 outPlate.setf(std::ios::scientific, std::ios::floatfield);
1127
1128 G4int i;
1129 G4double omega, varAngle, gamma;
1130 gamma = 10000.;
1131 varAngle = 1 / gamma / gamma;
1132 if(verboseLevel > 0)
1133 G4cout << "energy, keV" << "\t" << "Zmu for plate" << G4endl;
1134 for(i = 0; i < 100; ++i)
1135 {
1136 omega = (1.0 + i) * keV;
1137 if(verboseLevel > 1)
1138 G4cout << omega / keV << "\t"
1139 << GetPlateZmuProduct(omega, gamma, varAngle) << "\t";
1140 if(verboseLevel > 0)
1141 outPlate << omega / keV << "\t\t"
1142 << GetPlateZmuProduct(omega, gamma, varAngle) << G4endl;
1143 }
1144 return;
1145}

Referenced by GetPlateZmuProduct().

◆ GetPlateZmuProduct() [2/2]

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

Definition at line 1114 of file G4VXTRenergyLoss.cc.

1116{
1117 return GetPlateFormationZone(omega, gamma, varAngle) *
1119}

◆ GetProtonVector()

G4PhysicsLogVector * G4VXTRenergyLoss::GetProtonVector ( )
inline

Definition at line 193 of file G4VXTRenergyLoss.hh.

193{ return fProtonEnergyVector; };

◆ GetRandomAngle()

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

Definition at line 1457 of file G4VXTRenergyLoss.cc.

1458{
1459 G4int iTR, iAngle;
1460 G4double position, angle;
1461
1462 if(iTkin == fTotBin)
1463 --iTkin;
1464
1466
1467 for(iTR = 0; iTR < fBinTR; ++iTR)
1468 {
1469 if(energyXTR < fXTREnergyVector->GetLowEdgeEnergy(iTR))
1470 break;
1471 }
1472 if(iTR == fBinTR)
1473 --iTR;
1474
1475 position = (*(*fAngleForEnergyTable)(iTR))(0) * G4UniformRand();
1476 // position = (*(*fAngleForEnergyTable)(iTR))(1) * G4UniformRand(); // ATLAS TB
1477
1478 for(iAngle = 0;; ++iAngle)
1479 // for(iAngle = 1;; ++iAngle) // ATLAS TB
1480 {
1481 if(position >= (*(*fAngleForEnergyTable)(iTR))(iAngle))
1482 break;
1483 }
1484 angle = GetAngleXTR(iTR, position, iAngle);
1485 return angle;
1486}
G4double GetAngleXTR(G4int iTR, G4double position, G4int iAngle)

Referenced by PostStepDoIt().

◆ GetStackFactor()

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

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

Definition at line 1301 of file G4VXTRenergyLoss.cc.

1303{
1304 // return stack factor corresponding to one interface
1305 return std::real(OneInterfaceXTRdEdx(energy, gamma, varAngle));
1306}
G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma, G4double varAngle)

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

◆ GetTheMaxAngle()

G4double G4VXTRenergyLoss::GetTheMaxAngle ( )
inline

Definition at line 176 of file G4VXTRenergyLoss.hh.

176{ return fTheMaxAngle; };

◆ GetTheMaxEnergyTR()

G4double G4VXTRenergyLoss::GetTheMaxEnergyTR ( )
inline

Definition at line 166 of file G4VXTRenergyLoss.hh.

166{ return fTheMaxEnergyTR; };

◆ GetTheMinAngle()

G4double G4VXTRenergyLoss::GetTheMinAngle ( )
inline

Definition at line 174 of file G4VXTRenergyLoss.hh.

174{ return fTheMinAngle; };

◆ GetTheMinEnergyTR()

G4double G4VXTRenergyLoss::GetTheMinEnergyTR ( )
inline

Definition at line 164 of file G4VXTRenergyLoss.hh.

164{ return fTheMinEnergyTR; };

◆ GetTotBin()

G4int G4VXTRenergyLoss::GetTotBin ( )
inline

Definition at line 194 of file G4VXTRenergyLoss.hh.

194{ return fTotBin; };

◆ GetVarAngle()

G4double G4VXTRenergyLoss::GetVarAngle ( )
inline

Definition at line 150 of file G4VXTRenergyLoss.hh.

150{ return fVarAngle; };

◆ GetXTRenergy()

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

Definition at line 1423 of file G4VXTRenergyLoss.cc.

1424{
1425 G4double x1, x2, y1, y2, result;
1426
1427 if(iTransfer == 0)
1428 {
1429 result = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1430 }
1431 else
1432 {
1433 y1 = (*(*fEnergyDistrTable)(iPlace))(iTransfer - 1);
1434 y2 = (*(*fEnergyDistrTable)(iPlace))(iTransfer);
1435
1436 x1 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer - 1);
1437 x2 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1438
1439 if(x1 == x2)
1440 result = x2;
1441 else
1442 {
1443 if(y1 == y2)
1444 result = x1 + (x2 - x1) * G4UniformRand();
1445 else
1446 {
1447 result = x1 + (x2 - x1) * G4UniformRand();
1448 }
1449 }
1450 }
1451 return result;
1452}

Referenced by GetXTRrandomEnergy().

◆ GetXTRrandomEnergy()

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

Definition at line 1377 of file G4VXTRenergyLoss.cc.

1378{
1379 G4int iTransfer, iPlace;
1380 G4double transfer = 0.0, position, E1, E2, W1, W2, W;
1381
1382 iPlace = iTkin - 1;
1383
1384 if(iTkin == fTotBin) // relativistic plato, try from left
1385 {
1386 position = (*(*fEnergyDistrTable)(iPlace))(0) * G4UniformRand();
1387
1388 for(iTransfer = 0;; ++iTransfer)
1389 {
1390 if(position >= (*(*fEnergyDistrTable)(iPlace))(iTransfer))
1391 break;
1392 }
1393 transfer = GetXTRenergy(iPlace, position, iTransfer);
1394 }
1395 else
1396 {
1397 E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1);
1399 W = 1.0 / (E2 - E1);
1400 W1 = (E2 - scaledTkin) * W;
1401 W2 = (scaledTkin - E1) * W;
1402
1403 position = ((*(*fEnergyDistrTable)(iPlace))(0) * W1 +
1404 (*(*fEnergyDistrTable)(iPlace + 1))(0) * W2) *
1405 G4UniformRand();
1406
1407 for(iTransfer = 0;; ++iTransfer)
1408 {
1409 if(position >= ((*(*fEnergyDistrTable)(iPlace))(iTransfer) *W1 +
1410 (*(*fEnergyDistrTable)(iPlace + 1))(iTransfer) *W2))
1411 break;
1412 }
1413 transfer = GetXTRenergy(iPlace, position, iTransfer);
1414 }
1415 if(transfer < 0.0)
1416 transfer = 0.0;
1417 return transfer;
1418}
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 1288 of file G4VXTRenergyLoss.cc.

1291{
1292 G4double formationLength1, formationLength2;
1293 formationLength1 =
1294 1.0 / (1.0 / (gamma * gamma) + fSigma1 / (energy * energy) + varAngle);
1295 formationLength2 =
1296 1.0 / (1.0 / (gamma * gamma) + fSigma2 / (energy * energy) + varAngle);
1297 return (varAngle / energy) * (formationLength1 - formationLength2) *
1298 (formationLength1 - formationLength2);
1299}

Referenced by XTRNAngleSpectralDensity(), and XTRNSpectralAngleDensity().

◆ OneInterfaceXTRdEdx()

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

Definition at line 852 of file G4VXTRenergyLoss.cc.

854{
855 G4complex Z1 = GetPlateComplexFZ(energy, gamma, varAngle);
856 G4complex Z2 = GetGasComplexFZ(energy, gamma, varAngle);
857
858 G4complex zOut = (Z1 - Z2) * (Z1 - Z2) * (varAngle * energy / hbarc / hbarc);
859 return zOut;
860}
G4complex GetPlateComplexFZ(G4double, G4double, G4double)
G4complex GetGasComplexFZ(G4double, G4double, G4double)

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

◆ operator=()

G4VXTRenergyLoss & G4VXTRenergyLoss::operator= ( const G4VXTRenergyLoss & right)
delete

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 714 of file G4VXTRenergyLoss.cc.

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

◆ ProcessDescription()

void G4VXTRenergyLoss::ProcessDescription ( std::ostream & out) const
overridevirtual

Reimplemented from G4VProcess.

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

Definition at line 184 of file G4VXTRenergyLoss.cc.

185{
186 out << "Base class for 'fast' parameterisation model describing X-ray "
187 "transition\n"
188 "radiation. Angular distribution is very rough.\n";
189}

Referenced by DumpInfo().

◆ SetAlphaGas()

void G4VXTRenergyLoss::SetAlphaGas ( G4double ag)
inline

Definition at line 158 of file G4VXTRenergyLoss.hh.

158{ fAlphaGas = ag;};

◆ SetAlphaPlate()

void G4VXTRenergyLoss::SetAlphaPlate ( G4double ap)
inline

Definition at line 160 of file G4VXTRenergyLoss.hh.

◆ SetAngleRadDistr()

void G4VXTRenergyLoss::SetAngleRadDistr ( G4bool fatr)
inline

Definition at line 187 of file G4VXTRenergyLoss.hh.

187{ fAngleRadDistr = fatr;};

◆ SetCompton()

void G4VXTRenergyLoss::SetCompton ( G4bool pC)
inline

Definition at line 151 of file G4VXTRenergyLoss.hh.

151{ fCompton = pC; };

◆ SetEnergy()

void G4VXTRenergyLoss::SetEnergy ( G4double energy)
inline

Definition at line 147 of file G4VXTRenergyLoss.hh.

147{ fEnergy = energy; };

◆ SetFastAngle()

void G4VXTRenergyLoss::SetFastAngle ( G4bool fatr)
inline

Definition at line 185 of file G4VXTRenergyLoss.hh.

185{ fFastAngle = fatr;};

◆ SetGamma()

void G4VXTRenergyLoss::SetGamma ( G4double gamma)
inline

Definition at line 145 of file G4VXTRenergyLoss.hh.

145{ fGamma = gamma; };

◆ SetKrange()

void G4VXTRenergyLoss::SetKrange ( G4int kk)
inline

Definition at line 155 of file G4VXTRenergyLoss.hh.

155{ fKrange = kk;};

◆ SetMaxEnergyTR()

void G4VXTRenergyLoss::SetMaxEnergyTR ( G4double maxetr)
inline

Definition at line 170 of file G4VXTRenergyLoss.hh.

170{ fMaxEnergyTR = maxetr;};

◆ SetMaxThetaTR()

void G4VXTRenergyLoss::SetMaxThetaTR ( G4double maxatr)
inline

Definition at line 180 of file G4VXTRenergyLoss.hh.

180{ fMaxThetaTR = maxatr;};

◆ SetMinEnergyTR()

void G4VXTRenergyLoss::SetMinEnergyTR ( G4double minetr)
inline

Definition at line 168 of file G4VXTRenergyLoss.hh.

168{ fMinEnergyTR = minetr;};

◆ SetMinThetaTR()

void G4VXTRenergyLoss::SetMinThetaTR ( G4double minatr)
inline

Definition at line 178 of file G4VXTRenergyLoss.hh.

178{ fMinThetaTR = minatr;};

◆ SetTheMaxAngle()

void G4VXTRenergyLoss::SetTheMaxAngle ( G4double maxang)
inline

Definition at line 175 of file G4VXTRenergyLoss.hh.

175{ fTheMaxAngle = maxang;};

◆ SetTheMaxEnergyTR()

void G4VXTRenergyLoss::SetTheMaxEnergyTR ( G4double maxetr)
inline

Definition at line 165 of file G4VXTRenergyLoss.hh.

165{ fTheMaxEnergyTR = maxetr;};

◆ SetTheMinAngle()

void G4VXTRenergyLoss::SetTheMinAngle ( G4double minang)
inline

Definition at line 173 of file G4VXTRenergyLoss.hh.

173{ fTheMinAngle = minang;};

◆ SetTheMinEnergyTR()

void G4VXTRenergyLoss::SetTheMinEnergyTR ( G4double minetr)
inline

Definition at line 163 of file G4VXTRenergyLoss.hh.

163{ fTheMinEnergyTR = minetr;};

◆ SetVarAngle()

void G4VXTRenergyLoss::SetVarAngle ( G4double varAngle)
inline

Definition at line 149 of file G4VXTRenergyLoss.hh.

149{ fVarAngle = varAngle; };

◆ SpectralAngleXTRdEdx()

G4double G4VXTRenergyLoss::SpectralAngleXTRdEdx ( G4double varAngle)

Definition at line 864 of file G4VXTRenergyLoss.cc.

865{
866 G4double result = GetStackFactor(fEnergy, fGamma, varAngle);
867 if(result < 0.0)
868 result = 0.0;
869 return result;
870}

Referenced by BuildAngleForEnergyBank(), and SpectralXTRdEdx().

◆ SpectralXTRdEdx()

G4double G4VXTRenergyLoss::SpectralXTRdEdx ( G4double energy)
virtual

Reimplemented in G4GaussXTRadiator, G4RegularXTRadiator, G4TransparentRegXTRadiator, G4XTRRegularRadModel, and G4XTRTransparentRegRadModel.

Definition at line 874 of file G4VXTRenergyLoss.cc.

875{
876 G4int i;
877 static constexpr G4int iMax = 8;
878 G4double angleSum = 0.0;
879
880 G4double lim[iMax] = { 0.0, 0.01, 0.02, 0.05, 0.1, 0.2, 0.5, 1.0 };
881
882 for(i = 0; i < iMax; ++i)
883 lim[i] *= fMaxThetaTR;
884
886 integral;
887
888 fEnergy = energy;
889 {
890 for(i = 0; i < iMax - 1; ++i)
891 {
892 angleSum += integral.Legendre96(
893 this, &G4VXTRenergyLoss::SpectralAngleXTRdEdx, lim[i], lim[i + 1]);
894 }
895 }
896 return angleSum;
897}

Referenced by BuildEnergyTable().

◆ XTRNAngleDensity()

G4double G4VXTRenergyLoss::XTRNAngleDensity ( G4double varAngle)

Definition at line 1338 of file G4VXTRenergyLoss.cc.

1339{
1340 fVarAngle = varAngle;
1342 integral;
1345}
G4double Legendre96(T &typeT, F f, G4double a, G4double b)
G4double XTRNAngleSpectralDensity(G4double energy)

◆ XTRNAngleSpectralDensity()

G4double G4VXTRenergyLoss::XTRNAngleSpectralDensity ( G4double energy)

Definition at line 1331 of file G4VXTRenergyLoss.cc.

1332{
1333 return OneBoundaryXTRNdensity(energy, fGamma, fVarAngle) *
1335}
G4double OneBoundaryXTRNdensity(G4double energy, G4double gamma, G4double varAngle) const

Referenced by XTRNAngleDensity().

◆ XTRNSpectralAngleDensity()

G4double G4VXTRenergyLoss::XTRNSpectralAngleDensity ( G4double varAngle)

Definition at line 1310 of file G4VXTRenergyLoss.cc.

1311{
1312 return OneBoundaryXTRNdensity(fEnergy, fGamma, varAngle) *
1313 GetStackFactor(fEnergy, fGamma, varAngle);
1314}

Referenced by XTRNSpectralDensity().

◆ XTRNSpectralDensity()

G4double G4VXTRenergyLoss::XTRNSpectralDensity ( G4double energy)

Definition at line 1318 of file G4VXTRenergyLoss.cc.

1319{
1320 fEnergy = energy;
1322 integral;
1324 0.0, 0.2 * fMaxThetaTR) +
1325 integral.Legendre10(this, &G4VXTRenergyLoss::XTRNSpectralAngleDensity,
1326 0.2 * fMaxThetaTR, fMaxThetaTR);
1327}
G4double XTRNSpectralAngleDensity(G4double varAngle)

Member Data Documentation

◆ fAlphaGas

◆ fAlphaPlate

◆ fAngleBank

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

Definition at line 233 of file G4VXTRenergyLoss.hh.

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

◆ fAngleDistrTable

G4PhysicsTable* G4VXTRenergyLoss::fAngleDistrTable
protected

Definition at line 224 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 = CLHEP::fine_structure_const / CLHEP::pi
staticconstexprprotected

Definition at line 215 of file G4VXTRenergyLoss.hh.

Referenced by BuildEnergyTable(), and BuildGlobalAngleTable().

◆ fCompton

G4bool G4VXTRenergyLoss::fCompton
protected

◆ fEnergy

◆ fEnergyDistrTable

G4PhysicsTable* G4VXTRenergyLoss::fEnergyDistrTable
protected

Definition at line 225 of file G4VXTRenergyLoss.hh.

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

◆ fEnvelope

G4LogicalVolume* G4VXTRenergyLoss::fEnvelope
protected

Definition at line 223 of file G4VXTRenergyLoss.hh.

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

◆ fExitFlux

◆ fFastAngle

G4bool G4VXTRenergyLoss::fFastAngle
protected

◆ fGamma

◆ fGammaCutInKineticEnergy

G4double* G4VXTRenergyLoss::fGammaCutInKineticEnergy
protected

Definition at line 222 of file G4VXTRenergyLoss.hh.

◆ fGammaTkinCut

G4double G4VXTRenergyLoss::fGammaTkinCut
protected

◆ fGasPhotoAbsCof

G4SandiaTable* G4VXTRenergyLoss::fGasPhotoAbsCof
protected

Definition at line 230 of file G4VXTRenergyLoss.hh.

Referenced by ComputeGasPhotoAbsCof(), and GetGasLinearPhotoAbs().

◆ fGasThick

◆ fKrange

G4int G4VXTRenergyLoss::fKrange
protected

◆ fLambda

G4double G4VXTRenergyLoss::fLambda
protected

Definition at line 247 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss(), and GetMeanFreePath().

◆ fMatIndex1

G4int G4VXTRenergyLoss::fMatIndex1
protected

Definition at line 251 of file G4VXTRenergyLoss.hh.

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

◆ fMatIndex2

G4int G4VXTRenergyLoss::fMatIndex2
protected

Definition at line 252 of file G4VXTRenergyLoss.hh.

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

◆ fMaxEnergyTR

◆ fMaxProtonTkin

G4double G4VXTRenergyLoss::fMaxProtonTkin = 100. * CLHEP::TeV
staticconstexprprotected

Definition at line 210 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss().

◆ fMaxThetaTR

◆ fMinEnergyTR

◆ fMinProtonTkin

G4double G4VXTRenergyLoss::fMinProtonTkin = 100. * CLHEP::GeV
staticconstexprprotected

Definition at line 208 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss().

◆ fMinThetaTR

G4double G4VXTRenergyLoss::fMinThetaTR
protected

◆ fParticleChange

G4ParticleChange G4VXTRenergyLoss::fParticleChange
protected

Definition at line 232 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss(), and PostStepDoIt().

◆ fPlasmaCof

G4double G4VXTRenergyLoss::fPlasmaCof
staticconstexprprotected
Initial value:
=
4. * CLHEP::pi * CLHEP::fine_structure_const * CLHEP::hbarc * CLHEP::hbarc *
CLHEP::hbarc / CLHEP::electron_mass_c2

Definition at line 212 of file G4VXTRenergyLoss.hh.

Referenced by G4StrawTubeXTRadiator::G4StrawTubeXTRadiator(), and G4VXTRenergyLoss().

◆ fPlateNumber

◆ fPlatePhotoAbsCof

G4SandiaTable* G4VXTRenergyLoss::fPlatePhotoAbsCof
protected

Definition at line 229 of file G4VXTRenergyLoss.hh.

Referenced by ComputePlatePhotoAbsCof(), and GetPlateLinearPhotoAbs().

◆ fPlateThick

◆ fProtonEnergyVector

◆ fPtrGamma

G4ParticleDefinition* G4VXTRenergyLoss::fPtrGamma
protected

Definition at line 220 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss().

◆ fSigma1

◆ fSigma2

◆ fTheMaxAngle

G4double G4VXTRenergyLoss::fTheMaxAngle
protected

◆ fTheMaxEnergyTR

◆ fTheMinAngle

G4double G4VXTRenergyLoss::fTheMinAngle
protected

◆ fTheMinEnergyTR

◆ fTotalDist

G4double G4VXTRenergyLoss::fTotalDist
protected

Definition at line 239 of file G4VXTRenergyLoss.hh.

Referenced by BuildEnergyTable(), and G4VXTRenergyLoss().

◆ fTotBin

◆ fVarAngle

G4double G4VXTRenergyLoss::fVarAngle
protected

◆ fXTREnergyVector

G4PhysicsLogVector* G4VXTRenergyLoss::fXTREnergyVector
protected

Definition at line 228 of file G4VXTRenergyLoss.hh.

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

◆ secID

G4int G4VXTRenergyLoss::secID = -1
protected

Definition at line 259 of file G4VXTRenergyLoss.hh.

Referenced by G4VXTRenergyLoss().


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