ClosedOrbitFinder.h 32.1 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 "physics.h"
31

32
#include "AbsBeamline/Cyclotron.h"
frey_m's avatar
frey_m committed
33

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

37
/// Finds a closed orbit of a cyclotron for a given energy
38 39 40
template<typename Value_type, typename Size_type, class Stepper>
class ClosedOrbitFinder
{
41 42 43 44 45 46 47 48 49 50 51 52 53
    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
54
         * @param E0 is the potential energy (particle energy at rest) [MeV].
adelmann's avatar
adelmann committed
55
         * @param N specifies the number of splits (2pi/N), i.e number of integration steps
56
         * @param cycl is the cyclotron element
57 58
         * @param domain is a boolean (default: true). If "true" the closed orbit is computed over a single sector,
         * otherwise over 2*pi.
59
         */
60
        ClosedOrbitFinder(value_type E, value_type E0, size_type N,
61
                          Cyclotron* cycl, bool domain = true);
62 63

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

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

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

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

75 76 77 78 79 80 81
        /// 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°).
         */
82 83
        container_type getOrbit(value_type angle = 0);

84 85 86 87 88 89
        /// 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°).
90
         * @returns the momentum in \f$ \beta * \gamma \f$ units
91
         */
92
        container_type getMomentum(value_type angle = 0);
93 94 95 96 97 98 99

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

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

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

        /// 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
110
         * @param rguess initial radius guess in [mm]
111
         */
112
        bool findOrbit(value_type, size_type, value_type = -1.0);
113

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

117
    private:
118 119 120 121 122 123 124 125
        /// 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
126
        /// This function computes nzcross_ which is used to compute the tune in z-direction and the frequency error
127
        void computeVerticalOscillations();
128 129 130
        
        /// This function rotates the calculated closed orbit finder properties to the initial angle
        container_type rotate(value_type angle, container_type& orbitProperty);
131 132 133 134 135

        /// 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
136
        /// Stores current position in vertical direction for the solutions of the ODE with different initial values
137
        std::array<value_type,2> z_m; // z_m = [z1, z2]
frey_m's avatar
frey_m committed
138
        /// Stores current momenta in vertical direction for the solutions of the ODE with different initial values
139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158
        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;
        /// Stores the radial momentum
        container_type pr_m;
        /// 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;
159 160 161 162
        
        /// Is the potential energy [MeV]
        value_type E0_m;
        
163
        /// Is the nominal orbital frequency
164 165 166
        /* (see paper of Dr. C. Baumgarten: "Transverse-Longitudinal
         * Coupling by Space Charge in Cyclotrons" (2012), formula (1))
         */
167
        value_type wo_m;
adelmann's avatar
adelmann committed
168
        /// Number of integration steps
169 170 171 172 173 174 175 176 177 178 179 180 181
        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;

182
        /**
183 184 185 186
         * 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)
         */
187 188
        value_type lastOrbitVal_m;

189 190 191 192 193
        /**
         * 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)
         */
194
        value_type lastMomentumVal_m;
195 196

        /**
197 198 199
         * Boolean which is true if computeVerticalOscillations() executed, otherwise false. This is used for checking in
         * getTunes() and getFrequencyError().
         */
200 201
        bool vertOscDone_m;

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

208 209
        /// Defines the stepper for integration of the ODE's
        Stepper stepper_m;
210 211 212 213 214 215 216 217 218 219 220 221 222 223 224
        
