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

#include <G4XTRTransparentRegRadModel.hh>

Inheritance diagram for G4XTRTransparentRegRadModel:
Collaboration diagram for G4XTRTransparentRegRadModel:

Public Member Functions

 G4XTRTransparentRegRadModel (G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="XTRTransparentRegRadModel")
 
 ~G4XTRTransparentRegRadModel ()
 
G4double SpectralXTRdEdx (G4double energy) override
 
G4double GetStackFactor (G4double energy, G4double gamma, G4double varAngle) override
 
- Public Member Functions inherited from G4VXTRenergyLoss
 G4VXTRenergyLoss (G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="XTRenergyLoss", G4ProcessType type=fElectromagnetic)
 
virtual ~G4VXTRenergyLoss ()
 
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)
 
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 &)
 

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 ()
 
- Protected Attributes inherited from G4VXTRenergyLoss
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
 

Detailed Description

Definition at line 48 of file G4XTRTransparentRegRadModel.hh.

Constructor & Destructor Documentation

G4XTRTransparentRegRadModel::G4XTRTransparentRegRadModel ( G4LogicalVolume anEnvelope,
G4Material foilMat,
G4Material gasMat,
G4double  a,
G4double  b,
G4int  n,
const G4String processName = "XTRTransparentRegRadModel" 
)
explicit

Definition at line 40 of file G4XTRTransparentRegRadModel.cc.

43  :
44  G4VXTRenergyLoss(anEnvelope,foilMat,gasMat,a,b,n,processName)
45 {
46  G4cout<<"Regular transparent X-ray TR radiator EM process is called"<<G4endl;
47 
48  // Build energy and angular integral spectra of X-ray TR photons from
49  // a radiator
50  fExitFlux = true;
51  fAlphaPlate = 10000;
52  fAlphaGas = 1000;
53 
54  // BuildTable();
55 }
G4GLOB_DLL std::ostream G4cout
G4VXTRenergyLoss(G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="XTRenergyLoss", G4ProcessType type=fElectromagnetic)
#define G4endl
Definition: G4ios.hh:61
G4XTRTransparentRegRadModel::~G4XTRTransparentRegRadModel ( )

Definition at line 59 of file G4XTRTransparentRegRadModel.cc.

60 {
61  ;
62 }

Member Function Documentation

G4double G4XTRTransparentRegRadModel::GetStackFactor ( G4double  energy,
G4double  gamma,
G4double  varAngle 
)
overridevirtual

Reimplemented from G4VXTRenergyLoss.

Definition at line 145 of file G4XTRTransparentRegRadModel.cc.

147 {
148  /*
149  G4double result, Za, Zb, Ma, Mb, sigma;
150 
151  Za = GetPlateFormationZone(energy,gamma,varAngle);
152  Zb = GetGasFormationZone(energy,gamma,varAngle);
153  Ma = GetPlateLinearPhotoAbs(energy);
154  Mb = GetGasLinearPhotoAbs(energy);
155  sigma = Ma*fPlateThick + Mb*fGasThick;
156 
157  G4complex Ca(1.0+0.5*fPlateThick*Ma/fAlphaPlate,fPlateThick/Za/fAlphaPlate);
158  G4complex Cb(1.0+0.5*fGasThick*Mb/fAlphaGas,fGasThick/Zb/fAlphaGas);
159 
160  G4complex Ha = std::pow(Ca,-fAlphaPlate);
161  G4complex Hb = std::pow(Cb,-fAlphaGas);
162  G4complex H = Ha*Hb;
163  G4complex F1 = (1.0 - Ha)*(1.0 - Hb )/(1.0 - H)
164  * G4double(fPlateNumber) ;
165  G4complex F2 = (1.0-Ha)*(1.0-Ha)*Hb/(1.0-H)/(1.0-H)
166  * (1.0 - std::exp(-0.5*fPlateNumber*sigma)) ;
167  // *(1.0 - std::pow(H,fPlateNumber)) ;
168  G4complex R = (F1 + F2)*OneInterfaceXTRdEdx(energy,gamma,varAngle);
169  // G4complex R = F2*OneInterfaceXTRdEdx(energy,gamma,varAngle);
170  result = 2.0*std::real(R);
171  return result;
172  */
173  // numerically unstable result
174 
175  G4double result, Qa, Qb, Q, aZa, bZb, aMa, bMb, D, sigma;
176 
177  aZa = fPlateThick/GetPlateFormationZone(energy,gamma,varAngle);
178  bZb = fGasThick/GetGasFormationZone(energy,gamma,varAngle);
181  sigma = aMa*fPlateThick + bMb*fGasThick;
182  Qa = std::exp(-0.5*aMa);
183  Qb = std::exp(-0.5*bMb);
184  Q = Qa*Qb;
185 
186  G4complex Ha( Qa*std::cos(aZa), -Qa*std::sin(aZa) );
187  G4complex Hb( Qb*std::cos(bZb), -Qb*std::sin(bZb) );
188  G4complex H = Ha*Hb;
189  G4complex Hs = conj(H);
190  D = 1.0 /( (1. - Q)*(1. - Q) +
191  4.*Q*std::sin(0.5*(aZa + bZb))*std::sin(0.5*(aZa + bZb)) );
192  G4complex F1 = (1.0 - Ha)*(1.0 - Hb)*(1.0 - Hs)
194  G4complex F2 = (1.0 - Ha)*(1.0 - Ha)*Hb*(1.0 - Hs)*(1.0 - Hs)
195  // * (1.0 - std::pow(H,fPlateNumber)) * D*D;
196  * (1.0 - std::exp(-0.5*fPlateNumber*sigma)) * D*D;
197  G4complex R = (F1 + F2)*OneInterfaceXTRdEdx(energy,gamma,varAngle);
198  result = 2.0*std::real(R);
199  return result;
200 
201 }
G4double G4ParticleHPJENDLHEData::G4double result
G4double GetGasFormationZone(G4double, G4double, G4double)
G4double GetPlateLinearPhotoAbs(G4double)
G4double GetGasLinearPhotoAbs(G4double)
static double Q[]
G4double GetPlateFormationZone(G4double, G4double, G4double)
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4double energy(const ThreeVector &p, const G4double m)
double D(double temp)
G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma, G4double varAngle)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4XTRTransparentRegRadModel::SpectralXTRdEdx ( G4double  energy)
overridevirtual

