42#include <gsl/gsl_randist.h>
57 return deg_m->isInMaterial(
R(2));
68 return col_m->checkPoint(
R(0),
R(1));
75 explicit FlexCollimatorInsideTester(
ElementBase* el) {
79 return col_m->isStopped(
R);
85 constexpr long double operator"" _keV(
long double value) {
return value; }
86 constexpr long double operator"" _MeV(
long double value) {
return value * 1e3; }
87 constexpr long double operator"" _GeV(
long double value) {
return value * 1e6; }
92 std::string& material,
93 bool enableRutherford,
100 material_m(material),
117 stoppedPartStat_m(0),
119 totalPartsInMat_m(0),
123 enableRutherford_m(enableRutherford),
124 lowEnergyThr_m(lowEnergyThr)
128 rGen_m = gsl_rng_alloc(gsl_rng_default);
135 mySeed = tv.tv_sec + tv.tv_usec;
155 "Unsupported element type");
177 Z_m = material->getAtomicNumber();
178 A_m = material->getAtomicMass();
179 rho_m = material->getMassDensity();
180 X0_m = material->getRadiationLength();
181 I_m = material->getMeanExcitationEnergy();
195 const std::pair<Vector_t, double>& boundingSphere) {
216 "ScatteringPhysics::apply",
218 " is not supported for scattering interactions!");
243 if (onlyOneLoopOverParticles) {
262 onlyOneLoopOverParticles =
true;
264 }
while (onlyOneLoopOverParticles ==
false);
277 for (
size_t i = 0; i <
locParts_m.size(); ++i) {
329 bool includeFluctuations)
const {
340 const double gammaSqr =
std::pow(gamma, 2);
341 const double betaSqr = 1.0 - 1.0 / gammaSqr;
343 double Ekin_keV = (gamma - 1) * mass_keV;
345 double epsilon = 0.0;
347 const double deltas = deltat * beta *
Physics::c;
350 const double massRatio = massElectron_keV / mass_keV;
351 double Tmax = (2.0 * massElectron_keV * betaSqr * gammaSqr /
352 (
std::pow(gamma + massRatio, 2) - (gammaSqr - 1.0)));
356 if (Ekin_keV >= 0.6_MeV && Ekin_keV < 10.0_GeV) {
361 const double Ts = Ekin_keV / massProton_amu;
362 if (Ekin_keV > 10.0_keV) {
365 epsilon = (epsilon_low * epsilon_high) / (epsilon_low + epsilon_high);
366 }
else if (Ekin_keV > 1.0_keV) {
372 "for energy loss calculation." <<
endl);
375 if (Ekin_keV > 10.0_MeV && Ekin_keV < 1.0_GeV) {
378 }
else if (Ekin_keV > 1.0_keV && Ekin_keV <= 10.0_MeV) {
382 epsilon = (epsilon_low * epsilon_high) / (epsilon_low + epsilon_high);
386 "for energy loss calculation." <<
endl);
390 Ekin_keV += deltasrho * dEdx;
392 if (includeFluctuations) {
394 Ekin_keV += gsl_ran_gaussian(
rGen_m, sigma_E);
397 gamma = Ekin_keV / mass_keV + 1.0;
401 bool stopped = (Ekin_keV < lowEnergyThr_m * Units::MeV2keV || dEdx > 0);
428 double thetaru = 2.5 /
std::sqrt(gsl_rng_uniform(
rGen_m)) * 2.0 * theta0;
431 double normPtrans =
std::sqrt(P(0) * P(0) + P(1) * P(1));
455 constexpr double sqrtThreeInv = 0.57735026918962576451;
459 const double theta0 = (13.6e6 / (beta * normP *
mass_m) *
464 for (
unsigned int i = 0; i < 2; ++ i) {
469 double z1 = gsl_ran_gaussian(
rGen_m, 1.0);
470 double z2 = gsl_ran_gaussian(
rGen_m, 1.0);
473 z1 = gsl_ran_gaussian(
rGen_m, 1.0);
474 z2 = gsl_ran_gaussian(
rGen_m, 1.0);
477 double thetacou = z2 * theta0;
478 double shift = (z1 * sqrtThreeInv + z2) * deltas * theta0 * 0.5;
499 for (
size_t i = 0; i < nL; ++ i) {
522 const PART& particle,
bool pdead) {
524 unsigned int numLocalParticles = bunch->
getLocalNum();
529 bunch->
Bin[numLocalParticles] = 1;
531 bunch->
Bin[numLocalParticles] = -1;
534 bunch->
R[numLocalParticles] = particle.
Rincol;
535 bunch->
P[numLocalParticles] = particle.
Pincol;
536 bunch->
Q[numLocalParticles] = particle.
Qincol;
537 bunch->
M[numLocalParticles] = particle.
Mincol;
538 bunch->
Bf[numLocalParticles] = 0.0;
539 bunch->
Ef[numLocalParticles] = 0.0;
540 bunch->
dt[numLocalParticles] =
dT_m;
545 const std::pair<Vector_t, double>& boundingSphere)
551 double zmax = boundingSphere.first(2) + boundingSphere.second;
552 double zmin = boundingSphere.first(2) - boundingSphere.second;
559 std::set<size_t> partsToDel;
560 for (
size_t i = 0; i < nL; ++ i) {
561 if ((bunch->
Bin[i] == -1 || bunch->
Bin[i] == 1) &&
580 tau = bunch->
R[i](2) / stepWidth;
603 partsToDel.insert(i);
607 for (
auto it = partsToDel.begin(); it != partsToDel.end(); ++ it) {
627 <<
"--- ScatteringPhysics ---\n"
628 <<
"Name: " <<
name_m <<
" - "
631 <<
"Particle Statistics @ " << time.
time() <<
"\n"
662 if ((*inv).label == -1)
677 for (
size_t i = 0; i <
locParts_m.size(); ++i) {
695 for (
size_t i = 0; i <
locParts_m.size(); ++i) {
702 if (
R(2) < elementLength &&
703 R(2) + stepLength(2) > elementLength &&
707 double distance = elementLength -
R(2);
708 double tau = distance / stepLength(2);
719 for (
size_t i = 0; i <
locParts_m.size(); ++i) {
727 unsigned int locPartsInMat;
730 constexpr unsigned short numItems = 4;
731 unsigned int partStatistics[numItems] = {locPartsInMat,
736 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.
Tps< T > log(const Tps< T > &x)
Natural logarithm.
Tps< T > cos(const Tps< T > &x)
Cosine.
Tps< T > pow(const Tps< T > &x, int y)
Integer power.
Tps< T > sin(const Tps< T > &x)
Sine.
Tps< T > sqrt(const Tps< T > &x)
Square root.
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)
PETE_TUTree< FnArcTan, typename T::PETE_Expr_t > atan(const PETE_Expr< T > &l)
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
PETE_TBTree< FnFmod, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > fmod(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
PETE_TBTree< FnArcTan2, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > atan2(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
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.
std::string::iterator iterator
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()
virtual const std::string & getName() const
Get element name.
virtual double getElementLength() const
Get design length.
virtual ElementType getType() const =0
Get element type std::string.
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)
virtual bool checkHit(const Vector_t &R)=0
std::unique_ptr< InsideTester > hitTester_m
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