Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VXTRenergyLoss Class Reference

#include <G4VXTRenergyLoss.hh>

Inheritance diagram for G4VXTRenergyLoss:
Collaboration 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 &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

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
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Additional Inherited Members

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

Detailed Description

Definition at line 74 of file G4VXTRenergyLoss.hh.

Constructor & Destructor Documentation

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

68  :
69  G4VDiscreteProcess(processName, type),
70  fGammaCutInKineticEnergy(nullptr),
71  fGammaTkinCut(0.0),
72  fAngleDistrTable(nullptr),
73  fEnergyDistrTable(nullptr),
74  fPlatePhotoAbsCof(nullptr),
75  fGasPhotoAbsCof(nullptr),
76  fAngleForEnergyTable(nullptr)
77 {
78  verboseLevel = 1;
80 
81  fPtrGamma = nullptr;
84 
85  // Initialization of local constants
86  fTheMinEnergyTR = 1.0*keV;
87  fTheMaxEnergyTR = 100.0*keV;
88  fTheMaxAngle = 1.0e-2;
89  fTheMinAngle = 5.0e-6;
90  fBinTR = 50;
91 
92  fMinProtonTkin = 100.0*GeV;
93  fMaxProtonTkin = 100.0*TeV;
94  fTotBin = 50;
95 
96  // Proton energy vector initialization
99  fTotBin );
100 
103  fBinTR );
104 
106 
108 
109  fEnvelope = anEnvelope;
110 
111  fPlateNumber = n;
112  if(verboseLevel > 0)
113  G4cout<<"### G4VXTRenergyLoss: the number of TR radiator plates = "
114  <<fPlateNumber<<G4endl;
115  if(fPlateNumber == 0)
116  {
117  G4Exception("G4VXTRenergyLoss::G4VXTRenergyLoss()","VXTRELoss01",
118  FatalException,"No plates in X-ray TR radiator");
119  }
120  // default is XTR dEdx, not flux after radiator
121 
122  fExitFlux = false;
123  fAngleRadDistr = false;
124  fCompton = false;
125 
126  fLambda = DBL_MAX;
127  // Mean thicknesses of plates and gas gaps
128 
129  fPlateThick = a;
130  fGasThick = b;
132  if(verboseLevel > 0)
133  G4cout<<"total radiator thickness = "<<fTotalDist/cm<<" cm"<<G4endl;
134 
135  // index of plate material
136  fMatIndex1 = foilMat->GetIndex();
137  if(verboseLevel > 0)
138  G4cout<<"plate material = "<<foilMat->GetName()<<G4endl;
139 
140  // index of gas material
141  fMatIndex2 = gasMat->GetIndex();
142  if(verboseLevel > 0)
143  G4cout<<"gas material = "<<gasMat->GetName()<<G4endl;
144 
145  // plasma energy squared for plate material
146 
148  // fSigma1 = (20.9*eV)*(20.9*eV);
149  if(verboseLevel > 0)
150  G4cout<<"plate plasma energy = "<<std::sqrt(fSigma1)/eV<<" eV"<<G4endl;
151 
152  // plasma energy squared for gas material
153 
155  if(verboseLevel > 0)
156  G4cout<<"gas plasma energy = "<<std::sqrt(fSigma2)/eV<<" eV"<<G4endl;
157 
158  // Compute cofs for preparation of linear photo absorption
159 
162 
164 }
G4PhysicsTable * fEnergyDistrTable
G4int verboseLevel
Definition: G4VProcess.hh:368
G4LogicalVolume * fEnvelope
size_t GetIndex() const
Definition: G4Material.hh:262
static constexpr double hbarc
const G4String & GetName() const
Definition: G4Material.hh:178
static constexpr double electron_mass_c2
static constexpr double TeV
Definition: G4SIunits.hh:218
G4GLOB_DLL std::ostream G4cout
G4PhysicsLogVector * fProtonEnergyVector
G4double GetElectronDensity() const
Definition: G4Material.hh:217
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4SandiaTable * fGasPhotoAbsCof
G4double * fGammaCutInKineticEnergy
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
G4PhysicsLogVector * fXTREnergyVector
static constexpr double GeV
Definition: G4SIunits.hh:217
G4PhysicsTable * fAngleForEnergyTable
#define G4endl
Definition: G4ios.hh:61
G4PhysicsTable * fAngleDistrTable
static constexpr double pi
Definition: G4SIunits.hh:75
G4SandiaTable * fPlatePhotoAbsCof
G4ParticleDefinition * fPtrGamma
G4ParticleChange fParticleChange
static constexpr double fine_structure_const
#define DBL_MAX
Definition: templates.hh:83
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4VXTRenergyLoss::~G4VXTRenergyLoss ( )
virtual

Definition at line 168 of file G4VXTRenergyLoss.cc.

169 {
170  if(fEnvelope) delete fEnvelope;
171  delete fProtonEnergyVector;
172  delete fXTREnergyVector;
173  if(fEnergyDistrTable) {
175  delete fEnergyDistrTable;
176  }
177  if(fAngleRadDistr) {
179  delete fAngleDistrTable;
180  }
183  delete fAngleForEnergyTable;
184  }
185 }
G4PhysicsTable * fEnergyDistrTable
G4LogicalVolume * fEnvelope
G4PhysicsLogVector * fProtonEnergyVector
G4PhysicsLogVector * fXTREnergyVector
G4PhysicsTable * fAngleForEnergyTable
G4PhysicsTable * fAngleDistrTable
void clearAndDestroy()

Here is the call graph for this function:

Member Function Documentation

G4double G4VXTRenergyLoss::AngleSpectralXTRdEdx ( G4double  energy)

Definition at line 945 of file G4VXTRenergyLoss.cc.

946 {
948  if(result < 0) result = 0.0;
949  return result;
950 }
G4double G4ParticleHPJENDLHEData::G4double result
G4double energy(const ThreeVector &p, const G4double m)
double G4double
Definition: G4Types.hh:76
virtual G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle)

Here is the call graph for this function:

G4double G4VXTRenergyLoss::AngleXTRdEdx ( G4double  varAngle)

Definition at line 956 of file G4VXTRenergyLoss.cc.

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

Here is the caller graph for this function:

void G4VXTRenergyLoss::BuildAngleForEnergyBank ( )

Definition at line 391 of file G4VXTRenergyLoss.cc.

