ClosedOrbitFinder.h 33.9 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 "MagneticField.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].
55 56
         * @param wo is the nominal orbital frequency (see paper of Dr. C. Baumgarten: "Transverse-Longitudinal
         * Coupling by Space Charge in Cyclotrons" (2012), formula (1))
adelmann's avatar
adelmann committed
57
         * @param N specifies the number of splits (2pi/N), i.e number of integration steps
58 59 60 61
         * @param accuracy specifies the accuracy of the closed orbit
         * @param maxit is the maximal number of iterations done. Program stops if closed orbit not found within this time.
         * @param Emin is the minimum energy [MeV] needed in cyclotron
         * @param Emax is the maximum energy [MeV] reached in cyclotron
adelmann's avatar
adelmann committed
62
         * @param nSector is the number of sectors (--> symmetry) of cyclotron
63
         * @param fmapfn is the location of the file that specifies the magnetic field
frey_m's avatar
frey_m committed
64 65
	 * @param guess value of radius for closed orbit finder
         * @param type specifies the field format (e.g. CARBONCYCL)
66
         * @param scaleFactor for the magnetic field (default: 1.0)
67 68
         * @param domain is a boolean (default: true). If "true" the closed orbit is computed over a single sector,
         * otherwise over 2*pi.
69
         */
70 71 72
        ClosedOrbitFinder(value_type E, value_type E0, value_type wo, size_type N,
                          value_type accuracy, size_type maxit, value_type Emin, value_type Emax,
                          size_type nSector, const std::string& fmapfn, value_type guess,
frey_m's avatar
frey_m committed
73
                          const std::string& type, value_type scaleFactor = 1.0,
74
                          bool domain = true);
75 76

        /// Returns the inverse bending radius (size of container N+1)
77
        container_type getInverseBendingRadius(const value_type& angle = 0);
78 79

        /// Returns the step lengths of the path (size of container N+1)
80
        container_type getPathLength(const value_type& angle = 0);
81 82

        /// Returns the field index (size of container N+1)
83
        container_type getFieldIndex(const value_type& angle = 0);
84 85 86 87

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

88 89 90 91 92 93 94
        /// 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°).
         */
95 96
        container_type getOrbit(value_type angle = 0);

97 98 99 100 101 102
        /// 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°).
103
         * @returns the momentum in \f$ \beta * \gamma \f$ units
104
         */
105
        container_type getMomentum(value_type angle = 0);
106 107 108 109 110 111 112

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

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

adelmann's avatar
adelmann committed
113 114
        /// Returns the frequency error
        value_type getFrequencyError();
115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137

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

    private:
        /// 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
         */
        bool findOrbit(value_type, size_type);

        /// Fills in the values of h_m, ds_m, fidx_m. It gets called by in by constructor.
        void computeOrbitProperties();

        /// 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
138
        /// This function computes nzcross_ which is used to compute the tune in z-direction and the frequency error
139
        void computeVerticalOscillations();
140 141 142
        
        /// This function rotates the calculated closed orbit finder properties to the initial angle
        container_type rotate(value_type angle, container_type& orbitProperty);
143 144 145 146 147

        /// 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
148
        /// Stores current position in vertical direction for the solutions of the ODE with different initial values
149
        std::array<value_type,2> z_m; // z_m = [z1, z2]
frey_m's avatar
frey_m committed
150
        /// Stores current momenta in vertical direction for the solutions of the ODE with different initial values
151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170
        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;
171 172 173 174
        
        /// Is the potential energy [MeV]
        value_type E0_m;
        
175 176
        /// Is the nominal orbital frequency
        value_type wo_m;
adelmann's avatar
adelmann committed
177
        /// Number of integration steps
178 179 180 181 182 183 184 185 186 187 188 189 190
        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;

191 192 193
        /**
         * Boolean which tells if a closed orbit for this configuration could be found (get set by the function findOrbit)
         */
194 195 196 197 198 199 200
        bool converged_m;

        /// Minimum energy needed in cyclotron
        value_type Emin_m;

        /// Maximum energy reached in cyclotron
        value_type Emax_m;
201

adelmann's avatar
adelmann committed
202 203
        /// Number of sectors (symmetry)
        size_type nSector_m;
204 205

        /**
206 207 208 209
         * 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)
         */
210 211
        value_type lastOrbitVal_m;

212 213 214 215 216
        /**
         * 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)
         */
217
        value_type lastMomentumVal_m;
218 219

        /**
220 221 222
         * Boolean which is true if computeVerticalOscillations() executed, otherwise false. This is used for checking in
         * getTunes() and getFrequencyError().
         */
223 224
        bool vertOscDone_m;

225
        /**
226 227 228
         * Boolean which is true by default. "true": orbit integration over one sector only, "false": integration
         * over 2*pi
         */
adelmann's avatar
adelmann committed
229
        bool domain_m;
230

231 232 233
        /// Defines the stepper for integration of the ODE's
        Stepper stepper_m;

Andreas Adelmann's avatar
Andreas Adelmann committed
234 235
	/// a guesss for the clo finder
	value_type rguess_m;
236 237 238 239 240 241 242 243 244 245 246 247 248 249 250
        
        /*!
         * 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);
        };
251
        
252
        MagneticField bField_m;
253 254 255 256 257 258
};

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

259 260 261 262 263 264 265 266
template<typename Value_type, typename Size_type, class Stepper>
ClosedOrbitFinder<Value_type,
                  Size_type,
                  Stepper>::ClosedOrbitFinder(value_type E, value_type E0,
                                              value_type wo, size_type N,
                                              value_type accuracy, size_type maxit,
                                              value_type Emin, value_type Emax,
                                              size_type nSector, const std::string& fmapfn,
frey_m's avatar
frey_m committed
267 268
                                              value_type rguess, const std::string& type,
                                              value_type scaleFactor, bool domain)
269 270
: nxcross_m(0), nzcross_m(0), E_m(E), E0_m(E0), wo_m(wo), N_m(N), dtheta_m(Physics::two_pi/value_type(N)),
  gamma_m(E/E0+1.0), ravg_m(0), phase_m(0), converged_m(false), Emin_m(Emin), Emax_m(Emax), nSector_m(nSector),
271
  lastOrbitVal_m(0.0), lastMomentumVal_m(0.0),
frey_m's avatar
frey_m committed
272
  vertOscDone_m(false), domain_m(domain), stepper_m(), rguess_m(rguess), bField_m(fmapfn, nSector)
273
{
frey_m's avatar
frey_m committed
274 275 276 277 278 279 280 281
    
    if ( Emin_m > Emax_m )
        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.
//     if ( E_m < Emin_m )
//         throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()", "Kinetic energy smaller than minimum cyclotron energy");
282
     
frey_m's avatar
frey_m committed
283 284
    if ( E_m > Emax_m )
        throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()", "Kinetic energy exceeds cyclotron energy");
285

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

adelmann's avatar
adelmann committed
290 291 292 293
    // if domain_m = true --> integrate over a single sector
    if (domain_m) {
        N_m /=  nSector_m;
    }
294

295 296 297 298 299
    // 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
300 301
    r_m.reserve(N_m + 1);
    pr_m.reserve(N_m + 1);
302

303
    // reserve memory of N_m for the properties (--> size = 0, capacity = N_m)
adelmann's avatar
adelmann committed
304 305 306
    h_m.reserve(N_m);
    ds_m.reserve(N_m);
    fidx_m.reserve(N_m);
307 308
    
    // read in magnetic fieldmap
frey_m's avatar
frey_m committed
309
    bField_m.read(type, scaleFactor);
310

311
    // compute closed orbit
312
    converged_m = findOrbit(accuracy, maxit);
313

314 315 316 317 318
    // compute h, ds, fidx, rav (average radius)
    computeOrbitProperties();
}

template<typename Value_type, typename Size_type, class Stepper>
319
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
320
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getInverseBendingRadius(const value_type& angle)
321
{
322 323 324 325
    if (angle != 0.0)
        return rotate(angle, h_m);
    else
        return h_m;
326 327 328
}

template<typename Value_type, typename Size_type, class Stepper>
329
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
330
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getPathLength(const value_type& angle)
331
{
332 333 334 335
    if (angle != 0.0)
        return rotate(angle, ds_m);
    else
        return ds_m;
336 337 338
}

template<typename Value_type, typename Size_type, class Stepper>
339
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
340
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getFieldIndex(const value_type& angle)
341
{
342 343 344
    if (angle != 0.0)
        return rotate(angle, fidx_m);
    return ds_m;
345 346 347
}

template<typename Value_type, typename Size_type, class Stepper>
348 349 350
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);
351

352 353 354
    // compute nzcross_m
    if (!vertOscDone_m)
        computeVerticalOscillations();
355

356 357 358 359
    // compute vertical tune
    value_type nuz = computeTune(z_m,pz_m[1],nzcross_m);

    return std::make_pair(nur,nuz);
360 361 362
}

template<typename Value_type, typename Size_type, class Stepper>
363
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
364
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getOrbit(value_type angle)
365
{
366 367 368 369
    if (angle != 0.0)
        return rotate(angle, r_m);
    else
        return r_m;
370 371 372 373
}

template<typename Value_type, typename Size_type, class Stepper>
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
374
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getMomentum(value_type angle)
375 376
{
    container_type pr = pr_m;
377 378 379
    
    if (angle != 0.0)
        pr = rotate(angle, pr);
380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395
    
    // 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
     * 
396
     * The momentum in \beta * \gamma is obtained by dividing by "a"
397
     */
398
    value_type factor =  1.0 / acon_m(wo_m);
399
    std::for_each(pr.begin(), pr.end(), [factor](value_type& p) { p *= factor; });
400
    
401
    return pr;
402 403 404
}

