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

frey_m's avatar
frey_m committed
30 31
#include "AbstractObjects/OpalData.h"

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

34 35
// include headers for integration
#include <boost/numeric/odeint/integrate/integrate_n_steps.hpp>
frey_m's avatar
frey_m committed
36 37 38
#include <boost/filesystem.hpp>

extern Inform *gmsg;
39

40
/// Finds a closed orbit of a cyclotron for a given energy
41 42 43
template<typename Value_type, typename Size_type, class Stepper>
class ClosedOrbitFinder
{
44 45 46 47 48 49 50 51 52
    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;
frey_m's avatar
frey_m committed
53 54
        
        typedef std::function<void(const state_type&, state_type&, const double)> function_t;
55 56 57

        /// Sets the initial values for the integration and calls findOrbit().
        /*!
58
         * @param E0 is the potential energy (particle energy at rest) [MeV].
adelmann's avatar
adelmann committed
59
         * @param N specifies the number of splits (2pi/N), i.e number of integration steps
60
         * @param cycl is the cyclotron element
61 62
         * @param domain is a boolean (default: true). If "true" the closed orbit is computed over a single sector,
         * otherwise over 2*pi.
63
         */
frey_m's avatar
frey_m committed
64
        ClosedOrbitFinder(value_type E0, size_type N,
65
                          Cyclotron* cycl, bool domain = true);
66 67

        /// Returns the inverse bending radius (size of container N+1)
68
        container_type getInverseBendingRadius(const value_type& angle = 0);
69 70

        /// Returns the step lengths of the path (size of container N+1)
71
        container_type getPathLength(const value_type& angle = 0);
72 73

        /// Returns the field index (size of container N+1)
74
        container_type getFieldIndex(const value_type& angle = 0);
75 76 77 78

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

79 80 81 82 83 84 85
        /// 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°).
         */
86 87
        container_type getOrbit(value_type angle = 0);

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

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

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

        /// 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
111
         * @param ekin energy for which to find closed orbit (in tune mode: upper limit of range)
frey_m's avatar
frey_m committed
112
         * @param dE step increase [MeV]
113
         * @param rguess initial radius guess in [mm]
frey_m's avatar
frey_m committed
114
         * @param isTuneMode comptute tunes of all energies in one sweep
115
         */
116 117 118 119
        bool findOrbit(value_type accuracy, size_type maxit,
                       value_type ekin,
                       value_type dE = 1.0, value_type rguess = -1.0,
                       bool isTuneMode = false);
120

121
        /// Fills in the values of h_m, ds_m, fidx_m.
122
        void computeOrbitProperties(const value_type& E);
123

124
    private:
125 126 127 128 129 130 131 132
        /// 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);

frey_m's avatar
frey_m committed
133 134 135 136
        // Compute closed orbit for given energy
        bool findOrbitOfEnergy_m(const value_type&, container_type&, value_type&,
                                 const value_type&, size_type);
        
adelmann's avatar
adelmann committed
137
        /// This function computes nzcross_ which is used to compute the tune in z-direction and the frequency error
138
//         void computeVerticalOscillations();
139 140 141
        
        /// This function rotates the calculated closed orbit finder properties to the initial angle
        container_type rotate(value_type angle, container_type& orbitProperty);
142 143 144 145 146

        /// 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
147
        /// Stores current position in vertical direction for the solutions of the ODE with different initial values
148
        std::array<value_type,2> z_m; // z_m = [z1, z2]
frey_m's avatar
frey_m committed
149
        /// Stores current momenta in vertical direction for the solutions of the ODE with different initial values
150 151 152 153 154 155 156 157
        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;
158 159
        /// Stores the vertical oribt (size: N_m+1)
        container_type vz_m;
160 161
        /// Stores the radial momentum
        container_type pr_m;
162 163
        /// Stores the vertical momentum
        container_type vpz_m;
164 165 166 167 168 169 170 171
        /// 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

172
        /// Is the rest mass [MeV / c**2]
173 174
        value_type E0_m;
        
175
        /// Is the nominal orbital frequency
176 177 178
        /* (see paper of Dr. C. Baumgarten: "Transverse-Longitudinal
         * Coupling by Space Charge in Cyclotrons" (2012), formula (1))
         */
179
        value_type wo_m;
adelmann's avatar
adelmann committed
180
        /// Number of integration steps
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 average radius
        value_type ravg_m;

