60 G4int G4LivermorePolarizedComptonModel::maxZ = 99;
62 G4ShellData* G4LivermorePolarizedComptonModel::shellData = 0;
82 G4cout <<
"Livermore Polarized Compton is constructed " <<
G4endl;
88 fAtomDeexcitation = 0;
100 delete scatterFunctionData;
101 scatterFunctionData = 0;
102 for(
G4int i=0; i<maxZ; ++i) {
116 if (verboseLevel > 1)
117 G4cout <<
"Calling G4LivermorePolarizedComptonModel::Initialise()" <<
G4endl;
125 char* path = getenv(
"G4LEDATA");
132 for(
G4int i=0; i<numOfCouples; ++i) {
138 for (
G4int j=0; j<nelm; ++j) {
141 else if(Z > maxZ){ Z = maxZ; }
143 if( (!
data[Z]) ) { ReadData(Z, path); }
151 G4String file =
"/doppler/shell-doppler";
158 if(!scatterFunctionData)
162 G4String scatterFile =
"comp/ce-sf-";
164 scatterFunctionData->
LoadData(scatterFile);
170 if (verboseLevel > 2) {
174 if( verboseLevel>1 ) {
175 G4cout <<
"G4LivermoreComptonModel is initialized " << G4endl
182 if(isInitialised) {
return; }
186 isInitialised =
true;
198 void G4LivermorePolarizedComptonModel::ReadData(
size_t Z,
const char* path)
200 if (verboseLevel > 1)
202 G4cout <<
"G4LivermorePolarizedComptonModel::ReadData()"
205 if(
data[Z]) {
return; }
206 const char* datadir = path;
209 datadir = getenv(
"G4LEDATA");
212 G4Exception(
"G4LivermorePolarizedComptonModel::ReadData()",
214 "Environment variable G4LEDATA not defined");
222 data[
Z]->SetSpline(
false);
224 std::ostringstream ost;
225 ost << datadir <<
"/livermore/comp/ce-cs-" << Z <<
".dat";
226 std::ifstream fin(ost.str().c_str());
231 ed <<
"G4LivermorePolarizedComptonModel data file <" << ost.str().c_str()
232 <<
"> is not opened!" <<
G4endl;
235 ed,
"G4LEDATA version should be G4EMLOW6.34 or later");
238 if(verboseLevel > 3) {
239 G4cout <<
"File " << ost.str()
240 <<
" is opened by G4LivermorePolarizedComptonModel" <<
G4endl;
242 data[
Z]->Retrieve(fin,
true);
259 if (verboseLevel > 3)
260 G4cout <<
"Calling ComputeCrossSectionPerAtom() of G4LivermorePolarizedComptonModel" <<
G4endl;
268 if(intZ < 1 || intZ > maxZ) {
return cs; }
278 if(!pv) {
return cs; }
285 if(GammaEnergy <= e1) { cs = GammaEnergy/(e1*e1)*pv->
Value(e1); }
286 else if(GammaEnergy <= e2) { cs = pv->
Value(GammaEnergy)/GammaEnergy; }
287 else if(GammaEnergy > e2) { cs = pv->
Value(e2)/GammaEnergy; }
307 if (verboseLevel > 3)
308 G4cout <<
"Calling SampleSecondaries() of G4LivermorePolarizedComptonModel" <<
G4endl;
329 if(!(gammaPolarization0.
isOrthogonal(gammaDirection0, 1e-6))||(gammaPolarization0.
mag()==0))
331 gammaPolarization0 = GetRandomPolarization(gammaDirection0);
337 gammaPolarization0 = GetPerpendicularPolarization(gammaDirection0, gammaPolarization0);
355 G4double epsilon0Local = 1./(1. + 2*E0_m);
356 G4double epsilon0Sq = epsilon0Local*epsilon0Local;
357 G4double alpha1 = - std::log(epsilon0Local);
358 G4double alpha2 = 0.5*(1.- epsilon0Sq);
373 epsilon = std::sqrt(epsilonSq);
376 onecost = (1.-
epsilon)/(epsilon*E0_m);
377 sinThetaSqr = onecost*(2.-onecost);
380 if (sinThetaSqr > 1.)
383 <<
" -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
384 <<
"sin(theta)**2 = "
390 if (sinThetaSqr < 0.)
393 <<
" -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
394 <<
"sin(theta)**2 = "
402 G4double x = std::sqrt(onecost/2.) / (wlGamma/
cm);;
404 greject = (1. - epsilon*sinThetaSqr/(1.+ epsilonSq))*scatteringFunction;
413 G4double phi = SetPhi(epsilon,sinThetaSqr);
426 <<
" -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
436 <<
" -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
446 G4double sinTheta = std::sqrt (sinThetaSqr);
452 <<
" -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
462 <<
" -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
472 G4double dirx = sinTheta*std::cos(phi);
473 G4double diry = sinTheta*std::sin(phi);
486 static G4int maxDopplerIterations = 1000;
488 G4double photonEoriginal = epsilon * gammaEnergy0;
495 if (verboseLevel > 3) {
506 if (verboseLevel > 3) {
507 G4cout <<
"Shell ID= " << shellIdx
508 <<
" Ebind(keV)= " << bindingE/
keV <<
G4endl;
510 eMax = gammaEnergy0 - bindingE;
515 if (verboseLevel > 3) {
520 G4double pDoppler2 = pDoppler * pDoppler;
521 G4double var2 = 1. + onecost * E0_m;
522 G4double var3 = var2*var2 - pDoppler2;
523 G4double var4 = var2 - pDoppler2 * cosTheta;
524 G4double var = var4*var4 - var3 + pDoppler2 * var3;
528 G4double scale = gammaEnergy0 / var3;
530 if (
G4UniformRand() < 0.5) photonE = (var4 - varSqrt) * scale;
531 else photonE = (var4 + varSqrt) * scale;
537 }
while ( iteration <= maxDopplerIterations &&
538 (photonE < 0. || photonE > eMax || photonE < eMax*
G4UniformRand()) );
544 if (iteration >= maxDopplerIterations)
546 photonE = photonEoriginal;
550 gammaEnergy1 = photonE;
563 G4ThreeVector gammaPolarization1 = SetNewPolarization(epsilon,
570 gammaDirection1 = tmpDirection1;
574 SystemOfRefChange(gammaDirection0,gammaDirection1,
575 gammaPolarization0,gammaPolarization1);
577 if (gammaEnergy1 > 0.)
594 G4double ElecKineEnergy = gammaEnergy0 - gammaEnergy1 -bindingE;
598 if(ElecKineEnergy < 0.0) {
607 G4ThreeVector ElecDirection((gammaEnergy0 * gammaDirection0 -
608 gammaEnergy1 * gammaDirection1) * (1./ElecMomentum));
611 fvect->push_back(dp);
616 if (verboseLevel > 3) {
617 G4cout <<
"Started atomic de-excitation " << fAtomDeexcitation <<
G4endl;
620 if(fAtomDeexcitation && iteration < maxDopplerIterations) {
623 size_t nbefore = fvect->size();
627 size_t nafter = fvect->size();
628 if(nafter > nbefore) {
629 for (
size_t i=nbefore; i<nafter; ++i) {
631 if (bindingE >= ((*fvect)[i])->GetKineticEnergy())
634 bindingE -= ((*fvect)[i])->GetKineticEnergy();
649 G4Exception(
"G4LivermoreComptonModel::SampleSecondaries()",
675 b = energyRate + 1/energyRate;
677 phiProbability = 1 - (a/b)*(std::cos(phi)*std::cos(phi));
682 while ( rand2 > phiProbability );
698 return x < z ?
G4ThreeVector(-dy,dx,0) : G4ThreeVector(0,-dz,dy);
700 return y < z ?
G4ThreeVector(dz,0,-dx) : G4ThreeVector(-dy,dx,0);
719 c.
setX(std::cos(angle)*(a0.
x())+std::sin(angle)*b0.
x());
720 c.
setY(std::cos(angle)*(a0.
y())+std::sin(angle)*b0.
y());
721 c.
setZ(std::cos(angle)*(a0.
z())+std::sin(angle)*b0.
z());
731 G4ThreeVector G4LivermorePolarizedComptonModel::GetPerpendicularPolarization
745 return gammaPolarization - gammaPolarization.
dot(gammaDirection)/gammaDirection.
dot(gammaDirection) * gammaDirection;
759 G4double sinTheta = std::sqrt(sinSqrTh);
763 G4double normalisation = std::sqrt(1. - cosSqrPhi*sinSqrTh);
802 if (rand1<(epsilon+1.0/epsilon-2)/(2.0*(epsilon+1.0/epsilon)-4.0*sinSqrTh*cosSqrPhi))
817 G4double sinBeta = std::sqrt(1-cosBeta*cosBeta);
821 G4double xParallel = normalisation*cosBeta;
822 G4double yParallel = -(sinSqrTh*cosPhi*sinPhi)*cosBeta/normalisation;
823 G4double zParallel = -(costheta*sinTheta*cosPhi)*cosBeta/normalisation;
825 G4double yPerpendicular = (costheta)*sinBeta/normalisation;
826 G4double zPerpendicular = -(sinTheta*sinPhi)*sinBeta/normalisation;
828 G4double xTotal = (xParallel + xPerpendicular);
829 G4double yTotal = (yParallel + yPerpendicular);
830 G4double zTotal = (zParallel + zPerpendicular);
832 gammaPolarization1.
setX(xTotal);
833 gammaPolarization1.
setY(yTotal);
834 gammaPolarization1.
setZ(zTotal);
836 return gammaPolarization1;
842 void G4LivermorePolarizedComptonModel::SystemOfRefChange(
G4ThreeVector& direction0,
858 direction1 = (direction_x*Axis_X0 + direction_y*Axis_Y0 + direction_z*Axis_Z0).unit();
863 polarization1 = (polarization_x*Axis_X0 + polarization_y*Axis_Y0 + polarization_z*Axis_Z0).unit();
878 G4AutoLock l(&LivermorePolarizedComptonModelMutex);
881 if(!
data[Z]) { ReadData(Z); }
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4double LowEnergyLimit() const
static constexpr double h_Planck
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
static G4LossTableManager * Instance()
virtual ~G4LivermorePolarizedComptonModel()
std::vector< G4Element * > G4ElementVector
std::ostringstream G4ExceptionDescription
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
G4double HighEnergyLimit() const
double dot(const Hep3Vector &) const
static G4double angle[DIM]
virtual G4bool LoadData(const G4String &fileName)
G4ParticleDefinition * GetDefinition() const
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
size_t GetVectorLength() const
const G4ElementVector * GetElementVector() const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
#define G4MUTEX_INITIALIZER
G4int SelectRandomShell(G4int Z) const
double howOrthogonal(const Hep3Vector &v) const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double twopi
static constexpr double electron_mass_c2
void LoadData(const G4String &fileName)
const XML_Char const XML_Char * data
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4GLOB_DLL std::ostream G4cout
size_t GetTableSize() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
const G4ThreeVector & GetMomentumDirection() const
G4double BindingEnergy(G4int Z, G4int shellIndex) const
static constexpr double cm
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
void ProposePolarization(const G4ThreeVector &dir)
static constexpr double eV
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
G4double Energy(size_t index) const
std::vector< G4EmElementSelector * > * GetElementSelectors()
G4double Value(G4double theEnergy, size_t &lastidx) const
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
static G4ProductionCutsTable * GetProductionCutsTable()
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
virtual G4double FindValue(G4double x, G4int componentId=0) const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static constexpr double c_light
G4LivermorePolarizedComptonModel(const G4ParticleDefinition *p=0, const G4String &nam="LivermorePolarizedCompton")
const G4ThreeVector & GetPolarization() const
static constexpr double GeV
static G4Electron * Electron()
void SetProposedKineticEnergy(G4double proposedKinEnergy)
static constexpr double MeV
static constexpr double pi
size_t GetNumberOfElements() const
Hep3Vector cross(const Hep3Vector &) const
G4VAtomDeexcitation * AtomDeexcitation()
void ProposeTrackStatus(G4TrackStatus status)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetDeexcitationFlag(G4bool val)
static constexpr double fine_structure_const
static constexpr double barn
static constexpr double keV
double epsilon(double density, double temperature)
G4double RandomSelectMomentum(G4int Z, G4int shellIndex) const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
const G4Material * GetMaterial() const