ClosedOrbitFinder.h 29.5 KB
Newer Older
1 2 3 4
/**
 * @file ClosedOrbitFinder.h
 * The algorithm is based on the paper of M. M. Gordon: "Computation of closed orbits and basic focusing properties for
 * sector-focused cyclotrons and the design of 'cyclops'" (1983)
5 6
 * As template arguments one chooses the type of the variables and the integrator for the ODEs. The supported steppers can
 * be found on
7 8 9 10 11 12
 * http://www.boost.org/ where it is part of the library Odeint.
 *
 * @author Matthias Frey
 * @version 1.0
 */

13 14 15
#ifndef CLOSEDORBITFINDER_H
#define CLOSEDORBITFINDER_H

16
#include <algorithm>
17 18 19
#include <array>
#include <cmath>
#include <functional>
adelmann's avatar
adelmann committed
20
#include <limits>
21
#include <numeric>
adelmann's avatar
adelmann committed
22
#include <string>
23
#include <utility>
24 25
#include <vector>

26
#include "Utilities/Options.h"
27 28 29
#include "Utilities/Options.h"
#include "Utilities/OpalException.h"

30
#include "AbsBeamline/Cyclotron.h"
frey_m's avatar
frey_m committed
31

32 33 34
// include headers for integration
#include <boost/numeric/odeint/integrate/integrate_n_steps.hpp>

35
/// Finds a closed orbit of a cyclotron for a given energy
36 37 38
template<typename Value_type, typename Size_type, class Stepper>
class ClosedOrbitFinder
{
39 40 41 42 43 44 45 46 47 48 49 50 51
    public:
        /// Type of variables
        typedef Value_type value_type;
        /// Type for specifying sizes
        typedef Size_type size_type;
        /// Type of container for storing quantities (path length, orbit, etc.)
        typedef std::vector<value_type> container_type;
        /// Type for holding state of ODE values
        typedef std::vector<value_type> state_type;

        /// Sets the initial values for the integration and calls findOrbit().
        /*!
         * @param E is the energy [MeV] to which the closed orbit should be found
52
         * @param E0 is the potential energy (particle energy at rest) [MeV].
adelmann's avatar
adelmann committed
53
         * @param N specifies the number of splits (2pi/N), i.e number of integration steps
54
         * @param cycl is the cyclotron element
55 56
         * @param domain is a boolean (default: true). If "true" the closed orbit is computed over a single sector,
         * otherwise over 2*pi.
57
         */
58
        ClosedOrbitFinder(value_type E, value_type E0, size_type N,
59
                          Cyclotron* cycl, bool domain = true);
60 61

        /// Returns the inverse bending radius (size of container N+1)
62
        container_type getInverseBendingRadius(const value_type& angle = 0);
63 64

        /// Returns the step lengths of the path (size of container N+1)
65
        container_type getPathLength(const value_type& angle = 0);
66 67

        /// Returns the field index (size of container N+1)
68
        container_type getFieldIndex(const value_type& angle = 0);
69 70 71 72

        /// Returns the radial and vertical tunes (in that order)
        std::pair<value_type,value_type> getTunes();

73 74 75 76 77 78 79
        /// Returns the closed orbit (size of container N+1) starting at specific angle (only makes sense when computing
        /// the closed orbit for a whole turn) (default value: 0°).
        /// Attention: It computes the starting index of the array. If it's not an integer it just cuts the floating point
        /// part, i.e. it takes the next starting index below. There's no interpolation of the radius.
        /*!
         * @param angle is the start angle for the output. Has to be within [0°,360°[ (default: 0°).
         */
80 81
        container_type getOrbit(value_type angle = 0);

82 83 84 85 86 87
        /// Returns the momentum of the orbit (size of container N+1)starting at specific angle (only makes sense when
        /// computing the closed orbit for a whole turn) (default value: 0°), \f$ \left[ p_{r} \right] = \si{m}\f$.
        /// Attention: It computes the starting index of the array. If it's not an integer it just cuts the floating point
        /// part, i.e. it takes the next starting index below. There's no interpolation of the momentum.
        /*!
         * @param angle is the start angle for the output. Has to be within [0°,360°[ (default: 0°).
88
         * @returns the momentum in \f$ \beta * \gamma \f$ units
89
         */
90
        container_type getMomentum(value_type angle = 0);
91 92 93 94 95 96 97

