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

13 14 15
#ifndef CLOSEDORBITFINDER_H
#define CLOSEDORBITFINDER_H

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

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

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 101 102 103

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

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

adelmann's avatar
adelmann committed
104 105
        /// Returns the frequency error
        value_type getFrequencyError();
106 107 108 109 110 111 112 113

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

125
        /// Fills in the values of h_m, ds_m, fidx_m.
126 127
        void computeOrbitProperties();

128
    private:
129 130 131 132 133 134 135 136
        /// 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
137 138 139 140
        // 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
141
        /// This function computes nzcross_ which is used to compute the tune in z-direction and the frequency error
142
//         void computeVerticalOscillations();
143 144 145
        
        /// This function rotates the calculated closed orbit finder properties to the initial angle
        container_type rotate(value_type angle, container_type& orbitProperty);
146 147 148 149 150

        /// 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
151
        /// Stores current position in vertical direction for the solutions of the ODE with different initial values
152
        std::array<value_type,2> z_m; // z_m = [z1, z2]
frey_m's avatar
frey_m committed
153
        /// Stores current momenta in vertical direction for the solutions of the ODE with different initial values
154 155 156 157 158 159 160 161
        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;
162 163
        /// Stores the vertical oribt (size: N_m+1)
        container_type vz_m;
164 165
        /// Stores the radial momentum
        container_type pr_m;
166 167
        /// Stores the vertical momentum
        container_type vpz_m;
168 169 170 171 172 173 174 175 176 177
        /// 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;
178 179 180 181
        
        /// Is the potential energy [MeV]
        value_type E0_m;
        
182
        /// Is the nominal orbital frequency
183 184 185
        /* (see paper of Dr. C. Baumgarten: "Transverse-Longitudinal
         * Coupling by Space Charge in Cyclotrons" (2012), formula (1))
         */
186
        value_type wo_m;
adelmann's avatar
adelmann committed
187
        /// Number of integration steps
188 189 190 191 192 193 194 195 196 197 198 199 200
        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;

201
        /**
202 203 204 205
         * 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)
         */
206 207
        value_type lastOrbitVal_m;

208 209 210 211 212
        /**
         * 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)
         */
213
        value_type lastMomentumVal_m;
214 215

        /**
216 217 218
         * Boolean which is true by default. "true": orbit integration over one sector only, "false": integration
         * over 2*pi
         */
adelmann's avatar
adelmann committed
219
        bool domain_m;
220

221 222
        /// Defines the stepper for integration of the ODE's
        Stepper stepper_m;
223 224 225 226 227 228 229 230 231 232 233 234 235 236 237
        
        /*!
         * 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);
        };
238
        
239
        Cyclotron* cycl_m;
240 241 242 243 244 245
};

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

246 247 248
template<typename Value_type, typename Size_type, class Stepper>
ClosedOrbitFinder<Value_type,
                  Size_type,
frey_m's avatar
frey_m committed
249
                  Stepper>::ClosedOrbitFinder(value_type E0,
250
                                              size_type N, Cyclotron* cycl,
251
                                              bool domain)
252 253
    : nxcross_m(0)
    , nzcross_m(0)
frey_m's avatar
frey_m committed
254
    , E_m(0.0)
255
    , E0_m(E0)
256
    , wo_m(cycl->getRfFrequ()*1E6/cycl->getCyclHarm()*2.0*Physics::pi)
257 258
    , N_m(N)
    , dtheta_m(Physics::two_pi/value_type(N))
frey_m's avatar
frey_m committed
259
    , gamma_m(0.0)
260 261 262 263 264 265
    , ravg_m(0)
    , phase_m(0)
    , lastOrbitVal_m(0.0)
    , lastMomentumVal_m(0.0)
    , domain_m(domain)
    , stepper_m()
266
    , cycl_m(cycl)
267
{
frey_m's avatar
frey_m committed
268
    
269
    if ( cycl_m->getFMLowE() > cycl_m->getFMHighE() )
frey_m's avatar
frey_m committed
270 271 272
        throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()",
                            "Incorrect cyclotron energy (MeV) bounds: Maximum cyclotron energy smaller than minimum cyclotron energy.");
    
adelmann's avatar
adelmann committed
273 274
    // if domain_m = true --> integrate over a single sector
    if (domain_m) {
275
        N_m /=  cycl_m->getSymmetry();
adelmann's avatar
adelmann committed
276
    }
frey_m's avatar
frey_m committed
277 278 279
    
    cycl_m->read(cycl_m->getFieldFlag(cycl_m->getCyclotronType()),
                 cycl_m->getBScale());
280

281 282 283 284 285
    // 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
286 287
    r_m.reserve(N_m + 1);
    pr_m.reserve(N_m + 1);
288 289
    vz_m.reserve(N_m + 1);
    vpz_m.reserve(N_m + 1);
290

291
    // reserve memory of N_m for the properties (--> size = 0, capacity = N_m)
adelmann's avatar
adelmann committed
292 293 294
    h_m.reserve(N_m);
    ds_m.reserve(N_m);
    fidx_m.reserve(N_m);
295 296 297
}

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

template<typename Value_type, typename Size_type, class Stepper>
308
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
309
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getPathLength(const value_type& angle)
310
{
311 312 313 314
    if (angle != 0.0)
        return rotate(angle, ds_m);
    else
        return ds_m;
315 316 317
}

template<typename Value_type, typename Size_type, class Stepper>
318
inline typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::container_type
319
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getFieldIndex(const value_type& angle)
320
{
321 322
    if (angle != 0.0)
        return rotate(angle, fidx_m);
frey_m's avatar
frey_m committed
323
    return fidx_m;
324 325 326
}

template<typename Value_type, typename Size_type, class Stepper>
327 328 329
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);
330

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

    return std::make_pair(nur,nuz);
335 336 337
}

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

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

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

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

template<typename Value_type, typename Size_type, class Stepper>
394 395
typename ClosedOrbitFinder<Value_type, Size_type, Stepper>::value_type
    ClosedOrbitFinder<Value_type, Size_type, Stepper>::getFrequencyError()
396
{
397
    return phase_m;
398 399 400 401 402 403
}

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

frey_m's avatar
frey_m committed
404 405


406
template<typename Value_type, typename Size_type, class Stepper>
407 408
bool ClosedOrbitFinder<Value_type, Size_type, Stepper>::findOrbit(value_type accuracy,
                                                                  size_type maxit,
frey_m's avatar
frey_m committed
409
                                                                  value_type ekin,
frey_m's avatar
frey_m committed
410
                                                                  value_type dE,
frey_m's avatar
frey_m committed
411 412
                                                                  value_type rguess,
                                                                  bool isTuneMode)
413
{
414 415 416 417
    /* REMARK TO GORDON
     * q' = 1/b = 1/bcon
     * a' = a = acon
     */