392 {
393  if( this->GetProcessName() == "TranspRegXTRadiator" ||
394  this->GetProcessName() == "TranspRegXTRmodel" ||
395  this->GetProcessName() == "RegularXTRadiator" ||
396  this->GetProcessName() == "RegularXTRmodel" )
397  {
398  BuildAngleTable();
399  return;
400  }
401  G4int i, iTkin, iTR;
402  G4double angleSum = 0.0;
403 
404 
405  fGammaTkinCut = 0.0;
406 
407  // setting of min/max TR energies
408 
411 
414 
416  fMaxEnergyTR,
417  fBinTR );
418 
420 
421  G4cout.precision(4);
422  G4Timer timer;
423  timer.Start();
424 
425  for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
426  {
427 
428  fGamma = 1.0 + (fProtonEnergyVector->
429  GetLowEdgeEnergy(iTkin)/proton_mass_c2);
430 
431  fMaxThetaTR = 2500.0/(fGamma*fGamma) ; // theta^2
432 
433  fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
434 
437 
439 
440  for( iTR = 0; iTR < fBinTR; iTR++ )
441  {
442  angleSum = 0.0;
443  fEnergy = energyVector->GetLowEdgeEnergy(iTR);
444  G4PhysicsLinearVector* angleVector = new G4PhysicsLinearVector(0.0,
445  fMaxThetaTR,
446  fBinTR );
447 
448  angleVector ->PutValue(fBinTR - 1, angleSum);
449 
450  for( i = fBinTR - 2; i >= 0; i-- )
451  {
452  // Legendre96 or Legendre10
453 
454  angleSum += integral.Legendre10(
456  angleVector->GetLowEdgeEnergy(i),
457  angleVector->GetLowEdgeEnergy(i+1) );
458 
459  angleVector ->PutValue(i, angleSum);
460  }
461  fAngleForEnergyTable->insertAt(iTR, angleVector);
462  }
463  fAngleBank.push_back(fAngleForEnergyTable);
464  }
465  timer.Stop();
466  G4cout.precision(6);
467  if(verboseLevel > 0)
468  {
469  G4cout<<G4endl;
470  G4cout<<"total time for build X-ray TR angle for energy loss tables = "
471  <<timer.GetUserElapsed()<<" s"<<G4endl;
472  }
473  fGamma = 0.;
474  delete energyVector;
475 }
G4double Legendre10(T &typeT, F f, G4double a, G4double b)
G4int verboseLevel
Definition: G4VProcess.hh:368
static constexpr double proton_mass_c2
G4double GetLowEdgeEnergy(size_t binNumber) const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4PhysicsLogVector * fProtonEnergyVector
G4double GetUserElapsed() const
Definition: G4Timer.cc:130
void PutValue(size_t index, G4double theValue)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void Stop()
G4double SpectralAngleXTRdEdx(G4double varAngle)
std::vector< G4PhysicsTable * > fAngleBank
void insertAt(size_t, G4PhysicsVector *)
G4PhysicsTable * fAngleForEnergyTable
#define G4endl
Definition: G4ios.hh:61
void Start()
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4VXTRenergyLoss::BuildAngleTable ( )

Definition at line 482 of file G4VXTRenergyLoss.cc.

483 {
484  G4int iTkin, iTR;
486 
487  fGammaTkinCut = 0.0;
488 
489  // setting of min/max TR energies
490 
493 
496 
497  G4cout.precision(4);
498  G4Timer timer;
499  timer.Start();
500  if(verboseLevel > 0)
501  {
502  G4cout<<G4endl;
503  G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
504  G4cout<<G4endl;
505  }
506  for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
507  {
508 
509  fGamma = 1.0 + (fProtonEnergyVector->
510  GetLowEdgeEnergy(iTkin)/proton_mass_c2);
511 
512  fMaxThetaTR = 25.0/(fGamma*fGamma); // theta^2
513 
514  fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
515 
517  else
518  {
520  }
521 
523 
524  for( iTR = 0; iTR < fBinTR; iTR++ )
525  {
526  // energy = fMinEnergyTR*(iTR+1);
527 
528  energy = fXTREnergyVector->GetLowEdgeEnergy(iTR);
529 
530  G4PhysicsFreeVector* angleVector = GetAngleVector(energy,fBinTR);
531  // G4cout<<G4endl;
532 
533  fAngleForEnergyTable->insertAt(iTR,angleVector);
534  }
535 
536  fAngleBank.push_back(fAngleForEnergyTable);
537  }
538  timer.Stop();
539  G4cout.precision(6);
540  if(verboseLevel > 0)
541  {
542  G4cout<<G4endl;
543  G4cout<<"total time for build XTR angle for given energy tables = "
544  <<timer.GetUserElapsed()<<" s"<<G4endl;
545  }
546  fGamma = 0.;
547 
548  return;
549 }
G4int verboseLevel
Definition: G4VProcess.hh:368
static constexpr double proton_mass_c2
G4double GetLowEdgeEnergy(size_t binNumber) const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4PhysicsLogVector * fProtonEnergyVector
G4double GetUserElapsed() const
Definition: G4Timer.cc:130
G4PhysicsFreeVector * GetAngleVector(G4double energy, G4int n)
G4double energy(const ThreeVector &p, const G4double m)
void Stop()
G4PhysicsLogVector * fXTREnergyVector
std::vector< G4PhysicsTable * > fAngleBank
void insertAt(size_t, G4PhysicsVector *)
G4PhysicsTable * fAngleForEnergyTable
#define G4endl
Definition: G4ios.hh:61
void Start()
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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 
307 
308  fGammaTkinCut = 0.0;
309 
310  // setting of min/max TR energies
311 
314 
317 
319 
320  G4cout.precision(4);
321  G4Timer timer;
322  timer.Start();
323 
324  if(verboseLevel > 0)
325  {
326  G4cout<<G4endl;
327  G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
328  G4cout<<G4endl;
329  }
330  for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
331  {
333  fMaxEnergyTR,
334  fBinTR );
335 
336  fGamma = 1.0 + (fProtonEnergyVector->
337  GetLowEdgeEnergy(iTkin)/proton_mass_c2);
338 
339  fMaxThetaTR = 2500.0/(fGamma*fGamma) ; // theta^2
340 
341  fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
342 
345 
346  energySum = 0.0;
347 
348  energyVector->PutValue(fBinTR-1,energySum);
349 
350  for( iTR = fBinTR - 2; iTR >= 0; iTR-- )
351  {
352  // Legendre96 or Legendre10
353 
354  energySum += radiatorCof*fCofTR*integral.Legendre10(
356  energyVector->GetLowEdgeEnergy(iTR),
357  energyVector->GetLowEdgeEnergy(iTR+1) );
358 
359  energyVector->PutValue(iTR,energySum/fTotalDist);
360  }
361  iPlace = iTkin;
362  fEnergyDistrTable->insertAt(iPlace,energyVector);
363 
364  if(verboseLevel > 0)
365  {
366  G4cout
367  // <<iTkin<<"\t"
368  // <<"fGamma = "
369  <<fGamma<<"\t" // <<" fMaxThetaTR = "<<fMaxThetaTR
370  // <<"sumN = "
371  <<energySum // <<"; sumA = "<<angleSum
372  <<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 }
G4double Legendre10(T &typeT, F f, G4double a, G4double b)
G4PhysicsTable * fEnergyDistrTable
G4int verboseLevel
Definition: G4VProcess.hh:368
static constexpr double proton_mass_c2
G4double GetLowEdgeEnergy(size_t binNumber) const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4PhysicsLogVector * fProtonEnergyVector
G4double GetUserElapsed() const
Definition: G4Timer.cc:130
void PutValue(size_t index, G4double theValue)
void Stop()
virtual G4double SpectralXTRdEdx(G4double energy)
void insertAt(size_t, G4PhysicsVector *)
#define G4endl
Definition: G4ios.hh:61
G4PhysicsTable * fAngleDistrTable
void Start()
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4VXTRenergyLoss::BuildGlobalAngleTable ( )