        /// Returns the relativistic factor gamma
        value_type getGamma();

        /// Returns the average orbit radius
        value_type getAverageRadius();

adelmann's avatar
adelmann committed
98 99
        /// Returns the frequency error
        value_type getFrequencyError();
100 101 102 103 104 105 106 107

        /// Returns true if a closed orbit could be found
        bool isConverged();

        /// Computes the closed orbit
        /*!
         * @param accuracy specifies the accuracy of the closed orbit
         * @param maxit is the maximal number of iterations done for finding the closed orbit
frey_m's avatar
frey_m committed
108
         * @param energy step increase [MeV]
109
         * @param rguess initial radius guess in [mm]
110
         */
frey_m's avatar
frey_m committed
111
        bool findOrbit(value_type, size_type, value_type = 1.0, value_type = -1.0);
112

113
        /// Fills in the values of h_m, ds_m, fidx_m.
114 115
        void computeOrbitProperties();

116
    private:
117 118 119 120 121 122 123 124
        /// This function is called by the function getTunes().
        /*! Transfer matrix Y = [y11, y12; y21, y22] (see Gordon paper for more details).
         * @param y are the positions (elements y11 and y12 of Y)
         * @param py2 is the momentum of the second solution (element y22 of Y)
         * @param ncross is the number of sign changes (\#crossings of zero-line)
         */
        value_type computeTune(const std::array<value_type,2>&, value_type, size_type);

adelmann's avatar
adelmann committed
125
        /// This function computes nzcross_ which is used to compute the tune in z-direction and the frequency error
126
//         void computeVerticalOscillations();
127 128 129
        
        /// This function rotates the calculated closed orbit finder properties to the initial angle
        container_type rotate(value_type angle, container_type& orbitProperty);
130 131 132 133 134

        /// Stores current position in horizontal direction for the solutions of the ODE with different initial values
        std::array<value_type,2> x_m; // x_m = [x1, x2]
        /// Stores current momenta in horizontal direction for the solutions of the ODE with different initial values
        std::array<value_type,2> px_m; // px_m = [px1, px2]
frey_m's avatar
frey_m committed
135
        /// Stores current position in vertical direction for the solutions of the ODE with different initial values
136
        std::array<value_type,2> z_m; // z_m = [z1, z2]
frey_m's avatar
frey_m committed
137
        /// Stores current momenta in vertical direction for the solutions of the ODE with different initial values
138 139 140 141 142 143 144 145
        std::array<value_type,2> pz_m; // pz_m = [pz1, pz2]

        /// Stores the inverse bending radius
        container_type h_m;
        /// Stores the step length
        container_type ds_m;
        /// Stores the radial orbit (size: N_m+1)
        container_type r_m;
146 147
        /// Stores the vertical oribt (size: N_m+1)
        container_type vz_m;
148 149
        /// Stores the radial momentum
        container_type pr_m;
150 151
        /// Stores the vertical momentum
        container_type vpz_m;
152 153 154 155 156 157 158 159 160 161
        /// Stores the field index
        container_type fidx_m;

        /// Counts the number of zero-line crossings in horizontal direction (used for computing horizontal tune)
        size_type nxcross_m;
        /// Counts the number of zero-line crossings in vertical direction (used for computing vertical tune)
        size_type nzcross_m; //#crossings of zero-line in x- and z-direction

        /// Is the energy for which the closed orbit should be found
        value_type E_m;
162 163 164 165
        
        /// Is the potential energy [MeV]
        value_type E0_m;
        
166
        /// Is the nominal orbital frequency
167 168 169
        /* (see paper of Dr. C. Baumgarten: "Transverse-Longitudinal
         * Coupling by Space Charge in Cyclotrons" (2012), formula (1))
         */
170
        value_type wo_m;
adelmann's avatar
adelmann committed
171
        /// Number of integration steps
172 173 174 175 176 177 178 179 180 181 182 183 184
        size_type N_m;
        /// Is the angle step size
        value_type dtheta_m;

        /// Is the relativistic factor
        value_type gamma_m;

        /// Is the average radius
        value_type ravg_m;