template<typename Value_type, typename Size_type, class Stepper>
405 406 407
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::value_type
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getGamma()
{
408
    return gamma_m;
409 410 411
}

template<typename Value_type, typename Size_type, class Stepper>
412 413 414
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::value_type
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getAverageRadius()
{
415
    return ravg_m;
416 417 418
}

template<typename Value_type, typename Size_type, class Stepper>
419 420
typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::value_type
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getFrequencyError()
421
{
422 423 424
    // if the vertical oscillations aren't computed, we have to, since there we also compuote the frequency error.
    if(!vertOscDone_m)
        computeVerticalOscillations();
425

426
    return phase_m;
427 428 429 430
}

template<typename Value_type, typename Size_type, class Stepper>
inline bool ClosedOrbitFinder<Value_type, Size_type, Stepper>::isConverged() {
431
    return converged_m;
432
}
433 434 435 436 437 438 439

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

template<typename Value_type, typename Size_type, class Stepper>
bool ClosedOrbitFinder<Value_type, Size_type, Stepper>::findOrbit(value_type accuracy, size_type maxit) {
440 441 442 443 444
    /* REMARK TO GORDON
     * q' = 1/b = 1/bcon
     * a' = a = acon
     */

adelmann's avatar
adelmann committed
445
    // READ IN MAGNETIC FIELD: ONLY FOR STAND-ALONE PROGRAM
446
    
447
    value_type bint, brint, btint;
448

449 450 451
    // 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);
452

453
    // store acon and bcon locally
454 455
    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)