        /// Is the phase
        value_type phase_m;

191
        /**
192 193 194 195
         * 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)
         */
196 197
        value_type lastOrbitVal_m;

198 199 200 201 202
        /**
         * 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)
         */
203
        value_type lastMomentumVal_m;
204 205

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

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

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

236 237 238
template<typename Value_type, typename Size_type, class Stepper>
ClosedOrbitFinder<Value_type,
                  Size_type,
frey_m's avatar
frey_m committed
239
                  Stepper>::ClosedOrbitFinder(value_type E0,
240
                                              size_type N, Cyclotron* cycl,
241
                                              bool domain)
242 243 244
    : nxcross_m(0)
    , nzcross_m(0)
    , E0_m(E0)
245
    , wo_m(cycl->getRfFrequ()*1E6/cycl->getCyclHarm()*2.0*Physics::pi)
246 247 248 249 250 251 252 253
    , N_m(N)
    , dtheta_m(Physics::two_pi/value_type(N))
    , ravg_m(0)
    , phase_m(0)
    , lastOrbitVal_m(0.0)
    , lastMomentumVal_m(0.0)
    , 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
        throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()",
                            "Incorrect cyclotron energy (MeV) bounds: Maximum cyclotron energy smaller than minimum cyclotron energy.");
    
adelmann's avatar
adelmann committed
261 262
    // if domain_m = true --> integrate over a single sector
    if (domain_m) {
263
        N_m /=  cycl_m->getSymmetry();
adelmann's avatar
adelmann committed
264
    }
frey_m's avatar
frey_m committed
265 266 267
    
    cycl_m->read(cycl_m->getFieldFlag(cycl_m->getCyclotronType()),
                 cycl_m->getBScale());
268

269 270 271 272 273
    // 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
274 275
    r_m.reserve(N_m + 1);
    pr_m.reserve(N_m + 1);
276 277
    vz_m.reserve(N_m + 1);
    vpz_m.reserve(N_m + 1);
278

279
    // reserve memory of N_m for the properties (--> size = 0, capacity = N_m)
adelmann's avatar
adelmann committed
280 281 282
    h_m.reserve(N_m);
    ds_m.reserve(N_m);
    fidx_m.reserve(N_m);
283 284 285
}

template<typename Value_type, typename Size_type, class Stepper>
286
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
287
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getInverseBendingRadius(const value_type& angle)
288
{
289 290 291 292
    if (angle != 0.0)
        return rotate(angle, h_m);
    else
        return h_m;
293 294 295
}

template<typename Value_type, typename Size_type, class Stepper>
296
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
297
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getPathLength(const value_type& angle)
298
{
299 300 301 302
    if (angle != 0.0)
        return rotate(angle, ds_m);
    else
        return ds_m;
303 304 305
}

template<typename Value_type, typename Size_type, class Stepper>
306
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
307
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getFieldIndex(const value_type& angle)
308
{
309 310
    if (angle != 0.0)
        return rotate(angle, fidx_m);
frey_m's avatar
frey_m committed
311
    return fidx_m;
312 313 314
}

template<typename Value_type, typename Size_type, class Stepper>
315 316 317
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);
318

319 320 321 322
    // compute vertical tune
    value_type nuz = computeTune(z_m,pz_m[1],nzcross_m);

    return std::make_pair(nur,nuz);
323 324 325
}

template<typename Value_type, typename Size_type, class Stepper>
326
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
327
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getOrbit(value_type angle)
328
{
329 330 331 332
    if (angle != 0.0)
        return rotate(angle, r_m);
    else
        return r_m;
333 334 335 336
}

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

template<typename Value_type, typename Size_type, class Stepper>
368 369 370
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::value_type
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getAverageRadius()
{
371
    return ravg_m;
372 373 374
}

template<typename Value_type, typename Size_type, class Stepper>
375 376
typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::value_type
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getFrequencyError()
377
{
378
    return phase_m;
379 380 381 382 383 384
}

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

frey_m's avatar
frey_m committed
385 386


387
template<typename Value_type, typename Size_type, class Stepper>
388 389
bool ClosedOrbitFinder<Value_type, Size_type, Stepper>::findOrbit(value_type accuracy,
                                                                  size_type maxit,
frey_m's avatar
frey_m committed
390
                                                                  value_type ekin,
frey_m's avatar
frey_m committed
391
                                                                  value_type dE,
frey_m's avatar
frey_m committed
392 393
                                                                  value_type rguess,
                                                                  bool isTuneMode)