418
    
frey_m's avatar
frey_m committed
419 420
    E_m = ekin;
    gamma_m = E_m / E0_m + 1.0;
421

422 423 424
    // 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);
425 426
    vz_m.resize(N_m+1);
    vpz_m.resize(N_m+1);
427

428
    // store acon and bcon locally
429
    value_type acon = acon_m(wo_m);               // [acon] = m
frey_m's avatar
frey_m committed
430
    
431 432
    // amplitude of error; Gordon, formula (18) (a = a')
    value_type error = std::numeric_limits<value_type>::max();
433 434 435 436

    /*
     * Christian:
     * N = 1440 ---> N = 720 ---> dtheta = 2PI/720 --> nsteps = 721
437
     *
438
     * 0, 2, 4, ... ---> jeden zweiten berechnene: 1, 3, 5, ... interpolieren --> 1440 Werte
439
     *
440 441
     * Matthias:
     * N = 1440 --> dtheta = 2PI/1440 --> nsteps = 1441
442
     *
443
     * 0, 1, 2, 3, 4, 5, ... --> 1440 Werte
444
     *
445
     */
446

447
    // iterate until suggested energy (start with minimum energy)
448
    value_type E = cycl_m->getFMLowE();
449

450
    // energy dependent values
451
    value_type en = E / E0_m;                      // en = E/E0 = E/(mc^2) (E0 is potential energy)
452 453 454 455 456 457
    value_type gamma2 = (1.0 + en) * (1.0 + en);          // = 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)
Andreas Adelmann's avatar
Andreas Adelmann committed
458 459

    container_type init;
frey_m's avatar
frey_m committed
460 461 462 463 464 465
    //      r            pr   z    pz
    init = {beta * acon, 0.0, 0.0, 1.0};

    if (rguess >= 0.0) {
        init[0] = rguess * 0.001;
    }
frey_m's avatar
frey_m committed
466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481
    
    namespace fs = boost::filesystem;
    fs::path dir = OpalData::getInstance()->getInputBasename();
    dir = dir.parent_path() / "data";
    std::string tunefile = (dir / "tunes.dat").string();
    
    if ( isTuneMode ) {
        std::ofstream out(tunefile, std::ios::out);
        
        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;
    }
482