Definition at line 637 of file G4VXTRenergyLoss.cc.

638 {
639  G4int iTkin, iTR, iPlace;
640  G4double radiatorCof = 1.0; // for tuning of XTR yield
641  G4double angleSum;
643 
644  fGammaTkinCut = 0.0;
645 
646  // setting of min/max TR energies
647 
650 
653 
654  G4cout.precision(4);
655  G4Timer timer;
656  timer.Start();
657  if(verboseLevel > 0) {
658  G4cout<<G4endl;
659  G4cout<<"Lorentz Factor"<<"\t"<<"XTR photon number"<<G4endl;
660  G4cout<<G4endl;
661  }
662  for( iTkin = 0; iTkin < fTotBin; iTkin++ ) // Lorentz factor loop
663  {
664 
665  fGamma = 1.0 + (fProtonEnergyVector->
666  GetLowEdgeEnergy(iTkin)/proton_mass_c2);
667 
668  fMaxThetaTR = 25.0/(fGamma*fGamma); // theta^2
669 
670  fTheMinAngle = 1.0e-3; // was 5.e-6, e-6 !!!, e-5, e-4
671 
673  else
674  {
676  }
677  G4PhysicsLinearVector* angleVector = new G4PhysicsLinearVector(0.0,
678  fMaxThetaTR,
679  fBinTR );
680 
681  angleSum = 0.0;
682 
684 
685 
686  angleVector->PutValue(fBinTR-1,angleSum);
687 
688  for( iTR = fBinTR - 2; iTR >= 0; iTR-- )
689  {
690 
691  angleSum += radiatorCof*fCofTR*integral.Legendre96(
693  angleVector->GetLowEdgeEnergy(iTR),
694  angleVector->GetLowEdgeEnergy(iTR+1) );
695 
696  angleVector ->PutValue(iTR,angleSum);
697  }
698  if(verboseLevel > 1) {
699  G4cout
700  // <<iTkin<<"\t"
701  // <<"fGamma = "
702  <<fGamma<<"\t" // <<" fMaxThetaTR = "<<fMaxThetaTR
703  // <<"sumN = "<<energySum // <<"; sumA = "
704  <<angleSum
705  <<G4endl;
706  }
707  iPlace = iTkin;
708  fAngleDistrTable->insertAt(iPlace,angleVector);
709  }
710  timer.Stop();
711  G4cout.precision(6);
712  if(verboseLevel > 0) {
713  G4cout<<G4endl;
714  G4cout<<"total time for build X-ray TR angle tables = "
715  <<timer.GetUserElapsed()<<" s"<<G4endl;
716  }
717  fGamma = 0.;
718 
719  return;
720 }
G4double Legendre96(T &typeT, F f, G4double a, G4double b)
G4int verboseLevel
Definition: G4VProcess.hh:368
static constexpr double proton_mass_c2
G4double GetLowEdgeEnergy(size_t binNumber) const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4PhysicsLogVector * fProtonEnergyVector
G4double GetUserElapsed() const
Definition: G4Timer.cc:130
void PutValue(size_t index, G4double theValue)
void Stop()
void insertAt(size_t, G4PhysicsVector *)
#define G4endl
Definition: G4ios.hh:61
G4PhysicsTable * fAngleDistrTable
void Start()
double G4double
Definition: G4Types.hh:76
G4double AngleXTRdEdx(G4double varAngle)

Here is the call graph for this function:

void G4VXTRenergyLoss::BuildPhysicsTable ( const G4ParticleDefinition pd)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 274 of file G4VXTRenergyLoss.cc.

275 {
276  if(pd.GetPDGCharge() == 0.)
277  {
278  G4Exception("G4VXTRenergyLoss::BuildPhysicsTable", "Notification", JustWarning,
279  "XTR initialisation for neutral particle ?!" );
280  }
282 
283  if (fAngleRadDistr)
284  {
285  if(verboseLevel > 0)
286  {
287  G4cout<<"Build angle for energy distribution according the current radiator"
288  <<G4endl;
289  }
291  }
292 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4GLOB_DLL std::ostream G4cout
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
G4double GetPDGCharge() const

Here is the call graph for this function:

void G4VXTRenergyLoss::BuildTable ( )
inline

Definition at line 102 of file G4VXTRenergyLoss.hh.

102 {} ;
void G4VXTRenergyLoss::ComputeGasPhotoAbsCof ( )

Definition at line 1154 of file G4VXTRenergyLoss.cc.

1155 {
1156  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1157  const G4Material* mat = (*theMaterialTable)[fMatIndex2];
1159  return;
1160 }
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:229
G4SandiaTable * fGasPhotoAbsCof

Here is the call graph for this function:

Here is the caller graph for this function:

void G4VXTRenergyLoss::ComputePlatePhotoAbsCof ( )

Definition at line 1079 of file G4VXTRenergyLoss.cc.

1080 {
1081  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1082  const G4Material* mat = (*theMaterialTable)[fMatIndex1];
1084 
1085  return;
1086 }
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:229
G4SandiaTable * fPlatePhotoAbsCof

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 555 of file G4VXTRenergyLoss.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 1590 of file G4VXTRenergyLoss.cc.

1593 {
1594  G4double x1, x2, y1, y2, result;
1595 
1596  if(iTransfer == 0)
1597  {
1598  result = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1599  }
1600  else
1601  {
1602  y1 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer-1);
1603  y2 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer);
1604 
1605  x1 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
1606  x2 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1607 
1608  if ( x1 == x2 ) result = x2;
1609  else
1610  {
1611  if ( y1 == y2 ) result = x1 + (x2 - x1)*G4UniformRand();
1612  else
1613  {
1614  result = x1 + (position - y1)*(x2 - x1)/(y2 - y1);
1615  }
1616  }
1617  }
1618  return result;
1619 }
G4double G4ParticleHPJENDLHEData::G4double result
#define G4UniformRand()
Definition: Randomize.hh:97
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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