        /*!
         * 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);
        };
225
        
226
        Cyclotron* cycl_m;
227 228 229 230 231 232
};

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

233 234 235 236
template<typename Value_type, typename Size_type, class Stepper>
ClosedOrbitFinder<Value_type,
                  Size_type,
                  Stepper>::ClosedOrbitFinder(value_type E, value_type E0,
237
                                              size_type N, Cyclotron* cycl,
238
                                              bool domain)
239 240 241 242
    : nxcross_m(0)
    , nzcross_m(0)
    , E_m(E)
    , E0_m(E0)
243
    , wo_m(cycl->getRfFrequ()*1E6/cycl->getCyclHarm()*2.0*Physics::pi)
244 245 246 247 248 249 250 251 252 253
    , 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)
    , vertOscDone_m(false)
    , domain_m(domain)
    , stepper_m()
254
    , cycl_m(cycl)
255
{
frey_m's avatar
frey_m committed
256
    
257
    if ( cycl_m->getFMLowE() > cycl_m->getFMHighE() )
frey_m's avatar
frey_m committed
258 259 260 261
        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.
262
//     if ( E_m < cycl_m->getFMLowE() )
frey_m's avatar
frey_m committed
263
//         throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()", "Kinetic energy smaller than minimum cyclotron energy");
264
     
265
    if ( E_m > cycl_m->getFMHighE() )
frey_m's avatar
frey_m committed
266
        throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()", "Kinetic energy exceeds cyclotron energy");
267

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

adelmann's avatar
adelmann committed
272 273
    // if domain_m = true --> integrate over a single sector
    if (domain_m) {
274
        N_m /=  cycl_m->getSymmetry();
adelmann's avatar
adelmann committed
275
    }
276

277 278 279 280 281
    // 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
282 283
    r_m.reserve(N_m + 1);
    pr_m.reserve(N_m + 1);
284

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

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

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

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

template<typename Value_type, typename Size_type, class Stepper>
321 322 323
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);
324

325 326 327
    // compute nzcross_m
    if (!vertOscDone_m)
        computeVerticalOscillations();
328

329 330 331 332
    // compute vertical tune
    value_type nuz = computeTune(z_m,pz_m[1],nzcross_m);

    return std::make_pair(nur,nuz);
333 334 335
}

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

template<typename Value_type, typename Size_type, class Stepper>
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
347
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getMomentum(value_type angle)
348 349
{
    container_type pr = pr_m;
350 351 352
    
    if (angle != 0.0)
        pr = rotate(angle, pr);
353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368
    
    // 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
     * 
369
     * The momentum in \beta * \gamma is obtained by dividing by "a"
370
     */
371
    value_type factor =  1.0 / acon_m(wo_m);
372
    std::for_each(pr.begin(), pr.end(), [factor](value_type& p) { p *= factor; });
373
    
374
    return pr;
375 376 377
}

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

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

template<typename Value_type, typename Size_type, class Stepper>
392 393
typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::value_type
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getFrequencyError()
394
{
395 396 397
    // if the vertical oscillations aren't computed, we have to, since there we also compuote the frequency error.
    if(!vertOscDone_m)
        computeVerticalOscillations();
398

399
    return phase_m;
400 401 402 403 404 405 406
}

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

template<typename Value_type, typename Size_type, class Stepper>
407 408 409 410
bool ClosedOrbitFinder<Value_type, Size_type, Stepper>::findOrbit(value_type accuracy,
                                                                  size_type maxit,
                                                                  value_type rguess)
{
411 412 413 414 415
    /* REMARK TO GORDON
     * q' = 1/b = 1/bcon
     * a' = a = acon
     */

frey_m's avatar
frey_m committed
416 417
    cycl_m->read(cycl_m->getFieldFlag(cycl_m->getCyclotronType()),
                 cycl_m->getBScale());
418
    
419
    value_type bint, brint, btint;
420

421 422 423
    // 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);
424

425
    // store acon and bcon locally
426 427
    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)
428 429 430 431 432 433 434 435 436 437 438 439

    // 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

    // 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)
440
    {
441 442 443
        r_m[idx] = y[0];
        pr_m[idx] = y[1];

444
        // count number of crossings (excluding starting point --> idx>0)
445 446 447 448 449 450
        nxcross_m += (idx > 0) * (y[4] * xold < 0);
        xold = y[4];
        ++idx;
    };

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

459 460 461 462
        // Gordon, formula (5c)
        ptheta = std::sqrt(p2 - pr2);
        invptheta = 1.0 / ptheta;

463
        // interpolate values of magnetic field
464
        cycl_m->apply(y[0], 0.0, theta, brint, btint, bint);
465

466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481
        bint *= invbcon;
        brint *= invbcon;

        // 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];
        }
    };

    // define initial state container for integration: y = {r, pr, x1, px1, x2, px2}
    state_type y(6);
482

483 484 485 486 487 488 489 490
    // 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;
491 492 493 494

    /*
     * Christian:
     * N = 1440 ---> N = 720 ---> dtheta = 2PI/720 --> nsteps = 721
495
     *
496
     * 0, 2, 4, ... ---> jeden zweiten berechnene: 1, 3, 5, ... interpolieren --> 1440 Werte
497
     *
498 499
     * Matthias:
     * N = 1440 --> dtheta = 2PI/1440 --> nsteps = 1441
500
     *
501
     * 0, 1, 2, 3, 4, 5, ... --> 1440 Werte
502
     *
503
     */
504

Andreas Adelmann's avatar
Andreas Adelmann committed
505 506 507
    // step size of energy
    value_type dE; 

508
    if (cycl_m->getFMLowE() == cycl_m->getFMHighE())
Andreas Adelmann's avatar
Andreas Adelmann committed
509 510
      dE = 0.0;
    else
511
      dE = (E_m - cycl_m->getFMLowE()) / (cycl_m->getFMHighE() - cycl_m->getFMLowE());