394
{
395 396 397 398
    /* REMARK TO GORDON
     * q' = 1/b = 1/bcon
     * a' = a = acon
     */
399

400 401 402
    // 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);
403 404
    vz_m.resize(N_m+1);
    vpz_m.resize(N_m+1);
405

406 407
    // store acon locally
    value_type acon = acon_m(wo_m);            // [acon] = m
408 409
    // amplitude of error; Gordon, formula (18) (a = a')
    value_type error = std::numeric_limits<value_type>::max();
410 411 412 413

    /*
     * Christian:
     * N = 1440 ---> N = 720 ---> dtheta = 2PI/720 --> nsteps = 721
414
     *
415
     * 0, 2, 4, ... ---> jeden zweiten berechnene: 1, 3, 5, ... interpolieren --> 1440 Werte
416
     *
417 418
     * Matthias:
     * N = 1440 --> dtheta = 2PI/1440 --> nsteps = 1441
419
     *
420
     * 0, 1, 2, 3, 4, 5, ... --> 1440 Werte
421
     *
422
     */
423

424 425 426
    value_type E;     // starting energy
    value_type E_fin; // final    energy

frey_m's avatar
frey_m committed
427
    if ( isTuneMode ) {
428 429 430 431 432
        E     = cycl_m->getFMLowE();
        E_fin = cycl_m->getFMHighE();
    } else {
        E     = ekin;
        E_fin = ekin;
frey_m's avatar
frey_m committed
433
    }
Andreas Adelmann's avatar
Andreas Adelmann committed
434

frey_m's avatar
frey_m committed
435 436 437 438
    namespace fs = boost::filesystem;
    fs::path dir = OpalData::getInstance()->getInputBasename();
    dir = dir.parent_path() / "data";
    std::string tunefile = (dir / "tunes.dat").string();
439

frey_m's avatar
frey_m committed
440 441
    if ( isTuneMode ) {
        std::ofstream out(tunefile, std::ios::out);
442

frey_m's avatar
frey_m committed
443 444 445 446 447 448 449
        out << std::left
            << std::setw(15) << "energy [MeV]"
            << std::setw(15) << "radius [m]"
            << std::setw(15) << "nu_r"
            << std::setw(15) << "nu_z"
            << std::endl;
    }
450

451 452 453
    // iterate until suggested energy (start with minimum energy)
    // increase energy by dE
    for (; E <= E_fin ; E+=dE) {
frey_m's avatar
frey_m committed
454
        error = std::numeric_limits<value_type>::max();
455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478

        // energy dependent values
        value_type en     = E / E0_m;        // en = E/E0 = E/(mc^2) (E0 is potential energy)
        value_type gamma  = en + 1.0;
        value_type gamma2 = gamma * gamma;   // = gamma^2
        value_type beta   = std::sqrt(1.0 - 1.0 / gamma2);

        // set initial values for radius and radial momentum for lowest energy Emin
        // orbit, [r] = m;  Gordon, formula (20)
        // radial momentum; Gordon, formula (20)

        container_type init;
        //      r            pr   z    pz
        init = {beta * acon, 0.0, 0.0, 1.0};

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

        std::fill(  r_m.begin(),   r_m.end(), 0);
        std::fill( pr_m.begin(),  pr_m.end(), 0);
        std::fill( vz_m.begin(),  vz_m.end(), 0);
        std::fill(vpz_m.begin(), vpz_m.end(), 0);

479
        // (re-)set inital values for r and pr
480 481 482
        r_m[0]   = init[0];
        pr_m[0]  = init[1];
        vz_m[0]  = init[2];
483
        vpz_m[0] = init[3];
484

frey_m's avatar
frey_m committed
485
        if ( !this->findOrbitOfEnergy_m(E, init, error, accuracy, maxit) ) {
486 487 488 489
            /* throw OpalException("ClosedOrbitFinder::findOrbitOfEnergy()", */
            /*                     "Didn't converge for energy " + std::to_string(E) + " MeV."); */
            *gmsg << "Didn't converge for energy " + std::to_string(E) + " MeV." << endl;
            continue;
frey_m's avatar
frey_m committed
490
        }
491

frey_m's avatar
frey_m committed
492
        if ( isTuneMode ) {
493 494
            this->computeOrbitProperties(E);

frey_m's avatar
frey_m committed
495 496
            std::pair<value_type , value_type > tunes = this->getTunes();

497 498
            value_type reo = this->getOrbit(   cycl_m->getPHIinit())[0];
            value_type peo = this->getMomentum(cycl_m->getPHIinit())[0];
frey_m's avatar
frey_m committed
499 500 501 502 503 504 505 506 507 508 509 510

            *gmsg << "* ----------------------------" << endl
                  << "* Closed orbit info (Gordon units):" << endl
                  << "*" << endl
                  << "* kinetic energy:   " << E              << " [MeV]" << endl
                  << "* average radius:   " << ravg_m         << " [m]" << endl
                  << "* initial radius:   " << reo            << " [m]" << endl
                  << "* initial momentum: " << peo            << " [Beta Gamma]" << endl
                  << "* frequency error:  " << phase_m        << endl
                  << "* horizontal tune:  " << tunes.first    << endl
                  << "* vertical tune:    " << tunes.second   << endl
                  << "* ----------------------------" << endl << endl;
511

frey_m's avatar
frey_m committed
512 513 514 515 516 517 518 519
            std::ofstream out(tunefile, std::ios::app);
            out << std::left
                << std::setw(15) << E
                << std::setw(15) << reo
                << std::setw(15) << tunes.first
                << std::setw(15) << tunes.second << std::endl;
            out.close();
        }
520
    }
521

522 523 524 525 526
    /* 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()
527

528 529 530
    // 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();
531

532 533 534 535 536
    /* 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();
537

538 539 540
    // returns true if converged, otherwise false
    return error < accuracy;
}
541

frey_m's avatar
frey_m committed
542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560
template<typename Value_type, typename Size_type, class Stepper>
bool ClosedOrbitFinder<Value_type, Size_type, Stepper>::findOrbitOfEnergy_m(
    const value_type& E,
    container_type& init,
    value_type& error,
    const value_type& accuracy,
    size_type maxit)
{
    value_type bint, brint, btint;
    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)
    
    value_type xold = 0.0;                              // for counting nxcross
    value_type zold = 0.0;                              // for counting nzcross

    // 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)
    {
561 562 563
        r_m[idx]   = y[0];
        pr_m[idx]  = y[1];
        vz_m[idx]  = y[6];
frey_m's avatar
frey_m committed
564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589
        vpz_m[idx] = y[7];

        // count number of crossings (excluding starting point --> idx>0)
        nxcross_m += (idx > 0) * (y[4] * xold < 0);
        xold = y[4];
        
        // number of times z2 changes sign
        nzcross_m += (idx > 0) * (y[10] * zold < 0);
        zold = y[10];
        
        ++idx;
    };
    
    // define initial state container for integration: y = {r, pr, x1, px1, x2, px2,
    //                                                      z, pz, z1, pz1, z2, pz2,
    //                                                      phase}
    state_type y(11);
    
    // 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};
    // if niterations > maxit --> stop iteration
    size_type niterations = 0;
    
    // energy dependent values
590 591
    value_type en     = E / E0_m;                      // en = E/E0 = E/(mc^2) (E0 is potential energy)
    value_type gamma  = en + 1.0;
frey_m's avatar
frey_m committed
592
    value_type p = acon_m(wo_m) * std::sqrt(en * (2.0 + en));     // momentum [p] = m; Gordon, formula (3)
593 594
    value_type gamma2    = gamma * gamma;           // = gamma^2
    value_type invgamma4 = 1.0 / (gamma2 * gamma2); // = 1/gamma^4
frey_m's avatar
frey_m committed
595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638
    
    
    // helper constants
    value_type p2 = p * p;                                // p^2 = p*p
    value_type pr2;                                     // squared radial momentum (pr^2 = pr*pr)
    value_type ptheta, invptheta;                       // Gordon, formula (5c)

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

        // Gordon, formula (5c)
        ptheta = std::sqrt(p2 - pr2);
        invptheta = 1.0 / ptheta;
        // interpolate values of magnetic field
        cycl_m->apply(y[0], y[6], theta, brint, btint, bint);
        
        bint *= invbcon;
        brint *= invbcon;
        btint *= 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];
        }
        
        // 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
639
        dydt[12] = y[0] * invptheta * gamma - 1;
frey_m's avatar
frey_m committed
640 641 642 643 644
        
    };

    // integrate until error smaller than user-define accuracy
    do {
645
        // (re-)set initial values
frey_m's avatar
frey_m committed
646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666
        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)
        z_m[0]  = 1.0;
        pz_m[0] = 0.0;
        z_m[1]  = 0.0;
        pz_m[1] = 1.0;
        phase_m = 0.0;
        nxcross_m = 0;               // counts the number of crossings of x-axis (excluding first step)
        nzcross_m = 0;
        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],
             init[2], init[3],
             z_m[0], pz_m[0], z_m[1], pz_m[1],
             phase_m
        };

667 668 669 670 671 672 673
        try {
            // 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);
        } catch(OpalException &) {
            break;
        }

frey_m's avatar
frey_m committed
674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690
        // write new state
        x_m[0] = y[2];
        px_m[0] = y[3];
        x_m[1] = y[4];
        px_m[1] = y[5];
        
        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);
        
        // 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)

691
        // correct initial values of r and pr
frey_m's avatar
frey_m committed
692 693 694 695 696 697 698 699 700 701 702 703
        value_type 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)

        // improved initial values; Gordon, formula (17) (here it's used for higher energies)
        init[0] += delta[0];
        init[1] += delta[1];
        
        // 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);
    
704 705 706
    if (error > accuracy)
        *gmsg << "findOrbit not converged after " << maxit << " iterations with error: " << error << ". Needed accuracy " << accuracy << endl;

frey_m's avatar
frey_m committed
707 708 709
    return (error < accuracy);
}

710
template<typename Value_type, typename Size_type, class Stepper>
711 712 713
Value_type ClosedOrbitFinder<Value_type, Size_type, Stepper>::computeTune(const std::array<value_type,2>& y,
                                                                          value_type py2, size_type ncross)
{
714
    // Y = [y1, y2; py1, py2]
715

716 717
    // cos(mu)
    value_type cos = 0.5 * (y[0] + py2);
718
    
719
    value_type mu;
720

721 722
    // sign of sin(mu) has to be equal to y2
    bool negative = std::signbit(y[1]);
723

724
    bool uneven = (ncross % 2);
725

726 727 728
    if (std::fabs(cos) > 1.0) {
        // store the number of crossings
        value_type n = ncross;
729

730 731
        if (uneven)
            n = ncross - 1;
732

733 734
        // Gordon, formula (36b)
        value_type muPrime = -std::acosh(std::fabs(cos));      // mu'
735
        mu = n * Physics::pi + muPrime;
736

737 738 739 740 741 742 743
    } 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.
        */
744

745
        // Gordon, formula (36)
746
        mu = ncross * Physics::pi + muPrime;
747

748 749
        // if sign(y[1]) > 0 && sign(sin(mu)) < 0
        if (!negative && std::signbit(std::sin(mu))) {
750
            mu = ncross * Physics::pi - muPrime;
751
        } else if (negative && !std::signbit(std::sin(mu))) {
752
            mu = ncross * Physics::pi - muPrime + Physics::two_pi;
753 754
        }
    }