Definition at line 1310 of file G4VXTRenergyLoss.cc.

1311 {
1312  G4double CrossSection = 0.0;
1313  if ( Z < 0.9999 ) return CrossSection;
1314  if ( GammaEnergy < 0.1*keV ) return CrossSection;
1315  if ( GammaEnergy > (100.*GeV/Z) ) return CrossSection;
1316 
1317  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
1318 
1319  static const G4double
1320  d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527 *barn, d4=-1.9798e+1*barn,
1321  e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
1322  f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
1323 
1324  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
1325  p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
1326 
1327  G4double T0 = 15.0*keV;
1328  if (Z < 1.5) T0 = 40.0*keV;
1329 
1330  G4double X = std::max(GammaEnergy, T0) / electron_mass_c2;
1331  CrossSection = p1Z*std::log(1.+2.*X)/X
1332  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
1333 
1334  // modification for low energy. (special case for Hydrogen)
1335 
1336  if (GammaEnergy < T0)
1337  {
1338  G4double dT0 = 1.*keV;
1339  X = (T0+dT0) / electron_mass_c2;
1340  G4double sigma = p1Z*std::log(1.+2.*X)/X
1341  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
1342  G4double c1 = -T0*(sigma-CrossSection)/(CrossSection*dT0);
1343  G4double c2 = 0.150;
1344  if (Z > 1.5) c2 = 0.375-0.0556*std::log(Z);
1345  G4double y = std::log(GammaEnergy/T0);
1346  CrossSection *= std::exp(-y*(c1+c2*y));
1347  }
1348  // G4cout << "e= " << GammaEnergy << " Z= " << Z << " cross= " << CrossSection << G4endl;
1349  return CrossSection;
1350 }
static const G4double d2
static constexpr double electron_mass_c2
static const G4double d1
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static constexpr double GeV
Definition: G4SIunits.hh:217
static const G4double T0[78]
double G4double
Definition: G4Types.hh:76
static constexpr double barn
Definition: G4SIunits.hh:105
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4VXTRenergyLoss::GetEnergy ( )
inline

Definition at line 165 of file G4VXTRenergyLoss.hh.

165 {return fEnergy;};
G4double G4VXTRenergyLoss::GetGamma ( )
inline

Definition at line 164 of file G4VXTRenergyLoss.hh.

164 {return fGamma;};
G4complex G4VXTRenergyLoss::GetGasComplexFZ ( G4double  omega,
G4double  gamma,
G4double  varAngle 
)

Definition at line 1130 of file G4VXTRenergyLoss.cc.

1133 {
1134  G4double cof, length,delta, real_v, image_v;
1135 
1136  length = 0.5*GetGasFormationZone(omega,gamma,varAngle);
1137  delta = length*GetGasLinearPhotoAbs(omega);
1138  cof = 1.0/(1.0 + delta*delta);
1139 
1140  real_v = length*cof;
1141  image_v = real_v*delta;
1142 
1143  G4complex zone(real_v,image_v);
1144  return zone;
1145 }
G4double GetGasFormationZone(G4double, G4double, G4double)
G4double GetGasLinearPhotoAbs(G4double)
std::complex< G4double > G4complex
Definition: G4Types.hh:81
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4VXTRenergyLoss::GetGasCompton ( G4double  omega)

Definition at line 1286 of file G4VXTRenergyLoss.cc.

1287 {
1288  G4int i, numberOfElements;
1289  G4double xSection = 0., nowZ, sumZ = 0.;
1290 
1291  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1292  numberOfElements = (*theMaterialTable)[fMatIndex2]->GetNumberOfElements();
1293 
1294  for( i = 0; i < numberOfElements; i++ )
1295  {
1296  nowZ = (*theMaterialTable)[fMatIndex2]->GetElement(i)->GetZ();
1297  sumZ += nowZ;
1298  xSection += GetComptonPerAtom(omega,nowZ); // *nowZ;
1299  }
1300  xSection /= sumZ;
1301  xSection *= (*theMaterialTable)[fMatIndex2]->GetElectronDensity();
1302  return xSection;
1303 }
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
int G4int
Definition: G4Types.hh:78
G4double GetComptonPerAtom(G4double, G4double)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 1115 of file G4VXTRenergyLoss.cc.

1118 {
1119  G4double cof, lambda;
1120  lambda = 1.0/gamma/gamma + varAngle + fSigma2/omega/omega;
1121  cof = 2.0*hbarc/omega/lambda;
1122  return cof;
1123 }
static constexpr double hbarc
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

G4double G4VXTRenergyLoss::GetGasLinearPhotoAbs ( G4double  omega)

Definition at line 1167 of file G4VXTRenergyLoss.cc.

1168 {
1169  G4double omega2, omega3, omega4;
1170 
1171  omega2 = omega*omega;
1172  omega3 = omega2*omega;
1173  omega4 = omega2*omega2;
1174 
1175  const G4double* SandiaCof = fGasPhotoAbsCof->GetSandiaCofForMaterial(omega);
1176  G4double cross = SandiaCof[0]/omega + SandiaCof[1]/omega2 +
1177  SandiaCof[2]/omega3 + SandiaCof[3]/omega4;
1178  return cross;
1179 
1180 }
G4double GetSandiaCofForMaterial(G4int, G4int) const
G4SandiaTable * fGasPhotoAbsCof
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4VXTRenergyLoss::GetGasZmuProduct ( )

Definition at line 1237 of file G4VXTRenergyLoss.cc.

1238 {
1239  std::ofstream outGas("gasZmu.dat", std::ios::out );
1240  outGas.setf( std::ios::scientific, std::ios::floatfield );
1241  G4int i;
1242  G4double omega, varAngle, gamma;
1243  gamma = 10000.;
1244  varAngle = 1/gamma/gamma;
1245  if(verboseLevel > 0)
1246  G4cout<<"energy, keV"<<"\t"<<"Zmu for gas"<<G4endl;
1247  for(i=0;i<100;i++)
1248  {
1249  omega = (1.0 + i)*keV;
1250  if(verboseLevel > 1)
1251  G4cout<<omega/keV<<"\t"<<GetGasZmuProduct(omega,gamma,varAngle)<<"\t";
1252  if(verboseLevel > 0)
1253  outGas<<omega/keV<<"\t\t"<<GetGasZmuProduct(omega,gamma,varAngle)<<G4endl;
1254  }
1255  return;
1256 }
G4int verboseLevel
Definition: G4VProcess.hh:368
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static constexpr double keV
Definition: G4SIunits.hh:216
G4double G4VXTRenergyLoss::GetGasZmuProduct ( G4double  omega,
G4double  gamma,
G4double  varAngle 
)

