42#include <gsl/gsl_randist.h>
54 explicit DegraderInsideTester(ElementBase* el) {
55 deg_m =
static_cast<Degrader*
>(el);
57 bool checkHit(
const Vector_t& R)
override {
58 return deg_m->isInside(R);
65 explicit CollimatorInsideTester(ElementBase* el) {
66 col_m =
static_cast<CCollimator*
>(el);
68 bool checkHit(
const Vector_t& R)
override {
69 return col_m->checkPoint(
R(0),
R(1));
76 explicit FlexCollimatorInsideTester(ElementBase* el) {
77 col_m =
static_cast<FlexibleCollimator*
>(el);
79 bool checkHit(
const Vector_t& R)
override {
80 return col_m->isStopped(R);
83 FlexibleCollimator* col_m;
86 constexpr long double operator"" _keV(
long double value) {
return value; }
87 constexpr long double operator"" _MeV(
long double value) {
return value * 1e3; }
88 constexpr long double operator"" _GeV(
long double value) {
return value * 1e6; }
93 std::string& material,
94 bool enableRutherford,
129 rGen_m = gsl_rng_alloc(gsl_rng_default);
136 mySeed = tv.tv_sec + tv.tv_usec;
156 "Unsupported element type");
178 Z_m = material->getAtomicNumber();
179 A_m = material->getAtomicMass();
180 rho_m = material->getMassDensity();
181 X0_m = material->getRadiationLength();
182 I_m = material->getMeanExcitationEnergy();
196 const std::pair<Vector_t, double>& boundingSphere) {
217 "ScatteringPhysics::apply",
219 " is not supported for scattering interactions!");
244 if (onlyOneLoopOverParticles) {
263 onlyOneLoopOverParticles =
true;
265 }
while (onlyOneLoopOverParticles ==
false);
278 for (
size_t i = 0; i <
locParts_m.size(); ++i) {
330 bool includeFluctuations)
const {
341 const double gammaSqr = std::pow(gamma, 2);
342 const double betaSqr = 1.0 - 1.0 / gammaSqr;
343 double beta = std::sqrt(betaSqr);
344 double Ekin_keV = (gamma - 1) * mass_keV;
346 double epsilon = 0.0;
348 const double deltas = deltat * beta *
Physics::c;
351 const double massRatio = massElectron_keV / mass_keV;
352 double Tmax = (2.0 * massElectron_keV * betaSqr * gammaSqr /
353 (std::pow(gamma + massRatio, 2) - (gammaSqr - 1.0)));
357 if (Ekin_keV >= 0.6_MeV && Ekin_keV < 10.0_GeV) {
359 (0.5 * std::log(2 * massElectron_keV * betaSqr * gammaSqr * Tmax / std::pow(
I_m *
Units::eV2keV, 2)) - betaSqr));
362 const double Ts = Ekin_keV / massProton_amu;
363 if (Ekin_keV > 10.0_keV) {
364 const double epsilon_low =
A2_c * std::pow(Ts, 0.45);
365 const double epsilon_high = (
A3_c / Ts) * std::log(1 + (
A4_c / Ts) + (
A5_c * Ts));
366 epsilon = (epsilon_low * epsilon_high) / (epsilon_low + epsilon_high);
367 }
else if (Ekin_keV > 1.0_keV) {
368 epsilon =
A1_c * std::pow(Ts, 0.5);
373 "for energy loss calculation." <<
endl);
376 if (Ekin_keV > 10.0_MeV && Ekin_keV < 1.0_GeV) {
378 (0.5 * std::log(2 * massElectron_keV * betaSqr * gammaSqr * Tmax / std::pow(
I_m *
Units::eV2keV, 2)) - betaSqr));
379 }
else if (Ekin_keV > 1.0_keV && Ekin_keV <= 10.0_MeV) {
382 const double epsilon_high = (
B3_c / T) * std::log(1 + (
B4_c / T) + (
B5_c * T));
383 epsilon = (epsilon_low * epsilon_high) / (epsilon_low + epsilon_high);
387 "for energy loss calculation." <<
endl);
391 Ekin_keV += deltasrho * dEdx;
393 if (includeFluctuations) {
395 Ekin_keV += gsl_ran_gaussian(
rGen_m, sigma_E);
398 gamma = Ekin_keV / mass_keV + 1.0;
399 beta = std::sqrt(1.0 - 1.0 / std::pow(gamma, 2));
402 bool stopped = (Ekin_keV < lowEnergyThr_m * Units::MeV2keV || dEdx > 0);
418 R(0) = R(0) + shift * std::cos(Psixz);
419 R(2) = R(2) - shift * std::sin(Psixz);
423 P(0) = Px * std::cos(thetacou) + P(2) * std::sin(thetacou);
424 P(2) = -Px * std::sin(thetacou) + P(2) * std::cos(thetacou);
429 double thetaru = 2.5 / std::sqrt(gsl_rng_uniform(
rGen_m)) * 2.0 * theta0;
432 double normPtrans = std::sqrt(P(0) * P(0) + P(1) * P(1));
433 double Theta = std::atan(normPtrans / std::abs(P(2)));
434 double CosT = std::cos(Theta);
435 double SinT = std::sin(Theta);
437 Vector_t X(std::cos(phiru)*std::sin(thetaru),
438 std::sin(phiru)*std::sin(thetaru),
446 P =
X * CosT +
cross(W,
X) * SinT + W *
dot(W,
X) * (1.0 - CosT);
456 constexpr double sqrtThreeInv = 0.57735026918962576451;
460 const double theta0 = (13.6e6 / (beta * normP *
mass_m) *
462 (1.0 + 0.038 * std::log(deltas /
X0_m)));
465 for (
unsigned int i = 0; i < 2; ++ i) {
470 double z1 = gsl_ran_gaussian(
rGen_m, 1.0);
471 double z2 = gsl_ran_gaussian(
rGen_m, 1.0);
473 while(std::abs(z2) > 3.5) {
474 z1 = gsl_ran_gaussian(
rGen_m, 1.0);
475 z2 = gsl_ran_gaussian(
rGen_m, 1.0);
478 double thetacou = z2 * theta0;
479 double shift = (z1 * sqrtThreeInv + z2) * deltas * theta0 * 0.5;
498 const double elementLength =
element_ref_m->getElementLength();
500 for (
size_t i = 0; i < nL; ++ i) {
523 const PART& particle,
bool pdead) {
525 unsigned int numLocalParticles = bunch->
getLocalNum();
530 bunch->
Bin[numLocalParticles] = 1;
532 bunch->
Bin[numLocalParticles] = -1;
535 bunch->
R[numLocalParticles] = particle.Rincol;
536 bunch->
P[numLocalParticles] = particle.Pincol;
537 bunch->
Q[numLocalParticles] = particle.Qincol;
538 bunch->
M[numLocalParticles] = particle.Mincol;
539 bunch->
Bf[numLocalParticles] = 0.0;
540 bunch->
Ef[numLocalParticles] = 0.0;
541 bunch->
dt[numLocalParticles] =
dT_m;
546 const std::pair<Vector_t, double>& boundingSphere)
552 double zmax = boundingSphere.first(2) + boundingSphere.second;
553 double zmin = boundingSphere.first(2) - boundingSphere.second;
554 if (zmax < 0.0 || zmin >
element_ref_m->getElementLength()) {
560 std::set<size_t> partsToDel;
561 for (
size_t i = 0; i < nL; ++ i) {
562 if ((bunch->
Bin[i] == -1 || bunch->
Bin[i] == 1) &&
581 tau = bunch->
R[i](2) / stepWidth;
589 x.DTincol = bunch->
dt[i] * tau;
590 x.IDincol = bunch->
ID[i];
591 x.Binincol = bunch->
Bin[i];
592 x.Rincol = bunch->
R[i];
593 x.Pincol = bunch->
P[i];
594 x.Qincol = bunch->
Q[i];
595 x.Mincol = bunch->
M[i];
596 x.Bfincol = bunch->
Bf[i];
597 x.Efincol = bunch->
Ef[i];
604 partsToDel.insert(i);
608 for (
auto it = partsToDel.begin(); it != partsToDel.end(); ++ it) {
628 <<
"--- ScatteringPhysics ---\n"
629 <<
"Name: " <<
name_m <<
" - "
632 <<
"Particle Statistics @ " << time.
time() <<
"\n"
648 return x.label > y.label;
660 std::vector<PART>::iterator inv =
locParts_m.begin();
663 if ((*inv).label == -1)
678 for (
size_t i = 0; i <
locParts_m.size(); ++i) {
694 const double elementLength =
element_ref_m->getElementLength();
696 for (
size_t i = 0; i <
locParts_m.size(); ++i) {
703 if ( R(2) < elementLength &&
704 R(2) + stepLength(2) > elementLength &&
708 double distance = elementLength - R(2);
709 double tau = distance / stepLength(2);
720 for (
size_t i = 0; i <
locParts_m.size(); ++i) {
728 unsigned int locPartsInMat;
731 constexpr unsigned short numItems = 4;
732 unsigned int partStatistics[numItems] = {locPartsInMat,
737 allreduce(&(partStatistics[0]), numItems, std::plus<unsigned int>());
Vector3D cross(const Vector3D &lhs, const Vector3D &rhs)
Vector cross product.
double dot(const Vector3D &lhs, const Vector3D &rhs)
Vector dot product.
struct @313225227106163032040204343033076100017014020135 PART
Tps< T > cos(const Tps< T > &x)
Cosine.
Tps< T > sin(const Tps< T > &x)
Sine.
T euclidean_norm(const Vector< T > &)
Euclidean norm.
void allreduce(const T *input, T *output, int count, Op op)
PETE_TBTree< OpNE, Index::PETE_Expr_t, PETE_Scalar< double > > ne(const Index &idx, double x)
Inform & level2(Inform &inf)
Inform & level4(Inform &inf)
Inform & endl(Inform &inf)
constexpr double amu
The atomic mass unit energy equivalent in GeV.
constexpr double two_pi
The value of.
constexpr double Avo
The Avogadro's number.
constexpr double m_p
The proton rest mass in GeV.
constexpr double m_e
The electron rest mass in GeV.
constexpr double c
The velocity of light in m/s.
constexpr double pi
The value of.
constexpr double r_e
The classical electron radius in m.
int seed
The current random seed.
Vector_t getBeta(Vector_t p)
double getGamma(Vector_t p)
std::string toStringWithThousandSep(T value, char sep='\'')
ParticleAttrib< Vector_t > Ef
ParticleAttrib< int > Bin
double getQ() const
Access to reference data.
size_t getLocalNum() const
ParticleAttrib< double > M
ParticleAttrib< Vector_t > P
ParticleAttrib< double > Q
ParticleType getPType() const
void createWithID(unsigned id)
ParticleAttrib< double > dt
void performDestroy(bool updateLocalNum=false)
void destroy(size_t M, size_t I, bool doNow=false)
ParticleAttrib< Vector_t > Bf
static OpalData * getInstance()
Vector_t rotateFrom(const Vector_t &r) const
Vector_t transformFrom(const Vector_t &r) const
Vector_t rotateTo(const Vector_t &r) const
static std::shared_ptr< Material > getMaterial(const std::string &name)
static std::string getParticleTypeString(const ParticleType &type)
std::unique_ptr< InsideTester > hitTester_m
ParticleMatterInteractionHandler(std::string name, ElementBase *elref)
ElementBase * element_ref_m
bool allParticleInMat_m
if all particles are in matter stay inside the particle matter interaction
void addBackToBunch(PartBunchBase< double, 3 > *bunch)
void configureMaterialParameters()
The material of the collimator.
void setTimeStepForLeavingParticles()
void copyFromBunch(PartBunchBase< double, 3 > *bunch, const std::pair< Vector_t, double > &boundingSphere)
virtual void apply(PartBunchBase< double, 3 > *bunch, const std::pair< Vector_t, double > &boundingSphere) override
void addParticleBackToBunch(PartBunchBase< double, 3 > *bunch, const PART &particle, bool pdead=false)
unsigned int bunchToMatStat_m
void deleteParticleFromLocalVector()
IpplTimings::TimerRef DegraderLoopTimer_m
void applyRotation(Vector_t &P, Vector_t &R, double xplane, double thetacou)
virtual bool computeEnergyLoss(PartBunchBase< double, 3 > *bunch, Vector_t &P, const double deltat, bool includeFluctuations=true) const override
ScatteringPhysics(const std::string &name, ElementBase *element, std::string &mat, bool enableRutherford, double lowEnergyThr)
void computeInteraction(PartBunchBase< double, 3 > *bunch)
unsigned int totalPartsInMat_m
void applyRandomRotation(Vector_t &P, double theta0)
void computeCoulombScattering(Vector_t &R, Vector_t &P, double dt)
std::vector< PART > locParts_m
unsigned int rediffusedStat_m
IpplTimings::TimerRef DegraderApplyTimer_m
unsigned int stoppedPartStat_m
virtual bool stillActive() override
IpplTimings::TimerRef DegraderDestroyTimer_m
virtual std::string getName() override
std::unique_ptr< LossDataSink > lossDs_m
virtual void print(Inform &os) override
std::string time() const
Return time.
std::ios_base::fmtflags FmtFlags_t
static TimerRef getTimer(const char *nm)
static void stopTimer(TimerRef t)
static void startTimer(TimerRef t)
Vektor< double, 3 > Vector_t