OPAL (Object Oriented Parallel Accelerator Library) 2022.1
OPAL
MultipoleTCurvedVarRadius.cpp
Go to the documentation of this file.
1/*
2 * Copyright (c) 2017, Titus Dascalu
3 * Copyright (c) 2018, Martin Duy Tat
4 * All rights reserved.
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 * 1. Redistributions of source code must retain the above copyright notice,
8 * this list of conditions and the following disclaimer.
9 * 2. Redistributions in binary form must reproduce the above copyright notice,
10 * this list of conditions and the following disclaimer in the documentation
11 * and/or other materials provided with the distribution.
12 * 3. Neither the name of STFC nor the names of its contributors may be used to
13 * endorse or promote products derived from this software without specific
14 * prior written permission.
15 *
16 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
17 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
18 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
19 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE
20 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
21 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
22 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
23 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
24 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
25 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
26 * POSSIBILITY OF SUCH DAMAGE.
27 */
28
29
30#include <vector>
31#include "gsl/gsl_sf_pow_int.h"
34
35using namespace endfieldmodel;
36
39 maxOrderX_m(10),
40 varRadiusGeometry_m(1.0, 1.0, 1.0, 1.0, 1.0),
41 angle_m(0.0) {
42}
43
45 const MultipoleTCurvedVarRadius &right):
46 MultipoleTBase(right),
47 maxOrderX_m(right.maxOrderX_m),
48 recursion_m(right.recursion_m),
49 varRadiusGeometry_m(right.varRadiusGeometry_m),
50 angle_m(right.angle_m) {
52}
53
55}
56
58 return new MultipoleTCurvedVarRadius(*this);
59}
60
62 std::vector<double> fringeLength = getFringeLength();
64 getLength() / 2,
65 fringeLength[0],
66 fringeLength[1],
67 (getLength() / angle_m));
68 std::vector<double> r = t.getTransformation();
69 R[0] = r[0];
70 R[1] = r[1];
71 R[2] = r[2];
72}
73
75 const Vector_t &R) {
76 double length = getLength();
77 std::vector<double> fringeLength = getFringeLength();
78 double prefactor = (length / angle_m) *
79 (tanh((length / 2) / fringeLength[0])
80 + tanh((length / 2) / fringeLength[1]));
81 double theta = fringeLength[0] * log(cosh((R[2]
82 + (length / 2)) / fringeLength[0]))
83 - fringeLength[1] * log(cosh((R[2]
84 - (length / 2)) / fringeLength[1]));
85 theta /= prefactor;
86 double Bx = B[0], Bs = B[2];
87 B[0] = Bx * cos(theta) - Bs * sin(theta);
88 B[2] = Bx * sin(theta) + Bs * cos(theta);
89}
90
91void MultipoleTCurvedVarRadius::setMaxOrder(const std::size_t &maxOrder) {
93 std::size_t N = recursion_m.size();
94 while (maxOrder >= N) {
98 recursion_m.push_back(r);
99 N = recursion_m.size();
100 }
101}
102
104 if (getFringeDeriv(0, s) > 1.0e-12 && angle_m != 0.0) {
105 return getLength() * getFringeDeriv(0, 0)
106 / (getFringeDeriv(0, s) * angle_m);
107 } else {
108 return 1e300; // Return -1 if radius is infinite
109 }
110}
111
113 const double &s) {
114 return (1 + x / getRadius(s));
115}
116
117double MultipoleTCurvedVarRadius::getFn(const std::size_t &n,
118 const double &x,
119 const double &s) {
120 if (n == 0) {
121 return getTransDeriv(0, x) * getFringeDeriv(0, s);
122 }
123 double rho = getLength() / angle_m;
124 double S_0 = getFringeDeriv(0, 0);
125 double y = getFringeDeriv(0, s) / (S_0 * rho);
126 double func = 0.0;
127 std::vector<double> fringeDerivatives;
128 for (std::size_t j = 0;
129 j <= recursion_m.at(n).getMaxSDerivatives();
130 j++) {
131 fringeDerivatives.push_back(getFringeDeriv(j, s) / (S_0 * rho));
132 }
133 for (std::size_t i = 0;
134 i <= recursion_m.at(n).getMaxXDerivatives();
135 i++) {
136 double temp = 0.0;
137 for (std::size_t j = 0;
138 j <= recursion_m.at(n).getMaxSDerivatives();
139 j++) {
140 temp += recursion_m.at(n)
141 .evaluatePolynomial(x, y, i, j, fringeDerivatives)
142 * fringeDerivatives.at(j);
143 }
144 func += temp * getTransDeriv(i, x);
145 }
146 func *= gsl_sf_pow_int(-1.0, n) * S_0 * rho;
147 return func;
148}
Tps< T > log(const Tps< T > &x)
Natural logarithm.
Definition: TpsMath.h:182
Tps< T > cos(const Tps< T > &x)
Cosine.
Definition: TpsMath.h:129
Tps< T > cosh(const Tps< T > &x)
Hyperbolic cosine.
Definition: TpsMath.h:222
Tps< T > tanh(const Tps< T > &x)
Hyperbolic tangent.
Definition: TpsMath.h:240
Tps< T > sin(const Tps< T > &x)
Sine.
Definition: TpsMath.h:111
const std::string name
PartBunchBase< double, 3 > * RefPartBunch_m
Definition: Component.h:191
std::size_t getTransMaxOrder() const
std::vector< double > getFringeLength() const
virtual void setMaxOrder(const std::size_t &maxOrder)
double getLength() const
double getTransDeriv(const std::size_t &n, const double &x)
double getFringeDeriv(const std::size_t &n, const double &s)
std::vector< polynomial::RecursionRelationTwo > recursion_m
virtual double getScaleFactor(const double &x, const double &s) override
virtual void transformBField(Vector_t &B, const Vector_t &R) override
virtual double getRadius(const double &s) override
virtual void setMaxOrder(const std::size_t &maxOrder) override
virtual ElementBase * clone() const override
virtual void transformCoords(Vector_t &R) override
virtual double getFn(const std::size_t &n, const double &x, const double &s) override
MultipoleTCurvedVarRadius(const std::string &name)
std::vector< double > getTransformation() const
void truncate(std::size_t highestXorder)
void resizeX(const std::size_t &xDerivatives)