Definition at line 1226 of file G4VXTRenergyLoss.cc.

1229 {
1230  return GetGasFormationZone(omega,gamma,varAngle)*GetGasLinearPhotoAbs(omega);
1231 }
G4double GetGasFormationZone(G4double, G4double, G4double)
G4double GetGasLinearPhotoAbs(G4double)

Here is the call graph for this function:

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

Implements G4VDiscreteProcess.

Definition at line 201 of file G4VXTRenergyLoss.cc.

204 {
205  G4int iTkin, iPlace;
206  G4double lambda, sigma, kinEnergy, mass, gamma;
207  G4double charge, chargeSq, massRatio, TkinScaled;
208  G4double E1,E2,W,W1,W2;
209 
210  *condition = NotForced;
211 
212  if( aTrack.GetVolume()->GetLogicalVolume() != fEnvelope ) lambda = DBL_MAX;
213  else
214  {
215  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
216  kinEnergy = aParticle->GetKineticEnergy();
217  mass = aParticle->GetDefinition()->GetPDGMass();
218  gamma = 1.0 + kinEnergy/mass;
219  if(verboseLevel > 1)
220  {
221  G4cout<<" gamma = "<<gamma<<"; fGamma = "<<fGamma<<G4endl;
222  }
223 
224  if ( std::fabs( gamma - fGamma ) < 0.05*gamma ) 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)) break;
235  }
236  iPlace = iTkin - 1;
237 
238  if(iTkin == 0) lambda = DBL_MAX; // Tkin is too small, neglect of TR photon generation
239  else // general case: Tkin between two vectors of the material
240  {
241  if(iTkin == fTotBin)
242  {
243  sigma = (*(*fEnergyDistrTable)(iPlace))(0)*chargeSq;
244  }
245  else
246  {
247  E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1);
249  W = 1.0/(E2 - E1);
250  W1 = (E2 - TkinScaled)*W;
251  W2 = (TkinScaled - E1)*W;
252  sigma = ( (*(*fEnergyDistrTable)(iPlace ))(0)*W1 +
253  (*(*fEnergyDistrTable)(iPlace+1))(0)*W2 )*chargeSq;
254 
255  }
256  if (sigma < DBL_MIN) lambda = DBL_MAX;
257  else lambda = 1./sigma;
258  fLambda = lambda;
259  fGamma = gamma;
260  if(verboseLevel > 1)
261  {
262  G4cout<<" lambda = "<<lambda/mm<<" mm"<<G4endl;
263  }
264  }
265  }
266  }
267  return lambda;
268 }
G4double condition(const G4ErrorSymMatrix &m)
G4PhysicsTable * fEnergyDistrTable
static constexpr double mm
Definition: G4SIunits.hh:115
G4int verboseLevel
Definition: G4VProcess.hh:368
static constexpr double proton_mass_c2
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4LogicalVolume * fEnvelope
G4ParticleDefinition * GetDefinition() const
G4double GetLowEdgeEnergy(size_t binNumber) const
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4PhysicsLogVector * fProtonEnergyVector
G4LogicalVolume * GetLogicalVolume() const
G4double GetPDGMass() const
#define DBL_MIN
Definition: templates.hh:75
G4VPhysicalVolume * GetVolume() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

void G4VXTRenergyLoss::GetNumberOfPhotons ( )

Definition at line 1442 of file G4VXTRenergyLoss.cc.

1443 {
1444  G4int iTkin;
1445  G4double gamma, numberE;
1446 
1447  std::ofstream outEn("numberE.dat", std::ios::out );
1448  outEn.setf( std::ios::scientific, std::ios::floatfield );
1449 
1450  std::ofstream outAng("numberAng.dat", std::ios::out );
1451  outAng.setf( std::ios::scientific, std::ios::floatfield );
1452 
1453  for(iTkin=0;iTkin<fTotBin;iTkin++) // Lorentz factor loop
1454  {
1455  gamma = 1.0 + (fProtonEnergyVector->
1456  GetLowEdgeEnergy(iTkin)/proton_mass_c2);
1457  numberE = (*(*fEnergyDistrTable)(iTkin))(0);
1458  // numberA = (*(*fAngleDistrTable)(iTkin))(0);
1459  if(verboseLevel > 1)
1460  G4cout<<gamma<<"\t\t"<<numberE<<"\t" // <<numberA
1461  <<G4endl;
1462  if(verboseLevel > 0)
1463  outEn<<gamma<<"\t\t"<<numberE<<G4endl;
1464  }
1465  return;
1466 }
G4int verboseLevel
Definition: G4VProcess.hh:368
static constexpr double proton_mass_c2
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
G4PhysicsLogVector * fProtonEnergyVector
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4complex G4VXTRenergyLoss::GetPlateComplexFZ ( G4double  omega,
G4double  gamma,
G4double  varAngle 
)

Definition at line 1057 of file G4VXTRenergyLoss.cc.

1060 {
1061  G4double cof, length,delta, real_v, image_v;
1062 
1063  length = 0.5*GetPlateFormationZone(omega,gamma,varAngle);
1064  delta = length*GetPlateLinearPhotoAbs(omega);
1065  cof = 1.0/(1.0 + delta*delta);
1066 
1067  real_v = length*cof;
1068  image_v = real_v*delta;
1069 
1070  G4complex zone(real_v,image_v);
1071  return zone;
1072 }
G4double GetPlateLinearPhotoAbs(G4double)
G4double GetPlateFormationZone(G4double, G4double, G4double)
std::complex< G4double > G4complex
Definition: G4Types.hh:81
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4VXTRenergyLoss::GetPlateCompton ( G4double  omega)

Definition at line 1262 of file G4VXTRenergyLoss.cc.

1263 {
1264  G4int i, numberOfElements;
1265  G4double xSection = 0., nowZ, sumZ = 0.;
1266 
1267  const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
1268  numberOfElements = (*theMaterialTable)[fMatIndex1]->GetNumberOfElements();
1269 
1270  for( i = 0; i < numberOfElements; i++ )
1271  {
1272  nowZ = (*theMaterialTable)[fMatIndex1]->GetElement(i)->GetZ();
1273  sumZ += nowZ;
1274  xSection += GetComptonPerAtom(omega,nowZ); // *nowZ;
1275  }
1276  xSection /= sumZ;
1277  xSection *= (*theMaterialTable)[fMatIndex1]->GetElectronDensity();
1278  return xSection;
1279 }
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
int G4int
Definition: G4Types.hh:78
G4double GetComptonPerAtom(G4double, G4double)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 1043 of file G4VXTRenergyLoss.cc.