Reimplemented from G4VXTRenergyLoss.

Definition at line 68 of file G4XTRTransparentRegRadModel.cc.

69 {
70  G4double result, sum = 0., tmp, cof1, cof2, cofMin, cofPHC,aMa, bMb, sigma;
71  G4int k, kMax, kMin;
72 
75 
76  if(fCompton)
77  {
78  aMa += GetPlateCompton(energy);
79  bMb += GetGasCompton(energy);
80  }
81  aMa *= fPlateThick;
82  bMb *= fGasThick;
83 
84  sigma = aMa + bMb;
85 
86  cofPHC = 4.*pi*hbarc;
87  tmp = (fSigma1 - fSigma2)/cofPHC/energy;
88  cof1 = fPlateThick*tmp;
89  cof2 = fGasThick*tmp;
90 
93  cofMin /= cofPHC;
94 
95  // if (fGamma < 1200) kMin = G4int(cofMin); // 1200 ?
96  // else kMin = 1;
97 
98 
99  kMin = G4int(cofMin);
100  if (cofMin > kMin) kMin++;
101 
102  // tmp = (fPlateThick + fGasThick)*energy*fMaxThetaTR;
103  // tmp /= cofPHC;
104  // kMax = G4int(tmp);
105  // if(kMax < 0) kMax = 0;
106  // kMax += kMin;
107 
108 
109  kMax = kMin + 19; // 5; // 9; // kMin + G4int(tmp);
110 
111  // tmp /= fGamma;
112  // if( G4int(tmp) < kMin ) kMin = G4int(tmp);
113  // G4cout<<"kMin = "<<kMin<<"; kMax = "<<kMax<<G4endl;
114 
115  for( k = kMin; k <= kMax; k++ )
116  {
117  tmp = pi*fPlateThick*(k + cof2)/(fPlateThick + fGasThick);
118  result = (k - cof1)*(k - cof1)*(k + cof2)*(k + cof2);
119 
120  if( k == kMin && kMin == G4int(cofMin) )
121  {
122  sum += 0.5*std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
123  }
124  else
125  {
126  sum += std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
127  }
128  // G4cout<<"k = "<<k<<"; sum = "<<sum<<G4endl;
129  }
130  result = 4.*( cof1 + cof2 )*( cof1 + cof2 )*sum/energy;
131  result *= ( 1. - std::exp(-fPlateNumber*sigma) )/( 1. - std::exp(-sigma) );
132  return result;
133 }
G4double G4ParticleHPJENDLHEData::G4double result
G4double GetPlateLinearPhotoAbs(G4double)
G4double GetPlateCompton(G4double)
static constexpr double hbarc
G4double GetGasLinearPhotoAbs(G4double)
int G4int
Definition: G4Types.hh:78
G4double energy(const ThreeVector &p, const G4double m)
G4double GetGasCompton(G4double)
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:


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