Andreas Adelmann's avatar
Andreas Adelmann committed
512

513
    // iterate until suggested energy (start with minimum energy)
514
    value_type E = cycl_m->getFMLowE();
515

adelmann's avatar
adelmann committed
516 517
    // energy increase not more than 0.25
    dE = (dE > 0.25) ? 0.25 : dE;
518 519

    // energy dependent values
520
    value_type en = E / E0_m;                      // en = E/E0 = E/(mc^2) (E0 is potential energy)
521 522 523 524 525 526 527 528 529
    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
530 531

    container_type init;
532
    if (rguess < 0)
Andreas Adelmann's avatar
Andreas Adelmann committed
533 534
      init = {beta * acon, 0.0};
    else
535
      init = {rguess * 0.001, 0.0};
536

537 538
    do {
        
539
        // (re-)set inital values for r and pr
540
        r_m[0] = init[0];
541
        pr_m[0] = init[1];
542

543 544 545 546 547 548 549 550 551 552 553 554
        // 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)
            nxcross_m = 0;               // counts the number of crossings of x-axis (excluding first step)
            idx = 0;                     // index for looping over r and pr arrays

            // fill container with initial states
            y = {init[0],init[1], x_m[0], px_m[0], x_m[1], px_m[1]};
555

556 557
            // 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);
558

559 560 561 562 563
            // write new state
            x_m[0] = y[2];
            px_m[0] = y[3];
            x_m[1] = y[4];
            px_m[1] = y[5];
564

565 566 567 568
            // 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)
569

570 571 572 573
            // 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)
574

575 576 577
            // improved initial values; Gordon, formula (17) (here it's used for higher energies)
            init[0] += delta[0];
            init[1] += delta[1];
578

579 580 581
            // 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);
582

583 584
        // reset iteration counter
        niterations = 0;
585

586 587
        // reset correction term
        delta[0] = delta[1] = 0.0;
adelmann's avatar
adelmann committed
588 589 590 591 592 593

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

595
        // set constants for new energy E
596
        en = E / E0_m;                     // en = E/E0 = E/(mc^2) (E0 is potential energy)
597 598 599 600
        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);
601 602


603
    } while (E != E_m);
604

605 606 607 608 609
    /* 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()
610

611 612 613
    // 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();
614

615

616 617 618
    // returns true if converged, otherwise false
    return error < accuracy;
}
619 620

template<typename Value_type, typename Size_type, class Stepper>
621 622 623
Value_type ClosedOrbitFinder<Value_type, Size_type, Stepper>::computeTune(const std::array<value_type,2>& y,
                                                                          value_type py2, size_type ncross)
{
624
    // Y = [y1, y2; py1, py2]
625

626 627
    // cos(mu)
    value_type cos = 0.5 * (y[0] + py2);
628
    
629
    value_type mu;
630

631 632
    // sign of sin(mu) has to be equal to y2
    bool negative = std::signbit(y[1]);
633

634
    bool uneven = (ncross % 2);
635

636 637 638
    if (std::fabs(cos) > 1.0) {
        // store the number of crossings
        value_type n = ncross;
639

640 641
        if (uneven)
            n = ncross - 1;
642

643 644
        // Gordon, formula (36b)
        value_type muPrime = -std::acosh(std::fabs(cos));      // mu'
645
        mu = n * Physics::pi + muPrime;
646

647 648 649 650 651 652 653
    } 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.
        */
654

655
        // Gordon, formula (36)
656
        mu = ncross * Physics::pi + muPrime;
657

658 659
        // if sign(y[1]) > 0 && sign(sin(mu)) < 0
        if (!negative && std::signbit(std::sin(mu))) {
660
            mu = ncross * Physics::pi - muPrime;
661
        } else if (negative && !std::signbit(std::sin(mu))) {
662
            mu = ncross * Physics::pi - muPrime + Physics::two_pi;
663 664
        }
    }
665

666
    // nu = mu/theta, where theta = integration domain
667

668
    /* domain_m = true --> only integrated over a single sector --> multiply by cycl_m->getSymmetry() to
adelmann's avatar
adelmann committed
669 670 671
     * get correct tune.
     */
    if (domain_m)
672
        mu *= cycl_m->getSymmetry();
673

674
    return mu * Physics::u_two_pi;
675 676 677
}

template<typename Value_type, typename Size_type, class Stepper>
678
void ClosedOrbitFinder<Value_type, Size_type, Stepper>::computeOrbitProperties() {
679
    /*
680 681 682 683 684
     * 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
     */
685

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

689 690 691
    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)
692 693 694
    value_type p2 = p * p;
    value_type theta = 0.0;                                             // angle for interpolating
    value_type ptheta;
