OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
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 <functional>
34 #include <limits>
35 #include <numeric>
36 #include <string>
37 #include <utility>
38 #include <vector>
39 
40 #include "Utilities/Options.h"
42 #include "Physics/Physics.h"
43 
45 
46 #include "AbsBeamline/Cyclotron.h"
47 
48 // include headers for integration
49 #include <boost/numeric/odeint/integrate/integrate_n_steps.hpp>
50 #include <boost/filesystem.hpp>
51 
52 extern Inform *gmsg;
53 
54 template<typename Value_type, typename Size_type, class Stepper>
56 {
57  public:
59  typedef Value_type value_type;
61  typedef Size_type size_type;
63  typedef std::vector<value_type> container_type;
65  typedef std::vector<value_type> state_type;
66 
67  typedef std::function<void(const state_type&, state_type&, const double)> function_t;
68 
70 
81  Cyclotron* cycl, bool domain = true, int Nsectors = 1);
82 
85 
87  container_type getPathLength(const value_type& angle = 0);
88 
90  container_type getFieldIndex(const value_type& angle = 0);
91 
93  std::pair<value_type,value_type> getTunes();
94 
99 
103 
108 
113 
116 
119 
121  bool isConverged();
122 
124 
132  bool findOrbit(value_type accuracy, size_type maxit,
133  value_type ekin,
134  value_type dE = 1.0, value_type rguess = -1.0,
135  bool isTuneMode = false);
136 
138  void computeOrbitProperties(const value_type& E);
139 
140  private:
142 
147  value_type computeTune(const std::array<value_type,2>&, value_type, size_type);
148 
149  // Compute closed orbit for given energy
151  const value_type&, size_type);
152 
154 // void computeVerticalOscillations();
155 
157  container_type rotate(value_type angle, const container_type& orbitProperty);
158 
160  std::array<value_type,2> x_m; // x_m = [x1, x2]
162  std::array<value_type,2> px_m; // px_m = [px1, px2]
164  std::array<value_type,2> z_m; // z_m = [z1, z2]
166  std::array<value_type,2> pz_m; // pz_m = [pz1, pz2]
167 
182 
186  size_type nzcross_m; //#crossings of zero-line in x- and z-direction
187 
190 
191  // Is the particle charge [e]
193 
195  /* (see paper of Dr. C. Baumgarten: "Transverse-Longitudinal
196  * Coupling by Space Charge in Cyclotrons" (2012), formula (1))
197  */
203 
206 
209 
215  /* value_type lastOrbitVal_m; */
216 
222  /* value_type lastMomentumVal_m; */
223 
228  bool domain_m;
234 
237 
243  std::function<double(double)> acon_m = [](double wo) { return Physics::c / wo; };
244 
246 
249  std::function<double(double, double)> bcon_m = [this](double e0, double wo) {
250  return e0 * 1.0e7 / (q_m * Physics::c * Physics::c / wo);
251  };
252 
254 };
255 
256 // -----------------------------------------------------------------------------------------------------------------------
257 // PUBLIC MEMBER FUNCTIONS
258 // -----------------------------------------------------------------------------------------------------------------------
259 
260 template<typename Value_type, typename Size_type, class Stepper>
261 ClosedOrbitFinder<Value_type,
262  Size_type,
264  value_type q,
265  size_type N, Cyclotron* cycl,
266  bool domain, int nSectors)
267  : nxcross_m(0)
268  , nzcross_m(0)
269  , E0_m(E0)
270  , q_m(q)
271  , wo_m(cycl->getRfFrequ(0)*1E6/cycl->getCyclHarm()*2.0*Physics::pi)
272  , N_m(N)
273  , dtheta_m(Physics::two_pi/value_type(N))
274  , ravg_m(0)
275  , phase_m(0)
276  /* , lastOrbitVal_m(0.0) */
277  /* , lastMomentumVal_m(0.0) */
278  , domain_m(domain)
279  , nSectors_m(nSectors)
280  , stepper_m()
281  , cycl_m(cycl)
282 {
283 
284  if ( cycl_m->getFMLowE() > cycl_m->getFMHighE() )
285  throw OpalException("ClosedOrbitFinder::ClosedOrbitFinder()",
286  "Incorrect cyclotron energy (MeV) bounds: Maximum cyclotron energy smaller than minimum cyclotron energy.");
287 
288  // if domain_m = true --> integrate over a single sector
289  if (domain_m) {
290  N_m /= cycl_m->getSymmetry();
291  }
292 
294 
295  // reserve storage for the orbit and momentum (--> size = 0, capacity = N_m+1)
296  /*
297  * we need N+1 storage, since dtheta = 2pi/N (and not 2pi/(N-1)) that's why we need N+1 integration steps
298  * to return to the origin (but the return size is N_m)
299  */
300  r_m.reserve(N_m + 1);
301  pr_m.reserve(N_m + 1);
302  vz_m.reserve(N_m + 1);
303  vpz_m.reserve(N_m + 1);
304 
305  // reserve memory of N_m for the properties (--> size = 0, capacity = N_m)
306  h_m.reserve(N_m);
307  ds_m.reserve(N_m);
308  fidx_m.reserve(N_m);
309 }
310 
311 template<typename Value_type, typename Size_type, class Stepper>
314 {
315  if (angle != 0.0)
316  return rotate(angle, h_m);
317  else
318  return h_m;
319 }
320 
321 template<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 template<typename Value_type, typename Size_type, class Stepper>
334 {
335  if (angle != 0.0)
336  return rotate(angle, fidx_m);
337  return fidx_m;
338 }
339 
340 template<typename Value_type, typename Size_type, class Stepper>
342  // compute radial tune
343  value_type nur = computeTune(x_m,px_m[1],nxcross_m);
344 
345  // compute vertical tune
346  value_type nuz = computeTune(z_m,pz_m[1],nzcross_m);
347 
348  return std::make_pair(nur,nuz);
349 }
350 
351 template<typename Value_type, typename Size_type, class Stepper>
354 {
355  if (angle != 0.0)
356  return rotate(angle, r_m);
357  else
358  return r_m;
359 }
360 
361 template<typename Value_type, typename Size_type, class Stepper>
364 {
365  container_type pr = pr_m;
366 
367  if (angle != 0.0)
368  pr = rotate(angle, pr);
369 
370  // change units from meters to \beta * \gamma
371  /* in Gordon paper:
372  *
373  * p = \gamma * \beta * a
374  *
375  * where a = c / \omega_{0} with \omega_{0} = 2 * \pi * \nu_{0} = 2 * \pi * \nu_{rf} / h
376  *
377  * c: speed of light
378  * h: harmonic number
379  * v_{rf}: nomial rf frequency
380  *
381  * Units:
382  *
383  * [a] = m --> [p] = m
384  *
385  * The momentum in \beta * \gamma is obtained by dividing by "a"
386  */
387  value_type factor = 1.0 / acon_m(wo_m);
388  std::for_each(pr.begin(), pr.end(), [factor](value_type& p) { p *= factor; });
389 
390  return pr;
391 }
392 
393 template<typename Value_type, typename Size_type, class Stepper>
396 {
397  return ravg_m;
398 }
399 
400 template<typename Value_type, typename Size_type, class Stepper>
403 {
404  return phase_m;
405 }
406 
407 // -----------------------------------------------------------------------------------------------------------------------
408 // PRIVATE MEMBER FUNCTIONS
409 // -----------------------------------------------------------------------------------------------------------------------
410 
411 
412 
413 template<typename Value_type, typename Size_type, class Stepper>
415  size_type maxit,
416  value_type ekin,
417  value_type dE,
418  value_type rguess,
419  bool isTuneMode)
420 {
421  /* REMARK TO GORDON
422  * q' = 1/b = 1/bcon
423  * a' = a = acon
424  */
425 
426  // resize vectors (--> size = N_m+1, capacity = N_m+1), note: we do N_m+1 integration steps
427  r_m.resize(N_m+1);
428  pr_m.resize(N_m+1);
429  vz_m.resize(N_m+1);
430  vpz_m.resize(N_m+1);
431 
432  // store acon locally
433  value_type acon = acon_m(wo_m); // [acon] = m
434  // amplitude of error; Gordon, formula (18) (a = a')
436 
437  /*
438  * Christian:
439  * N = 1440 ---> N = 720 ---> dtheta = 2PI/720 --> nsteps = 721
440  *
441  * 0, 2, 4, ... ---> jeden zweiten berechnene: 1, 3, 5, ... interpolieren --> 1440 Werte
442  *
443  * Matthias:
444  * N = 1440 --> dtheta = 2PI/1440 --> nsteps = 1441
445  *
446  * 0, 1, 2, 3, 4, 5, ... --> 1440 Werte
447  *
448  */
449 
450  value_type E = cycl_m->getFMLowE(); // starting energy
451  value_type E_fin = ekin; // final energy
452  const value_type eps = dE * 1.0e-1; // articial constant for stopping criteria
453 
454  if (isTuneMode) {
455  E_fin = cycl_m->getFMHighE();
456  }
457 
458  namespace fs = boost::filesystem;
459  fs::path dir = OpalData::getInstance()->getInputBasename();
460  dir = dir.parent_path() / OpalData::getInstance()->getAuxiliaryOutputDirectory();
461  std::string tunefile = (dir / "tunes.dat").string();
462 
463  if ( isTuneMode ) {
464  std::ofstream out(tunefile, std::ios::out);
465 
466  out << std::left
467  << std::setw(15) << "# energy[MeV]"
468  << std::setw(15) << "radius_ini[m]"
469  << std::setw(15) << "momentum_ini[Beta Gamma]"
470  << std::setw(15) << "radius_avg[m]"
471  << std::setw(15) << "nu_r"
472  << std::setw(15) << "nu_z"
473  << std::endl;
474  }
475  // initial guess
476  container_type init;
477  enum Guess {NONE, FIRST, SECOND};
478  Guess guess = NONE;
479  value_type rn1 = 0.0, pn1 = 0.0; // normalised r, pr value of previous closed orbit
480  value_type rn2 = 0.0, pn2 = 0.0; // normalised r, pr value of second to previous closed orbit
481 
482  // iterate until suggested energy (start with minimum energy)
483  // increase energy by dE
484  *gmsg << level3 << "Start iteration to find closed orbit of energy " << E_fin << " MeV "
485  << "in steps of " << dE << " MeV." << endl;
486 
487  for (; E <= E_fin + eps; E += dE) {
488 
490 
491  // energy dependent values
492  value_type en = E / E0_m; // en = E/E0 = E/(mc^2) (E0 is potential energy)
493  value_type gamma = en + 1.0;
494  value_type gamma2 = gamma * gamma; // = gamma^2
495  value_type beta = std::sqrt(1.0 - 1.0 / gamma2);
496  value_type p = acon_m(wo_m) * std::sqrt(en * (2.0 + en)); // momentum [p] = m; Gordon, formula (3)
497 
498  if (guess == NONE) {
499  // set initial values for radius and radial momentum for lowest energy Emin
500  // orbit, [r] = m; Gordon, formula (20)
501  // radial momentum; Gordon, formula (20)
502  // r pr z pz
503  init = {beta * acon, 0.0, 0.0, 1.0};
504  if (rguess >= 0.0) {
505  init[0] = rguess * 0.001;
506  }
507  guess = FIRST;
508  } else if (guess == FIRST) {
509  // new initial values based on previous one, formula (21)
510  init[0] = (beta*acon) * rn1;
511  init[1] = p*pn1;
512  guess = SECOND;
513  } else if (guess == SECOND) {
514  // second extrapolation, formula (21)
515  init[0] = (beta*acon) * (rn1 + (rn1-rn2));
516  init[1] = p*(pn1 + (pn1-pn2));
517  }
518 
519  std::fill( r_m.begin(), r_m.end(), 0);
520  std::fill( pr_m.begin(), pr_m.end(), 0);
521  std::fill( vz_m.begin(), vz_m.end(), 0);
522  std::fill(vpz_m.begin(), vpz_m.end(), 0);
523 
524  // (re-)set inital values for r and pr
525  r_m[0] = init[0];
526  pr_m[0] = init[1];
527  vz_m[0] = init[2];
528  vpz_m[0] = init[3];
529 
530  *gmsg << level3 << " Try to find orbit for " << E << " MeV ... ";
531 
532  if ( !this->findOrbitOfEnergy_m(E, init, error, accuracy, maxit) ) {
533  *gmsg << endl << "ClosedOrbitFinder didn't converge for energy " + std::to_string(E) + " MeV." << endl;
534  guess = NONE;
535  continue;
536  }
537 
538  *gmsg << level3 << "Successfully found." << endl;
539 
540  // store for next initial guess
541  rn2 = rn1;
542  pn2 = pn1;
543  rn1 = r_m[0] / (acon * beta);
544  pn1 = pr_m[0] / p;
545 
546  if ( isTuneMode ) {
547 
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  *gmsg << level3 << "Finished closed orbit finder successfully." << endl;
579 
580  /* store last entry, since it is needed in computeVerticalOscillations(), because we have to do the same
581  * number of integrations steps there.
582  */
583  /* lastOrbitVal_m = r_m[N_m]; // needed in computeVerticalOscillations() */
584  /* lastMomentumVal_m = pr_m[N_m]; // needed in computeVerticalOscillations() */
585 
586  // remove last entry (since we don't have to store [0,2pi], but [0,2pi[) --> size = N_m, capacity = N_m+1
587  r_m.pop_back();
588  pr_m.pop_back();
589 
590  /* domain_m = true --> only integrated over a single sector
591  * --> multiply by cycl_m->getSymmetry() to get correct phase_m
592  */
593  if (domain_m)
594  phase_m *= cycl_m->getSymmetry();
595 
596  // returns true if converged, otherwise false
597  return error < accuracy;
598 }
599 
600 template<typename Value_type, typename Size_type, class Stepper>
602  const value_type& E,
603  container_type& init,
604  value_type& error,
605  const value_type& accuracy,
606  size_type maxit)
607 {
608  /* *gmsg << "rguess : " << init[0] << endl; */
609  /* *gmsg << "prguess: " << init[1] << endl; */
610 
611  value_type bint, brint, btint;
612  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)
613 
614  value_type xold = 0.0; // for counting nxcross
615  value_type zold = 0.0; // for counting nzcross
616 
617  // index for reaching next element of the arrays r and pr (no nicer way found yet)
618  size_type idx = 0;
619  // observer for storing the current value after each ODE step (e.g. Runge-Kutta step) into the containers of r and pr
620  auto store = [&](state_type& y, const value_type /*t*/)
621  {
622  r_m[idx] = y[0];
623  pr_m[idx] = y[1];
624  vz_m[idx] = y[6];
625  vpz_m[idx] = y[7];
626 
627  // count number of crossings (excluding starting point --> idx>0)
628  nxcross_m += (idx > 0) * (y[4] * xold < 0);
629  xold = y[4];
630 
631  // number of times z2 changes sign
632  nzcross_m += (idx > 0) * (y[10] * zold < 0);
633  zold = y[10];
634 
635  ++idx;
636  };
637 
638  // define initial state container for integration: y = {r, pr, x1, px1, x2, px2,
639  // z, pz, z1, pz1, z2, pz2,
640  // phase}
641  state_type y(11);
642 
643  // difference of last and first value of r (1. element) and pr (2. element)
644  container_type err(2);
645  // correction term for initial values: r = r + dr, pr = pr + dpr; Gordon, formula (17)
646  container_type delta = {0.0, 0.0};
647  // if niterations > maxit --> stop iteration
648  size_type niterations = 0;
649 
650  // energy dependent values
651  value_type en = E / E0_m; // en = E/E0 = E/(mc^2) (E0 is potential energy)
652  value_type gamma = en + 1.0;
653  value_type p = acon_m(wo_m) * std::sqrt(en * (2.0 + en)); // momentum [p] = m; Gordon, formula (3)
654  value_type gamma2 = gamma * gamma; // = gamma^2
655  value_type invgamma4 = 1.0 / (gamma2 * gamma2); // = 1/gamma^4
656  value_type p2 = p * p; // p^2 = p*p
657 
658  // helper constants
659  value_type pr2; // squared radial momentum (pr^2 = pr*pr)
660  value_type ptheta, invptheta; // azimuthal momentum
661 
662  // define the six ODEs (using lambda function)
663  function_t orbit_integration = [&](const state_type &y,
664  state_type &dydt,
665  const double theta)
666  {
667  pr2 = y[1] * y[1];
668  if (p2 < pr2) {
669  //*gmsg << theta << " " << p2 << " " << pr2 << endl;
670  throw OpalException("ClosedOrbitFinder::findOrbitOfEnergy_m()",
671  "p_{r}^2 > p^{2} (defined in Gordon paper) --> Square root of negative number.");
672  }
673  // Gordon, formula (5c)
674  ptheta = std::sqrt(p2 - pr2);
675  invptheta = 1.0 / ptheta;
676  // average field over the number of sectors
677  brint=0.0, btint=0.0, bint=0.0;
678  for (int i = 0; i<nSectors_m; i++) {
679  double angle = theta + i * Physics::two_pi / nSectors_m;
680  double tmpbr, tmpbt, tmpb;
681  // interpolate values of magnetic field
682  cycl_m->apply(y[0], y[6], angle, tmpbr, tmpbt, tmpb);
683  brint += tmpbr;
684  btint += tmpbt;
685  bint += tmpb;
686  }
687  brint /= nSectors_m;
688  btint /= nSectors_m;
689  bint /= nSectors_m;
690  // multiply by 1 / b
691  bint *= invbcon;
692  brint *= invbcon;
693  btint *= invbcon;
694 
695  // Gordon, formula (5a)
696  dydt[0] = y[0] * y[1] * invptheta;
697  // Gordon, formula (5b) (typo in paper! second equal sign is a minus)
698  dydt[1] = ptheta - y[0] * bint;
699  // Gordon, formulas (9a) and (9b)
700  for (size_type i = 2; i < 5; i += 2) {
701  dydt[i] = (y[1] * y[i] + y[0] * p2 * y[i+1] * invptheta * invptheta) * invptheta;
702  dydt[i+1] = - y[1] * y[i+1] * invptheta - (bint + y[0] * brint) * y[i];
703  }
704 
705  // Gordon, formulas (22a) and (22b)
706  for (size_type i = 6; i < 12; i += 2) {
707  dydt[i] = y[0] * y[i+1] * invptheta;
708  dydt[i+1] = (y[0] * brint - y[1] * invptheta * btint) * y[i];
709  }
710 
711  // integrate phase
712  dydt[12] = y[0] * invptheta * gamma - 1;
713 
714  };
715 
716  // integrate until error smaller than user-defined accuracy
717  do {
718  // (re-)set initial values
719  x_m[0] = 1.0; // x1; Gordon, formula (10)
720  px_m[0] = 0.0; // px1; Gordon, formula (10)
721  x_m[1] = 0.0; // x2; Gordon, formula (10)
722  px_m[1] = 1.0; // px2; Gordon, formula (10)
723  z_m[0] = 1.0;
724  pz_m[0] = 0.0;
725  z_m[1] = 0.0;
726  pz_m[1] = 1.0;
727  phase_m = 0.0;
728  nxcross_m = 0; // counts the number of crossings of x-axis (excluding first step)
729  nzcross_m = 0;
730  idx = 0; // index for looping over r and pr arrays
731 
732  // fill container with initial states
733  y = {init[0],init[1],
734  x_m[0], px_m[0], x_m[1], px_m[1],
735  init[2], init[3],
736  z_m[0], pz_m[0], z_m[1], pz_m[1],
737  phase_m
738  };
739 
740  try {
741  // integrate from 0 to 2*pi / nSectors (one has to get back to the "origin")
742  boost::numeric::odeint::integrate_n_steps(stepper_m, orbit_integration,y,0.0,dtheta_m,N_m,store);
743  } catch(OpalException & ex) {
744  *gmsg << ex.where() << " " << ex.what() << endl;
745  break;
746  }
747 
748  // write new state
749  x_m[0] = y[2];
750  px_m[0] = y[3];
751  x_m[1] = y[4];
752  px_m[1] = y[5];
753 
754  z_m[0] = y[8];
755  pz_m[0] = y[9];
756  z_m[1] = y[10];
757  pz_m[1] = y[11];
758  phase_m = y[12] * Physics::u_two_pi; // / (2.0 * Physics::pi);
759 
760  // compute error (compare values of orbit and momentum for theta = 0 and theta = 2*pi)
761  // (Note: size = N_m+1 --> last entry is N_m)
762  err[0] = r_m[N_m] - r_m[0]; // Gordon, formula (14)
763  err[1] = pr_m[N_m] - pr_m[0]; // Gordon, formula (14)
764 
765  // correct initial values of r and pr
766  value_type invdenom = 1.0 / (x_m[0] + px_m[1] - 2.0);
767  delta[0] = ((px_m[1] - 1.0) * err[0] - x_m[1] * err[1]) * invdenom; // dr; Gordon, formula (16a)
768  delta[1] = (( x_m[0] - 1.0) * err[1] - px_m[0] * err[0]) * invdenom; // dpr; Gordon, formula (16b)
769 
770  // improved initial values; Gordon, formula (17) (here it's used for higher energies)
771  init[0] += delta[0];
772  init[1] += delta[1];
773 
774  // compute amplitude of the error (Gordon, formula (18)
775  error = std::sqrt(delta[0] * delta[0] + delta[1] * delta[1] * invgamma4) / r_m[0];
776 
777  // *gmsg << "iteration " << niterations << " error: " << error << endl;
778 
779  } while ((error > accuracy) && (niterations++ < maxit));
780 
781  if (error > accuracy)
782  *gmsg << "findOrbit not converged after " << niterations << " iterations with error: " << error << ". Needed accuracy " << accuracy << endl;
783 
784  return (error < accuracy);
785 }
786 
787 template<typename Value_type, typename Size_type, class Stepper>
788 Value_type ClosedOrbitFinder<Value_type, Size_type, Stepper>::computeTune(const std::array<value_type,2>& y,
789  value_type py2, size_type ncross)
790 {
791  // Y = [y1, y2; py1, py2]
792 
793  // cos(mu)
794  value_type cos = 0.5 * (y[0] + py2);
795 
796  value_type mu;
797 
798  // sign of sin(mu) has to be equal to y2
799  bool negative = std::signbit(y[1]);
800 
801  bool uneven = (ncross % 2);
802 
803  value_type abscos = std::abs(cos);
804  value_type muPrime;
805  if (abscos > 1.0) {
806  // store the number of crossings
807  if (uneven)
808  ncross = ncross + 1;
809 
810  // Gordon, formula (36b)
811  muPrime = -std::acosh(abscos); // mu'
812 
813  } else {
814  muPrime = (uneven) ? std::acos(-cos) : std::acos(cos); // mu'
815  /* It has to be fulfilled: 0<= mu' <= pi
816  * But since |cos(mu)| <= 1, we have
817  * -1 <= cos(mu) <= 1 --> 0 <= mu <= pi (using above programmed line), such
818  * that condition is already fulfilled.
819  */
820  }
821 
822  // Gordon, formula (36)
823  mu = ncross * Physics::pi + muPrime;
824 
825  if (abscos < 1.0) {
826  // if sign(y[1]) > 0 && sign(sin(mu)) < 0
827  if (!negative && std::signbit(std::sin(mu))) {
828  mu = ncross * Physics::pi - muPrime;
829  } else if (negative && !std::signbit(std::sin(mu))) {
830  mu = ncross * Physics::pi - muPrime + Physics::two_pi;
831  }
832  }
833 
834  // nu = mu/theta, where theta = integration domain
835 
836  /* domain_m = true --> only integrated over a single sector --> multiply by cycl_m->getSymmetry() to
837  * get correct tune.
838  */
839  if (domain_m)
840  mu *= cycl_m->getSymmetry();
841 
842  return mu * Physics::u_two_pi;
843 }
844 
845 template<typename Value_type, typename Size_type, class Stepper>
847  /*
848  * The formulas for h, fidx and ds are from the paper:
849  * "Tranverse-Longitudinal Coupling by Space Charge in Cyclotrons"
850  * written by Dr. Christian Baumgarten (2012, PSI)
851  * p. 6
852  */
853 
854  // READ IN MAGNETIC FIELD: ONLY FOR STAND-ALONE PROGRAM
855  value_type bint, brint, btint; // B, dB/dr, dB/dtheta
856 
857  value_type invbcon = 1.0 / bcon_m(E0_m, wo_m);
858  value_type en = E / E0_m; // en = E/E0 = E/(mc^2) (E0 is potential energy)
859  value_type p = acon_m(wo_m) * std::sqrt(en * (2.0 + en)); // momentum [p] = m; Gordon, formula (3)
860  value_type p2 = p * p;
861  value_type theta = 0.0; // angle for interpolating
862  value_type ptheta;
863 
864  // resize of container (--> size = N_m, capacity = N_m)
865  h_m.resize(N_m);
866  fidx_m.resize(N_m);
867  ds_m.resize(N_m);
868 
869  for (size_type i = 0; i < N_m; ++i) {
870  // interpolate magnetic field
871  cycl_m->apply(r_m[i], vz_m[i], theta, brint, btint, bint);
872  bint *= invbcon;
873  brint *= invbcon;
874  btint *= invbcon;
875 
876  // inverse bending radius
877  h_m[i] = bint / p;
878 
879  // local field index
880  if (p < pr_m[i])
881  throw OpalException("ClosedOrbitFinder::computeOrbitProperties()",
882  "p_{r}^2 > p^{2} " + std::to_string(p) + " " + std::to_string(pr_m[i]) + " (defined in Gordon paper) --> Square root of negative number.");
883 
884  ptheta = std::sqrt(p2 - pr_m[i] * pr_m[i]);
885 
886  fidx_m[i] = (brint * ptheta - btint * pr_m[i] / r_m[i]) / p2; //(bint*bint);
887 
888  // path length element
889  ds_m[i] = std::hypot(r_m[i] * pr_m[i] / ptheta,r_m[i]) * dtheta_m; // C++11 function
890 
891  // increase angle
892  theta += dtheta_m;
893  }
894 
895  // compute average radius
896  ravg_m = std::accumulate(r_m.begin(),r_m.end(),0.0) / value_type(r_m.size());
897 }
898 
899 template<typename Value_type, typename Size_type, class Stepper>
902 
903  container_type orbitPropertyCopy = orbitProperty;
904 
905  // compute the number of steps per degree
906  value_type deg_step = N_m / 360.0;
907 
908  if (angle < 0.0) {
909  angle = 360.0 + angle;
910  }
911 
912  // compute starting point
913  unsigned int start = deg_step * angle;
914 
915  start %= orbitProperty.size();
916 
917  // copy end to start
918  std::copy(orbitProperty.begin() + start, orbitProperty.end(), orbitPropertyCopy.begin());
919 
920  // copy start to end
921  std::copy_n(orbitProperty.begin(), start, orbitPropertyCopy.end() - start);
922 
923  return orbitPropertyCopy;
924 
925 }
926 
927 #endif
T * value_type(const SliceIterator< T > &)
Tps< T > cos(const Tps< T > &x)
Cosine.
Definition: TpsMath.h:129
Tps< T > sin(const Tps< T > &x)
Sine.
Definition: TpsMath.h:111
Tps< T > sqrt(const Tps< T > &x)
Square root.
Definition: TpsMath.h:91
Inform * gmsg
Definition: Main.cpp:62
const double pi
Definition: fftpack.cpp:894
bool for_each(const BareFieldIterator< T, D > &p, SameFieldID s, C)
Definition: AssignDefs.h:30
T::PETE_Expr_t::PETE_Return_t max(const PETE_Expr< T > &expr, NDIndex< D > &loc)
Definition: ReductionLoc.h:84
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
PETE_TUTree< FnArcCos, typename T::PETE_Expr_t > acos(const PETE_Expr< T > &l)
Definition: PETE.h:725
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:51
constexpr double pi
The value of.
Definition: Physics.h:30
FRONT * fs
Definition: hypervolume.cpp:59
std::string getInputBasename()
get input file name without extension
Definition: OpalData.cpp:658
static OpalData * getInstance()
Definition: OpalData.cpp:195
std::string getAuxiliaryOutputDirectory() const
get the name of the the additional data directory
Definition: OpalData.cpp:650
void read(const double &scaleFactor)
Definition: Cyclotron.cpp:765
virtual double getFMHighE() const
Definition: Cyclotron.cpp:360
virtual double getBScale() const
Definition: Cyclotron.cpp:271
virtual double getSymmetry() const
Definition: Cyclotron.cpp:247
virtual double getFMLowE() const
Definition: Cyclotron.cpp:352
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.
Definition: OpalException.h:28
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.
Definition: Inform.h:42