OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
ClosedOrbitFinder.h
Go to the documentation of this file.
1//
2// Class ClosedOrbitFinder
3// This class finds a closed orbit of a cyclotron for a given energy.
4// The algorithm is based on the paper of M. M. Gordon: "Computation of
5// closed orbits and basic focusing properties for sector-focused cyclotrons
6// and the design of 'cyclops'" (1983)
7// As template arguments one chooses the type of the variables and the
8// integrator for the ODEs. The supported steppers can be found on
9// http://www.boost.org/ where it is part of the library Odeint.
10//
11// Copyright (c) 2014, Matthias Frey, ETH Zürich
12// All rights reserved
13//
14// Implemented as part of the Semester thesis by Matthias Frey
15// "Matched Distributions in Cyclotrons"
16//
17// This file is part of OPAL.
18//
19// OPAL is free software: you can redistribute it and/or modify
20// it under the terms of the GNU General Public License as published by
21// the Free Software Foundation, either version 3 of the License, or
22// (at your option) any later version.
23//
24// You should have received a copy of the GNU General Public License
25// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
26//
27#ifndef CLOSEDORBITFINDER_H
28#define CLOSEDORBITFINDER_H
29
30#include <algorithm>
31#include <array>
32#include <cmath>
33#include <filesystem>
34#include <functional>
35#include <limits>
36#include <numeric>
37#include <string>
38#include <utility>
39#include <vector>
40
44#include "Utilities/Options.h"
45#include "Physics/Physics.h"
46#include "Physics/Units.h"
47
48// include headers for integration
49#include <boost/numeric/odeint/integrate/integrate_n_steps.hpp>
50
51extern Inform *gmsg;
52
53template<typename Value_type, typename Size_type, class Stepper>
55
56public:
58 typedef Value_type value_type;
60 typedef Size_type size_type;
62 typedef std::vector<value_type> container_type;
64 typedef std::vector<value_type> state_type;
65
66 typedef std::function<void(const state_type&, state_type&, const double)> function_t;
67
69
80 Cyclotron* cycl, bool domain = true, int Nsectors = 1);
81
84
86 container_type getPathLength(const value_type& angle = 0);
87
89 container_type getFieldIndex(const value_type& angle = 0);
90
92 std::pair<value_type,value_type> getTunes();
93
98
102
107
112
115
118
121
123
131 bool findOrbit(value_type accuracy, size_type maxit,
132 value_type ekin,
133 value_type dE = 1.0, value_type rguess = -1.0,
134 bool isTuneMode = false);
135
137 void computeOrbitProperties(const value_type& E);
138
139private:
141
146 value_type computeTune(const std::array<value_type,2>&, value_type, size_type);
147
148 // Compute closed orbit for given energy
150 const value_type&, size_type);
151
153// void computeVerticalOscillations();
154
156 container_type rotate(value_type angle, const container_type& orbitProperty);
157
159 std::array<value_type,2> x_m; // x_m = [x1, x2]
161 std::array<value_type,2> px_m; // px_m = [px1, px2]
163 std::array<value_type,2> z_m; // z_m = [z1, z2]
165 std::array<value_type,2> pz_m; // pz_m = [pz1, pz2]
166
181
185 size_type nzcross_m; //#crossings of zero-line in x- and z-direction
186
189
190 // Is the particle charge [e]
192
194 /* (see paper of Dr. C. Baumgarten: "Transverse-Longitudinal
195 * Coupling by Space Charge in Cyclotrons" (2012), formula (1))
196 */
202
205
208
214 /* value_type lastOrbitVal_m; */
215
221 /* value_type lastMomentumVal_m; */
222
233
236
242 std::function<double(double)> acon_m = [](double wo) { return Physics::c / wo; };
243
245
248 std::function<double(double, double)> bcon_m = [this](double e0, double wo) {
249 return e0 * 1.0e7 / (q_m * Physics::c * Physics::c / wo);
250 };
251
253};
254
255// -----------------------------------------------------------------------------------------------------------------------
256// PUBLIC MEMBER FUNCTIONS
257// -----------------------------------------------------------------------------------------------------------------------
258
259template<typename Value_type, typename Size_type, class Stepper>
260ClosedOrbitFinder<Value_type,
261 Size_type,
263 value_type q,
264 size_type N, Cyclotron* cycl,
265 bool domain, int nSectors)
266 : nxcross_m(0)
267 , nzcross_m(0)
268 , E0_m(E0)
269 , q_m(q)
270 , wo_m(cycl->getRfFrequ()[0]*Units::MHz2Hz/cycl->getCyclHarm()*2.0*Physics::pi)
271 , N_m(N)
272 , dtheta_m(Physics::two_pi/value_type(N))
273 , ravg_m(0)
274 , phase_m(0)
275 /* , lastOrbitVal_m(0.0) */
276 /* , lastMomentumVal_m(0.0) */
277 , domain_m(domain)
278 , nSectors_m(nSectors)
279 , stepper_m()
280 , cycl_m(cycl)
281{
282
283 if ( cycl_m->getFMLowE() > cycl_m->getFMHighE() )
284 throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()",
285 "Incorrect cyclotron energy (MeV) bounds: Maximum cyclotron energy smaller than minimum cyclotron energy.");
286
287 // if domain_m = true --> integrate over a single sector
288 if (domain_m) {
289 N_m /= cycl_m->getSymmetry();
290 }
291
292 cycl_m->read(cycl_m->getBScale());
293
294 // reserve storage for the orbit and momentum (--> size = 0, capacity = N_m+1)
295 /*
296 * we need N+1 storage, since dtheta = 2pi/N (and not 2pi/(N-1)) that's why we need N+1 integration steps
297 * to return to the origin (but the return size is N_m)
298 */
299 r_m.reserve(N_m + 1);
300 pr_m.reserve(N_m + 1);
301 vz_m.reserve(N_m + 1);
302 vpz_m.reserve(N_m + 1);
303
304 // reserve memory of N_m for the properties (--> size = 0, capacity = N_m)
305 h_m.reserve(N_m);
306 ds_m.reserve(N_m);
307 fidx_m.reserve(N_m);
308}
309
310template<typename Value_type, typename Size_type, class Stepper>
313{
314 if (angle != 0.0) {
315 return rotate(angle, h_m);
316 } else {
317 return h_m;
318 }
319}
320
321template<typename Value_type, typename Size_type, class Stepper>
324{
325 if (angle != 0.0) {
326 return rotate(angle, ds_m);
327 } else {
328 return ds_m;
329 }
330}
331
332template<typename Value_type, typename Size_type, class Stepper>
335{
336 if (angle != 0.0) {
337 return rotate(angle, fidx_m);
338 }
339 return fidx_m;
340}
341
342template<typename Value_type, typename Size_type, class Stepper>
344 // compute radial tune
346
347 // compute vertical tune
349
350 return std::make_pair(nur,nuz);
351}
352
353template<typename Value_type, typename Size_type, class Stepper>
356{
357 if (angle != 0.0) {
358 return rotate(angle, r_m);
359 } else {
360 return r_m;
361 }
362}
363
364template<typename Value_type, typename Size_type, class Stepper>
367{
368 container_type pr = pr_m;
369
370 if (angle != 0.0) {
371 pr = rotate(angle, pr);
372 }
373 // change units from meters to \beta * \gamma
374 /* in Gordon paper:
375 *
376 * p = \gamma * \beta * a
377 *
378 * where a = c / \omega_{0} with \omega_{0} = 2 * \pi * \nu_{0} = 2 * \pi * \nu_{rf} / h
379 *
380 * c: speed of light
381 * h: harmonic number
382 * v_{rf}: nomial rf frequency
383 *
384 * Units:
385 *
386 * [a] = m --> [p] = m
387 *
388 * The momentum in \beta * \gamma is obtained by dividing by "a"
389 */
390 value_type factor = 1.0 / acon_m(wo_m);
391 std::for_each(pr.begin(), pr.end(), [factor](value_type& p) { p *= factor; });
392
393 return pr;
394}
395
396template<typename Value_type, typename Size_type, class Stepper>
402
403template<typename Value_type, typename Size_type, class Stepper>
409
410// -----------------------------------------------------------------------------------------------------------------------
411// PRIVATE MEMBER FUNCTIONS
412// -----------------------------------------------------------------------------------------------------------------------
413
414template<typename Value_type, typename Size_type, class Stepper>
416 size_type maxit,
417 value_type ekin,
418 value_type dE,
419 value_type rguess,
420 bool isTuneMode)
421{
422 /* REMARK TO GORDON
423 * q' = 1/b = 1/bcon
424 * a' = a = acon
425 */
426
427 // resize vectors (--> size = N_m+1, capacity = N_m+1), note: we do N_m+1 integration steps
428 r_m.resize(N_m+1);
429 pr_m.resize(N_m+1);
430 vz_m.resize(N_m+1);
431 vpz_m.resize(N_m+1);
432
433 // store acon locally
434 value_type acon = acon_m(wo_m); // [acon] = m
435 // amplitude of error; Gordon, formula (18) (a = a')
436 value_type error = std::numeric_limits<value_type>::max();
437
438 /*
439 * Christian:
440 * N = 1440 ---> N = 720 ---> dtheta = 2PI/720 --> nsteps = 721
441 *
442 * 0, 2, 4, ... ---> jeden zweiten berechnene: 1, 3, 5, ... interpolieren --> 1440 Werte
443 *
444 * Matthias:
445 * N = 1440 --> dtheta = 2PI/1440 --> nsteps = 1441
446 *
447 * 0, 1, 2, 3, 4, 5, ... --> 1440 Werte
448 *
449 */
450
451 value_type E = cycl_m->getFMLowE(); // starting energy
452 value_type E_fin = ekin; // final energy
453 const value_type eps = dE * 1.0e-1; // articial constant for stopping criteria
454
455 if (isTuneMode) {
456 E_fin = cycl_m->getFMHighE();
457 }
458
459 namespace fs = std::filesystem;
460 fs::path dir = OpalData::getInstance()->getInputBasename();
461 dir = dir.parent_path() / OpalData::getInstance()->getAuxiliaryOutputDirectory();
462 std::string tunefile = (dir / "tunes.dat").string();
463
464 if ( isTuneMode ) {
465 std::ofstream out(tunefile, std::ios::out);
466
467 out << std::left
468 << std::setw(15) << "# energy[MeV]"
469 << std::setw(15) << "radius_ini[m]"
470 << std::setw(15) << "momentum_ini[Beta Gamma]"
471 << std::setw(15) << "radius_avg[m]"
472 << std::setw(15) << "nu_r"
473 << std::setw(15) << "nu_z"
474 << std::endl;
475 }
476 // initial guess
477 container_type init;
478 enum Guess {NONE, FIRST, SECOND};
479 Guess guess = NONE;
480 value_type rn1 = 0.0, pn1 = 0.0; // normalised r, pr value of previous closed orbit
481 value_type rn2 = 0.0, pn2 = 0.0; // normalised r, pr value of second to previous closed orbit
482
483 // iterate until suggested energy (start with minimum energy)
484 // increase energy by dE
485 *gmsg << level3 << "Start iteration to find closed orbit of energy "
486 << E_fin << " MeV " << "in steps of " << dE << " MeV." << endl;
487
488 for (; E <= E_fin + eps; E += dE) {
489
490 error = std::numeric_limits<value_type>::max();
491
492 // energy dependent values
493 value_type en = E / E0_m; // en = E/E0 = E/(mc^2) (E0 is potential energy)
494 value_type gamma = en + 1.0;
495 value_type gamma2 = gamma * gamma; // = gamma^2
496 value_type beta = std::sqrt(1.0 - 1.0 / gamma2);
497 value_type p = acon_m(wo_m) * std::sqrt(en * (2.0 + en)); // momentum [p] = m; Gordon, formula (3)
498
499 if (guess == NONE) {
500 // set initial values for radius and radial momentum for lowest energy Emin
501 // orbit, [r] = m; Gordon, formula (20)
502 // radial momentum; Gordon, formula (20)
503 // r pr z pz
504 init = {beta * acon, 0.0, 0.0, 1.0};
505 if (rguess >= 0.0) {
506 init[0] = rguess;
507 }
508 guess = FIRST;
509 } else if (guess == FIRST) {
510 // new initial values based on previous one, formula (21)
511 init[0] = (beta*acon) * rn1;
512 init[1] = p*pn1;
513 guess = SECOND;
514 } else if (guess == SECOND) {
515 // second extrapolation, formula (21)
516 init[0] = (beta*acon) * (rn1 + (rn1-rn2));
517 init[1] = p*(pn1 + (pn1-pn2));
518 }
519
520 std::fill( r_m.begin(), r_m.end(), 0);
521 std::fill( pr_m.begin(), pr_m.end(), 0);
522 std::fill( vz_m.begin(), vz_m.end(), 0);
523 std::fill(vpz_m.begin(), vpz_m.end(), 0);
524
525 // (re-)set inital values for r and pr
526 r_m[0] = init[0];
527 pr_m[0] = init[1];
528 vz_m[0] = init[2];
529 vpz_m[0] = init[3];
530
531 *gmsg << level3 << " Try to find orbit for " << E << " MeV ... ";
532
533 if ( !this->findOrbitOfEnergy_m(E, init, error, accuracy, maxit) ) {
534 *gmsg << endl << "ClosedOrbitFinder didn't converge for energy " + std::to_string(E) + " MeV." << endl;
535 guess = NONE;
536 continue;
537 }
538
539 *gmsg << level3 << "Successfully found." << endl;
540
541 // store for next initial guess
542 rn2 = rn1;
543 pn2 = pn1;
544 rn1 = r_m[0] / (acon * beta);
545 pn1 = pr_m[0] / p;
546
547 if ( isTuneMode ) {
548 this->computeOrbitProperties(E);
549 value_type reo = this->getOrbit( cycl_m->getPHIinit())[0];
550 value_type peo = this->getMomentum(cycl_m->getPHIinit())[0];
551 std::pair<value_type, value_type> tunes = this->getTunes();
552
553 *gmsg << std::left
554 << "* ----------------------------" << endl
555 << "* Closed orbit info (Gordon units):" << endl
556 << "*" << endl
557 << "* kinetic energy: " << std::setw(12) << E << " [MeV]" << endl
558 << "* average radius: " << std::setw(12) << ravg_m << " [m]" << endl
559 << "* initial radius: " << std::setw(12) << reo << " [m]" << endl
560 << "* initial momentum: " << std::setw(12) << peo << " [Beta Gamma]" << endl
561 << "* frequency error: " << phase_m << endl
562 << "* horizontal tune: " << tunes.first << endl
563 << "* vertical tune: " << tunes.second << endl
564 << "* ----------------------------" << endl << endl;
565
566 std::ofstream out(tunefile, std::ios::app);
567 out << std::left
568 << std::setw(15) << E
569 << std::setw(15) << reo
570 << std::setw(15) << peo
571 << std::setw(15) << ravg_m
572 << std::setw(15) << tunes.first
573 << std::setw(15) << tunes.second << std::endl;
574 out.close();
575 }
576 }
577
578 bool isConvergent = error < accuracy;
579 if (isConvergent) {
580 *gmsg << level3 << "Finished closed orbit finder successfully." << endl;
581 }
582
583 /* store last entry, since it is needed in computeVerticalOscillations(), because we have to do the same
584 * number of integrations steps there.
585 */
586 /* lastOrbitVal_m = r_m[N_m]; // needed in computeVerticalOscillations() */
587 /* lastMomentumVal_m = pr_m[N_m]; // needed in computeVerticalOscillations() */
588
589 // remove last entry (since we don't have to store [0,2pi], but [0,2pi[) --> size = N_m, capacity = N_m+1
590 r_m.pop_back();
591 pr_m.pop_back();
592
593 /* domain_m = true --> only integrated over a single sector
594 * --> multiply by cycl_m->getSymmetry() to get correct phase_m
595 */
596 if (domain_m) {
597 phase_m *= cycl_m->getSymmetry();
598 }
599
600 // returns true if converged, otherwise false
601 return isConvergent;
602}
603
604template<typename Value_type, typename Size_type, class Stepper>
606 const value_type& E,
607 container_type& init,
608 value_type& error,
609 const value_type& accuracy,
610 size_type maxit)
611{
612 /* *gmsg << "rguess : " << init[0] << endl; */
613 /* *gmsg << "prguess: " << init[1] << endl; */
614
615 value_type bint, brint, btint;
616 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)
617
618 value_type xold = 0.0; // for counting nxcross
619 value_type zold = 0.0; // for counting nzcross
620
621 // index for reaching next element of the arrays r and pr (no nicer way found yet)
622 size_type idx = 0;
623 // observer for storing the current value after each ODE step (e.g. Runge-Kutta step) into the containers of r and pr
624 auto store = [&](state_type& y, const value_type /*t*/) {
625 r_m[idx] = y[0];
626 pr_m[idx] = y[1];
627 vz_m[idx] = y[6];
628 vpz_m[idx] = y[7];
629
630 // count number of crossings (excluding starting point --> idx>0)
631 nxcross_m += (idx > 0) * (y[4] * xold < 0);
632 xold = y[4];
633
634 // number of times z2 changes sign
635 nzcross_m += (idx > 0) * (y[10] * zold < 0);
636 zold = y[10];
637
638 ++idx;
639 };
640
641 // define initial state container for integration: y = {r, pr, x1, px1, x2, px2,
642 // z, pz, z1, pz1, z2, pz2,
643 // phase}
644 state_type y(11);
645
646 // difference of last and first value of r (1. element) and pr (2. element)
647 container_type err(2);
648 // correction term for initial values: r = r + dr, pr = pr + dpr; Gordon, formula (17)
649 container_type delta = {0.0, 0.0};
650 // if niterations > maxit --> stop iteration
651 size_type niterations = 0;
652
653 // energy dependent values
654 value_type en = E / E0_m; // en = E/E0 = E/(mc^2) (E0 is potential energy)
655 value_type gamma = en + 1.0;
656 value_type p = acon_m(wo_m) * std::sqrt(en * (2.0 + en)); // momentum [p] = m; Gordon, formula (3)
657 value_type gamma2 = gamma * gamma; // = gamma^2
658 value_type invgamma4 = 1.0 / (gamma2 * gamma2); // = 1/gamma^4
659 value_type p2 = p * p; // p^2 = p*p
660
661 // helper constants
662 value_type pr2; // squared radial momentum (pr^2 = pr*pr)
663 value_type ptheta, invptheta; // azimuthal momentum
664
665 // define the six ODEs (using lambda function)
666 function_t orbit_integration = [&](const state_type &y,
667 state_type &dydt,
668 const double theta)
669 {
670 pr2 = y[1] * y[1];
671 if (p2 < pr2) {
672 //*gmsg << theta << " " << p2 << " " << pr2 << endl;
673 throw OpalException("ClosedOrbitFinder::findOrbitOfEnergy_m()",
674 "p_{r}^2 > p^{2} (defined in Gordon paper) --> Square root of negative number.");
675 }
676 // Gordon, formula (5c)
677 ptheta = std::sqrt(p2 - pr2);
678 invptheta = 1.0 / ptheta;
679 // average field over the number of sectors
680 brint=0.0, btint=0.0, bint=0.0;
681 for (int i = 0; i<nSectors_m; i++) {
682 double angle = theta + i * Physics::two_pi / nSectors_m;
683 double tmpbr, tmpbt, tmpb;
684 // interpolate values of magnetic field
685 cycl_m->apply(y[0], y[6], angle, tmpbr, tmpbt, tmpb);
686 brint += tmpbr;
687 btint += tmpbt;
688 bint += tmpb;
689 }
690 brint /= nSectors_m;
691 btint /= nSectors_m;
692 bint /= nSectors_m;
693 // multiply by 1 / b
694 bint *= invbcon;
695 brint *= invbcon;
696 btint *= invbcon;
697
698 // Gordon, formula (5a)
699 dydt[0] = y[0] * y[1] * invptheta;
700 // Gordon, formula (5b) (typo in paper! second equal sign is a minus)
701 dydt[1] = ptheta - y[0] * bint;
702 // Gordon, formulas (9a) and (9b)
703 for (size_type i = 2; i < 5; i += 2) {
704 dydt[i] = (y[1] * y[i] + y[0] * p2 * y[i+1] * invptheta * invptheta) * invptheta;
705 dydt[i+1] = - y[1] * y[i+1] * invptheta - (bint + y[0] * brint) * y[i];
706 }
707
708 // Gordon, formulas (22a) and (22b)
709 for (size_type i = 6; i < 12; i += 2) {
710 dydt[i] = y[0] * y[i+1] * invptheta;
711 dydt[i+1] = (y[0] * brint - y[1] * invptheta * btint) * y[i];
712 }
713
714 // integrate phase
715 dydt[12] = y[0] * invptheta * gamma - 1;
716
717 };
718
719 // integrate until error smaller than user-defined accuracy
720 do {
721 // (re-)set initial values
722 x_m[0] = 1.0; // x1; Gordon, formula (10)
723 px_m[0] = 0.0; // px1; Gordon, formula (10)
724 x_m[1] = 0.0; // x2; Gordon, formula (10)
725 px_m[1] = 1.0; // px2; Gordon, formula (10)
726 z_m[0] = 1.0;
727 pz_m[0] = 0.0;
728 z_m[1] = 0.0;
729 pz_m[1] = 1.0;
730 phase_m = 0.0;
731 nxcross_m = 0; // counts the number of crossings of x-axis (excluding first step)
732 nzcross_m = 0;
733 idx = 0; // index for looping over r and pr arrays
734
735 // fill container with initial states
736 y = {init[0],init[1],
737 x_m[0], px_m[0], x_m[1], px_m[1],
738 init[2], init[3],
739 z_m[0], pz_m[0], z_m[1], pz_m[1],
740 phase_m
741 };
742
743 try {
744 // integrate from 0 to 2*pi / nSectors (one has to get back to the "origin")
745 boost::numeric::odeint::integrate_n_steps(stepper_m, orbit_integration,y,0.0,dtheta_m,N_m,store);
746 } catch(OpalException & ex) {
747 *gmsg << ex.where() << " " << ex.what() << endl;
748 break;
749 }
750
751 // write new state
752 x_m[0] = y[2];
753 px_m[0] = y[3];
754 x_m[1] = y[4];
755 px_m[1] = y[5];
756
757 z_m[0] = y[8];
758 pz_m[0] = y[9];
759 z_m[1] = y[10];
760 pz_m[1] = y[11];
761 phase_m = y[12] * Physics::u_two_pi; // / (2.0 * Physics::pi);
762
763 // compute error (compare values of orbit and momentum for theta = 0 and theta = 2*pi)
764 // (Note: size = N_m+1 --> last entry is N_m)
765 err[0] = r_m[N_m] - r_m[0]; // Gordon, formula (14)
766 err[1] = pr_m[N_m] - pr_m[0]; // Gordon, formula (14)
767
768 // correct initial values of r and pr
769 value_type invdenom = 1.0 / (x_m[0] + px_m[1] - 2.0);
770 delta[0] = ((px_m[1] - 1.0) * err[0] - x_m[1] * err[1]) * invdenom; // dr; Gordon, formula (16a)
771 delta[1] = (( x_m[0] - 1.0) * err[1] - px_m[0] * err[0]) * invdenom; // dpr; Gordon, formula (16b)
772
773 // improved initial values; Gordon, formula (17) (here it's used for higher energies)
774 init[0] += delta[0];
775 init[1] += delta[1];
776
777 // compute amplitude of the error (Gordon, formula (18)
778 error = std::sqrt(delta[0] * delta[0] + delta[1] * delta[1] * invgamma4) / r_m[0];
779
780 // *gmsg << "iteration " << niterations << " error: " << error << endl;
781
782 } while ((error > accuracy) && (niterations++ < maxit));
783
784 if (error > accuracy) {
785 *gmsg << "findOrbit not converged after " << niterations
786 << " iterations with error: " << error
787 << ". Needed accuracy " << accuracy << endl;
788 }
789 return (error < accuracy);
790}
791
792template<typename Value_type, typename Size_type, class Stepper>
793Value_type ClosedOrbitFinder<Value_type, Size_type, Stepper>::computeTune(const std::array<value_type,2>& y,
794 value_type py2, size_type ncross)
795{
796 // Y = [y1, y2; py1, py2]
797
798 // cos(mu)
799 value_type cos = 0.5 * (y[0] + py2);
800
801 value_type mu;
802
803 // sign of sin(mu) has to be equal to y2
804 bool negative = std::signbit(y[1]);
805
806 bool uneven = (ncross % 2);
807
808 value_type abscos = std::abs(cos);
809 value_type muPrime;
810 if (abscos > 1.0) {
811 // store the number of crossings
812 if (uneven) {
813 ncross = ncross + 1;
814 }
815 // Gordon, formula (36b)
816 muPrime = -std::acosh(abscos); // mu'
817
818 } else {
819 muPrime = (uneven) ? std::acos(-cos) : std::acos(cos); // mu'
820 /* It has to be fulfilled: 0<= mu' <= pi
821 * But since |cos(mu)| <= 1, we have
822 * -1 <= cos(mu) <= 1 --> 0 <= mu <= pi (using above programmed line), such
823 * that condition is already fulfilled.
824 */
825 }
826
827 // Gordon, formula (36)
828 mu = ncross * Physics::pi + muPrime;
829
830 if (abscos < 1.0) {
831 // if sign(y[1]) > 0 && sign(sin(mu)) < 0
832 if (!negative && std::signbit(std::sin(mu))) {
833 mu = ncross * Physics::pi - muPrime;
834 } else if (negative && !std::signbit(std::sin(mu))) {
835 mu = ncross * Physics::pi - muPrime + Physics::two_pi;
836 }
837 }
838
839 // nu = mu/theta, where theta = integration domain
840
841 /* domain_m = true --> only integrated over a single sector --> multiply by cycl_m->getSymmetry() to
842 * get correct tune.
843 */
844 if (domain_m) {
845 mu *= cycl_m->getSymmetry();
846 }
847 return mu * Physics::u_two_pi;
848}
849
850template<typename Value_type, typename Size_type, class Stepper>
852 /*
853 * The formulas for h, fidx and ds are from the paper:
854 * "Tranverse-Longitudinal Coupling by Space Charge in Cyclotrons"
855 * written by Dr. Christian Baumgarten (2012, PSI)
856 * p. 6
857 */
858
859 // READ IN MAGNETIC FIELD: ONLY FOR STAND-ALONE PROGRAM
860 value_type bint, brint, btint; // B, dB/dr, dB/dtheta
861
862 value_type invbcon = 1.0 / bcon_m(E0_m, wo_m);
863 value_type en = E / E0_m; // en = E/E0 = E/(mc^2) (E0 is potential energy)
864 value_type p = acon_m(wo_m) * std::sqrt(en * (2.0 + en)); // momentum [p] = m; Gordon, formula (3)
865 value_type p2 = p * p;
866 value_type theta = 0.0; // angle for interpolating
867 value_type ptheta;
868
869 // resize of container (--> size = N_m, capacity = N_m)
870 h_m.resize(N_m);
871 fidx_m.resize(N_m);
872 ds_m.resize(N_m);
873
874 for (size_type i = 0; i < N_m; ++i) {
875 // interpolate magnetic field
876 cycl_m->apply(r_m[i], vz_m[i], theta, brint, btint, bint);
877 bint *= invbcon;
878 brint *= invbcon;
879 btint *= invbcon;
880
881 // inverse bending radius
882 h_m[i] = bint / p;
883
884 // local field index
885 if (p < pr_m[i])
886 throw OpalException("ClosedOrbitFinder::computeOrbitProperties()",
887 "p_{r}^2 > p^{2} " + std::to_string(p) + " " + std::to_string(pr_m[i]) + " (defined in Gordon paper) --> Square root of negative number.");
888
889 ptheta = std::sqrt(p2 - pr_m[i] * pr_m[i]);
890
891 fidx_m[i] = (brint * ptheta - btint * pr_m[i] / r_m[i]) / p2; //(bint*bint);
892
893 // path length element
894 ds_m[i] = std::hypot(r_m[i] * pr_m[i] / ptheta,r_m[i]) * dtheta_m; // C++11 function
895
896 // increase angle
897 theta += dtheta_m;
898 }
899
900 // compute average radius
901 ravg_m = std::accumulate(r_m.begin(),r_m.end(),0.0) / value_type(r_m.size());
902}
903
904template<typename Value_type, typename Size_type, class Stepper>
907
908 container_type orbitPropertyCopy = orbitProperty;
909
910 if (angle < 0.0) {
911 angle = Physics::two_pi + angle;
912 }
913
914 // compute starting point
915 unsigned int start = angle / dtheta_m;
916
917 start %= orbitProperty.size();
918
919 // copy end to start
920 std::copy(orbitProperty.begin() + start, orbitProperty.end(), orbitPropertyCopy.begin());
921
922 // copy start to end
923 std::copy_n(orbitProperty.begin(), start, orbitPropertyCopy.end() - start);
924
925 return orbitPropertyCopy;
926}
927
928#endif
Tps< T > cos(const Tps< T > &x)
Cosine.
Definition TpsMath.h:129
T * value_type(const SliceIterator< T > &)
Inform * gmsg
Definition Main.cpp:70
Inform * gmsg
Definition Main.cpp:70
const double pi
Definition fftpack.cpp:894
Inform & endl(Inform &inf)
Definition Inform.cpp:42
Inform & level3(Inform &inf)
Definition Inform.cpp:47
Definition Air.h:27
constexpr double u_two_pi
The value of.
Definition Physics.h:36
constexpr double two_pi
The value of.
Definition Physics.h:33
constexpr double c
The velocity of light in m/s.
Definition Physics.h:45
constexpr double pi
The value of.
Definition Physics.h:30
Definition Units.h:23
std::string getInputBasename()
get input file name without extension
Definition OpalData.cpp:674
static OpalData * getInstance()
Definition OpalData.cpp:196
std::string getAuxiliaryOutputDirectory() const
get the name of the the additional data directory
Definition OpalData.cpp:666
std::vector< value_type > container_type
Type of container for storing quantities (path length, orbit, etc.).
std::pair< value_type, value_type > getTunes()
Returns the radial and vertical tunes (in that order).
std::array< value_type, 2 > px_m
Stores current momenta in horizontal direction for the solutions of the ODE with different initial va...
size_type N_m
Number of integration steps.
ClosedOrbitFinder(value_type E0, value_type q, size_type N, Cyclotron *cycl, bool domain=true, int Nsectors=1)
Sets the initial values for the integration and calls findOrbit().
container_type vz_m
Stores the vertical oribt (size: N_m+1).
size_type nzcross_m
Counts the number of zero-line crossings in vertical direction (used for computing vertical tune).
value_type phase_m
Is the phase.
bool findOrbitOfEnergy_m(const value_type &, container_type &, value_type &, const value_type &, size_type)
std::vector< value_type > state_type
Type for holding state of ODE values.
void computeOrbitProperties(const value_type &E)
Fills in the values of h_m, ds_m, fidx_m.
std::array< value_type, 2 > z_m
Stores current position in vertical direction for the solutions of the ODE with different initial val...
container_type fidx_m
Stores the field index.
value_type wo_m
Is the nominal orbital frequency.
container_type ds_m
Stores the step length.
bool findOrbit(value_type accuracy, size_type maxit, value_type ekin, value_type dE=1.0, value_type rguess=-1.0, bool isTuneMode=false)
Computes the closed orbit.
container_type h_m
Stores the inverse bending radius.
std::function< double(double, double)> bcon_m
Cyclotron unit (Tesla).
container_type r_m
Stores the radial orbit (size: N_m+1).
Size_type size_type
Type for specifying sizes.
container_type vpz_m
Stores the vertical momentum.
value_type getAverageRadius()
Returns the average orbit radius.
std::function< double(double)> acon_m
value_type dtheta_m
Is the angle step size.
std::function< void(const state_type &, state_type &, const double)> function_t
container_type getFieldIndex(const value_type &angle=0)
Returns the field index (size of container N+1).
Value_type value_type
Type of variables.
container_type getMomentum(value_type angle=0)
container_type getOrbit(value_type angle=0)
container_type rotate(value_type angle, const container_type &orbitProperty)
This function computes nzcross_ which is used to compute the tune in z-direction and the frequency er...
bool isConverged()
Returns true if a closed orbit could be found.
container_type getPathLength(const value_type &angle=0)
Returns the step lengths of the path (size of container N+1).
std::array< value_type, 2 > pz_m
Stores current momenta in vertical direction for the solutions of the ODE with different initial valu...
std::array< value_type, 2 > x_m
Stores current position in horizontal direction for the solutions of the ODE with different initial v...
value_type ravg_m
Is the average radius.
container_type pr_m
Stores the radial momentum.
container_type getInverseBendingRadius(const value_type &angle=0)
Returns the inverse bending radius (size of container N+1).
value_type getFrequencyError()
Returns the frequency error.
value_type computeTune(const std::array< value_type, 2 > &, value_type, size_type)
This function is called by the function getTunes().
value_type E0_m
Is the rest mass [MeV / c**2].
size_type nxcross_m
Counts the number of zero-line crossings in horizontal direction (used for computing horizontal tune)...
Stepper stepper_m
Defines the stepper for integration of the ODE's.
The base class for all OPAL exceptions.
virtual const std::string & what() const
Return the message string for the exception.
virtual const std::string & where() const
Return the name of the method or function which detected the exception.