695

696 697 698 699 700 701 702
    // 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
703
        cycl_m->apply(r_m[i], 0.0, theta, brint, btint, bint);
704 705 706
        bint *= invbcon;
        brint *= invbcon;
        btint *= invbcon;
frey_m's avatar
frey_m committed
707
        
708 709 710 711 712 713 714 715 716 717 718 719
        // 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;
720
    }
721 722 723

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

template<typename Value_type, typename Size_type, class Stepper>
727
void ClosedOrbitFinder<Value_type, Size_type, Stepper>::computeVerticalOscillations() {
728

729 730 731 732 733
    vertOscDone_m = true;

    // READ IN MAGNETIC FIELD: ONLY FOR STAND-ALONE PROGRAM
    value_type bint, brint, btint; // B, dB/dr, dB/dtheta

734 735
    value_type en = E_m / E0_m;                                  // en = E/E0 = E/(mc^2) with potential energy E0
    value_type p = acon_m(wo_m) * std::sqrt(en *(en + 2.0));     // Gordon, formula (3)
736 737 738 739 740 741 742
    value_type p2 = p * p;                                              // p^2 = p*p
    size_type idx = 0;                                                  // index for going through container
    value_type pr2;                                                     // pr^2 = pr*pr
    value_type ptheta, invptheta;                                       // Gordon, formula (5c)
    value_type zold = 0.0;                                              // for counting nzcross

    // store bcon locally
743
    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)
744 745

    // define the ODEs (using lambda function)
746 747 748 749
    std::function<void(const state_type&, state_type&, const double)> vertical = [&](const state_type &y,
                                                                                     state_type &dydt,
                                                                                     const double theta)
    {
750
        pr2 = y[1] * y[1];
751
        if (p2 < pr2) {
752 753
            throw OpalException("ClosedOrbitFinder::computeVerticalOscillations()",
                                "p_{r}^2 > p^{2} (defined in Gordon paper) --> Square root of negative number.");
754
        }
755

756 757 758 759 760
        // Gordon, formula (5c)
        ptheta = std::sqrt(p2 - pr2);
        invptheta = 1.0 / ptheta;

        // intepolate values of magnetic field
frey_m's avatar
frey_m committed
761
        cycl_m->apply(y[0], y[2], theta, brint, btint, bint);
frey_m's avatar
frey_m committed
762
        
763 764 765
        bint *= invbcon;
        brint *= invbcon;
        btint *= invbcon;
766

767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796
        // We have to integrate r and pr again, otherwise we don't have the Runge-Kutta of the B-field
        // Gordon, formula (5a)
        dydt[0] = y[0] * y[1] * invptheta;
        // Gordon, formula (5b)
        dydt[1] = ptheta - y[0] * bint;

        // Gordon, formulas (22a) and (22b)
        for (size_type i = 2; i < 5; 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[6] = y[0] * invptheta * gamma_m - 1;
    };

    // to get next index for r and pr (to iterate over container)
    auto next = [&](state_type& y, const value_type t) {
        // number of times z2 changes sign
        nzcross_m += (idx > 0) * (y[4] * zold < 0);
        zold = y[4];
        ++idx;
    };

    // set initial state container for integration: y = {r, pr, z1, pz1, z2, pz2, phase}
    state_type y = {r_m[0], pr_m[0], 1.0, 0.0, 0.0, 1.0, 0.0};

    // add last element for integration (since we have to return to the initial point (--> size = N_m+1, capacity = N_m+1)
    r_m.push_back(lastOrbitVal_m);
    pr_m.push_back(lastMomentumVal_m);
797

798 799
    // integrate: assume no imperfections --> only integrate over a single sector (dtheta_m = 2pi/N_m)
    boost::numeric::odeint::integrate_n_steps(stepper_m,vertical,y,0.0,dtheta_m,N_m,next);
800

801 802 803
    // remove last element again (--> size = N_m, capacity = N_m+1)
    r_m.pop_back();
    pr_m.pop_back();
804

805 806 807 808 809
    // write new state
    z_m[0] = y[2];
    pz_m[0] = y[3];
    z_m[1] = y[4];
    pz_m[1] = y[5];
810
    phase_m = y[6] * Physics::u_two_pi; // / (2.0 * Physics::pi);
811

adelmann's avatar
adelmann committed
812
    /* domain_m = true --> only integrated over a single sector
813
     * --> multiply by cycl_m->getSymmetry() to get correct phase_m
adelmann's avatar
adelmann committed
814 815
     */
    if (domain_m)
816
        phase_m *= cycl_m->getSymmetry();
817 818
}

819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840
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
841
#endif