OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
SigmaGenerator.cpp
Go to the documentation of this file.
1//
2// Class: SigmaGenerator
3// The SigmaGenerator class uses the class <b>ClosedOrbitFinder</b> to get the parameters(inverse bending radius, path length
4// field index and tunes) to initialize the sigma matrix.
5// The main function of this class is <b>match(double, unsigned int)</b>, where it iteratively tries to find a matched
6// distribution for given
7// emittances, energy and current. The computation stops when the L2-norm is smaller than a user-defined tolerance. \n
8// In default mode it prints all space charge maps, cyclotron maps and second moment matrices. The orbit properties, i.e.
9// tunes, average radius, orbit radius, inverse bending radius, path length, field index and frequency error, are printed
10// as well.
11//
12// Copyright (c) 2014, 2018, Matthias Frey, Cristopher Cortes, ETH Zürich
13// All rights reserved
14//
15// Implemented as part of the Semester thesis by Matthias Frey
16// "Matched Distributions in Cyclotrons"
17//
18// Some adaptations done as part of the Bachelor thesis by Cristopher Cortes
19// "Limitations of a linear transfer map method for finding matched distributions in high intensity cyclotrons"
20//
21// This file is part of OPAL.
22//
23// OPAL is free software: you can redistribute it and/or modify
24// it under the terms of the GNU General Public License as published by
25// the Free Software Foundation, either version 3 of the License, or
26// (at your option) any later version.
27//
28// You should have received a copy of the GNU General Public License
29// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
30//
31#include "SigmaGenerator.h"
32
36#include "Utilities/Util.h"
37
39#include "ClosedOrbitFinder.h"
40#include "MapGenerator.h"
41
42#include <cmath>
43#include <filesystem>
44#include <iomanip>
45#include <limits>
46#include <numeric>
47#include <sstream>
48#include <utility>
49
50#include <boost/numeric/odeint/stepper/runge_kutta4.hpp>
51
52extern Inform *gmsg;
53
55 double ex,
56 double ey,
57 double ez,
58 double E,
59 double m,
60 double q,
61 const Cyclotron* cycl,
62 unsigned int N,
63 unsigned int Nsectors,
64 unsigned int truncOrder,
65 bool write)
66 : I_m(I)
67 , wo_m(cycl->getRfFrequ()[0]*Units::MHz2Hz/cycl->getCyclHarm()*2.0*Physics::pi)
68 , E_m(E)
69 , gamma_m(E/m+1.0)
71 , nh_m(cycl->getCyclHarm())
72 , beta_m(std::sqrt(1.0-1.0/gamma2_m))
73 , mass_m(m)
74 , q_m(q)
75 , niterations_m(0)
76 , converged_m(false)
77 , Emin_m(cycl->getFMLowE())
78 , Emax_m(cycl->getFMHighE())
79 , N_m(N)
80 , nSectors_m(Nsectors)
81 , nStepsPerSector_m(N/cycl->getSymmetry())
82 , nSteps_m(N)
83 , error_m(std::numeric_limits<double>::max())
85 , write_m(write)
87 , rinit_m(0.0)
88 , prinit_m(0.0)
89{
90 // minimum beta*gamma
91 double bgam = Util::getBetaGamma(Emin_m, mass_m);
92
93 // set emittances (initialization like that due to old compiler version)
94 // [ex] = [ey] = [ez] = pi*mm*mrad --> [emittance] = m rad
95 // normalized emittance (--> multiply with beta*gamma)
99
100 // Define the Hamiltonian
102
103 // infinitesimal elements
110
111 H_m = [&](double h, double kx, double ky) {
112 return 0.5*px_m*px_m + 0.5*kx*x_m*x_m - h*x_m*delta_m +
113 0.5*pz_m*pz_m + 0.5*ky*z_m*z_m +
115 };
116
117 Hsc_m = [&](double sigx, double sigy, double sigz) {
118 double m = mass_m * Units::MeV2eV;
119
120 // formula (57)
121 double lam = Physics::two_pi*Physics::c / (wo_m * nh_m); // wavelength, [lam] = m
122 // [K3] = m
123 double K3 = 3.0 * std::abs(q_m) * I_m * lam
124 / (20.0 * std::sqrt(5.0) * Physics::pi * Physics::epsilon_0 * m
126
127 // formula (30), (31)
128 // [sigma(0,0)] = m^{2} --> [sx] = [sy] = [sz] = m
129 // In the cyclotron community z is the vertical and y the longitudinal
130 // direction.
131 // x: horizontal
132 // y/l: longitudinal
133 // z: vertical
134 double sx = std::sqrt(std::abs(sigx));
135 double sy = std::sqrt(std::abs(sigy));
136 double sz = std::sqrt(std::abs(sigz));
137
138 double tmp = sx * sy; // [tmp] = m^{2}
139
140 double f = std::sqrt(tmp) / (3.0 * gamma_m * sz); // [f] = 1
141 double kxy = K3 * std::abs(1.0 - f) / ((sx + sy) * sz); // [kxy] = 1/m
142
143 double Kx = kxy / sx;
144 double Ky = kxy / sy;
145 double Kz = K3 * f / (tmp * sz);
146
147 return -0.5 * Kx * x_m * x_m
148 -0.5 * Kz * z_m * z_m
149 -0.5 * Ky * l_m * l_m * gamma2_m;
150 };
151}
152
153
154bool SigmaGenerator::match(double accuracy,
155 unsigned int maxit,
156 unsigned int maxitOrbit,
157 Cyclotron* cycl,
158 double denergy,
159 double rguess,
160 bool full)
161{
162 /* compute the equilibrium orbit for energy E_
163 * and get the the following properties:
164 * - inverse bending radius h
165 * - step sizes of path ds
166 * - tune nuz
167 */
168
169 try {
170 if ( !full )
172
173 // object for space charge map and cyclotron map
174 MapGenerator<double,
175 unsigned int,
176 Series,
177 Map,
179 SpaceCharge> mapgen(nSteps_m);
180
181 // compute cyclotron map and space charge map for each angle and store them into a vector
182 std::vector<matrix_t> Mcycs(nSteps_m), Mscs(nSteps_m);
183
184 ClosedOrbitFinder<double, unsigned int,
185 boost::numeric::odeint::runge_kutta4<container_t> > cof(mass_m, q_m, N_m, cycl, false, nSectors_m);
186
187 if ( !cof.findOrbit(accuracy, maxitOrbit, E_m, denergy, rguess) ) {
188 throw OpalException("SigmaGenerator::match()",
189 "Closed orbit finder didn't converge.");
190 }
191
192 cof.computeOrbitProperties(E_m);
193
194 double angle = cycl->getPHIinit();
195 container_t h = cof.getInverseBendingRadius(angle);
196 container_t r = cof.getOrbit(angle);
197 container_t peo = cof.getMomentum(angle);
198 container_t ds = cof.getPathLength(angle);
199 container_t fidx = cof.getFieldIndex(angle);
200
201 // write properties to file
202 writeOrbitOutput_m(r, peo, h, fidx, ds);
203
204 rinit_m = r[0];
205 prinit_m = peo[0];
206
207 std::string fpath = Util::combineFilePath({
209 "maps"
210 });
211 if (!std::filesystem::exists(fpath)) {
212 std::filesystem::create_directory(fpath);
213 }
214
215 std::pair<double,double> tunes = cof.getTunes();
216 double ravg = cof.getAverageRadius();
217
218 // write to terminal
219 *gmsg << "* ----------------------------" << endl
220 << "* Closed orbit info:" << endl
221 << "*" << endl
222 << "* average radius: " << ravg << " [m]" << endl
223 << "* initial radius: " << r[0] << " [m]" << endl
224 << "* initial momentum: " << peo[0] << " [Beta Gamma]" << endl
225 << "* frequency error: " << cof.getFrequencyError()*100 <<" [ % ] "<< endl
226 << "* horizontal tune: " << tunes.first << endl
227 << "* vertical tune: " << tunes.second << endl
228 << "* ----------------------------" << endl << endl;
229
230 // initialize sigma matrices (for each angle one) (first guess)
231 initialize(tunes.second,ravg);
232
233 // for writing
234 std::ofstream writeMturn, writeMcyc, writeMsc;
235
236 if (write_m) {
237
238 std::string energy = float2string(E_m);
239
240 std::string fname = Util::combineFilePath({
242 "maps",
243 "OneTurnMapsForEnergy" + energy + "MeV.dat"
244 });
245
246 writeMturn.open(fname, std::ios::out);
247
248 fname = Util::combineFilePath({
250 "maps",
251 "SpaceChargeMapPerAngleForEnergy" + energy + "MeV_iteration_0.dat"
252 });
253
254 writeMsc.open(fname, std::ios::out);
255
256 fname = Util::combineFilePath({
258 "maps",
259 "CyclotronMapPerAngleForEnergy" + energy + "MeV.dat"
260 });
261
262 writeMcyc.open(fname, std::ios::out);
263 }
264
265 // calculate only for a single sector (a nSector_-th) of the whole cyclotron
266 for (unsigned int i = 0; i < nSteps_m; ++i) {
267 Mcycs[i] = mapgen.generateMap(H_m(h[i],
268 h[i]*h[i]+fidx[i],
269 -fidx[i]),
270 ds[i],truncOrder_m);
271
272
273 Mscs[i] = mapgen.generateMap(Hsc_m(sigmas_m[i](0,0),
274 sigmas_m[i](2,2),
275 sigmas_m[i](4,4)),
276 ds[i],truncOrder_m);
277
278 writeMatrix(writeMcyc, Mcycs[i]);
279 writeMatrix(writeMsc, Mscs[i]);
280 }
281
282 if (write_m)
283 writeMsc.close();
284
285 // one turn matrix
286 mapgen.combine(Mscs,Mcycs);
287 matrix_t Mturn = mapgen.getMap();
288
289 writeMatrix(writeMturn, Mturn);
290
291 // (inverse) transformation matrix
292 sparse_matrix_t R(6, 6), invR(6, 6);
293
294 // new initial sigma matrix
295 matrix_t newSigma(6,6);
296
297 // for exiting loop
298 bool stop = false;
299
300 double weight = 0.5;
301
302 while (error_m > accuracy && !stop) {
303 // decouple transfer matrix and compute (inverse) tranformation matrix
304 vector_t eigen = decouple(Mturn, R,invR);
305
306 // construct new initial sigma-matrix
307 newSigma = updateInitialSigma(Mturn, eigen, R, invR);
308
309 // compute new sigma matrices for all angles (except for initial sigma)
310 updateSigma(Mscs,Mcycs);
311
312 // compute error with mm^2 and (mrad)^2
313 error_m = L2ErrorNorm(sigmas_m[0] * 1e6, newSigma * 1e6);
314
315 // write new initial sigma-matrix into vector
316 sigmas_m[0] = weight*newSigma + (1.0-weight)*sigmas_m[0];
317
318 if (write_m) {
319
320 std::string energy = float2string(E_m);
321
322 std::string fname = Util::combineFilePath({
324 "maps",
325 "SpaceChargeMapPerAngleForEnergy" + energy + "MeV_iteration_"
326 + std::to_string(niterations_m + 1)
327 + ".dat"
328 });
329
330 writeMsc.open(fname, std::ios::out);
331 }
332
333 // compute new space charge maps
334 for (unsigned int i = 0; i < nSteps_m; ++i) {
335 Mscs[i] = mapgen.generateMap(Hsc_m(sigmas_m[i](0,0),
336 sigmas_m[i](2,2),
337 sigmas_m[i](4,4)),
338 ds[i],truncOrder_m);
339
340 writeMatrix(writeMsc, Mscs[i]);
341 }
342
343 if (write_m) {
344 writeMsc.close();
345 }
346
347 // construct new one turn transfer matrix M
348 mapgen.combine(Mscs,Mcycs);
349 Mturn = mapgen.getMap();
350
351 writeMatrix(writeMturn, Mturn);
352
353 // check if number of iterations has maxit exceeded.
354 stop = (niterations_m++ > maxit);
355 }
356
357 // store converged sigma-matrix
358 sigma_m.resize(6,6,false);
359 sigma_m.swap(newSigma);
360
361 // returns if the sigma matrix has converged
362 converged_m = error_m < accuracy;
363
364 // Close files
365 if (write_m) {
366 writeMturn.close();
367 writeMcyc.close();
368 }
369
370 } catch(const std::exception& e) {
371 throw OpalException("SigmaGenerator::match()", e.what());
372 }
373
374 if ( converged_m && write_m ) {
375 // write tunes
376 std::string fname = Util::combineFilePath({
378 "MatchedDistributions.dat"
379 });
380
381 std::ofstream writeSigmaMatched(fname, std::ios::out);
382
383 std::array<double,3> emit = this->getEmittances();
384
385 writeSigmaMatched << "* Converged (Ex, Ey, Ez) = (" << emit[0]
386 << ", " << emit[1] << ", " << emit[2]
387 << ") pi mm mrad for E= " << E_m << " (MeV)"
388 << std::endl << "* Sigma-Matrix " << std::endl;
389
390 for(unsigned int i = 0; i < sigma_m.size1(); ++ i) {
391 writeSigmaMatched << std::setprecision(4) << std::setw(11)
392 << sigma_m(i,0);
393 for(unsigned int j = 1; j < sigma_m.size2(); ++ j) {
394 writeSigmaMatched << " & " << std::setprecision(4)
395 << std::setw(11) << sigma_m(i,j);
396 }
397 writeSigmaMatched << " \\\\" << std::endl;
398 }
399
400 writeSigmaMatched << std::endl;
401 writeSigmaMatched.close();
402 }
403
404 return converged_m;
405}
406
407
411 sparse_matrix_t& invR)
412{
413 // copy one turn matrix
414 matrix_t M(Mturn);
415
416 // reduce 6x6 matrix to 4x4 matrix
418
419 // compute symplex part
420 matrix_t Ms = rdm_m.symplex(M);
421
422 // diagonalize and compute transformation matrices
423 rdm_m.diagonalize(Ms,R,invR);
424
425 /*
426 * formula (38) in paper of Dr. Christian Baumgarten:
427 * Geometrical method of decoupling
428 *
429 * [0, alpha, 0, 0;
430 * F_{d} = -beta, 0, 0, 0;
431 * 0, 0, 0, gamma;
432 * 0, 0, -delta, 0]
433 *
434 *
435 */
436 vector_t eigen(4);
437 eigen(0) = Ms(0,1); // alpha
438 eigen(1) = - Ms(1,0); // beta
439 eigen(2) = Ms(2,3); // gamma
440 eigen(3) = - Ms(3,2); // delta
441 return eigen;
442}
443
444
446 const matrix_t& sigma)
447{
448 // compute sigma matrix after one turn
449 matrix_t newSigma = matt_boost::gemmm<matrix_t>(Mturn,
450 sigma,
451 boost::numeric::ublas::trans(Mturn));
452
453 // return L2 error
454 return L2ErrorNorm(sigma,newSigma);
455}
456
457
458void SigmaGenerator::initialize(double nuz, double ravg)
459{
460 /*
461 * The initialization is based on the analytical solution of
462 * using a spherical symmetric beam in the paper:
463 * Transverse-longitudinal coupling by space charge in cyclotrons
464 * by Dr. Christian Baumgarten
465 * (formulas: (46), (56), (57))
466 */
467
468
469 /* Units:
470 * ----------------------------------------------
471 * [wo] = 1/s
472 * [nh] = 1
473 * [q0] = 1 e
474 * [I] = A
475 * [eps0] = (A*s)^{2}/(N*m^{2})
476 * [E0] = MeV/(c^{2}) (with speed of light c)
477 * [beta] = 1
478 * [gamma] = 1
479 * [m] = eV/c^2
480 *
481 * [lam] = m
482 * [K3] = m
483 * [alpha] = 1
484 * ----------------------------------------------
485 */
486
487 // helper constants
488 double invbg = 1.0 / (beta_m * gamma_m);
489
490 double mass = mass_m * Units::MeV2eV;
491
492 // emittance [ex] = [ey] = [ez] = m rad
493 double ex = emittance_m[0] * invbg; // [ex] = m rad
494 double ey = emittance_m[1] * invbg; // [ey] = m rad
495 double ez = emittance_m[2] * invbg; // [ez] = m rad
496
497 // initial guess of emittance, [e] = m rad
498 double guessedEmittance = std::cbrt(ex * ey * ez); // cbrt computes cubic root (C++11) <cmath>
499
500 // cyclotron radius [rcyc] = m
501 double rcyc = ravg / beta_m;
502
503 // "average" inverse bending radius
504 double avgInverseBendingRadius = 1.0 / ravg;
505
506 // formula (57)
507 double lam = Physics::two_pi * Physics::c / (wo_m * nh_m); // wavelength, [lam] = m
508
509 // m * c^3 --> c^2 in [m] = eV / c^2 cancel --> m * c in denominator
510 double K3 = 3.0 * std::abs(q_m) * I_m * lam
511 / (20.0 * std::sqrt(5.0) * Physics::pi * Physics::epsilon_0 * mass
512 * Physics::c * beta_m * beta_m * gamma2_m * gamma_m); // [K3] = m
513
514 // c in denominator cancels with unit of [m] = eV / c^2 --> we need to multiply
515 // with c in order to get dimensionless quantity
516 double alpha = (std::abs(q_m) * Physics::mu_0 * I_m * Physics::c
517 / (5.0 * std::sqrt(10.0) * mass * gamma_m * nh_m)
518 * std::sqrt(rcyc * rcyc * rcyc / (std::pow(guessedEmittance, 3)))); // [alpha] = 1/rad --> [alpha] = 1
519
520 double sig0 = std::sqrt(2.0 * rcyc * guessedEmittance) / gamma_m; // [sig0] = m*sqrt(rad) --> [sig0] = m
521
522 // formula (56)
523 double sig; // [sig] = m
524 if (alpha >= 2.5) {
525 sig = sig0 * std::cbrt(1.0 + alpha); // cbrt computes cubic root (C++11) <cmath>
526 } else if (alpha >= 0) {
527 sig = sig0 * (1 + alpha * (0.25 - 0.03125 * alpha));
528 } else {
529 throw OpalException("SigmaGenerator::initialize()",
530 "Negative alpha value: " + std::to_string(alpha) + " < 0");
531 }
532
533 // K = Kx = Ky = Kz
534 double K = K3 * gamma_m / (3.0 * sig * sig * sig); // formula (46), [K] = 1/m^{2}
535 double kx = std::pow(avgInverseBendingRadius, 2) * gamma2_m;// formula (46) (assumption of an isochronous cyclotron), [kx] = 1/m^{2}
536
537 double a = 0.5 * kx - K; // formula (22) (with K = Kx = Kz), [a] = 1/m^{2}
538 double b = K * K; // formula (22) (with K = Kx = Kz and kx = h^2*gamma^2), [b] = 1/m^{4}
539
540
541 // b must be positive, otherwise no real-valued frequency
542 if (b < 0)
543 throw OpalException("SigmaGenerator::initialize()",
544 "Negative value --> No real-valued frequency.");
545
546 double tmp = a * a - b; // [tmp] = 1/m^{4}
547 if (tmp < 0)
548 throw OpalException("SigmaGenerator::initialize()",
549 "Square root of negative number.");
550
551 tmp = std::sqrt(tmp); // [tmp] = 1/m^{2}
552
553 if (a < tmp)
554 throw OpalException("SigmaGenerator::initialize()",
555 "Square root of negative number.");
556
557 if (std::pow(avgInverseBendingRadius, 2) * nuz * nuz <= K)
558 throw OpalException("SigmaGenerator::initialize()",
559 "h^{2} * nu_{z}^{2} <= K (i.e. square root of negative number)");
560
561 double Omega = std::sqrt(a + tmp); // formula (22), [Omega] = 1/m
562 double omega = std::sqrt(a - tmp); // formula (22), [omega] = 1/m
563
564 double A = avgInverseBendingRadius / (Omega * Omega + K); // formula (26), [A] = m
565 double B = avgInverseBendingRadius / (omega * omega + K); // formula (26), [B] = m
566 double invAB = 1.0 / (B - A); // [invAB] = 1/m
567
568 // construct initial sigma-matrix (formula (29, 30, 31)
569 matrix_t sigma = boost::numeric::ublas::zero_matrix<double>(6);
570
571 // formula (30), [sigma(0,0)] = m^2 rad
572 sigma(0,0) = invAB * (B * ex / Omega + A * ez / omega);
573
574 // [sigma(0,5)] = [sigma(5,0)] = m rad
575 sigma(0,5) = sigma(5,0) = invAB * (ex / Omega + ez / omega);
576
577 // [sigma(1,1)] = rad
578 sigma(1,1) = invAB * (B * ex * Omega + A * ez * omega);
579
580 // [sigma(1,4)] = [sigma(4,1)] = m rad
581 sigma(1,4) = sigma(4,1) = invAB * (ex * Omega+ez * omega) / (K * gamma2_m);
582
583 // formula (31), [sigma(2,2)] = m rad
584 sigma(2,2) = ey / (std::sqrt(std::pow(avgInverseBendingRadius,2) * nuz * nuz - K));
585
586 sigma(3,3) = (ey * ey) / sigma(2,2);
587
588 // [sigma(4,4)] = m^{2} rad
589 sigma(4,4) = invAB * (A * ex * Omega + B * ez * omega) / (K * gamma2_m);
590
591 // formula (30), [sigma(5,5)] = rad
592 sigma(5,5) = invAB * (ex / (B * Omega) + ez / (A * omega));
593
594 // fill in initial guess of the sigma matrix (for each angle the same guess)
595 sigmas_m.resize(nSteps_m);
596 for (typename std::vector<matrix_t>::iterator it = sigmas_m.begin(); it != sigmas_m.end(); ++it) {
597 *it = sigma;
598 }
599
600 if (write_m) {
601 std::string energy = float2string(E_m);
602
603 std::string fname = Util::combineFilePath({
605 "maps",
606 "InitialSigmaPerAngleForEnergy" + energy + "MeV.dat"
607 });
608
609 std::ofstream writeInit(fname, std::ios::out);
610 writeMatrix(writeInit, sigma);
611 writeInit.close();
612 }
613}
614
615
618 const vector_t& eigen,
620 sparse_matrix_t& invR)
621{
622 /*
623 * This function is based on the paper of Dr. Christian Baumgarten:
624 * Transverse-Longitudinal Coupling by Space Charge in Cyclotrons (2012)
625 */
626
627 /*
628 * Function input:
629 * - M: one turn transfer matrix
630 * - eigen = {alpha, beta, gamma, delta}
631 * - R: transformation matrix (in paper: E)
632 * - invR: inverse transformation matrix (in paper: E^{-1}
633 */
634
635 // normalize emittances
636 double invbg = 1.0 / (beta_m * gamma_m);
637 double ex = emittance_m[0] * invbg;
638 double ey = emittance_m[1] * invbg;
639 double ez = emittance_m[2] * invbg;
640
641 // alpha^2-beta*gamma = 1
642
643 /* 0 eigen(0) 0 0
644 * eigen(1) 0 0 0
645 * 0 0 0 eigen(2)
646 * 0 0 eigen(3) 0
647 *
648 * M = cos(mux)*[1, 0; 0, 1] + sin(mux)*[alpha, beta; -gamma, -alpha], Book, p. 242
649 *
650 * -----------------------------------------------------------------------------------
651 * X-DIRECTION and L-DIRECTION
652 * -----------------------------------------------------------------------------------
653 * --> eigen(0) = sin(mux)*betax
654 * --> eigen(1) = -gammax*sin(mux)
655 *
656 * What is sin(mux)? --> alphax = 0 --> -alphax^2+betax*gammax = betax*gammax = 1
657 *
658 * Convention: betax > 0
659 *
660 * 1) betax = 1/gammax
661 * 2) eig0 = sin(mux)*betax
662 * 3) eig1 = -gammax*sin(mux)
663 *
664 * eig0 = sin(mux)/gammax
665 * eig1 = -gammax*sin(mux) <--> 1/gammax = -sin(mux)/eig1
666 *
667 * eig0 = -sin(mux)^2/eig1 --> -sin(mux)^2 = eig0*eig1 --> sin(mux) = sqrt(-eig0*eig1)
668 * --> gammax = -eig1/sin(mux)
669 * --> betax = eig0/sin(mux)
670 */
671
672
673 // x-direction
674 //double alphax = 0.0;
675 double betax = std::sqrt(std::fabs(eigen(0) / eigen(1)));
676 double gammax = 1.0 / betax;
677
678 // l-direction
679 //double alphal = 0.0;
680 double betal = std::sqrt(std::fabs(eigen(2) / eigen(3)));
681 double gammal = 1.0 / betal;
682
683 /*
684 * -----------------------------------------------------------------------------------
685 * Z-DIRECTION
686 * -----------------------------------------------------------------------------------
687 *
688 * m22 m23
689 * m32 m33
690 *
691 * m22 = cos(muz) + alpha*sin(muz)
692 * m33 = cos(muz) - alpha*sin(muz)
693 *
694 * --> cos(muz) = 0.5*(m22 + m33)
695 * sin(muz) = sign(m32)*sqrt(1-cos(muz)^2)
696 *
697 * m22-m33 = 2*alpha*sin(muz) --> alpha = 0.5*(m22-m33)/sin(muz)
698 *
699 * m23 = betaz*sin(muz) --> betaz = m23/sin(muz)
700 * m32 = -gammaz*sin(muz) --> gammaz = -m32/sin(muz)
701 */
702
703 double cosz = 0.5 * (M(2,2) + M(3,3));
704
705 double sign = (std::signbit(M(2,3))) ? double(-1) : double(1);
706
707 double invsinz = sign / std::sqrt(std::fabs( 1.0 - cosz * cosz));
708
709 double alphaz = 0.5 * (M(2,2) - M(3,3)) * invsinz;
710 double betaz = M(2,3) * invsinz;
711 double gammaz = - M(3,2) * invsinz;
712
713 // Convention beta>0
714 if (std::signbit(betaz)) // singbit = true if beta<0, else false
715 betaz *= -1.0;
716
717 // diagonal matrix with eigenvalues
718 matrix_t D = boost::numeric::ublas::zero_matrix<double>(6,6);
719 // x-direction
720 D(0,1) = betax * ex;
721 D(1,0) = - gammax * ex;
722 // z-direction
723 D(2,2) = alphaz * ey;
724 D(3,3) = - alphaz * ey;
725 D(2,3) = betaz * ey;
726 D(3,2) = - gammaz * ey;
727 // l-direction
728 D(4,5) = betal * ez;
729 D(5,4) = - gammal * ez;
730
731 // expand 4x4 transformation matrices to 6x6
734
735 // symplectic matrix
736 sparse_matrix_t S(6,6,6);
737 S(0,1) = S(2,3) = S(4,5) = 1;
738 S(1,0) = S(3,2) = S(5,4) = -1;
739
740 matrix_t sigma = matt_boost::gemmm<matrix_t>(-invR,D,R);
741 sigma = boost::numeric::ublas::prod(sigma,S);
742
743 if (write_m) {
744 std::string energy = float2string(E_m);
745
746 std::string fname = Util::combineFilePath({
748 "maps",
749 "InitialSigmaPerAngleForEnergy" + energy + "MeV.dat"
750 });
751
752 std::ofstream writeInit(fname, std::ios::app);
753 writeMatrix(writeInit, sigma);
754 writeInit.close();
755 }
756
757 return sigma;
758}
759
760
761void SigmaGenerator::updateSigma(const std::vector<matrix_t>& Mscs,
762 const std::vector<matrix_t>& Mcycs)
763{
764 std::ofstream writeSigma;
765
766 if (write_m) {
767 std::string energy = float2string(E_m);
768
769 std::string fname = Util::combineFilePath({
771 "maps",
772 "SigmaPerAngleForEnergy" + energy + "MeV_iteration_"
773 + std::to_string(niterations_m + 1)
774 + ".dat"
775 });
776
777 writeSigma.open(fname,std::ios::out);
778 }
779
780 // initial sigma is already computed
781 writeMatrix(writeSigma, sigmas_m[0]);
782
783 for (unsigned int i = 1; i < nSteps_m; ++i) {
784 // transfer matrix for one angle
785 matrix_t M = boost::numeric::ublas::prod(Mscs[i - 1],Mcycs[i - 1]);
786 // transfer the matrix sigma
788 boost::numeric::ublas::trans(M));
789
790 writeMatrix(writeSigma, sigmas_m[i]);
791 }
792
793 if (write_m) {
794 writeSigma.close();
795 }
796}
797
798
800 const matrix_t& newS)
801{
802 // compute difference
803 matrix_t diff = newS - oldS;
804
805 // sum squared error up and take square root
806 return std::sqrt(std::inner_product(diff.data().begin(),
807 diff.data().end(),
808 diff.data().begin(),
809 0.0));
810}
811
812
813std::string SigmaGenerator::float2string(double val) {
814 std::ostringstream out;
815 out << std::setprecision(6) << val;
816 return out.str();
817}
818
819
821 const container_t& r,
822 const container_t& peo,
823 const container_t& h,
824 const container_t& fidx,
825 const container_t& ds)
826{
827 if (!write_m)
828 return;
829
830 std::string energy = float2string(E_m);
831 std::string fname = Util::combineFilePath({
833 "PropertiesForEnergy" + energy + "MeV.dat"
834 });
835
836 std::ofstream writeProperties(fname, std::ios::out);
837
838 writeProperties << std::left
839 << std::setw(25) << "orbit radius"
840 << std::setw(25) << "orbit momentum"
841 << std::setw(25) << "inverse bending radius"
842 << std::setw(25) << "field index"
843 << std::setw(25) << "path length" << std::endl;
844
845 for (unsigned int i = 0; i < r.size(); ++i) {
846 writeProperties << std::setprecision(10) << std::left
847 << std::setw(25) << r[i]
848 << std::setw(25) << peo[i]
849 << std::setw(25) << h[i]
850 << std::setw(25) << fidx[i]
851 << std::setw(25) << ds[i] << std::endl;
852 }
853 writeProperties.close();
854}
855
856
857void SigmaGenerator::writeMatrix(std::ofstream& os, const matrix_t& m) {
858 if (!write_m)
859 return;
860
861 for (unsigned int i = 0; i < 6; ++i) {
862 for (unsigned int j = 0; j < 6; ++j)
863 os << m(i, j) << " ";
864 }
865 os << std::endl;
866}
Tps< T > sqrt(const Tps< T > &x)
Square root.
Definition TpsMath.h:91
int Tps< T >::truncOrder
Definition Tps.hpp:217
Inform * gmsg
Definition Main.cpp:70
T::PETE_Expr_t::PETE_Return_t max(const PETE_Expr< T > &expr, NDIndex< D > &loc)
const double pi
Definition fftpack.cpp:894
std::complex< double > a
FnArg FnReal sign(a))
Inform & endl(Inform &inf)
Definition Inform.cpp:42
STL namespace.
Definition Air.h:27
constexpr double two_pi
The value of.
Definition Physics.h:33
constexpr double epsilon_0
The permittivity of vacuum in As/Vm.
Definition Physics.h:51
constexpr double mu_0
The permeability of vacuum in Vs/Am.
Definition Physics.h:48
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
constexpr double mrad2rad
Definition Units.h:140
constexpr double mm2m
Definition Units.h:29
constexpr double MeV2eV
Definition Units.h:74
std::string combineFilePath(std::initializer_list< std::string > ilist)
Definition Util.cpp:197
double getBetaGamma(double Ekin, double mass)
Definition Util.h:61
BOOST_UBLAS_INLINE M gemmm(const ublas::matrix_expression< E1 > &e1, const ublas::matrix_expression< E2 > &e2, const ublas::matrix_expression< E3 > &e3)
Generalized matrix-matrix-matrix multiplication .
static OpalData * getInstance()
Definition OpalData.cpp:196
std::string getAuxiliaryOutputDirectory() const
get the name of the the additional data directory
Definition OpalData.cpp:666
virtual double getPHIinit() const
static FTps makeVariable(int var)
static void setGlobalTruncOrder(int order)
This class generates the matrices for the one turn matrix of a cyclotron.
unsigned int niterations_m
Is the number of iterations needed for convergence.
double E_m
Stores the user-define energy, .
bool write_m
Decides for writing output (default: true).
void writeMatrix(std::ofstream &, const matrix_t &)
double isEigenEllipse(const matrix_t &Mturn, const matrix_t &sigma)
Checks if the sigma-matrix is an eigenellipse and returns the L2 error.
std::vector< double > container_t
Container for storing the properties for each angle.
std::array< double, 3 > getEmittances() const
double gamma2_m
Relativistic factor squared.
FTps< double, 2 *3 > Series
Type of the truncated power series.
void expand(matrix &)
Hamiltonian H_m
Stores the Hamiltonian of the cyclotron.
double error_m
Error of computation.
matrix_t updateInitialSigma(const matrix_t &, const vector_t &, sparse_matrix_t &, sparse_matrix_t &)
Computes the new initial sigma matrix.
unsigned int nSteps_m
Number of integration steps in total.
RealDiracMatrix::sparse_matrix_t sparse_matrix_t
Type for storing the sparse maps.
unsigned int nSectors_m
Number of (symmetric) sectors the field is averaged over.
unsigned int nStepsPerSector_m
Number of integration steps per sector (--> also: number of maps per sector).
double I_m
Stores the value of the current, .
std::function< Series(double, double, double)> SpaceCharge
Type of the Hamiltonian for the space charge.
double Emin_m
Minimum energy needed in cyclotron, .
unsigned int truncOrder_m
Truncation order of the power series for the Hamiltonian (cyclotron and space charge).
SigmaGenerator(double I, double ex, double ey, double ez, double E, double m, double q, const Cyclotron *cycl, unsigned int N, unsigned int Nsectors, unsigned int truncOrder, bool write=true)
Constructs an object of type SigmaGenerator.
double L2ErrorNorm(const matrix_t &, const matrix_t &)
Returns the L2-error norm between the old and new sigma-matrix.
RealDiracMatrix rdm_m
RealDiracMatrix-class member used for decoupling
FVps< double, 2 *3 > Map
Type of a map.
double mass_m
Is the mass of the particles, .
Series x_m
All variables x, px, z, pz, l, delta.
void reduce(matrix &)
SpaceCharge Hsc_m
Stores the Hamiltonian for the space charge.
double q_m
Is the particle charge [e].
RealDiracMatrix::matrix_t matrix_t
Type for storing maps.
std::vector< matrix_t > sigmas_m
Vector storing the second moment matrix for each angle.
double Emax_m
Maximum energy reached in cyclotron, .
matrix_t sigma_m
Stores the stationary distribution (sigma matrix).
double nh_m
Harmonic number, .
vector_t decouple(const matrix_t &Mturn, sparse_matrix_t &R, sparse_matrix_t &invR)
Block diagonalizes the symplex part of the one turn transfer matrix.
double wo_m
Is the orbital frequency, .
std::function< Series(double, double, double)> Hamiltonian
Type of the Hamiltonian for the cyclotron.
unsigned int N_m
Number of integration steps for closed orbit computation.
void initialize(double, double)
bool converged_m
Is true if converged, false otherwise.
void writeOrbitOutput_m(const container_t &r_turn, const container_t &peo, const container_t &h_turn, const container_t &fidx_turn, const container_t &ds_turn)
Called within SigmaGenerator::match().
void updateSigma(const std::vector< matrix_t > &, const std::vector< matrix_t > &)
Computes new sigma matrices (one for each angle).
std::string float2string(double val)
Transforms a floating point value to a string.
RealDiracMatrix::vector_t vector_t
Type for storing vectors.
double beta_m
Velocity (c/v), .
std::array< double, 3 > emittance_m
bool match(double accuracy, unsigned int maxit, unsigned int maxitOrbit, Cyclotron *cycl, double denergy, double rguess, bool full)
Searches for a matched distribution.
The base class for all OPAL exceptions.