        /// Is the phase
        value_type phase_m;

185
        /**
186 187 188 189
         * Stores the last orbit value (since we have to return to the beginning to check the convergence in the
         * findOrbit() function. This last value is then deleted from the array but is stored in lastOrbitVal_m to
         * compute the vertical oscillations)
         */
190 191
        value_type lastOrbitVal_m;

192 193 194 195 196
        /**
         * Stores the last momentum value (since we have to return to the beginning to check the convergence in the
         * findOrbit() function. This last value is then deleted from the array but is stored in lastMomentumVal_m to
         * compute the vertical oscillations)
         */
197
        value_type lastMomentumVal_m;
198 199

        /**
200 201 202
         * Boolean which is true by default. "true": orbit integration over one sector only, "false": integration
         * over 2*pi
         */
adelmann's avatar
adelmann committed
203
        bool domain_m;
204

205 206
        /// Defines the stepper for integration of the ODE's
        Stepper stepper_m;
207 208 209 210 211 212 213 214 215 216 217 218 219 220 221
        
        /*!
         * This quantity is defined in the paper "Transverse-Longitudinal Coupling by Space Charge in Cyclotrons" 
         * of Dr. Christian Baumgarten (2012)
         * The lambda function takes the orbital frequency \f$ \omega_{o} \f$ (also defined in paper) as input argument.
         */
        std::function<double(double)> acon_m = [](double wo) { return Physics::c / wo; };
        
        /// Cyclotron unit \f$ \left[T\right] \f$ (Tesla)
        /*!
         * The lambda function takes the orbital frequency \f$ \omega_{o} \f$ as input argument.
         */
        std::function<double(double, double)> bcon_m = [](double e0, double wo) {
            return e0 * 1.0e7 / (/* physics::q0 */ 1.0 * Physics::c * Physics::c / wo);
        };
222
        
223
        Cyclotron* cycl_m;
224 225 226 227 228 229
};

// -----------------------------------------------------------------------------------------------------------------------
// PUBLIC MEMBER FUNCTIONS
// -----------------------------------------------------------------------------------------------------------------------

230 231 232 233
template<typename Value_type, typename Size_type, class Stepper>
ClosedOrbitFinder<Value_type,
                  Size_type,
                  Stepper>::ClosedOrbitFinder(value_type E, value_type E0,
234
                                              size_type N, Cyclotron* cycl,
235
                                              bool domain)
236 237 238 239
    : nxcross_m(0)
    , nzcross_m(0)
    , E_m(E)
    , E0_m(E0)
240
    , wo_m(cycl->getRfFrequ()*1E6/cycl->getCyclHarm()*2.0*Physics::pi)
241 242 243 244 245 246 247 248 249
    , N_m(N)
    , dtheta_m(Physics::two_pi/value_type(N))
    , gamma_m(E/E0+1.0)
    , ravg_m(0)
    , phase_m(0)
    , lastOrbitVal_m(0.0)
    , lastMomentumVal_m(0.0)
    , domain_m(domain)
    , stepper_m()
250
    , cycl_m(cycl)
251
{
frey_m's avatar
frey_m committed
252
    
253
    if ( cycl_m->getFMLowE() > cycl_m->getFMHighE() )
frey_m's avatar
frey_m committed
254 255 256 257
        throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()",
                            "Incorrect cyclotron energy (MeV) bounds: Maximum cyclotron energy smaller than minimum cyclotron energy.");
    
//     // Even if the numbers are equal --> if statement is true.
258
//     if ( E_m < cycl_m->getFMLowE() )
frey_m's avatar
frey_m committed
259
//         throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()", "Kinetic energy smaller than minimum cyclotron energy");
260
     
261
    if ( E_m > cycl_m->getFMHighE() )
frey_m's avatar
frey_m committed
262
        throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()", "Kinetic energy exceeds cyclotron energy");
263

adelmann's avatar
adelmann committed
264 265
    // velocity: beta = v/c = sqrt(1-1/(gamma*gamma))
    if (gamma_m == 0)
266
        throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()", "Relativistic factor equal zero.");
267

adelmann's avatar
adelmann committed
268 269
    // if domain_m = true --> integrate over a single sector
    if (domain_m) {
270
        N_m /=  cycl_m->getSymmetry();
adelmann's avatar
adelmann committed
271
    }
frey_m's avatar
frey_m committed
272 273 274
    
    cycl_m->read(cycl_m->getFieldFlag(cycl_m->getCyclotronType()),
                 cycl_m->getBScale());