755

756
    // nu = mu/theta, where theta = integration domain
757

758
    /* domain_m = true --> only integrated over a single sector --> multiply by cycl_m->getSymmetry() to
adelmann's avatar
adelmann committed
759 760 761
     * get correct tune.
     */
    if (domain_m)
762
        mu *= cycl_m->getSymmetry();
763

764
    return mu * Physics::u_two_pi;
765 766 767
}

template<typename Value_type, typename Size_type, class Stepper>
768
void ClosedOrbitFinder<Value_type, Size_type, Stepper>::computeOrbitProperties(const value_type& E) {
769
    /*
770 771 772 773 774
     * 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
     */
775

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

779
    value_type invbcon = 1.0 / bcon_m(E0_m, wo_m);
780
    value_type en = E / E0_m;                                  // en = E/E0 = E/(mc^2) (E0 is potential energy)
781
    value_type p = acon_m(wo_m) * std::sqrt(en * (2.0 + en));    // momentum [p] = m; Gordon, formula (3)
782 783 784
    value_type p2 = p * p;
    value_type theta = 0.0;                                             // angle for interpolating
    value_type ptheta;
785

786 787 788 789 790 791 792
    // 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
793
        cycl_m->apply(r_m[i], vz_m[i], theta, brint, btint, bint);
794 795 796
        bint *= invbcon;
        brint *= invbcon;
        btint *= invbcon;
frey_m's avatar
frey_m committed
797
        
798 799 800 801
        // inverse bending radius
        h_m[i] = bint / p;

        // local field index
802 803 804 805
        if (p < pr_m[i])
            throw OpalException("ClosedOrbitFinder::computeTune()",
                                "p_{r}^2 > p^{2} " + std::to_string(p) + " " + std::to_string(pr_m[i]) + " (defined in Gordon paper) --> Square root of negative number.");

806
        ptheta = std::sqrt(p2 - pr_m[i] * pr_m[i]);
807

808 809 810 811 812 813 814
        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;
815
    }
816 817 818

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

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