31#include "gsl/gsl_poly.h"
47 messageHeader_m(
" * "),
48 pusher_m(right.pusher_m),
49 designRadius_m(right.designRadius_m),
50 exitAngle_m(right.exitAngle_m),
51 fieldIndex_m(right.fieldIndex_m),
52 startField_m(right.startField_m),
53 endField_m(right.endField_m),
54 widthEntranceFringe_m(right.widthEntranceFringe_m),
55 widthExitFringe_m(right.widthExitFringe_m),
56 reinitialize_m(right.reinitialize_m),
57 entranceParameter1_m(right.entranceParameter1_m),
58 entranceParameter2_m(right.entranceParameter2_m),
59 entranceParameter3_m(right.entranceParameter3_m),
60 exitParameter1_m(right.exitParameter1_m),
61 exitParameter2_m(right.exitParameter2_m),
62 exitParameter3_m(right.exitParameter3_m),
63 engeCoeffsEntry_m(right.engeCoeffsEntry_m),
64 engeCoeffsExit_m(right.engeCoeffsExit_m),
65 entryFieldValues_m(nullptr),
66 exitFieldValues_m(nullptr),
67 entryFieldAccel_m(nullptr),
68 exitFieldAccel_m(nullptr),
69 deltaBeginEntry_m(right.deltaBeginEntry_m),
70 deltaEndEntry_m(right.deltaEndEntry_m),
71 polyOrderEntry_m(right.polyOrderEntry_m),
72 deltaBeginExit_m(right.deltaBeginExit_m),
73 deltaEndExit_m(right.deltaEndExit_m),
74 polyOrderExit_m(right.polyOrderExit_m),
75 cosEntranceAngle_m(right.cosEntranceAngle_m),
76 sinEntranceAngle_m(right.sinEntranceAngle_m),
77 tanEntranceAngle_m(right.tanEntranceAngle_m),
78 tanExitAngle_m(right.tanExitAngle_m),
79 beginToEnd_m(right.beginToEnd_m),
80 beginToEnd_lcs_m(right.beginToEnd_lcs_m),
81 toEntranceRegion_m(right.toEntranceRegion_m),
82 toExitRegion_m(right.toExitRegion_m),
83 computeAngleTrafo_m(right.computeAngleTrafo_m),
84 maxAngle_m(right.maxAngle_m),
85 nSlices_m(right.nSlices_m){
90 messageHeader_m(
" * "),
97 widthEntranceFringe_m(0.0),
98 widthExitFringe_m(0.0),
99 reinitialize_m(false),
100 entranceParameter1_m(0.0),
101 entranceParameter2_m(0.0),
102 entranceParameter3_m(0.0),
103 exitParameter1_m(0.0),
104 exitParameter2_m(0.0),
105 exitParameter3_m(0.0),
106 entryFieldValues_m(nullptr),
107 exitFieldValues_m(nullptr),
108 entryFieldAccel_m(nullptr),
109 exitFieldAccel_m(nullptr),
110 deltaBeginEntry_m(0.0),
111 deltaEndEntry_m(0.0),
113 deltaBeginExit_m(0.0),
116 cosEntranceAngle_m(1.0),
117 sinEntranceAngle_m(0.0),
118 tanEntranceAngle_m(0.0),
126 for (
unsigned int i = 0; i < 3u; ++ i) {
184 return apply(
R, P, t, E, B);
200 <<
"======================================================================\n"
201 <<
"===== " << std::left << std::setw(64) << std::setfill(
'=') <<
name <<
"\n"
202 <<
"======================================================================\n";
208 double bendAngleX = 0.0;
209 double bendAngleY = 0.0;
211 print(msg, bendAngleX, bendAngleY);
223 <<
"======================================================================\n"
224 <<
"======================================================================\n"
227 ERRORMSG(
"There is something wrong with your field map '"
229 endField = startField - 1.0e-3;
265 while (deltaS < bendLength) {
291 const std::vector<double> &engeCoeff,
294 double &engeFuncDeriv,
295 double &engeFuncSecDerivNorm) {
298 double expSumDeriv = 0.0;
299 double expSumSecDeriv = 0.0;
301 if (polyOrder >= 2) {
303 expSum = engeCoeff.at(0)
304 + engeCoeff.at(1) * zNormalized;
305 expSumDeriv = engeCoeff.at(1);
307 for (
int index = 2; index <= polyOrder; index++) {
308 expSum += engeCoeff.at(index) *
std::pow(zNormalized, index);
309 expSumDeriv += index * engeCoeff.at(index)
311 expSumSecDeriv += index * (index - 1) * engeCoeff.at(index)
315 }
else if (polyOrder == 1) {
317 expSum = engeCoeff.at(0)
318 + engeCoeff.at(1) * zNormalized;
319 expSumDeriv = engeCoeff.at(1);
322 expSum = engeCoeff.at(0);
325 engeFunc = 1.0 / (1.0 + engeExp);
329 expSumDeriv /=
gap_m;
331 double engeExpDeriv = expSumDeriv * engeExp;
332 double engeExpSecDeriv = (expSumSecDeriv +
std::pow(expSumDeriv, 2.0))
334 double engeFuncSq =
std::pow(engeFunc, 2.0);
336 engeFuncDeriv = -engeExpDeriv * engeFuncSq;
340 engeFuncSecDerivNorm = -engeExpSecDeriv * engeFunc
344 engeFuncSecDerivNorm = 0.0;
349 engeFuncSecDerivNorm = 0.0;
402 B(1) = (engeFunc - 0.5 * engeFuncSecDeriv *
std::pow(Rprime(1), 2.0));
404 B(2) = engeFuncDeriv * Rprime(1);
439 B(1) = (engeFunc - 0.5 * engeFuncSecDeriv *
std::pow(Rprime(1), 2.0));
441 B(2) = engeFuncDeriv * Rprime(1);
450 bool horizontallyInside =
false;
453 if (verticallyInside) {
458 if (inEntranceRegion && inExitRegion) {
463 inExitRegion =
false;
465 inEntranceRegion =
false;
468 if (inEntranceRegion) {
470 }
else if (inExitRegion) {
481 Vector_t BEntrance(0.0), BExit(0.0);
487 if (inEntranceRegion) {
488 horizontallyInside =
true;
489 if (verticallyInside) {
495 horizontallyInside =
true;
496 if (verticallyInside) {
503 B = BEntrance + BExit;
505 bool hitMaterial = (horizontallyInside && (!verticallyInside));
511 std::ofstream trajectoryOutput;
517 trajectoryOutput.open(fname);
518 trajectoryOutput.precision(12);
519 trajectoryOutput <<
"# " << std::setw(18) <<
"s"
520 << std::setw(20) <<
"x"
521 << std::setw(20) <<
"z"
522 << std::setw(20) <<
"By"
539 const double stepSize = betaGamma / gamma *
Physics::c * dt;
542 while (deltaS < bendLength) {
555 trajectoryOutput << std::setw(20) << deltaS + 0.5 * stepSize
556 << std::setw(20) <<
X(0)
557 << std::setw(20) <<
X(2)
558 << std::setw(20) << bField(1)
582 double fieldStep = (
angle_m - actualBendAngle) * effectiveFieldAmplitude;
583 if (ratioSquared < 1.0) {
584 fieldStep = (
angle_m - actualBendAngle) * effectiveFieldAmplitude *
std::sqrt(1.0 - ratioSquared);
606 if (error > 1.0e-6) {
608 double deltaStep = 0.0;
610 deltaStep = -
gap_m / 2.0;
612 deltaStep =
gap_m / 2.0;
615 double bendAngle1 = actualBendAngle;
617 double delta2 = deltaStep;
642 unsigned int iterations = 1;
644 while (error > 1.0e-6 && iterations < 100) {
646 double delta = (delta1 + delta2) / 2.0;
654 if (error > 1.0e-6) {
656 if (bendAngle1 -
angle_m < 0.0) {
658 if (newBendAngle -
angle_m < 0.0) {
659 bendAngle1 = newBendAngle;
668 if (newBendAngle -
angle_m < 0.0) {
672 bendAngle1 = newBendAngle;
689 const double tolerance = 1
e-7;
692 if (error < tolerance)
697 double amplitude2 = amplitude1;
698 double bendAngle1 = actualBendAngle;
699 double bendAngle2 = bendAngle1;
703 amplitude1 = amplitude2;
704 bendAngle1 = bendAngle2;
706 amplitude2 += stepSign * fieldStep;
715 unsigned int iterations = 1;
716 while (error > tolerance && iterations < 100) {
723 if (error > tolerance) {
725 if (bendAngle1 -
angle_m < 0.0) {
727 if (newBendAngle -
angle_m < 0.0) {
728 bendAngle1 = newBendAngle;
737 if (newBendAngle -
angle_m < 0.0) {
741 bendAngle1 = newBendAngle;
752 bool reinitialize =
false;
770 if (refCharge < 0.0) {
777 reinitialize =
false;
807 effectiveAngle >= 0.0 && effectiveAngle <
maxAngle_m) {
808 if (effectiveAngle < 0.5 *
maxAngle_m)
return R(2) >= 0.0;
809 return Rprime(2) < 0.0;
826 return (Rprime(2) >= 0 &&
848 <<
"Start of field map: "
850 <<
" m (in s coordinates)"
852 msg <<
"End of field map: "
854 <<
" m (in s coordinates)"
856 msg <<
"Entrance edge of magnet: "
858 <<
" m (in s coordinates)"
861 <<
"Reference Trajectory Properties"
863 <<
"======================================================================"
865 msg <<
"Bend angle magnitude: "
871 msg <<
"Entrance edge angle: "
877 msg <<
"Exit edge angle: "
883 msg <<
"Bend design radius: "
887 msg <<
"Bend design energy: "
892 <<
"Bend Field and Rotation Properties"
894 <<
"======================================================================"
896 msg <<
"Field amplitude: "
900 msg <<
"Field index: "
903 msg <<
"Rotation about z axis: "
910 <<
"Reference Trajectory Properties Through Bend Magnet with Fringe Fields"
912 <<
"======================================================================"
914 msg <<
"Reference particle is bent: "
918 <<
" degrees) in x plane"
920 msg <<
"Reference particle is bent: "
924 <<
" degrees) in y plane"
929 <<
"Effective Field Map\n"
930 <<
"======================================================================\n"
966 unsigned int numStepsEntry =
std::ceil(entryLength / stepSize) + 3;
967 double stepSizeEntry = entryLength / (numStepsEntry - 3);
968 std::vector<double> zvalsEntry(numStepsEntry);
969 std::vector<std::vector<double> > fieldValuesEntry(3);
972 unsigned int numStepsExit =
std::ceil(exitLength / stepSize) + 3;
973 double stepSizeExit = exitLength / (numStepsExit - 3);
974 std::vector<double> zvalsExit(numStepsExit);
975 std::vector<std::vector<double> > fieldValuesExit(3);
982 for (
unsigned int i = 0; i < 3u; ++ i) {
983 fieldValuesEntry[i].resize(numStepsEntry);
984 fieldValuesExit[i].resize(numStepsExit);
990 for (
unsigned int j = 0; j < numStepsEntry; ++ j) {
994 zvalsEntry[j] = zvalsEntry[j - 1] + stepSizeEntry;
999 fieldValuesEntry[0][j],
1000 fieldValuesEntry[1][j],
1001 fieldValuesEntry[2][j]);
1002 fieldValuesEntry[2][j] *= fieldValuesEntry[0][j];
1005 for (
unsigned int j = 0; j < numStepsExit; ++ j) {
1009 zvalsExit[j] = zvalsExit[j - 1] + stepSizeExit;
1014 fieldValuesExit[0][j],
1015 fieldValuesExit[1][j],
1016 fieldValuesExit[2][j]);
1017 fieldValuesExit[2][j] *= fieldValuesExit[0][j];
1020 for (
unsigned int i = 0; i < 3u; ++ i) {
1021 gsl_spline_init(
entryFieldValues_m[i], &(zvalsEntry[0]), &(fieldValuesEntry[i][0]), numStepsEntry);
1022 gsl_spline_init(
exitFieldValues_m[i], &(zvalsExit[0]), &(fieldValuesExit[i][0]), numStepsExit);
1086 std::sin(0.5 * (0.5 * bendAngle - entranceAngle)) * rotationAxis);
1102 if (rotationCenter(2) < 0.0) {
1105 gsl_poly_solve_quadratic(
dot(P,P),
1111 maxAngleEntranceAperture = tau * P;
1115 gsl_poly_solve_quadratic(
dot(P,P),
1121 maxAngleEntranceAperture = tau * P;
1127 gsl_poly_solve_quadratic(
dot(P,P),
1137 gsl_poly_solve_quadratic(
dot(P,P),
1146 maxAngleEntranceAperture -= rotationCenter;
1154 maxAngleExitAperture -= rotationCenter;
1221 ERRORMSG(
"If using \"1DPROFILE1-DEFAULT\" field map you must set the "
1222 "chord length of the bend using the L attribute in the OPAL "
1250 WARNMSG(
"Warning: bend design energy is zero. Treating as drift."
1261 WARNMSG(
"Warning: bend strength and defined reference trajectory "
1262 <<
"chord length are not consistent. Treating bend as drift."
1271 WARNMSG(
"Warning bend angle/strength is zero. Treating as drift."
1281 std::vector<Vector_t> outline;
1283 unsigned int numSteps = 2;
1294 gsl_poly_solve_quadratic(
dot(P,P),
1300 Vector_t upperCornerAtEntry = tau1 * P;
1303 gsl_poly_solve_quadratic(
dot(P,P),
1315 double dAngle = 0.5 * totalAngle / (1.0 * numSteps - 1.0);
1317 outline.push_back(upperCornerAtEntry);
1319 for (
unsigned int i = 0; i < numSteps; ++ i, angle += dAngle) {
1321 outline.push_back(rot.
rotate(upperCornerAtEntry - rotationCenter) + rotationCenter);
1323 outline.push_back(upperCornerAtExit);
1335 gsl_poly_solve_quadratic(
dot(P,P),
1341 Vector_t lowerCornerAtEntry = tau1 * P;
1344 gsl_poly_solve_quadratic(
dot(P,P),
1355 double dAngle = 0.5 * totalAngle / (1.0 * numSteps - 1.0);
1357 outline.push_back(lowerCornerAtExit);
1358 double angle = 0.5 * totalAngle;
1359 for (
unsigned int i = 0; i < numSteps; ++ i, angle -= dAngle) {
1361 outline.push_back(rot.
rotate(lowerCornerAtEntry - rotationCenter) + rotationCenter);
1363 outline.push_back(lowerCornerAtEntry);
1366 std::ofstream outlineOutput;
1372 outlineOutput.open(fname);
1373 outlineOutput.precision(8);
1375 for (
auto a: outline) {
1376 outlineOutput << std::setw(18) <<
a(2)
1377 << std::setw(18) <<
a(0)
1380 outlineOutput << std::setw(18) << outline.front()(2)
1381 << std::setw(18) << outline.front()(0)
1390 std::vector<Vector_t> outline =
getOutline();
1392 unsigned int size = outline.size();
1393 unsigned int last = size - 1;
1394 unsigned int numSteps = (size - 14) / 2;
1395 unsigned int midIdx = numSteps + 7;
1397 for (
unsigned int i = 0; i < 6; ++ i) {
1398 mesh.
vertices_m.push_back(outline[i] + 2 * hgap);
1400 for (
unsigned int i = 6; i < 6 + numSteps; ++ i) {
1401 mesh.
vertices_m.push_back(outline[i] + hgap);
1403 for (
unsigned int i = 6 + numSteps; i < 13 + numSteps; ++ i) {
1404 mesh.
vertices_m.push_back(outline[i] + 2 * hgap);
1406 for (
unsigned int i = 13 + numSteps; i < 13 + 2 * numSteps; ++ i) {
1407 mesh.
vertices_m.push_back(outline[i] + hgap);
1409 mesh.
vertices_m.push_back(outline[last] + 2 * hgap);
1411 for (
unsigned int i = 0; i < 6; ++ i) {
1412 mesh.
vertices_m.push_back(outline[i] - 2 * hgap);
1414 for (
unsigned int i = 6; i < 6 + numSteps; ++ i) {
1415 mesh.
vertices_m.push_back(outline[i] - hgap);
1417 for (
unsigned int i = 6 + numSteps; i < 13 + numSteps; ++ i) {
1418 mesh.
vertices_m.push_back(outline[i] - 2 * hgap);
1420 for (
unsigned int i = 13 + numSteps; i < 13 + 2 * numSteps; ++ i) {
1421 mesh.
vertices_m.push_back(outline[i] - hgap);
1423 mesh.
vertices_m.push_back(outline[last] - 2 * hgap);
1428 mesh.
vertices_m.push_back(outline[midIdx + 4]);
1466 for (
unsigned int i = 6; i < 5 + numSteps; ++ i) {
1475 for (
unsigned int i = 0; i < size - 2; ++ i) {
1476 if (i == 4u || i == 5u ||
1477 i == 5 + numSteps || i == 6 + numSteps ||
1478 i == 11 + numSteps || i == 12 + numSteps)
continue;
1480 unsigned int next = (i + 1) % size;
1519 mesh.
decorations_m.push_back(std::make_pair(0.5 * (P1 + P2) - 0.25 * dir1,
1520 0.5 * (P1 + P2) + 0.25 * dir1));
1524 double det = -dir1[0] * dir2[2] + dir1[2] * dir2[0];
1525 inv(0, 0) = -dir2[2] /
det;
1526 inv(0, 2) = dir2[0] /
det;
1528 inv(2, 0) = -dir1[2] /
det;
1529 inv(2, 2) = dir1[0] /
det;
1531 Vector_t crossPoint = P1 + Tau[0] * dir1;
1535 Vector_t R = crossPoint - rotationCenter;
1539 gsl_poly_solve_quadratic(
dot(P,P),
1545 std::swap(tau1, tau2);
1547 Vector_t lowerCornerFringeLimit = crossPoint + tau1 * P;
1550 gsl_poly_solve_quadratic(
dot(P,P),
1556 std::swap(tau1, tau2);
1558 Vector_t upperCornerFringeLimit = crossPoint + tau1 * P;
1560 mesh.
decorations_m.push_back(std::make_pair(lowerCornerFringeLimit,
1561 upperCornerFringeLimit));
1568 gsl_poly_solve_quadratic(
dot(P,P),
1576 gsl_poly_solve_quadratic(
dot(P,P),
1586 gsl_poly_solve_quadratic(
dot(P,P),
1594 gsl_poly_solve_quadratic(
dot(P,P),
1602 mesh.
decorations_m.push_back(std::make_pair(lowerCornerFringeLimitEntrance, upperCornerFringeLimitEntrance));
1603 mesh.
decorations_m.push_back(std::make_pair(lowerCornerFringeLimitExit, upperCornerFringeLimitExit));
1648 double entrParam1, entrParam2, entrParam3;
1654 std::array<double,2> entFFL;
1655 entFFL[0]=entrParam2-entrParam1;
1656 entFFL[1]=entrParam3-entrParam2;
1668 std::array<double,2> extFFL;
1670 double exitParam1, exitParam2, exitParam3;
1676 extFFL[0]=exitParam3-exitParam2;
1677 extFFL[1]=exitParam2-exitParam1;
1689 std::vector<Vector_t> outline =
getOutline();
1693 for (
int i : {-1, 1}) {
Vector3D cross(const Vector3D &lhs, const Vector3D &rhs)
Vector cross product.
double dot(const Vector3D &lhs, const Vector3D &rhs)
Vector dot product.
Quaternion getQuaternion(Vector_t u, Vector_t ref)
Tps< T > cos(const Tps< T > &x)
Cosine.
Tps< T > pow(const Tps< T > &x, int y)
Integer power.
Tps< T > exp(const Tps< T > &x)
Exponential.
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.
T::PETE_Expr_t::PETE_Return_t max(const PETE_Expr< T > &expr, NDIndex< D > &loc)
T::PETE_Expr_t::PETE_Return_t min(const PETE_Expr< T > &expr, NDIndex< D > &loc)
PETE_TUTree< FnCeil, typename T::PETE_Expr_t > ceil(const PETE_Expr< T > &l)
PETE_TUTree< FnArcSin, typename T::PETE_Expr_t > asin(const PETE_Expr< T > &l)
PETE_TBTree< FnCopysign, PETE_Scalar< Vektor< T1, Dim > >, typename T2::PETE_Expr_t > copysign(const Vektor< T1, Dim > &l, const PETE_Expr< T2 > &r)
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)
T det(const AntiSymTenzor< T, 3 > &)
Inform & level2(Inform &inf)
Inform & endl(Inform &inf)
constexpr double two_pi
The value of.
constexpr double e
The value of.
constexpr double c
The velocity of light in m/s.
constexpr double pi
The value of.
T isinf(T x)
isinf function with adjusted return type
T isnan(T x)
isnan function with adjusted return type
bool writeBendTrajectories
std::string combineFilePath(std::initializer_list< std::string > ilist)
std::string toUpper(const std::string &str)
const PartData * getReference() const
double getQ() const
Access to reference data.
std::string getInputBasename()
get input file name without extension
static OpalData * getInstance()
std::string getAuxiliaryOutputDirectory() const
get the name of the the additional data directory
void calculateRefTrajectory(double &angleX, double &angleY)
virtual bool apply(const size_t &i, const double &t, Vector_t &E, Vector_t &B) override
Apply field to particles with coordinates in magnet frame.
Vector_t transformToExitRegion(const Vector_t &R) const
void readFieldMap(Inform &msg)
int polyOrderExit_m
function origin that Enge function ends.
Vector_t calcExitFringeField(const Vector_t &R, double deltaX)
virtual ElementType getType() const override=0
Get element type std::string.
virtual bool findChordLength(double &chordLength)=0
gsl_interp_accel * entryFieldAccel_m
bool inMagnetCentralRegion(const Vector_t &R) const
std::vector< double > engeCoeffsEntry_m
Enge coefficients for map entry and exit regions.
std::array< double, 2 > getEntranceFringeFieldLength() const
Get entrance fringe field length.
virtual void initialise(PartBunchBase< double, 3 > *bunch, double &startField, double &endField) override
double deltaEndEntry_m
function origin where Enge function starts.
void setExitAngle(double exitAngle)
bool inMagnetExitRegion(const Vector_t &R) const
CoordinateSystemTrafo beginToEnd_m
CoordinateSystemTrafo getBeginToEnd_local() const
double deltaBeginExit_m
Enge function order for entry region.
Vector_t calcCentralField(const Vector_t &R, double deltaX)
void setEngeOriginDelta(double delta)
double entranceParameter3_m
bool calculateMapField(const Vector_t &R, Vector_t &B)
double widthEntranceFringe_m
End of magnet field map in s coordinates (m).
gsl_spline ** exitFieldValues_m
std::vector< double > engeCoeffsExit_m
double entranceParameter1_m
double designRadius_m
through the bend.
double sinEntranceAngle_m
void setBendEffectiveLength(double startField, double endField)
void setNSlices(const std::size_t &nSlices)
Vector_t calcEntranceFringeField(const Vector_t &R, double deltaX)
bool setupBendGeometry(double &startField, double &endField)
bool findIdealBendParameters(double chordLength)
std::string messageHeader_m
double tanEntranceAngle_m
virtual void setEntranceAngle(double entranceAngle) override
double startField_m
Dipole field index.
virtual void goOnline(const double &kineticEnergy) override
double endField_m
Start of magnet field map in s coordinates (m).
CoordinateSystemTrafo toEntranceRegion_m
bool isPositionInExitField(const Vector_t &R) const
void setFieldBoundaries(double startField, double endField)
gsl_interp_accel * exitFieldAccel_m
double deltaEndExit_m
function origin that Enge function starts.
Vector_t transformToEntranceRegion(const Vector_t &R) const
bool inMagnetEntranceRegion(const Vector_t &R) const
void adjustFringeFields(double ratio)
double exitAngle_m
Bend design radius (m).
double fieldIndex_m
and the exit face of the magnet (radians).
MeshData getSurfaceMesh() const
std::array< double, 2 > getExitFringeFieldLength() const
Get exit fringe field length.
CoordinateSystemTrafo beginToEnd_lcs_m
double entranceParameter2_m
int polyOrderEntry_m
function origin that Enge function ends.
double estimateFieldAdjustmentStep(double actualBendAngle)
virtual bool isInside(const Vector_t &r) const override
virtual CoordinateSystemTrafo getEdgeToEnd() const override
bool setupDefaultFieldMap()
BoundingBox getBoundingBoxInLabCoords() const override
std::vector< Vector_t > getOutline() const
CoordinateSystemTrafo toExitRegion_m
gsl_spline ** entryFieldValues_m
bool initializeFieldMap()
CoordinateSystemTrafo computeAngleTrafo_m
void calcEngeFunction(double zNormalized, const std::vector< double > &engeCoeff, int polyOrder, double &engeFunc, double &engeFuncDeriv, double &engeFuncSecDerivNorm)
std::size_t getNSlices() const
void setupPusher(PartBunchBase< double, 3 > *bunch)
void print(Inform &msg, double bendAngleX, double bendAngle)
virtual bool applyToReferenceParticle(const Vector_t &R, const Vector_t &P, const double &t, Vector_t &E, Vector_t &B) override
void findBendEffectiveLength(double startField, double endField)
void setGapFromFieldMap()
bool isPositionInEntranceField(const Vector_t &R) const
double cosEntranceAngle_m
Enge function order for entry region.
double calculateBendAngle()
double calcGamma() const
Calculate gamma from design energy.
std::vector< Vector_t > refTrajMap_m
Map of reference particle trajectory.
const bool fast_m
Flag to turn on fast field calculation.
double calcBendAngle(double chordLength, double radius) const
Calculate bend angle from chord length and design radius.
double fieldAmplitude_m
Field amplitude.
double calcBetaGamma() const
Calculate beta*gamma from design energy.
double getChordLength() const
double designEnergy_m
Bend design energy (eV).
double getBendAngle() const
double gap_m
Full vertical gap of the magnets.
Fieldmap * fieldmap_m
Magnet field map.
double getEntranceAngle() const
double calcFieldAmplitude(double radius) const
Calculate field amplitude from design energy and radius.
double calcDesignRadius(double fieldAmplitude) const
Calculate design radius from design energy and field amplitude.
double getFullGap() const
double angle_m
Bend angle.
PartBunchBase< double, 3 > * RefPartBunch_m
virtual const std::string & getName() const
Get element name.
bool getFlagDeleteOnTransverseExit() const
virtual double getElementLength() const
Get design length.
virtual void setElementLength(double length)
Set design length.
std::pair< ApertureType, std::vector< double > > aperture_m
virtual CoordinateSystemTrafo getEdgeToBegin() const
double getRotationAboutZ() const
CoordinateSystemTrafo csTrafoGlobal2Local_m
Vector_t rotateFrom(const Vector_t &r) const
Vector_t transformFrom(const Vector_t &r) const
Vector_t transformTo(const Vector_t &r) const
Vector_t rotate(const Vector_t &) const
Quaternion conjugate() const
virtual void getInfo(Inform *msg)=0
static Fieldmap * getFieldmap(std::string Filename, bool fast=false)
virtual void get1DProfile1EntranceParam(double &entranceParameter1, double &entranceParameter2, double &entranceParameter3)
virtual double getFieldGap()
virtual void get1DProfile1EngeCoeffs(std::vector< double > &engeCoeffsEntry, std::vector< double > &engeCoeffsExit)
virtual void get1DProfile1ExitParam(double &exitParameter1, double &exitParameter2, double &exitParameter3)
void enlargeToContainPosition(const Vector_t &position)
std::vector< std::pair< Vector_t, Vector_t > > decorations_m
std::vector< Vektor< unsigned int, 3 > > triangles_m
std::vector< Vector_t > vertices_m
void kick(const Vector_t &R, Vector_t &P, const Vector_t &Ef, const Vector_t &Bf, const double &dt) const
void initialise(const PartData *ref)
void push(Vector_t &R, const Vector_t &P, const double &dt) const
Vektor< double, 3 > Vector_t