69 fGammaCutInKineticEnergy(0),
75 fAngleForEnergyTable(0)
112 G4cout<<
"### G4VXTRenergyLoss: the number of TR radiator plates = "
116 G4Exception(
"G4VXTRenergyLoss::G4VXTRenergyLoss()",
"VXTRELoss01",
197 G4double charge, chargeSq, massRatio, TkinScaled;
208 gamma = 1.0 + kinEnergy/mass;
214 if ( std::fabs( gamma -
fGamma ) < 0.05*gamma ) lambda =
fLambda;
218 chargeSq = charge*charge;
220 TkinScaled = kinEnergy*massRatio;
222 for(iTkin = 0; iTkin <
fTotBin; iTkin++)
224 if( TkinScaled < fProtonEnergyVector->GetLowEdgeEnergy(iTkin))
break;
228 if(iTkin == 0) lambda =
DBL_MAX;
233 sigma = (*(*fEnergyDistrTable)(iPlace))(0)*chargeSq;
240 W1 = (E2 - TkinScaled)*W;
241 W2 = (TkinScaled - E1)*W;
242 sigma = ( (*(*fEnergyDistrTable)(iPlace ))(0)*W1 +
247 else lambda = 1./sigma;
269 "XTR initialisation for neutral particle ?!" );
277 G4cout<<
"Build angle for energy distribution according the current radiator"
291 G4int iTkin, iTR, iPlace;
317 G4cout<<
"Lorentz Factor"<<
"\t"<<
"XTR photon number"<<
G4endl;
320 for( iTkin = 0; iTkin <
fTotBin; iTkin++ )
340 for( iTR =
fBinTR - 2; iTR >= 0; iTR-- )
370 G4cout<<
"total time for build X-ray TR energy loss tables = "
415 for( iTkin = 0; iTkin <
fTotBin; iTkin++ )
430 for( iTR = 0; iTR <
fBinTR; iTR++ )
438 angleVector ->
PutValue(fBinTR - 1, angleSum);
440 for( i = fBinTR - 2; i >= 0; i-- )
449 angleVector ->
PutValue(i, angleSum);
460 G4cout<<
"total time for build X-ray TR angle for energy loss tables = "
493 G4cout<<
"Lorentz Factor"<<
"\t"<<
"XTR photon number"<<
G4endl;
496 for( iTkin = 0; iTkin <
fTotBin; iTkin++ )
514 for( iTR = 0; iTR <
fBinTR; iTR++ )
533 G4cout<<
"total time for build XTR angle for given energy tables = "
547 G4double theta=0., result,
tmp=0., cof1, cof2, cofMin, cofPHC, angleSum = 0.;
548 G4int iTheta, k, kMin;
561 kMin =
G4int(cofMin);
562 if (cofMin > kMin) kMin++;
568 G4cout<<
"n-1 = "<<n-1<<
"; theta = "
571 <<
"; angleSum = "<<angleSum<<
G4endl;
575 for( iTheta = n - 1; iTheta >= 1; iTheta-- )
578 k = iTheta- 1 + kMin;
582 result = (k - cof1)*(k - cof1)*(k + cof2)*(k + cof2);
584 tmp = std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
586 if( k == kMin && kMin ==
G4int(cofMin) )
590 else if(iTheta == n-1);
595 theta = std::abs(k-cofMin)*cofPHC/energy/(fPlateThick +
fGasThick);
599 G4cout<<
"iTheta = "<<iTheta<<
"; k = "<<k<<
"; theta = "
600 <<std::sqrt(theta)*
fGamma<<
"; tmp = "
602 <<
"; angleSum = "<<angleSum<<
G4endl;
604 angleVector->
PutValue( iTheta, theta, angleSum );
613 G4cout<<
"iTheta = "<<iTheta<<
"; theta = "
614 <<std::sqrt(theta)*
fGamma<<
"; tmp = "
616 <<
"; angleSum = "<<angleSum<<
G4endl;
618 angleVector->
PutValue( iTheta, theta, angleSum );
629 G4int iTkin, iTR, iPlace;
649 G4cout<<
"Lorentz Factor"<<
"\t"<<
"XTR photon number"<<
G4endl;
652 for( iTkin = 0; iTkin <
fTotBin; iTkin++ )
678 for( iTR =
fBinTR - 2; iTR >= 0; iTR-- )
686 angleVector ->
PutValue(iTR,angleSum);
704 G4cout<<
"total time for build X-ray TR angle tables = "
722 G4double energyTR, theta,theta2, phi, dirX, dirY, dirZ;
729 G4cout<<
"Start of G4VXTRenergyLoss::PostStepDoIt "<<
G4endl;
730 G4cout<<
"name of current material = "
737 G4cout<<
"Go out from G4VXTRenergyLoss::PostStepDoIt: wrong volume "<<
G4endl;
750 G4double gamma = 1.0 + kinEnergy/mass;
757 G4double TkinScaled = kinEnergy*massRatio;
762 for( iTkin = 0; iTkin <
fTotBin; iTkin++ )
764 if(TkinScaled < fProtonEnergyVector->GetLowEdgeEnergy(iTkin))
break;
772 G4cout<<
"Go out from G4VXTRenergyLoss::PostStepDoIt:iTkin = "<<iTkin<<
G4endl;
790 if(theta2 > 0.) theta = std::sqrt(theta2);
793 else theta = std::fabs(G4RandGauss::shoot(0.0,
pi/gamma));
797 if( theta >= 0.1 ) theta = 0.1;
803 dirX = std::sin(theta)*std::cos(phi);
804 dirY = std::sin(theta)*std::sin(phi);
805 dirZ = std::cos(theta);
812 directionTR, energyTR);
832 position += distance*directionTR;
836 startTime, position );
880 if(result < 0.0) result = 0.0;
893 G4double lim[8] = { 0.0, 0.01, 0.02, 0.05, 0.1, 0.2, 0.5, 1.0 };
919 for( i = 0; i < iMax-1; i++ )
938 if(result < 0) result = 0.0;
951 G4double sum = 0., tmp1, tmp2,
tmp=0., cof1, cof2, cofMin, cofPHC, energy1, energy2;
952 G4int k, kMax, kMin, i;
961 cofMin = std::sqrt(cof1*cof2);
964 kMin =
G4int(cofMin);
965 if (cofMin > kMin) kMin++;
971 for( k = kMin; k <= kMax; k++ )
974 tmp2 = std::sqrt(tmp1*tmp1-cof1*cof2);
975 energy1 = (tmp1+tmp2)/cof1;
976 energy2 = (tmp1-tmp2)/cof1;
978 for(i = 0; i < 2; i++)
985 tmp2 = std::sin(tmp1);
986 tmp = energy1*tmp2*tmp2;
988 tmp1 = hbarc*energy1/( energy1*energy1*(1./
fGamma/
fGamma + varAngle) + fSigma2 );
989 tmp *= (tmp1-tmp2)*(tmp1-tmp2);
990 tmp1 = cof1/(4*
hbarc) - cof2/(4*hbarc*energy1*energy1);
991 tmp2 = std::abs(tmp1);
992 if(tmp2 > 0.) tmp /= tmp2;
1000 tmp2 = std::sin(tmp1);
1001 tmp = energy2*tmp2*tmp2;
1003 tmp1 = hbarc*energy2/( energy2*energy2*(1./
fGamma/
fGamma + varAngle) + fSigma2 );
1004 tmp *= (tmp1-tmp2)*(tmp1-tmp2);
1005 tmp1 = cof1/(4*
hbarc) - cof2/(4*hbarc*energy2*energy2);
1006 tmp2 = std::abs(tmp1);
1007 if(tmp2 > 0.) tmp /= tmp2;
1016 result /= hbarc*
hbarc;
1038 lambda = 1.0/gamma/gamma + varAngle +
fSigma1/omega/omega;
1051 G4double cof, length,delta, real_v, image_v;
1055 cof = 1.0/(1.0 + delta*delta);
1057 real_v = length*cof;
1058 image_v = real_v*delta;
1090 omega2 = omega*omega;
1091 omega3 = omega2*omega;
1092 omega4 = omega2*omega2;
1095 G4double cross = SandiaCof[0]/omega + SandiaCof[1]/omega2 +
1096 SandiaCof[2]/omega3 + SandiaCof[3]/omega4;
1110 lambda = 1.0/gamma/gamma + varAngle +
fSigma2/omega/omega;
1124 G4double cof, length,delta, real_v, image_v;
1128 cof = 1.0/(1.0 + delta*delta);
1130 real_v = length*cof;
1131 image_v = real_v*delta;
1161 omega2 = omega*omega;
1162 omega3 = omega2*omega;
1163 omega4 = omega2*omega2;
1166 G4double cross = SandiaCof[0]/omega + SandiaCof[1]/omega2 +
1167 SandiaCof[2]/omega3 + SandiaCof[3]/omega4;
1191 std::ofstream outPlate(
"plateZmu.dat", std::ios::out );
1192 outPlate.setf( std::ios::scientific, std::ios::floatfield );
1197 varAngle = 1/gamma/gamma;
1202 omega = (1.0 + i)*
keV;
1229 std::ofstream outGas(
"gasZmu.dat", std::ios::out );
1230 outGas.setf( std::ios::scientific, std::ios::floatfield );
1234 varAngle = 1/gamma/gamma;
1239 omega = (1.0 + i)*
keV;
1254 G4int i, numberOfElements;
1255 G4double xSection = 0., nowZ, sumZ = 0.;
1258 numberOfElements = (*theMaterialTable)[
fMatIndex1]->GetNumberOfElements();
1260 for( i = 0; i < numberOfElements; i++ )
1262 nowZ = (*theMaterialTable)[
fMatIndex1]->GetElement(i)->GetZ();
1267 xSection *= (*theMaterialTable)[
fMatIndex1]->GetElectronDensity();
1278 G4int i, numberOfElements;
1279 G4double xSection = 0., nowZ, sumZ = 0.;
1282 numberOfElements = (*theMaterialTable)[
fMatIndex2]->GetNumberOfElements();
1284 for( i = 0; i < numberOfElements; i++ )
1286 nowZ = (*theMaterialTable)[
fMatIndex2]->GetElement(i)->GetZ();
1291 xSection *= (*theMaterialTable)[
fMatIndex2]->GetElectronDensity();
1303 if ( Z < 0.9999 )
return CrossSection;
1304 if ( GammaEnergy < 0.1*
keV )
return CrossSection;
1305 if ( GammaEnergy > (100.*
GeV/Z) )
return CrossSection;
1307 static const G4double a = 20.0 ,
b = 230.0 ,
c = 440.0;
1310 d1= 2.7965e-1*
barn, d2=-1.8300e-1*
barn, d3= 6.7527 *
barn, d4=-1.9798e+1*
barn,
1311 e1= 1.9756e-5*
barn, e2=-1.0205e-2*
barn, e3=-7.3913e-2*
barn, e4= 2.7079e-2*
barn,
1314 G4double p1Z = Z*(d1 + e1*Z + f1*Z*
Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
1315 p3Z = Z*(d3 + e3*Z + f3*Z*
Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
1318 if (Z < 1.5) T0 = 40.0*
keV;
1321 CrossSection = p1Z*std::log(1.+2.*X)/X
1322 + (p2Z + p3Z*X + p4Z*X*
X)/(1. + a*X +
b*X*X +
c*X*X*X);
1326 if (GammaEnergy < T0)
1330 G4double sigma = p1Z*std::log(1.+2*X)/X
1331 + (p2Z + p3Z*X + p4Z*X*
X)/(1. + a*X +
b*X*X +
c*X*X*X);
1332 G4double c1 = -T0*(sigma-CrossSection)/(CrossSection*dT0);
1334 if (Z > 1.5) c2 = 0.375-0.0556*std::log(Z);
1336 CrossSection *= std::exp(-y*(c1+c2*y));
1339 return CrossSection;
1357 G4double formationLength1, formationLength2;
1358 formationLength1 = 1.0/
1362 formationLength2 = 1.0/
1366 return (varAngle/energy)*(formationLength1 - formationLength2)
1367 *(formationLength1 - formationLength2);
1437 std::ofstream outEn(
"numberE.dat", std::ios::out );
1438 outEn.setf( std::ios::scientific, std::ios::floatfield );
1440 std::ofstream outAng(
"numberAng.dat", std::ios::out );
1441 outAng.setf( std::ios::scientific, std::ios::floatfield );
1443 for(iTkin=0;iTkin<
fTotBin;iTkin++)
1447 numberE = (*(*fEnergyDistrTable)(iTkin))(0);
1450 G4cout<<gamma<<
"\t\t"<<numberE<<
"\t"
1453 outEn<<gamma<<
"\t\t"<<numberE<<
G4endl;
1465 G4int iTransfer, iPlace;
1476 for(iTransfer=0;;iTransfer++)
1487 W1 = (E2 - scaledTkin)*W;
1488 W2 = (scaledTkin - E1)*W;
1490 position =( (*(*fEnergyDistrTable)(iPlace))(0)*W1 +
1495 for(iTransfer=0;;iTransfer++)
1505 if(transfer < 0.0 ) transfer = 0.0;
1522 result = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1526 y1 = (*(*fEnergyDistrTable)(iPlace))(iTransfer-1);
1527 y2 = (*(*fEnergyDistrTable)(iPlace))(iTransfer);
1529 x1 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
1530 x2 = (*fEnergyDistrTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1532 if ( x1 == x2 ) result =
x2;
1555 if (iTkin ==
fTotBin) iTkin--;
1559 for( iTR = 0; iTR <
fBinTR; iTR++ )
1561 if( energyXTR < fXTREnergyVector->GetLowEdgeEnergy(iTR) )
break;
1563 if (iTR == fBinTR) iTR--;
1565 position = (*(*fAngleForEnergyTable)(iTR))(0)*
G4UniformRand();
1567 for(iAngle = 0;;iAngle++)
1588 result = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1592 y1 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer-1);
1593 y2 = (*(*fAngleForEnergyTable)(iPlace))(iTransfer);
1595 x1 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
1596 x2 = (*fAngleForEnergyTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
1598 if ( x1 == x2 ) result =
x2;
1604 result = x1 + (position -
y1)*(x2 - x1)/(y2 -
y1);