275

276 277 278 279 280
    // reserve storage for the orbit and momentum (--> size = 0, capacity = N_m+1)
    /*
     * we need N+1 storage, since dtheta = 2pi/N (and not 2pi/(N-1)) that's why we need N+1 integration steps
     * to return to the origin (but the return size is N_m)
     */
adelmann's avatar
adelmann committed
281 282
    r_m.reserve(N_m + 1);
    pr_m.reserve(N_m + 1);
283 284
    vz_m.reserve(N_m + 1);
    vpz_m.reserve(N_m + 1);
285

286
    // reserve memory of N_m for the properties (--> size = 0, capacity = N_m)
adelmann's avatar
adelmann committed
287 288 289
    h_m.reserve(N_m);
    ds_m.reserve(N_m);
    fidx_m.reserve(N_m);
290 291 292
}

template<typename Value_type, typename Size_type, class Stepper>
293
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
294
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getInverseBendingRadius(const value_type& angle)
295
{
296 297 298 299
    if (angle != 0.0)
        return rotate(angle, h_m);
    else
        return h_m;
300 301 302
}

template<typename Value_type, typename Size_type, class Stepper>
303
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
304
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getPathLength(const value_type& angle)
305
{
306 307 308 309
    if (angle != 0.0)
        return rotate(angle, ds_m);
    else
        return ds_m;
310 311 312
}

template<typename Value_type, typename Size_type, class Stepper>
313
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
314
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getFieldIndex(const value_type& angle)
315
{
316 317
    if (angle != 0.0)
        return rotate(angle, fidx_m);
frey_m's avatar
frey_m committed
318
    return fidx_m;
319 320 321
}

template<typename Value_type, typename Size_type, class Stepper>
322 323 324
std::pair<Value_type,Value_type> ClosedOrbitFinder<Value_type, Size_type, Stepper>::getTunes() {
    // compute radial tune
    value_type nur = computeTune(x_m,px_m[1],nxcross_m);
325

326 327 328 329
    // compute vertical tune
    value_type nuz = computeTune(z_m,pz_m[1],nzcross_m);

    return std::make_pair(nur,nuz);
330 331 332
}

template<typename Value_type, typename Size_type, class Stepper>
333
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
334
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getOrbit(value_type angle)
335
{
336 337 338 339
    if (angle != 0.0)
        return rotate(angle, r_m);
    else
        return r_m;
340 341 342 343
}

template<typename Value_type, typename Size_type, class Stepper>
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
344
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getMomentum(value_type angle)
345 346
{
    container_type pr = pr_m;
347 348 349
    
    if (angle != 0.0)
        pr = rotate(angle, pr);
350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365
    
    // change units from meters to \beta * \gamma
    /* in Gordon paper:
     * 
     * p = \gamma * \beta * a
     * 
     * where a = c / \omega_{0} with \omega_{0} = 2 * \pi * \nu_{0} = 2 * \pi * \nu_{rf} / h
     * 
     * c: speed of light
     * h: harmonic number
     * v_{rf}: nomial rf frequency
     * 
     * Units:
     * 
     * [a] = m --> [p] = m
     * 
366
     * The momentum in \beta * \gamma is obtained by dividing by "a"
367
     */
368
    value_type factor =  1.0 / acon_m(wo_m);
369
    std::for_each(pr.begin(), pr.end(), [factor](value_type& p) { p *= factor; });
370
    
371
    return pr;
372 373 374
}

template<typename Value_type, typename Size_type, class Stepper>
375 376 377
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::value_type
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getGamma()
{
378
    return gamma_m;
379 380 381
}

template<typename Value_type, typename Size_type, class Stepper>
382 383 384
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::value_type
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getAverageRadius()
{
385
    return ravg_m;
386 387 388
}

template<typename Value_type, typename Size_type, class Stepper>
389 390
typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::value_type
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getFrequencyError()
391
{
392
    return phase_m;
393 394 395 396 397 398 399
}

// -----------------------------------------------------------------------------------------------------------------------
// PRIVATE MEMBER FUNCTIONS
// -----------------------------------------------------------------------------------------------------------------------