456 457 458 459 460 461 462 463 464 465 466 467

    // 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)
468
    {
469 470 471
        r_m[idx] = y[0];
        pr_m[idx] = y[1];

472
        // count number of crossings (excluding starting point --> idx>0)
473 474 475 476 477 478
        nxcross_m += (idx > 0) * (y[4] * xold < 0);
        xold = y[4];
        ++idx;
    };

    // define the six ODEs (using lambda function)
479 480 481 482
    std::function<void(const state_type&, state_type&, const double)> orbit_integration = [&](const state_type &y,
                                                                                              state_type &dydt,
                                                                                              const double theta)
    {
483 484
        pr2 = y[1] * y[1];
        if (p2 < pr2)
485
            throw OpalException("ClosedOrbitFinder::findOrbit()", "p_{r}^2 > p^{2} (defined in Gordon paper) --> Square root of negative number.");
486

487 488 489 490
        // Gordon, formula (5c)
        ptheta = std::sqrt(p2 - pr2);
        invptheta = 1.0 / ptheta;

491
        // interpolate values of magnetic field
492
        bField_m.interpolate(bint, brint, btint, y[0], theta);
493

494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509
        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);
510

511 512 513 514 515 516 517 518
    // 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;
519 520 521 522

    /*
     * Christian:
     * N = 1440 ---> N = 720 ---> dtheta = 2PI/720 --> nsteps = 721
523
     *
524
     * 0, 2, 4, ... ---> jeden zweiten berechnene: 1, 3, 5, ... interpolieren --> 1440 Werte
525
     *
526 527
     * Matthias:
     * N = 1440 --> dtheta = 2PI/1440 --> nsteps = 1441
528
     *
529
     * 0, 1, 2, 3, 4, 5, ... --> 1440 Werte
530
     *
531
     */
532