1046 {
1047  G4double cof, lambda;
1048  lambda = 1.0/gamma/gamma + varAngle + fSigma1/omega/omega;
1049  cof = 2.0*hbarc/omega/lambda;
1050  return cof;
1051 }
static constexpr double hbarc
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

G4double G4VXTRenergyLoss::GetPlateLinearPhotoAbs ( G4double  omega)

Definition at line 1095 of file G4VXTRenergyLoss.cc.

1096 {
1097  // G4int i;
1098  G4double omega2, omega3, omega4;
1099 
1100  omega2 = omega*omega;
1101  omega3 = omega2*omega;
1102  omega4 = omega2*omega2;
1103 
1104  const G4double* SandiaCof = fPlatePhotoAbsCof->GetSandiaCofForMaterial(omega);
1105  G4double cross = SandiaCof[0]/omega + SandiaCof[1]/omega2 +
1106  SandiaCof[2]/omega3 + SandiaCof[3]/omega4;
1107  return cross;
1108 }
G4double GetSandiaCofForMaterial(G4int, G4int) const
G4SandiaTable * fPlatePhotoAbsCof
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4VXTRenergyLoss::GetPlateZmuProduct ( )

Definition at line 1199 of file G4VXTRenergyLoss.cc.

1200 {
1201  std::ofstream outPlate("plateZmu.dat", std::ios::out );
1202  outPlate.setf( std::ios::scientific, std::ios::floatfield );
1203 
1204  G4int i;
1205  G4double omega, varAngle, gamma;
1206  gamma = 10000.;
1207  varAngle = 1/gamma/gamma;
1208  if(verboseLevel > 0)
1209  G4cout<<"energy, keV"<<"\t"<<"Zmu for plate"<<G4endl;
1210  for(i=0;i<100;i++)
1211  {
1212  omega = (1.0 + i)*keV;
1213  if(verboseLevel > 1)
1214  G4cout<<omega/keV<<"\t"<<GetPlateZmuProduct(omega,gamma,varAngle)<<"\t";
1215  if(verboseLevel > 0)
1216  outPlate<<omega/keV<<"\t\t"<<GetPlateZmuProduct(omega,gamma,varAngle)<<G4endl;
1217  }
1218  return;
1219 }
G4int verboseLevel
Definition: G4VProcess.hh:368
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static constexpr double keV
Definition: G4SIunits.hh:216
G4double G4VXTRenergyLoss::GetPlateZmuProduct ( G4double  omega,
G4double  gamma,
G4double  varAngle 
)

Definition at line 1187 of file G4VXTRenergyLoss.cc.

1190 {
1191  return GetPlateFormationZone(omega,gamma,varAngle)
1192  * GetPlateLinearPhotoAbs(omega);
1193 }
G4double GetPlateLinearPhotoAbs(G4double)
G4double GetPlateFormationZone(G4double, G4double, G4double)

Here is the call graph for this function:

G4PhysicsLogVector* G4VXTRenergyLoss::GetProtonVector ( )
inline

Definition at line 174 of file G4VXTRenergyLoss.hh.

174 { return fProtonEnergyVector;};
G4PhysicsLogVector * fProtonEnergyVector
G4double G4VXTRenergyLoss::GetRandomAngle ( G4double  energyXTR,
G4int  iTkin 
)

Definition at line 1560 of file G4VXTRenergyLoss.cc.

1561 {
1562  G4int iTR, iAngle;
1564 
1565  if (iTkin == fTotBin) iTkin--;
1566 
1568 
1569  for( iTR = 0; iTR < fBinTR; iTR++ )
1570  {
1571  if( energyXTR < fXTREnergyVector->GetLowEdgeEnergy(iTR) ) break;
1572  }
1573  if (iTR == fBinTR) iTR--;
1574 
1575  position = (*(*fAngleForEnergyTable)(iTR))(0)*G4UniformRand();
1576 
1577  for(iAngle = 0;;iAngle++)
1578  {
1579  if(position >= (*(*fAngleForEnergyTable)(iTR))(iAngle)) break;
1580  }
1581  angle = GetAngleXTR(iTR,position,iAngle);
1582  return angle;
1583 }
static G4double angle[DIM]
int G4int
Definition: G4Types.hh:78
#define position
Definition: xmlparse.cc:622
#define G4UniformRand()
Definition: Randomize.hh:97
G4double GetAngleXTR(G4int iTR, G4double position, G4int iAngle)
std::vector< G4PhysicsTable * > fAngleBank
G4PhysicsTable * fAngleForEnergyTable
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

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

Definition at line 1381 of file G4VXTRenergyLoss.cc.

1383 {
1384  // return stack factor corresponding to one interface
1385 
1386  return std::real( OneInterfaceXTRdEdx(energy,gamma,varAngle) );
1387 }
G4double energy(const ThreeVector &p, const G4double m)
G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma, G4double varAngle)

Here is the call graph for this function:

Here is the caller graph for this function:

G4int G4VXTRenergyLoss::GetTotBin ( )
inline

Definition at line 175 of file G4VXTRenergyLoss.hh.

175 {return fTotBin;};
G4double G4VXTRenergyLoss::GetVarAngle ( )
inline

Definition at line 166 of file G4VXTRenergyLoss.hh.

166 {return fVarAngle;};
G4double G4VXTRenergyLoss::GetXTRenergy ( G4int  iPlace,
G4double  position,
G4int  iTransfer 
)

Definition at line 1524 of file G4VXTRenergyLoss.cc.

1527 {
1528  G4double x1, x2, y1, y2, result;
1529 
1530  if(iTransfer == 0)
1531  {
1532  result = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1533  }
1534  else
1535  {
1536  y1 = (*(*fEnergyDistrTable)(iPlace))(iTransfer-1);
1537  y2 = (*(*fEnergyDistrTable)(iPlace))(iTransfer);
1538 
1539  x1 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
1540  x2 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1541 
1542  if ( x1 == x2 ) result = x2;
1543  else
1544  {
1545  if ( y1 == y2 ) result = x1 + (x2 - x1)*G4UniformRand();
1546  else
1547  {
1548  // result = x1 + (position - y1)*(x2 - x1)/(y2 - y1);
1549  result = x1 + (x2 - x1)*G4UniformRand();
1550  }
1551  }
1552  }
1553  return result;
1554 }
G4double G4ParticleHPJENDLHEData::G4double result
#define G4UniformRand()
Definition: Randomize.hh:97
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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