template<typename Value_type, typename Size_type, class Stepper>
400 401
bool ClosedOrbitFinder<Value_type, Size_type, Stepper>::findOrbit(value_type accuracy,
                                                                  size_type maxit,
frey_m's avatar
frey_m committed
402
                                                                  value_type dE,
403 404
                                                                  value_type rguess)
{
405 406 407 408
    /* REMARK TO GORDON
     * q' = 1/b = 1/bcon
     * a' = a = acon
     */
409
    
410
    value_type bint, brint, btint;
411

412 413 414
    // resize vectors (--> size = N_m+1, capacity = N_m+1), note: we do N_m+1 integration steps
    r_m.resize(N_m+1);
    pr_m.resize(N_m+1);
415 416
    vz_m.resize(N_m+1);
    vpz_m.resize(N_m+1);
417

418
    // store acon and bcon locally
419 420
    value_type acon = acon_m(wo_m);               // [acon] = m
    value_type invbcon = 1.0 / bcon_m(E0_m, wo_m);        // [bcon] = MeV*s/(C*m^2) = 10^6 T = 10^7 kG (kilo Gauss)
421 422 423 424 425 426 427

    // helper constants
    value_type p2;                                      // p^2 = p*p
    value_type pr2;                                     // squared radial momentum (pr^2 = pr*pr)
    value_type ptheta, invptheta;                       // Gordon, formula (5c)
    value_type invdenom;                                // denominator for computing dr,dpr
    value_type xold = 0.0;                              // for counting nxcross
428
    value_type zold = 0.0;                              // for counting nzcross
429 430 431 432 433

    // index for reaching next element of the arrays r and pr (no nicer way found yet)
    size_type idx = 0;
    // observer for storing the current value after each ODE step (e.g. Runge-Kutta step) into the containers of r and pr
    auto store = [&](state_type& y, const value_type t)
434
    {
435 436
        r_m[idx] = y[0];
        pr_m[idx] = y[1];
437 438
        vz_m[idx] = y[6];
        vpz_m[idx] = y[7];
439

440
        // count number of crossings (excluding starting point --> idx>0)
441 442
        nxcross_m += (idx > 0) * (y[4] * xold < 0);
        xold = y[4];
443 444 445 446 447
        
        // number of times z2 changes sign
        nzcross_m += (idx > 0) * (y[10] * zold < 0);
        zold = y[10];
        
448 449 450 451
        ++idx;
    };

    // define the six ODEs (using lambda function)
452 453 454 455
    std::function<void(const state_type&, state_type&, const double)> orbit_integration = [&](const state_type &y,
                                                                                              state_type &dydt,
                                                                                              const double theta)
    {
456 457
        pr2 = y[1] * y[1];
        if (p2 < pr2)
458 459
            throw OpalException("ClosedOrbitFinder::findOrbit()",
                                "p_{r}^2 > p^{2} (defined in Gordon paper) --> Square root of negative number.");
460

461 462 463
        // Gordon, formula (5c)
        ptheta = std::sqrt(p2 - pr2);
        invptheta = 1.0 / ptheta;
464
        // interpolate values of magnetic field
465 466
        cycl_m->apply(y[0], y[6], theta, brint, btint, bint);
        
467 468
        bint *= invbcon;
        brint *= invbcon;
469
        btint *= invbcon;
470 471 472 473 474 475 476 477 478 479

        // Gordon, formula (5a)
        dydt[0] = y[0] * y[1] * invptheta;
        // Gordon, formula (5b)
        dydt[1] = ptheta - y[0] * bint;
        // Gordon, formulas (9a) and (9b)
        for (size_type i = 2; i < 5; i += 2) {
            dydt[i] = (y[1] * y[i] + y[0] * p2 * y[i+1] * invptheta * invptheta) * invptheta;
            dydt[i+1] = - y[1] * y[i+1] * invptheta - (bint + y[0] * brint) * y[i];
        }
480 481 482 483 484 485 486 487 488 489
        
        // Gordon, formulas (22a) and (22b)
        for (size_type i = 6; i < 12; i += 2) {
            dydt[i] = y[0] * y[i+1] * invptheta;
            dydt[i+1] = (y[0] * brint - y[1] * invptheta * btint) * y[i];
        }

        // integrate phase
        dydt[12] = y[0] * invptheta * gamma_m - 1;
        
490 491
    };

492 493 494 495
    // define initial state container for integration: y = {r, pr, x1, px1, x2, px2,
    //                                                      z, pz, z1, pz1, z2, pz2,
    //                                                      phase}
    state_type y(11);
496

497 498 499 500 501 502 503 504
    // difference of last and first value of r (1. element) and pr (2. element)
    container_type err(2);
    // correction term for initial values: r = r + dr, pr = pr + dpr; Gordon, formula (17)
    container_type delta = {0.0, 0.0};
    // amplitude of error; Gordon, formula (18) (a = a')
    value_type error = std::numeric_limits<value_type>::max();
    // if niterations > maxit --> stop iteration
    size_type niterations = 0;
505 506 507 508

    /*
     * Christian:
     * N = 1440 ---> N = 720 ---> dtheta = 2PI/720 --> nsteps = 721
509
     *
510
     * 0, 2, 4, ... ---> jeden zweiten berechnene: 1, 3, 5, ... interpolieren --> 1440 Werte
511
     *
512 513
     * Matthias:
     * N = 1440 --> dtheta = 2PI/1440 --> nsteps = 1441
514
     *
515
     * 0, 1, 2, 3, 4, 5, ... --> 1440 Werte
516
     *
517
     */
518

519
    // iterate until suggested energy (start with minimum energy)
520
    value_type E = cycl_m->getFMLowE();
521

522
    // energy dependent values
523
    value_type en = E / E0_m;                      // en = E/E0 = E/(mc^2) (E0 is potential energy)
524 525 526 527 528 529 530 531 532
    value_type p = acon * std::sqrt(en * (2.0 + en));     // momentum [p] = m; Gordon, formula (3)
    value_type gamma2 = (1.0 + en) * (1.0 + en);          // = gamma^2
    value_type beta = std::sqrt(1.0 - 1.0 / gamma2);
    p2 = p * p;                                           // p^2 = p*p
    value_type invgamma4 = 1.0 / (gamma2 * gamma2);       // = 1/gamma^4

    // set initial values for radius and radial momentum for lowest energy Emin
    // orbit, [r] = m; Gordon, formula (20)
    // radial momentum; Gordon, formula (20)
Andreas Adelmann's avatar
Andreas Adelmann committed
533 534

    container_type init;
frey_m's avatar
frey_m committed
535 536 537 538 539 540
    //      r            pr   z    pz
    init = {beta * acon, 0.0, 0.0, 1.0};

    if (rguess >= 0.0) {
        init[0] = rguess * 0.001;
    }
541

542 543
    do {
        
544
        // (re-)set inital values for r and pr
545
        r_m[0] = init[0];
546
        pr_m[0] = init[1];
547 548
        vz_m[0] = init[2];
        vpz_m[0] = init[3];
549

550 551 552 553 554 555 556
        // integrate until error smaller than user-define accuracy
        do {
            // (re-)set inital values
            x_m[0]  = 1.0;               // x1; Gordon, formula (10)
            px_m[0] = 0.0;               // px1; Gordon, formula (10)
            x_m[1]  = 0.0;               // x2; Gordon, formula (10)
            px_m[1] = 1.0;               // px2; Gordon, formula (10)
frey_m's avatar
frey_m committed
557 558 559 560 561
            z_m[0]  = 1.0;
            pz_m[0] = 0.0;
            z_m[1]  = 0.0;
            pz_m[1] = 1.0;
            phase_m = 0.0;
562
            nxcross_m = 0;               // counts the number of crossings of x-axis (excluding first step)
563
            nzcross_m = 0;
564 565 566
            idx = 0;                     // index for looping over r and pr arrays

            // fill container with initial states
567 568 569
            y = {init[0],init[1],
                 x_m[0], px_m[0], x_m[1], px_m[1],
                 init[2], init[3],
frey_m's avatar
frey_m committed
570 571
                 z_m[0], pz_m[0], z_m[1], pz_m[1],
                 phase_m
572
            };
573

574 575
            // integrate from 0 to 2*pi (one has to get back to the "origin")
            boost::numeric::odeint::integrate_n_steps(stepper_m,orbit_integration,y,0.0,dtheta_m,N_m,store);
576

577 578 579 580 581
            // write new state
            x_m[0] = y[2];
            px_m[0] = y[3];
            x_m[1] = y[4];
            px_m[1] = y[5];
582 583 584 585 586 587
            
            z_m[0] = y[8];
            pz_m[0] = y[9];
            z_m[1] = y[10];
            pz_m[1] = y[11];
            phase_m = y[12] * Physics::u_two_pi; // / (2.0 * Physics::pi);
588

589 590 591 592
            // compute error (compare values of orbit and momentum for theta = 0 and theta = 2*pi)
            // (Note: size = N_m+1 --> last entry is N_m)
            err[0] = r_m[N_m] - r_m[0];      // Gordon, formula (14)
            err[1] = pr_m[N_m] - pr_m[0];    // Gordon, formula (14)
593

594 595 596 597
            // correct inital values of r and pr
            invdenom = 1.0 / (x_m[0] + px_m[1] - 2.0);
            delta[0] = ((px_m[1] - 1.0) * err[0] - x_m[1] * err[1]) * invdenom; // dr; Gordon, formula (16a)
            delta[1] = ((x_m[0] - 1.0) * err[1] - px_m[0] * err[0]) * invdenom; // dpr; Gordon, formula (16b)
598

599 600 601
            // improved initial values; Gordon, formula (17) (here it's used for higher energies)
            init[0] += delta[0];
            init[1] += delta[1];
602

603 604 605
            // compute amplitude of the error
            error = std::sqrt(delta[0] * delta[0] + delta[1] * delta[1] * invgamma4) / r_m[0];
        } while (error > accuracy && niterations++ < maxit);
606

607 608
        // reset iteration counter
        niterations = 0;
609

610 611
        // reset correction term
        delta[0] = delta[1] = 0.0;
adelmann's avatar
adelmann committed
612 613 614 615 616 617

        // increase energy by dE
        if (E_m <= E + dE)
            E = E_m;
        else
            E += dE;
618

619
        // set constants for new energy E
620
        en = E / E0_m;                     // en = E/E0 = E/(mc^2) (E0 is potential energy)
621 622 623 624
        p = acon * std::sqrt(en * (2.0 + en));    // momentum [p] = m; Gordon, formula (3)
        p2 = p * p;                               // p^2 = p*p
        gamma2 = (1.0 + en) * (1.0 + en);
        invgamma4 = 1.0 / (gamma2 * gamma2);
625 626


627
    } while (E != E_m);
628

629 630 631 632 633
    /* store last entry, since it is needed in computeVerticalOscillations(), because we have to do the same
     * number of integrations steps there.
     */
    lastOrbitVal_m = r_m[N_m];           // needed in computeVerticalOscillations()
    lastMomentumVal_m = pr_m[N_m];       // needed in computeVerticalOscillations()
634

635 636 637
    // remove last entry (since we don't have to store [0,2pi], but [0,2pi[)  --> size = N_m, capacity = N_m+1
    r_m.pop_back();
    pr_m.pop_back();
638 639 640 641 642 643 644
    
    
    /* domain_m = true --> only integrated over a single sector
     * --> multiply by cycl_m->getSymmetry() to get correct phase_m
     */
    if (domain_m)
        phase_m *= cycl_m->getSymmetry();
645

646

647 648 649
    // returns true if converged, otherwise false
    return error < accuracy;
}
650 651