Andreas Adelmann's avatar
Andreas Adelmann committed
533 534 535 536 537 538 539 540
    // step size of energy
    value_type dE; 

    if (Emin_m == Emax_m)
      dE = 0.0;
    else
      dE = (E_m - Emin_m) / (Emax_m - Emin_m);

541 542
    // iterate until suggested energy (start with minimum energy)
    value_type E = Emin_m;
543

adelmann's avatar
adelmann committed
544 545
    // energy increase not more than 0.25
    dE = (dE > 0.25) ? 0.25 : dE;
546 547

    // energy dependent values
548
    value_type en = E / E0_m;                      // en = E/E0 = E/(mc^2) (E0 is potential energy)
549 550 551 552 553 554 555 556 557
    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
558 559 560 561 562

    container_type init;
    if (rguess_m < 0)
      init = {beta * acon, 0.0};
    else
Andreas Adelmann's avatar
Andreas Adelmann committed
563
      init = {rguess_m/1000.0, 0.0};
564 565 566

    // store initial values for updating values for higher energies
    container_type previous_init = {0.0, 0.0};
567

568 569
    do {
        
570
        // (re-)set inital values for r and pr
571
        r_m[0] = init[0];
572
        pr_m[0] = init[1];
573

574 575 576 577 578 579 580 581 582 583 584 585
        // 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]};
586

587 588
            // 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);
589

590 591 592 593 594
            // write new state
            x_m[0] = y[2];
            px_m[0] = y[3];
            x_m[1] = y[4];
            px_m[1] = y[5];
595

596 597 598 599
            // 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)
600

601 602 603 604
            // 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)
605

606 607 608
            // improved initial values; Gordon, formula (17) (here it's used for higher energies)
            init[0] += delta[0];
            init[1] += delta[1];
609

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

614 615
        // reset iteration counter
        niterations = 0;
616

617 618
        // reset correction term
        delta[0] = delta[1] = 0.0;
adelmann's avatar
adelmann committed
619 620 621 622 623 624

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

626
        // set constants for new energy E
627
        en = E / E0_m;                     // en = E/E0 = E/(mc^2) (E0 is potential energy)
628 629 630 631
        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);
632 633


634
    } while (E != E_m);
635

636 637 638 639 640
    /* 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()
641

642 643 644
    // 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();
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

adelmann's avatar
adelmann committed
699 700 701 702 703
    /* domain_m = true --> only integrated over a single sector --> multiply by nSector_m to
     * get correct tune.
     */
    if (domain_m)
        mu *= nSector_m;
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
        bField_m.interpolate(bint, brint, btint, r_m[i], theta);
735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750
        bint *= invbcon;
        brint *= invbcon;
        btint *= invbcon;

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

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

760 761 762 763 764
    vertOscDone_m = true;

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

765 766
    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)
767 768 769 770 771 772 773
    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
774
    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)
775 776

    // define the ODEs (using lambda function)
777 778 779 780
    std::function<void(const state_type&, state_type&, const double)> vertical = [&](const state_type &y,
                                                                                     state_type &dydt,
                                                                                     const double theta)
    {
781
        pr2 = y[1] * y[1];
782
        if (p2 < pr2) {
783 784
            throw OpalException("ClosedOrbitFinder::computeVerticalOscillations()",
                                "p_{r}^2 > p^{2} (defined in Gordon paper) --> Square root of negative number.");
785
        }
786

787 788 789 790 791
        // Gordon, formula (5c)
        ptheta = std::sqrt(p2 - pr2);
        invptheta = 1.0 / ptheta;

        // intepolate values of magnetic field
792
        bField_m.interpolate(bint, brint, btint, y[0], theta);
frey_m's avatar
frey_m committed
793
        
794 795 796
        bint *= invbcon;
        brint *= invbcon;
        btint *= invbcon;
797

798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827
        // 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);
828

829 830
    // 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);
831

832 833 834
    // remove last element again (--> size = N_m, capacity = N_m+1)
    r_m.pop_back();
    pr_m.pop_back();
835

836 837 838 839 840
    // write new state
    z_m[0] = y[2];
    pz_m[0] = y[3];
    z_m[1] = y[4];
    pz_m[1] = y[5];
841
    phase_m = y[6] * Physics::u_two_pi; // / (2.0 * Physics::pi);
842

adelmann's avatar
adelmann committed
843 844 845 846 847
    /* domain_m = true --> only integrated over a single sector
     * --> multiply by nSector_m to get correct phase_m
     */
    if (domain_m)
        phase_m *= nSector_m;
848 849
}

850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871
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
872
#endif