Definition at line 1473 of file G4VXTRenergyLoss.cc.

1474 {
1475  G4int iTransfer, iPlace;
1476  G4double transfer = 0.0, position, E1, E2, W1, W2, W;
1477 
1478  iPlace = iTkin - 1;
1479 
1480  // G4cout<<"iPlace = "<<iPlace<<endl;
1481 
1482  if(iTkin == fTotBin) // relativistic plato, try from left
1483  {
1484  position = (*(*fEnergyDistrTable)(iPlace))(0)*G4UniformRand();
1485 
1486  for(iTransfer=0;;iTransfer++)
1487  {
1488  if(position >= (*(*fEnergyDistrTable)(iPlace))(iTransfer)) break;
1489  }
1490  transfer = GetXTRenergy(iPlace,position,iTransfer);
1491  }
1492  else
1493  {
1494  E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1);
1496  W = 1.0/(E2 - E1);
1497  W1 = (E2 - scaledTkin)*W;
1498  W2 = (scaledTkin - E1)*W;
1499 
1500  position =( (*(*fEnergyDistrTable)(iPlace))(0)*W1 +
1501  (*(*fEnergyDistrTable)(iPlace+1))(0)*W2 )*G4UniformRand();
1502 
1503  // G4cout<<position<<"\t";
1504 
1505  for(iTransfer=0;;iTransfer++)
1506  {
1507  if( position >=
1508  ( (*(*fEnergyDistrTable)(iPlace))(iTransfer)*W1 +
1509  (*(*fEnergyDistrTable)(iPlace+1))(iTransfer)*W2) ) break;
1510  }
1511  transfer = GetXTRenergy(iPlace,position,iTransfer);
1512 
1513  }
1514  // G4cout<<"XTR transfer = "<<transfer/keV<<" keV"<<endl;
1515  if(transfer < 0.0 ) transfer = 0.0;
1516  return transfer;
1517 }
G4PhysicsTable * fEnergyDistrTable
G4double GetLowEdgeEnergy(size_t binNumber) const
int G4int
Definition: G4Types.hh:78
#define position
Definition: xmlparse.cc:622
#define G4UniformRand()
Definition: Randomize.hh:97
G4PhysicsLogVector * fProtonEnergyVector
G4double GetXTRenergy(G4int iPlace, G4double position, G4int iTransfer)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4VXTRenergyLoss::IsApplicable ( const G4ParticleDefinition particle)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 192 of file G4VXTRenergyLoss.cc.

193 {
194  return ( particle.GetPDGCharge() != 0.0 );
195 }
G4double GetPDGCharge() const

Here is the call graph for this function:

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

Definition at line 1364 of file G4VXTRenergyLoss.cc.

1366 {
1367  G4double formationLength1, formationLength2;
1368  formationLength1 = 1.0/
1369  (1.0/(gamma*gamma)
1370  + fSigma1/(energy*energy)
1371  + varAngle);
1372  formationLength2 = 1.0/
1373  (1.0/(gamma*gamma)
1374  + fSigma2/(energy*energy)
1375  + varAngle);
1376  return (varAngle/energy)*(formationLength1 - formationLength2)
1377  *(formationLength1 - formationLength2);
1378 
1379 }
G4double energy(const ThreeVector &p, const G4double m)
double G4double
Definition: G4Types.hh:76

Here is the caller graph for this function:

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

Definition at line 868 of file G4VXTRenergyLoss.cc.

871 {
872  G4complex Z1 = GetPlateComplexFZ(energy,gamma,varAngle);
873  G4complex Z2 = GetGasComplexFZ(energy,gamma,varAngle);
874 
875  G4complex zOut = (Z1 - Z2)*(Z1 - Z2)
876  * (varAngle*energy/hbarc/hbarc);
877  return zOut;
878 
879 }
static constexpr double hbarc
G4complex GetPlateComplexFZ(G4double, G4double, G4double)
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double energy(const ThreeVector &p, const G4double m)
G4complex GetGasComplexFZ(G4double, G4double, G4double)

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VDiscreteProcess.

Definition at line 728 of file G4VXTRenergyLoss.cc.

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

Here is the call graph for this function:

void G4VXTRenergyLoss::SetAngleRadDistr ( G4bool  pAngleRadDistr)
inline

Definition at line 171 of file G4VXTRenergyLoss.hh.

171 {fAngleRadDistr=pAngleRadDistr;};
void G4VXTRenergyLoss::SetCompton ( G4bool  pC)
inline

Definition at line 172 of file G4VXTRenergyLoss.hh.

172 {fCompton=pC;};
void G4VXTRenergyLoss::SetEnergy ( G4double  energy)
inline

Definition at line 169 of file G4VXTRenergyLoss.hh.

169 {fEnergy = energy;};
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

void G4VXTRenergyLoss::SetGamma ( G4double  gamma)
inline

Definition at line 168 of file G4VXTRenergyLoss.hh.

168 {fGamma = gamma;};
void G4VXTRenergyLoss::SetVarAngle ( G4double  varAngle)
inline

Definition at line 170 of file G4VXTRenergyLoss.hh.

170 {fVarAngle = varAngle;};
G4double G4VXTRenergyLoss::SpectralAngleXTRdEdx ( G4double  varAngle)

Definition at line 887 of file G4VXTRenergyLoss.cc.

888 {
890  if(result < 0.0) result = 0.0;
891  return result;
892 }
G4double G4ParticleHPJENDLHEData::G4double result
double G4double
Definition: G4Types.hh:76
virtual G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4VXTRenergyLoss::SpectralXTRdEdx ( G4double  energy)
virtual

Reimplemented in G4TransparentRegXTRadiator, G4RegularXTRadiator, G4XTRTransparentRegRadModel, and G4XTRRegularRadModel.

Definition at line 898 of file G4VXTRenergyLoss.cc.

899 {
900  G4int i, iMax = 8;
901  G4double angleSum = 0.0;
902 
903  G4double lim[8] = { 0.0, 0.01, 0.02, 0.05, 0.1, 0.2, 0.5, 1.0 };
904 
905  for( i = 0; i < iMax; i++ ) lim[i] *= fMaxThetaTR;
906 
908 
909  fEnergy = energy;
910  /*
911  if( fAngleRadDistr && ( fEnergy == fEnergyForAngle ) )
912  {
913  fAngleVector ->PutValue(fBinTR - 1, angleSum);
914 
915  for( i = fBinTR - 2; i >= 0; i-- )
916  {
917 
918  angleSum += integral.Legendre10(
919  this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
920  fAngleVector->GetLowEdgeEnergy(i),
921  fAngleVector->GetLowEdgeEnergy(i+1) );
922 
923  fAngleVector ->PutValue(i, angleSum);
924  }
925  }
926  else
927  */
928  {
929  for( i = 0; i < iMax-1; i++ )
930  {
931  angleSum += integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
932  lim[i],lim[i+1]);
933  // result += integral.Legendre10(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
934  // lim[i],lim[i+1]);
935  }
936  }
937  return angleSum;
938 }
G4double Legendre96(T &typeT, F f, G4double a, G4double b)
int G4int
Definition: G4Types.hh:78
G4double energy(const ThreeVector &p, const G4double m)
G4double SpectralAngleXTRdEdx(G4double varAngle)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4VXTRenergyLoss::XTRNAngleDensity ( G4double  varAngle)