template<typename Value_type, typename Size_type, class Stepper>
652 653 654
Value_type ClosedOrbitFinder<Value_type, Size_type, Stepper>::computeTune(const std::array<value_type,2>& y,
                                                                          value_type py2, size_type ncross)
{
655
    // Y = [y1, y2; py1, py2]
656

657 658
    // cos(mu)
    value_type cos = 0.5 * (y[0] + py2);
659
    
660
    value_type mu;
661

662 663
    // sign of sin(mu) has to be equal to y2
    bool negative = std::signbit(y[1]);
664

665
    bool uneven = (ncross % 2);
666

667 668 669
    if (std::fabs(cos) > 1.0) {
        // store the number of crossings
        value_type n = ncross;
670

671 672
        if (uneven)
            n = ncross - 1;
673

674 675
        // Gordon, formula (36b)
        value_type muPrime = -std::acosh(std::fabs(cos));      // mu'
676
        mu = n * Physics::pi + muPrime;
677

678 679 680 681 682 683 684
    } else {
        value_type muPrime = (uneven) ? std::acos(-cos) : std::acos(cos);    // mu'
        /* It has to be fulfilled: 0<= mu' <= pi
        * But since |cos(mu)| <= 1, we have
        * -1 <= cos(mu) <= 1 --> 0 <= mu <= pi (using above programmed line), such
        * that condition is already fulfilled.
        */
685

686
        // Gordon, formula (36)
687
        mu = ncross * Physics::pi + muPrime;
688

689 690
        // if sign(y[1]) > 0 && sign(sin(mu)) < 0
        if (!negative && std::signbit(std::sin(mu))) {
691
            mu = ncross * Physics::pi - muPrime;
692
        } else if (negative && !std::signbit(std::sin(mu))) {
693
            mu = ncross * Physics::pi - muPrime + Physics::two_pi;
694 695
        }
    }
696

697
    // nu = mu/theta, where theta = integration domain
698

699
    /* domain_m = true --> only integrated over a single sector --> multiply by cycl_m->getSymmetry() to
adelmann's avatar
adelmann committed
700 701 702
     * get correct tune.
     */
    if (domain_m)
703
        mu *= cycl_m->getSymmetry();
704

705
    return mu * Physics::u_two_pi;
706 707 708
}