483 484
    do {
        
frey_m's avatar
frey_m committed
485 486
        error = std::numeric_limits<value_type>::max();
        
487
        // (re-)set inital values for r and pr
488
        r_m[0] = init[0];
489
        pr_m[0] = init[1];
490 491
        vz_m[0] = init[2];
        vpz_m[0] = init[3];
frey_m's avatar
frey_m committed
492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517
        
        if ( !this->findOrbitOfEnergy_m(E, init, error, accuracy, maxit) ) {
            throw OpalException("ClosedOrbitFinder::findOrbitOfEnergy()",
                                "Didn't converge for energy " + std::to_string(E) + " MeV.");
        }
        
        if ( isTuneMode ) {
            this->computeOrbitProperties();
    
            std::pair<value_type , value_type > tunes = this->getTunes();
    
            value_type reo = this->getOrbit(cycl_m->getPHIinit())[0];
            value_type  peo = this->getMomentum(cycl_m->getPHIinit())[0];


            *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;
518
            
frey_m's avatar
frey_m committed
519 520 521 522 523 524 525 526 527
            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();
        }
        
adelmann's avatar
adelmann committed
528 529 530 531 532
        // increase energy by dE
        if (E_m <= E + dE)
            E = E_m;
        else
            E += dE;
frey_m's avatar
frey_m committed
533
    
534
    } while (E != E_m);
535

536 537 538 539 540
    /* 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()
541

542 543 544
    // 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();
545 546 547 548 549 550 551
    
    
    /* 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();
552

553

554 555 556
    // returns true if converged, otherwise false
    return error < accuracy;
}
557

frey_m's avatar
frey_m committed
558 559 560 561 562 563 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 590 591 592 593 594 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 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717
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)
    {
        r_m[idx] = y[0];
        pr_m[idx] = y[1];
        vz_m[idx] = y[6];
        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
    value_type en = E / 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)
    value_type gamma2 = (1.0 + en) * (1.0 + en);          // = gamma^2
    value_type invgamma4 = 1.0 / (gamma2 * gamma2);       // = 1/gamma^4
    
    
    // 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
        dydt[12] = y[0] * invptheta * gamma_m - 1;
        
    };

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

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

        // correct inital values of r and pr
        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);
    
    return (error < accuracy);
}

718
template<typename Value_type, typename Size_type, class Stepper>
719 720 721
Value_type ClosedOrbitFinder<Value_type, Size_type, Stepper>::computeTune(const std::array<value_type,2>& y,
                                                                          value_type py2, size_type ncross)
{
722
    // Y = [y1, y2; py1, py2]
723

724 725
    // cos(mu)
    value_type cos = 0.5 * (y[0] + py2);
726
    
727
    value_type mu;
728

729 730
    // sign of sin(mu) has to be equal to y2
    bool negative = std::signbit(y[1]);
731

732
    bool uneven = (ncross % 2);
733

734 735 736
    if (std::fabs(cos) > 1.0) {
        // store the number of crossings
        value_type n = ncross;
737

738 739
        if (uneven)
            n = ncross - 1;
740

741 742
        // Gordon, formula (36b)
        value_type muPrime = -std::acosh(std::fabs(cos));      // mu'
743
        mu = n * Physics::pi + muPrime;
744

745 746 747 748 749 750 751
    } 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.
        */
752

753
        // Gordon, formula (36)
754
        mu = ncross * Physics::pi + muPrime;
755

756 757
        // if sign(y[1]) > 0 && sign(sin(mu)) < 0
        if (!negative && std::signbit(std::sin(mu))) {
758
            mu = ncross * Physics::pi - muPrime;
759
        } else if (negative && !std::signbit(std::sin(mu))) {
760
            mu = ncross * Physics::pi - muPrime + Physics::two_pi;
761 762
        }
    }
763

764
    // nu = mu/theta, where theta = integration domain
765

766
    /* domain_m = true --> only integrated over a single sector --> multiply by cycl_m->getSymmetry() to
adelmann's avatar
adelmann committed
767 768 769
     * get correct tune.
     */
    if (domain_m)
770
        mu *= cycl_m->getSymmetry();
771

772
    return mu * Physics::u_two_pi;
773 774 775
}

template<typename Value_type, typename Size_type, class Stepper>
776
void ClosedOrbitFinder<Value_type, Size_type, Stepper>::computeOrbitProperties() {
777
    /*
778 779 780 781 782
     * 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
     */
783

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

787 788 789
    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)
790 791 792
    value_type p2 = p * p;
    value_type theta = 0.0;                                             // angle for interpolating
    value_type ptheta;
793

794 795 796 797 798 799 800
    // 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
801
        cycl_m->apply(r_m[i], vz_m[i], theta, brint, btint, bint);
802 803 804
        bint *= invbcon;
        brint *= invbcon;
        btint *= invbcon;
frey_m's avatar
frey_m committed
805
        
806 807 808 809 810 811 812 813 814 815 816 817
        // 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;
818
    }
819 820 821

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

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