OPAL (Object Oriented Parallel Accelerator Library)  2024.1
OPAL
MultipoleT.cpp
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2017, Titus Dascalu
3  * Copyright (c) 2018, Martin Duy Tat
4  * Copyright (c) 2019-2023, Chris Rogers
5  * All rights reserved.
6  * Redistribution and use in source and binary forms, with or without
7  * modification, are permitted provided that the following conditions are met:
8  * 1. Redistributions of source code must retain the above copyright notice,
9  * this list of conditions and the following disclaimer.
10  * 2. Redistributions in binary form must reproduce the above copyright notice,
11  * this list of conditions and the following disclaimer in the documentation
12  * and/or other materials provided with the distribution.
13  * 3. Neither the name of STFC nor the names of its contributors may be used to
14  * endorse or promote products derived from this software without specific
15  * prior written permission.
16  *
17  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
18  * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
19  * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
20  * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
21  * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
22  * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
23  * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
24  * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
25  * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
26  * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
27  * POSSIBILITY OF SUCH DAMAGE.
28  */
29 
30 
31 #include "Algorithms/PartBunch.h"
32 #include "MultipoleT.h"
33 #include <gsl/gsl_math.h>
38 
39 using namespace endfieldmodel;
40 
41 MultipoleT::MultipoleT(const std::string &name):
42  Component(name),
43  fringeField_l(endfieldmodel::Tanh()),
44  fringeField_r(endfieldmodel::Tanh()),
45  maxOrder_m(5),
46  maxOrderX_m(10),
47  transMaxOrder_m(0),
48  planarArcGeometry_m(1., 1.),
49  length_m(1.0),
50  angle_m(0.0),
51  entranceAngle_m(0.0),
52  rotation_m(0.0),
53  variableRadius_m(false),
54  boundingBoxLength_m(0.0),
55  verticalApert_m(0.5),
56  horizApert_m(0.5) {
57 }
58 
60  Component(right),
61  fringeField_l(right.fringeField_l),
62  fringeField_r(right.fringeField_r),
63  maxOrder_m(right.maxOrder_m),
64  maxOrderX_m(right.maxOrderX_m),
65  recursion_VarRadius_m(right.recursion_VarRadius_m),
66  recursion_ConstRadius_m(right.recursion_ConstRadius_m),
67  transMaxOrder_m(right.transMaxOrder_m),
68  transProfile_m(right.transProfile_m),
69  planarArcGeometry_m(right.planarArcGeometry_m),
70  length_m(right.length_m),
71  angle_m(right.angle_m),
72  entranceAngle_m(right.entranceAngle_m),
73  rotation_m(right.rotation_m),
74  variableRadius_m(right.variableRadius_m),
75  boundingBoxLength_m(right.boundingBoxLength_m),
76  verticalApert_m(right.verticalApert_m),
77  horizApert_m(right.horizApert_m),
78  dummy() {
80 }
81 
83 }
84 
86  MultipoleT* newMultipole = new MultipoleT(*this);
87  newMultipole->initialise();
88  return newMultipole;
89 }
90 
92  RefPartBunch_m = nullptr;
93 }
94 
95 bool MultipoleT::apply(const Vector_t &R, const Vector_t &/*P*/,
96  const double &/*t*/,Vector_t &/*E*/, Vector_t &B) {
98  Vector_t R_prime = rotateFrame(R);
100  //R_prime[2] *= -1; //OPAL uses different sign convention...
101  Vector_t X = R_prime;
102  R_prime = transformCoords(X);
103  if (insideAperture(R_prime)) {
105  double theta;
106  double prefactor = (length_m / angle_m) * (tanh((fringeField_l.getX0())
107  / fringeField_l.getLambda()) +
109  / fringeField_r.getLambda()));
110  if (angle_m == 0.0) {
111  theta = 0.0;
112  } else if (!variableRadius_m) {
113  theta = R_prime[2] * angle_m / length_m;
114  } else { // variableRadius_m == true
115  theta = fringeField_l.getLambda() * log(cosh((R_prime[2] +
117  fringeField_r.getLambda() * log(cosh((R_prime[2] -
119  theta /= prefactor;
120  }
121  double Bx = getBx(R_prime);
122  double Bs = getBs(R_prime);
123  B[0] = Bx * cos(theta) - Bs * sin(theta);
124  B[2] = Bx * sin(theta) + Bs * cos(theta);
125  B[1] = getBz(R_prime);
126  // B[2] *= -1; //OPAL uses different sign convention
127  return false;
128  } else {
129  for(int i = 0; i < 3; i++) {
130  B[i] = 0.0;
131  }
133  }
134 }
135 
136 bool MultipoleT::apply(const size_t &i, const double &t,
137  Vector_t &E, Vector_t &B) {
138  return apply(RefPartBunch_m->R[i], RefPartBunch_m->P[i], t, E, B);
139 }
140 
142  Vector_t R_prime(3), R_pprime(3);
147  // 1st rotation
148  R_prime[0] = R[0] * cos(rotation_m) + R[1] * sin(rotation_m);
149  R_prime[1] = - R[0] * sin(rotation_m) + R[1] * cos(rotation_m);
150  R_prime[2] = R[2];
151  // 2nd rotation
152  R_pprime[0] = R_prime[2] * sin(entranceAngle_m) +
153  R_prime[0] * cos(entranceAngle_m);
154  R_pprime[1] = R_prime[1];
155  R_pprime[2] = R_prime[2] * cos(entranceAngle_m) -
156  R_prime[0] * sin(entranceAngle_m);
157  return R_pprime;
158 
159 }
160 
166  Vector_t B_prime(3);
167  B_prime[0] = B[0] * cos(rotation_m) + B[1] * sin(rotation_m);
168  B_prime[1] = - B[0] * sin(rotation_m) + B[1] * cos(rotation_m);
169  B_prime[2] = B[2];
170  return B_prime;
171 }
172 
174  if (std::abs(R[1]) > verticalApert_m / 2.) {
175  return false;
176  } else if (std::abs(R[0]) > horizApert_m / 2.) {
177  return false;
178  } else if (R[2] < -boundingBoxLength_m/2.0) {
179  return false;
180  } else if (R[2] > boundingBoxLength_m/2.0) {
181  return false;
182  }
183  return true;
184 }
185 
187  Vector_t X(3);
188  // If radius is constant
189  if (!variableRadius_m) {
190  double radius = length_m / angle_m;
191  // Transform from Cartesian to Frenet-Serret along the magnet
192  double alpha = atan(R[2] / (R[0] + radius ));
193  if (alpha != 0.0 && angle_m != 0.0) {
194  X[0] = R[2] / sin(alpha) - radius;
195  X[1] = R[1];
196  X[2] = radius * alpha;// + boundingBoxLength_m;
197  } else {
198  X[0] = R[0];
199  X[1] = R[1];
200  X[2] = R[2];// + boundingBoxLength_m;
201  }
202  } else {
203  // If radius is variable
204  coordinatetransform::CoordinateTransform t(R[0], R[1], R[2],
208  (length_m / angle_m));
209  std::vector<double> r = t.getTransformation();
210  X[0] = r[0];
211  X[1] = r[1];
212  X[2] = r[2];
213  }
214  return X;
215 }
216 
217 void MultipoleT::setMaxOrder(std::size_t maxOrder) {
218  if (variableRadius_m && angle_m != 0.0) {
219  std::size_t N = recursion_VarRadius_m.size();
220  while (maxOrder >= N) {
221  polynomial::RecursionRelationTwo r(N, 2 * (N + maxOrderX_m + 1));
224  recursion_VarRadius_m.push_back(r);
225  N = recursion_VarRadius_m.size();
226  }
227  } else if (!variableRadius_m && angle_m != 0.0) {
228  std::size_t N = recursion_ConstRadius_m.size();
229  while (maxOrder >= N) {
230  polynomial::RecursionRelation r(N, 2 * (N + maxOrderX_m + 1));
233  recursion_ConstRadius_m.push_back(r);
234  N = recursion_ConstRadius_m.size();
235  }
236  }
237  maxOrder_m = maxOrder;
238 }
239 
240 void MultipoleT::setTransProfile(std::size_t n, double dTn) {
241  if (n > transMaxOrder_m) {
242  transMaxOrder_m = n;
243  transProfile_m.resize(n+1, 0.0);
244  }
245  transProfile_m[n] = dTn;
246 }
247 
248 bool MultipoleT::setFringeField(double s0, double lambda_l, double lambda_r) {
249  fringeField_l.Tanh::setLambda(lambda_l);
250  fringeField_l.Tanh::setX0(s0);
251  fringeField_l.Tanh::setTanhDiffIndices(2 * maxOrder_m + 1);
252 
253  fringeField_r.Tanh::setLambda(lambda_r);
254  fringeField_r.Tanh::setX0(s0);
255  fringeField_r.Tanh::setTanhDiffIndices(2 * maxOrder_m + 1);
256 
257  return true;
258 }
259 
260 double MultipoleT::getBz(const Vector_t &R) {
264  double Bz = 0.0;
265  if (angle_m == 0.0) {
266  // Straight geometry -> Use corresponding field expansion directly
267  for(std::size_t n = 0; n <= maxOrder_m; n++) {
268  double f_n = 0.0;
269  for(std::size_t i = 0; i <= n; i++) {
270  f_n += gsl_sf_choose(n, i) * getTransDeriv(2 * i, R[0]) *
271  getFringeDeriv(2 * n - 2 * i, R[2]);
272  }
273  f_n *= gsl_sf_pow_int(-1.0, n);
274  Bz += gsl_sf_pow_int(R[1], 2 * n) / gsl_sf_fact(2 * n) * f_n;
275  }
276  } else {
277  if (variableRadius_m == true and getFringeDeriv(0, R[2]) < 1.0e-12) {
278  // Return 0 if end of fringe field is reached
279  // This is to avoid functions being called at infinite radius
280  return 0.0;
281  }
282  // Curved geometry -> Use full machinery of differential operators
283  for(std::size_t n = 0; n <= maxOrder_m; n++) {
284  double f_n = getFn(n, R[0], R[2]);
285  Bz += gsl_sf_pow_int(R[1], 2 * n) / gsl_sf_fact(2 * n) * f_n;
286  }
287  }
288  return Bz;
289 }
290 
291 double MultipoleT::getBx(const Vector_t &R) {
295  double Bx = 0.0;
296  if (angle_m == 0.0) {
297  // Straight geometry -> Use corresponding field expansion directly
298  for(std::size_t n = 0; n <= maxOrder_m; n++) {
299  double f_n = 0.0;
300  for(std::size_t i = 0; i <= n; i++) {
301  f_n += gsl_sf_choose(n, i) * getTransDeriv(2 * i + 1, R[0]) *
302  getFringeDeriv(2 * n - 2 * i, R[2]);
303  }
304  f_n *= gsl_sf_pow_int(-1.0, n);
305  Bx += gsl_sf_pow_int(R[1], 2 * n + 1) /
306  gsl_sf_fact(2 * n + 1) * f_n;
307  }
308  } else {
309  if (variableRadius_m == true and getFringeDeriv(0, R[2]) < 1.0e-12) {
310  // Return 0 if end of fringe field is reached
311  // This is to avoid functions being called at infinite radius
312  return 0.0;
313  }
314  // Curved geometry -> Use full machinery of differential operators
315  for(std::size_t n = 0; n <= maxOrder_m; n++) {
316  double partialX_fn = getFnDerivX(n, R[0], R[2]);
317  Bx += partialX_fn * gsl_sf_pow_int(R[1], 2 * n + 1)
318  / gsl_sf_fact(2 * n + 1);
319  }
320  }
321  return Bx;
322 }
323 
324 double MultipoleT::getBs(const Vector_t &R) {
328  double Bs = 0.0;
329  if (angle_m == 0.0) {
330  // Straight geometry -> Use corresponding field expansion directly
331  for(std::size_t n = 0; n <= maxOrder_m; n++) {
332  double f_n = 0.0;
333  for(std::size_t i = 0; i <= n; i++) {
334  f_n += gsl_sf_choose(n, i) * getTransDeriv(2 * i, R[0]) *
335  getFringeDeriv(2 * n - 2 * i + 1, R[2]);
336  }
337  f_n *= gsl_sf_pow_int(-1.0, n);
338  Bs += gsl_sf_pow_int(R[1], 2 * n + 1) /
339  gsl_sf_fact(2 * n + 1) * f_n;
340  }
341  } else {
342  if (variableRadius_m == true and getFringeDeriv(0, R[2]) < 1.0e-12) {
343  // Return 0 if end of fringe field is reached
344  // This is to avoid functions being called at infinite radius
345  return 0.0;
346  }
347  // Curved geometry -> Use full machinery of differential operators
348  for(std::size_t n = 0; n <= maxOrder_m; n++) {
349  double partialS_fn = getFnDerivS(n, R[0], R[2]);
350  Bs += partialS_fn * gsl_sf_pow_int(R[1], 2 * n + 1)
351  / gsl_sf_fact(2 * n + 1);
352  }
353  Bs /= getScaleFactor(R[0], R[2]);
354  }
355  return Bs;
356 }
357 
358 double MultipoleT::getFringeDeriv(int n, double s) {
359  if (n <= 10) {
360  return (fringeField_l.getTanh(s, n) - fringeField_r.getNegTanh(s, n))
361  / 2;
362  } else {
363  return tanhderiv::integrate(s,
364  fringeField_l.Tanh::getX0(),
365  fringeField_l.Tanh::getLambda(),
366  fringeField_r.Tanh::getLambda(),
367  n);
368  }
369 }
370 
371 double MultipoleT::getTransDeriv(std::size_t n, double x) {
377  double func = 0;
378  std::vector<double> temp = transProfile_m;
379  if (n <= transMaxOrder_m) {
380  if (n != 0) {
381  for(std::size_t i = 1; i <= n; i++) {
382  for(std::size_t j = 0; j <= transMaxOrder_m; j++) {
383  if (j <= transMaxOrder_m - i) {
384  // Move terms to the left and multiply by power
385  temp[j] = temp[j + 1] * (j + 1);
386  } else {
387  // Put 0 at the end for missing higher powers
388  temp[j] = 0.0;
389  }
390  }
391  }
392  }
393  // Now use the vector to calculate value of the function
394  for(std::size_t k = 0; k <= transMaxOrder_m; k++) {
395  func += temp[k] * gsl_sf_pow_int(x, k);
396  }
397  }
398  return func;
399 }
400 
402  if (transMaxOrder_m < 1) {
403  transProfile_m.resize(1, 0.);
404  }
405  transProfile_m[0] = B0;
406 }
407 
408 void MultipoleT::accept(BeamlineVisitor& visitor) const {
409  visitor.visitMultipoleT(*this);
410 }
411 
412 void MultipoleT::getDimensions(double &/*zBegin*/, double &/*zEnd*/) const {
413 }
414 
415 void MultipoleT::setAperture(double vertAp, double horizAp) {
416  verticalApert_m = vertAp;
417  horizApert_m = horizAp;
418 }
419 
420 std::vector<double> MultipoleT::getAperture() const {
421  std::vector<double> temp(2, 0.0);
422  temp[0] = verticalApert_m;
423  temp[1] = horizApert_m;
424  return temp;
425 }
426 
427 std::vector<double> MultipoleT::getFringeLength() const {
428  std::vector<double> temp(2, 0.0);
429  temp[0] = fringeField_l.getLambda();
430  temp[1] = fringeField_r.getLambda();
431  return temp;
432 }
433 
434 double MultipoleT::getRadius(double s) {
435  double centralRadius = length_m / angle_m;
436  if (!variableRadius_m) {
437  return centralRadius;
438  }
439  if (getFringeDeriv(0, s) != 0.0) {
440  return centralRadius * getFringeDeriv(0, 0) / getFringeDeriv(0, s);
441  } else {
442  return -1; // Return -1 if radius is infinite
443  }
444 }
445 
446 double MultipoleT::getScaleFactor(double x, double s) {
447  if (!variableRadius_m) {
448  return 1 + x * angle_m/ length_m;
449  } else {
450  return 1 + x / getRadius(s);
451  }
452 }
453 
454 double MultipoleT::getFnDerivX(std::size_t n, double x, double s) {
455  if (n == 0) {
456  return getTransDeriv(1, x) * getFringeDeriv(0, s);
457  }
458  double deriv = 0.0;
459  double stepSize = 1e-3;
460  deriv += 1. * getFn(n, x - 2. * stepSize, s);
461  deriv += -8. * getFn(n, x - stepSize, s);
462  deriv += 8. * getFn(n, x + stepSize, s);
463  deriv += -1. * getFn(n, x + 2. * stepSize, s);
464  deriv /= 12 * stepSize;
465  return deriv;
466 }
467 
468 double MultipoleT::getFnDerivS(std::size_t n, double x, double s) {
469  if (n == 0) {
470  return getTransDeriv(0, x) * getFringeDeriv(1, s);
471  }
472  double deriv = 0.0;
473  double stepSize = 1e-3;
474  deriv += 1. * getFn(n, x, s - 2. * stepSize);
475  deriv += -8. * getFn(n, x, s - stepSize);
476  deriv += 8. * getFn(n, x, s + stepSize);
477  deriv += -1. * getFn(n, x, s + 2. * stepSize);
478  deriv /= 12 * stepSize;
479  return deriv;
480 }
481 
482 double MultipoleT::getFn(std::size_t n, double x, double s) {
483  if (n == 0) {
484  return getTransDeriv(0, x) * getFringeDeriv(0, s);
485  }
486  if (!variableRadius_m) {
487  double rho = length_m / angle_m;
488  double func = 0.0;
489 
490  for (std::size_t j = 0;
491  j <= recursion_ConstRadius_m.at(n).getMaxSDerivatives();
492  j++) {
493  double FringeDerivj = getFringeDeriv(2 * j, s);
494  for (std::size_t i = 0;
495  i <= recursion_ConstRadius_m.at(n).getMaxXDerivatives();
496  i++) {
497  if (recursion_ConstRadius_m.at(n).isPolynomialZero(i, j)) {
498  continue;
499  }
500  func += (recursion_ConstRadius_m.at(n)
501  .evaluatePolynomial(x / rho, i, j) *
502  getTransDeriv(i, x) * FringeDerivj) /
503  gsl_sf_pow_int(rho, 2 * n - i - 2 * j);
504  }
505  }
506  func *= gsl_sf_pow_int(-1.0, n);
507  return func;
508  } else {
509  double rho = length_m / angle_m;
510  double S_0 = getFringeDeriv(0, 0);
511  double y = getFringeDeriv(0, s) / (S_0 * rho);
512  double func = 0.0;
513  std::vector<double> fringeDerivatives;
514  for (std::size_t j = 0;
515  j <= recursion_VarRadius_m.at(n).getMaxSDerivatives();
516  j++) {
517  fringeDerivatives.push_back(getFringeDeriv(j, s) / (S_0 * rho));
518  }
519  for (std::size_t i = 0;
520  i <= recursion_VarRadius_m.at(n).getMaxXDerivatives();
521  i++) {
522  double temp = 0.0;
523  for (std::size_t j = 0;
524  j <= recursion_VarRadius_m.at(n).getMaxSDerivatives();
525  j++) {
526  temp += recursion_VarRadius_m.at(n)
527  .evaluatePolynomial(x, y, i, j, fringeDerivatives)
528  * fringeDerivatives.at(j);
529  }
530  func += temp * getTransDeriv(i, x);
531  }
532  func *= gsl_sf_pow_int(-1.0, n) * S_0 * rho;
533  return func;
534  }
535 }
536 
540 }
541 
543  double &/*startField*/,
544  double &/*endField*/) {
545  RefPartBunch_m = bunch;
546  initialise();
547 }
548 
549 bool MultipoleT::bends() const {
550  return (transProfile_m[0] != 0);
551 }
552 
554  return planarArcGeometry_m;
555 }
556 
558  return planarArcGeometry_m;
559 }
560 
562  return dummy;
563 }
564 
566  return dummy;
567 }
void truncate(std::size_t highestXorder)
std::vector< double > getFringeLength() const
Definition: MultipoleT.cpp:427
std::vector< double > transProfile_m
Definition: MultipoleT.h:265
BMultipoleField dummy
Definition: MultipoleT.h:299
double length_m
Definition: MultipoleT.h:282
void finalise() override
Definition: MultipoleT.cpp:91
double rotation_m
Definition: MultipoleT.h:285
virtual void setElementLength(double)
Set length.
double getFn(std::size_t n, double x, double s)
Definition: MultipoleT.cpp:482
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
Abstract base class for electromagnetic fields.
Definition: EMField.h:188
double getTransDeriv(std::size_t n, double x)
Definition: MultipoleT.cpp:371
ParticleAttrib< Vector_t > P
void setDipoleConstant(double B0)
Definition: MultipoleT.cpp:401
and that you know you can do these things To protect your we need to make restrictions that forbid anyone to deny you these rights or to ask you to surrender the rights These restrictions translate to certain responsibilities for you if you distribute copies of the or if you modify it For if you distribute copies of such a whether gratis or for a you must give the recipients all the rights that you have You must make sure that receive or can get the source code And you must show them these terms so they know their rights We protect your rights with two and(2) offer you this license which gives you legal permission to copy
double horizApert_m
Definition: MultipoleT.h:297
ElementBase * clone() const override
Definition: MultipoleT.cpp:85
double getTanh(double x, int n) const
Definition: Tanh.cpp:50
double getFringeDeriv(int n, double s)
Definition: MultipoleT.cpp:358
A simple arc in the XZ plane.
PartBunchBase< double, 3 > * RefPartBunch_m
Definition: Component.h:191
void resizeX(const std::size_t &xDerivatives)
void setCurvature(double)
Set curvature.
Vector_t rotateFrame(const Vector_t &R)
Definition: MultipoleT.cpp:141
bool variableRadius_m
Definition: MultipoleT.h:287
double getScaleFactor(double x, double s)
Definition: MultipoleT.cpp:446
PETE_TUTree< FnArcTan, typename T::PETE_Expr_t > atan(const PETE_Expr< T > &l)
Definition: PETE.h:727
endfieldmodel::Tanh fringeField_l
Definition: MultipoleT.h:251
void setMaxOrder(std::size_t maxOrder)
Definition: MultipoleT.cpp:217
Tps< T > tanh(const Tps< T > &x)
Hyperbolic tangent.
Definition: TpsMath.h:240
bool bends() const override
Definition: MultipoleT.cpp:549
double getLambda() const
Definition: Tanh.h:101
void truncate(std::size_t highestXorder)
bool getFlagDeleteOnTransverseExit() const
Definition: ElementBase.h:614
bool apply(const Vector_t &R, const Vector_t &P, const double &t, Vector_t &E, Vector_t &B) override
Definition: MultipoleT.cpp:95
void setAperture(double vertAp, double horizAp)
Definition: MultipoleT.cpp:415
double getBz(const Vector_t &R)
Definition: MultipoleT.cpp:260
std::size_t maxOrder_m
Definition: MultipoleT.h:256
double getNegTanh(double x, int n) const
Definition: Tanh.cpp:61
EMField & getField() override
Definition: MultipoleT.cpp:561
double getFnDerivS(std::size_t n, double x, double s)
Definition: MultipoleT.cpp:468
endfieldmodel::Tanh fringeField_r
Definition: MultipoleT.h:252
std::vector< double > getTransformation() const
double getX0() const
Definition: Tanh.h:104
#define X(arg)
Definition: fftpack.cpp:112
virtual void visitMultipoleT(const MultipoleT &)=0
Apply the algorithm to an arbitrary multipole.
Vector_t transformCoords(const Vector_t &R)
Definition: MultipoleT.cpp:186
double getBx(const Vector_t &R)
Definition: MultipoleT.cpp:291
double angle_m
Definition: MultipoleT.h:283
double getRadius(double s)
Definition: MultipoleT.cpp:434
Tps< T > cos(const Tps< T > &x)
Cosine.
Definition: TpsMath.h:129
std::size_t maxOrderX_m
Definition: MultipoleT.h:258
constexpr double alpha
The fine structure constant, no dimension.
Definition: Physics.h:72
double entranceAngle_m
Definition: MultipoleT.h:284
void initialise()
Definition: MultipoleT.cpp:537
Tps< T > log(const Tps< T > &x)
Natural logarithm.
Definition: TpsMath.h:182
std::vector< double > getAperture() const
Definition: MultipoleT.cpp:420
const std::string name
ParticlePos_t & R
double getBs(const Vector_t &R)
Definition: MultipoleT.cpp:324
MultipoleT(const std::string &name)
Definition: MultipoleT.cpp:41
double boundingBoxLength_m
Definition: MultipoleT.h:289
bool insideAperture(const Vector_t &R)
Definition: MultipoleT.cpp:173
std::vector< polynomial::RecursionRelation > recursion_ConstRadius_m
Definition: MultipoleT.h:261
double verticalApert_m
Definition: MultipoleT.h:296
bool setFringeField(double s0, double lambda_left, double lambda_right)
Definition: MultipoleT.cpp:248
constexpr double e
The value of .
Definition: Physics.h:39
std::size_t transMaxOrder_m
Definition: MultipoleT.h:263
double getFnDerivX(std::size_t n, double x, double s)
Definition: MultipoleT.cpp:454
void getDimensions(double &zBegin, double &zEnd) const override
Definition: MultipoleT.cpp:412
Interface for a single beam element.
Definition: Component.h:50
Tps< T > sin(const Tps< T > &x)
Sine.
Definition: TpsMath.h:111
Tps< T > cosh(const Tps< T > &x)
Hyperbolic cosine.
Definition: TpsMath.h:222
void setTransProfile(std::size_t n, double Bn)
Definition: MultipoleT.cpp:240
std::vector< polynomial::RecursionRelationTwo > recursion_VarRadius_m
Definition: MultipoleT.h:260
void resizeX(const std::size_t &xDerivatives)
Vector_t rotateFrameInverse(Vector_t &B)
Definition: MultipoleT.cpp:161
PlanarArcGeometry planarArcGeometry_m
Definition: MultipoleT.h:267
double integrate(const double &a, const double &s0, const double &lambdaleft, const double &lambdaright, const int &n)
Definition: tanhDeriv.cpp:66
void initialise(PartBunchBase< double, 3 > *, double &startField, double &endField) override
Definition: MultipoleT.cpp:542
PlanarArcGeometry & getGeometry() override
Definition: MultipoleT.cpp:553
void accept(BeamlineVisitor &visitor) const override
Definition: MultipoleT.cpp:408