template<typename Value_type, typename Size_type, class Stepper>
709
void ClosedOrbitFinder<Value_type, Size_type, Stepper>::computeOrbitProperties() {
710
    /*
711 712 713 714 715
     * The formulas for h, fidx and ds are from the paper:
     * "Tranverse-Longitudinal Coupling by Space Charge in Cyclotrons"
     * written by Dr. Christian Baumgarten (2012, PSI)
     * p. 6
     */
716

adelmann's avatar
adelmann committed
717
    // READ IN MAGNETIC FIELD: ONLY FOR STAND-ALONE PROGRAM
718
    value_type bint, brint, btint; // B, dB/dr, dB/dtheta
719

720 721 722
    value_type invbcon = 1.0 / bcon_m(E0_m, wo_m);
    value_type en = E_m / E0_m;                                  // en = E/E0 = E/(mc^2) (E0 is potential energy)
    value_type p = acon_m(wo_m) * std::sqrt(en * (2.0 + en));    // momentum [p] = m; Gordon, formula (3)
723 724 725
    value_type p2 = p * p;
    value_type theta = 0.0;                                             // angle for interpolating
    value_type ptheta;
726

727 728 729 730 731 732 733
    // resize of container (--> size = N_m, capacity = N_m)
    h_m.resize(N_m);
    fidx_m.resize(N_m);
    ds_m.resize(N_m);

    for (size_type i = 0; i < N_m; ++i) {
        // interpolate magnetic field
734
        cycl_m->apply(r_m[i], vz_m[i], theta, brint, btint, bint);
735 736 737
        bint *= invbcon;
        brint *= invbcon;
        btint *= invbcon;
frey_m's avatar
frey_m committed
738
        
739 740 741 742 743 744 745 746 747 748 749 750
        // inverse bending radius
        h_m[i] = bint / p;

        // local field index
        ptheta = std::sqrt(p2 - pr_m[i] * pr_m[i]);
        fidx_m[i] = (brint * ptheta - btint * pr_m[i] / r_m[i]) / p2; //(bint*bint);

        // path length element
        ds_m[i] = std::hypot(r_m[i] * pr_m[i] / ptheta,r_m[i]) * dtheta_m; // C++11 function

        // increase angle
        theta += dtheta_m;
751
    }
752 753 754

    // compute average radius
    ravg_m = std::accumulate(r_m.begin(),r_m.end(),0.0) / value_type(r_m.size());
755 756
}

757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778
template<typename Value_type, typename Size_type, class Stepper> 
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
ClosedOrbitFinder<Value_type, Size_type, Stepper>::rotate(value_type angle, container_type &orbitProperty) {

    container_type orbitPropertyCopy = orbitProperty;
    
    // compute the number of steps per degree
    value_type deg_step = N_m / 360.0;

    // compute starting point
    size_type start = deg_step * angle;

    // copy end to start
    std::copy(orbitProperty.begin() + start, orbitProperty.end(), orbitPropertyCopy.begin());
    
    // copy start to end
    std::copy_n(orbitProperty.begin(), start, orbitPropertyCopy.end() - start);

    return orbitPropertyCopy;

}

Andreas Adelmann's avatar
Andreas Adelmann committed
779
#endif