Definition at line 1429 of file G4VXTRenergyLoss.cc.

1430 {
1431  fVarAngle = varAngle;
1435 }
G4double Legendre96(T &typeT, F f, G4double a, G4double b)
G4double XTRNAngleSpectralDensity(G4double energy)

Here is the call graph for this function:

G4double G4VXTRenergyLoss::XTRNAngleSpectralDensity ( G4double  energy)

Definition at line 1419 of file G4VXTRenergyLoss.cc.

1420 {
1423 }
G4double OneBoundaryXTRNdensity(G4double energy, G4double gamma, G4double varAngle) const
G4double energy(const ThreeVector &p, const G4double m)
virtual G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4VXTRenergyLoss::XTRNSpectralAngleDensity ( G4double  varAngle)

Definition at line 1394 of file G4VXTRenergyLoss.cc.

1395 {
1396  return OneBoundaryXTRNdensity(fEnergy,fGamma,varAngle)*
1397  GetStackFactor(fEnergy,fGamma,varAngle);
1398 }
G4double OneBoundaryXTRNdensity(G4double energy, G4double gamma, G4double varAngle) const
virtual G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle)

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4VXTRenergyLoss::XTRNSpectralDensity ( G4double  energy)

Definition at line 1404 of file G4VXTRenergyLoss.cc.

1405 {
1406  fEnergy = energy;
1409  0.0,0.2*fMaxThetaTR) +
1411  0.2*fMaxThetaTR,fMaxThetaTR);
1412 }
G4double Legendre10(T &typeT, F f, G4double a, G4double b)
G4double Legendre96(T &typeT, F f, G4double a, G4double b)
G4double XTRNSpectralAngleDensity(G4double varAngle)
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

Member Data Documentation

G4double G4VXTRenergyLoss::fAlphaGas
protected

Definition at line 226 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fAlphaPlate
protected

Definition at line 225 of file G4VXTRenergyLoss.hh.

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

Definition at line 235 of file G4VXTRenergyLoss.hh.

G4PhysicsTable* G4VXTRenergyLoss::fAngleDistrTable
protected

Definition at line 186 of file G4VXTRenergyLoss.hh.

G4PhysicsTable* G4VXTRenergyLoss::fAngleForEnergyTable
protected

Definition at line 234 of file G4VXTRenergyLoss.hh.

G4bool G4VXTRenergyLoss::fAngleRadDistr
protected

Definition at line 213 of file G4VXTRenergyLoss.hh.

G4int G4VXTRenergyLoss::fBinTR
protected

Definition at line 199 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fCofTR
protected

Definition at line 210 of file G4VXTRenergyLoss.hh.

G4bool G4VXTRenergyLoss::fCompton
protected

Definition at line 214 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fEnergy
protected

Definition at line 205 of file G4VXTRenergyLoss.hh.

G4PhysicsTable* G4VXTRenergyLoss::fEnergyDistrTable
protected

Definition at line 187 of file G4VXTRenergyLoss.hh.

G4LogicalVolume* G4VXTRenergyLoss::fEnvelope
protected

Definition at line 185 of file G4VXTRenergyLoss.hh.

G4bool G4VXTRenergyLoss::fExitFlux
protected

Definition at line 212 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fGamma
protected

Definition at line 204 of file G4VXTRenergyLoss.hh.

G4double* G4VXTRenergyLoss::fGammaCutInKineticEnergy
protected

Definition at line 182 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fGammaTkinCut
protected

Definition at line 184 of file G4VXTRenergyLoss.hh.

G4SandiaTable* G4VXTRenergyLoss::fGasPhotoAbsCof
protected

Definition at line 230 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fGasThick
protected

Definition at line 224 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fLambda
protected

Definition at line 207 of file G4VXTRenergyLoss.hh.

G4int G4VXTRenergyLoss::fMatIndex1
protected

Definition at line 218 of file G4VXTRenergyLoss.hh.

G4int G4VXTRenergyLoss::fMatIndex2
protected

Definition at line 219 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fMaxEnergyTR
protected

Definition at line 195 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fMaxProtonTkin
protected

Definition at line 202 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fMaxThetaTR
protected

Definition at line 198 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fMinEnergyTR
protected

Definition at line 194 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fMinProtonTkin
protected

Definition at line 201 of file G4VXTRenergyLoss.hh.

G4ParticleChange G4VXTRenergyLoss::fParticleChange
protected

Definition at line 232 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fPlasmaCof
protected

Definition at line 209 of file G4VXTRenergyLoss.hh.

G4int G4VXTRenergyLoss::fPlateNumber
protected

Definition at line 220 of file G4VXTRenergyLoss.hh.

G4SandiaTable* G4VXTRenergyLoss::fPlatePhotoAbsCof
protected

Definition at line 228 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fPlateThick
protected

Definition at line 223 of file G4VXTRenergyLoss.hh.

G4PhysicsLogVector* G4VXTRenergyLoss::fProtonEnergyVector
protected

Definition at line 189 of file G4VXTRenergyLoss.hh.

G4ParticleDefinition* G4VXTRenergyLoss::fPtrGamma
protected

Definition at line 180 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fSigma1
protected

Definition at line 215 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fSigma2
protected

Definition at line 216 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fTheMaxAngle
protected

Definition at line 196 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fTheMaxEnergyTR
protected

Definition at line 193 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fTheMinAngle
protected

Definition at line 197 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fTheMinEnergyTR
protected

Definition at line 192 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fTotalDist
protected

Definition at line 222 of file G4VXTRenergyLoss.hh.

G4int G4VXTRenergyLoss::fTotBin
protected

Definition at line 203 of file G4VXTRenergyLoss.hh.

G4double G4VXTRenergyLoss::fVarAngle
protected

Definition at line 206 of file G4VXTRenergyLoss.hh.

G4PhysicsLogVector* G4VXTRenergyLoss::fXTREnergyVector
protected

Definition at line 190 of file G4VXTRenergyLoss.hh.


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