OPAL (Object Oriented Parallel Accelerator Library)  2021.1.99
OPAL
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
SigmaGenerator Class Reference

#include <SigmaGenerator.h>

Collaboration diagram for SigmaGenerator:
Collaboration graph
[legend]

Public Types

typedef RealDiracMatrix::matrix_t matrix_t
 Type for storing maps. More...
 
typedef RealDiracMatrix::sparse_matrix_t sparse_matrix_t
 Type for storing the sparse maps. More...
 
typedef RealDiracMatrix::vector_t vector_t
 Type for storing vectors. More...
 
typedef std::vector< double > container_t
 Container for storing the properties for each angle. More...
 
typedef FTps< double, 2 *3 > Series
 Type of the truncated power series. More...
 
typedef FVps< double, 2 *3 > Map
 Type of a map. More...
 
typedef std::function< Series(double, double, double)> Hamiltonian
 Type of the Hamiltonian for the cyclotron. More...
 
typedef std::function< Series(double, double, double)> SpaceCharge
 Type of the Hamiltonian for the space charge. More...
 

Public Member Functions

 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. More...
 
bool match (double accuracy, unsigned int maxit, unsigned int maxitOrbit, Cyclotron *cycl, double denergy, double rguess, bool full)
 Searches for a matched distribution. More...
 
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. More...
 
double isEigenEllipse (const matrix_t &Mturn, const matrix_t &sigma)
 Checks if the sigma-matrix is an eigenellipse and returns the L2 error. More...
 
matrix_tgetSigma ()
 Returns the converged stationary distribution. More...
 
unsigned int getIterations () const
 Returns the number of iterations needed for convergence (if not converged, it returns zero) More...
 
double getError () const
 Returns the error (if the program didn't converged, one can call this function to check the error) More...
 
std::array< double, 3 > getEmittances () const
 
const double & getInjectionRadius () const
 
const double & getInjectionMomentum () const
 

Private Member Functions

void initialize (double, double)
 
matrix_t updateInitialSigma (const matrix_t &, const vector_t &, sparse_matrix_t &, sparse_matrix_t &)
 Computes the new initial sigma matrix. More...
 
void updateSigma (const std::vector< matrix_t > &, const std::vector< matrix_t > &)
 Computes new sigma matrices (one for each angle) More...
 
double L2ErrorNorm (const matrix_t &, const matrix_t &)
 Returns the L2-error norm between the old and new sigma-matrix. More...
 
std::string float2string (double val)
 Transforms a floating point value to a string. More...
 
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(). More...
 
void writeMatrix (std::ofstream &, const matrix_t &)
 
template<class matrix >
void reduce (matrix &)
 
template<class matrix >
void expand (matrix &)
 

Private Attributes

double I_m
 Stores the value of the current, \( \left[I\right] = A \). More...
 
std::array< double, 3 > emittance_m
 
double wo_m
 Is the orbital frequency, \( \left[\omega_{o}\right] = \frac{1}{s} \). More...
 
double E_m
 Stores the user-define energy, \( \left[E\right] = MeV \). More...
 
double gamma_m
 
double gamma2_m
 Relativistic factor squared. More...
 
double nh_m
 Harmonic number, \( \left[N_{h}\right] = 1 \). More...
 
double beta_m
 Velocity (c/v), \( \left[\beta\right] = 1 \). More...
 
double m_m
 Is the mass of the particles, \( \left[m\right] = \frac{MeV}{c^{2}} \). More...
 
double q_m
 Is the particle charge [e]. More...
 
unsigned int niterations_m
 Is the number of iterations needed for convergence. More...
 
bool converged_m
 Is true if converged, false otherwise. More...
 
double Emin_m
 Minimum energy needed in cyclotron, \( \left[E_{min}\right] = MeV \). More...
 
double Emax_m
 Maximum energy reached in cyclotron, \( \left[E_{max}\right] = MeV \). More...
 
unsigned int N_m
 Number of integration steps for closed orbit computation. More...
 
unsigned int nSectors_m
 Number of (symmetric) sectors the field is averaged over. More...
 
unsigned int nStepsPerSector_m
 Number of integration steps per sector (--> also: number of maps per sector) More...
 
unsigned int nSteps_m
 Number of integration steps in total. More...
 
double error_m
 Error of computation. More...
 
unsigned int truncOrder_m
 Truncation order of the power series for the Hamiltonian (cyclotron and space charge) More...
 
bool write_m
 Decides for writing output (default: true) More...
 
matrix_t sigma_m
 Stores the stationary distribution (sigma matrix) More...
 
std::vector< matrix_tsigmas_m
 Vector storing the second moment matrix for each angle. More...
 
Hamiltonian H_m
 Stores the Hamiltonian of the cyclotron. More...
 
SpaceCharge Hsc_m
 Stores the Hamiltonian for the space charge. More...
 
Series x_m
 All variables x, px, z, pz, l, delta. More...
 
Series px_m
 
Series z_m
 
Series pz_m
 
Series l_m
 
Series delta_m
 
double rinit_m
 
double prinit_m
 
RealDiracMatrix rdm_m
 RealDiracMatrix-class member used for decoupling More...
 

Detailed Description

Definition at line 50 of file SigmaGenerator.h.

Member Typedef Documentation

◆ container_t

typedef std::vector<double> SigmaGenerator::container_t

Container for storing the properties for each angle.

Definition at line 60 of file SigmaGenerator.h.

◆ Hamiltonian

typedef std::function<Series(double,double,double)> SigmaGenerator::Hamiltonian

Type of the Hamiltonian for the cyclotron.

Definition at line 66 of file SigmaGenerator.h.

◆ Map

typedef FVps<double,2*3> SigmaGenerator::Map

Type of a map.

Definition at line 64 of file SigmaGenerator.h.

◆ matrix_t

Type for storing maps.

Definition at line 54 of file SigmaGenerator.h.

◆ Series

typedef FTps<double,2*3> SigmaGenerator::Series

Type of the truncated power series.

Definition at line 62 of file SigmaGenerator.h.

◆ SpaceCharge

typedef std::function<Series(double,double,double)> SigmaGenerator::SpaceCharge

Type of the Hamiltonian for the space charge.

Definition at line 68 of file SigmaGenerator.h.

◆ sparse_matrix_t

Type for storing the sparse maps.

Definition at line 56 of file SigmaGenerator.h.

◆ vector_t

Type for storing vectors.

Definition at line 58 of file SigmaGenerator.h.

Constructor & Destructor Documentation

◆ SigmaGenerator()

SigmaGenerator::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.

Parameters
Ispecifies the current for which a matched distribution should be found, \( [I] = A \)
exis the emittance in x-direction (horizontal), \( \left[\varepsilon_{x}\right] = \pi\ mm\ mrad \)
eyis the emittance in y-direction (longitudinal), \( \left[\varepsilon_{y}\right] = \pi\ mm\ mrad \)
ezis the emittance in z-direction (vertical), \( \left[\varepsilon_{z}\right] = \pi\ mm\ mrad \)
Eis the energy, \( \left[E\right] = MeV \)
mis the mass of the particles \( \left[m\right] = \frac{MeV}{c^{2}} \)
qis the particle charge [e]
cyclis the cyclotron element
Nis the number of integration steps (closed orbit computation). That's why its also the number of maps (for each integration step a map)
Nsectorsis the number of sectors that the field map is averaged over (closed orbit computation)
truncOrderis the truncation order for power series of the Hamiltonian
writeis a boolean (default: true). If true all maps of all iterations are stored, otherwise not.

Definition at line 55 of file SigmaGenerator.cpp.

References abs(), beta_m, Physics::c, delta_m, Emin_m, emittance_m, Physics::epsilon_0, gamma2_m, gamma_m, Util::getBetaGamma(), H_m, Hsc_m, I_m, l_m, m_m, FTps< double, 2 *3 >::makeVariable(), nh_m, Physics::pi, px_m, pz_m, q_m, FTps< double, 2 *3 >::setGlobalTruncOrder(), sqrt(), truncOrder_m, Physics::two_pi, wo_m, x_m, and z_m.

Here is the call graph for this function:

Member Function Documentation

◆ decouple()

SigmaGenerator::vector_t SigmaGenerator::decouple ( const matrix_t Mturn,
sparse_matrix_t R,
sparse_matrix_t invR 
)

Block diagonalizes the symplex part of the one turn transfer matrix.

It computes the transformation matrix R and its inverse invR.

Parameters
Mturnis a 6x6 dimensional one turn transfer matrix
Ris the 6x6 dimensional transformation matrix (gets computed)
invRis the 6x6 dimensional inverse transformation (gets computed)

Definition at line 411 of file SigmaGenerator.cpp.

References RealDiracMatrix::diagonalize(), Attrib::Distribution::R, rdm_m, and RealDiracMatrix::symplex().

Referenced by match().

Here is the call graph for this function:

◆ expand()

template<class matrix >
void SigmaGenerator::expand ( matrix &  M)
private

Definition at line 352 of file SigmaGenerator.h.

◆ float2string()

std::string SigmaGenerator::float2string ( double  val)
private

Transforms a floating point value to a string.

Parameters
valis the floating point value which is transformed to a string

Definition at line 815 of file SigmaGenerator.cpp.

Referenced by initialize(), match(), updateInitialSigma(), updateSigma(), and writeOrbitOutput_m().

◆ getEmittances()

std::array< double, 3 > SigmaGenerator::getEmittances ( ) const
inline

Returns the emittances (ex,ey,ez) in \( \pi\ mm\ mrad \) for which the code converged (since the whole simulation is done on normalized emittances)

Definition at line 307 of file SigmaGenerator.h.

References beta_m, emittance_m, gamma_m, and Physics::pi.

Referenced by match().

◆ getError()

double SigmaGenerator::getError ( ) const
inline

Returns the error (if the program didn't converged, one can call this function to check the error)

Definition at line 300 of file SigmaGenerator.h.

References error_m.

◆ getInjectionMomentum()

const double& SigmaGenerator::getInjectionMomentum ( ) const
inline

Definition at line 141 of file SigmaGenerator.h.

References prinit_m.

◆ getInjectionRadius()

const double& SigmaGenerator::getInjectionRadius ( ) const
inline

Definition at line 137 of file SigmaGenerator.h.

References rinit_m.

◆ getIterations()

unsigned int SigmaGenerator::getIterations ( ) const
inline

Returns the number of iterations needed for convergence (if not converged, it returns zero)

Definition at line 293 of file SigmaGenerator.h.

References converged_m, and niterations_m.

◆ getSigma()

SigmaGenerator::matrix_t & SigmaGenerator::getSigma ( )
inline

Returns the converged stationary distribution.

Definition at line 286 of file SigmaGenerator.h.

References sigma_m.

◆ initialize()

void SigmaGenerator::initialize ( double  nuz,
double  ravg 
)
private

Initializes a first guess of the sigma matrix with the assumption of a spherical symmetric beam (ex = ey = ez). For each angle split the same initial guess is taken.

Parameters
nuzis the vertical tune
ravgis the average radius of the closed orbit

Definition at line 460 of file SigmaGenerator.cpp.

References a, abs(), Physics::alpha, beta_m, Physics::c, Util::combineFilePath(), Physics::e, E_m, emittance_m, Physics::epsilon_0, float2string(), gamma2_m, gamma_m, OpalData::getAuxiliaryOutputDirectory(), OpalData::getInstance(), I_m, m_m, Physics::mu_0, nh_m, nSteps_m, Physics::pi, q_m, sigmas_m, sqrt(), Physics::two_pi, wo_m, write_m, and writeMatrix().

Referenced by match().

Here is the call graph for this function:

◆ isEigenEllipse()

double SigmaGenerator::isEigenEllipse ( const matrix_t Mturn,
const matrix_t sigma 
)

Checks if the sigma-matrix is an eigenellipse and returns the L2 error.

The idea of this function is taken from Dr. Christian Baumgarten's program.

Parameters
Mturnis the one turn transfer matrix
sigmais the sigma matrix to be tested

Definition at line 447 of file SigmaGenerator.cpp.

References L2ErrorNorm().

Here is the call graph for this function:

◆ L2ErrorNorm()

double SigmaGenerator::L2ErrorNorm ( const matrix_t oldS,
const matrix_t newS 
)
private

Returns the L2-error norm between the old and new sigma-matrix.

Parameters
oldSis the old sigma matrix
newSis the new sigma matrix

Definition at line 801 of file SigmaGenerator.cpp.

References sqrt().

Referenced by isEigenEllipse(), and match().

Here is the call graph for this function:

◆ match()

bool SigmaGenerator::match ( double  accuracy,
unsigned int  maxit,
unsigned int  maxitOrbit,
Cyclotron cycl,
double  denergy,
double  rguess,
bool  full 
)

Searches for a matched distribution.

Returns "true" if a matched distribution within the accuracy could be found, returns "false" otherwise.

Parameters
accuracyis used for computing the equilibrium orbit (to a certain accuracy)
maxitis the maximum number of iterations (the program stops if within this value no stationary distribution was found)
maxitOrbitis the maximum number of iterations for finding closed orbit
angledefines the start of the sector (one can choose any angle between 0° and 360°)
denergyenergy step size for closed orbit finder [MeV]
rguessvalue of radius for closed orbit finder
typespecifies the magnetic field format (e.g. CARBONCYCL)
fullmatch over full turn not just single sector

Definition at line 156 of file SigmaGenerator.cpp.

References Util::combineFilePath(), converged_m, decouple(), Physics::e, E_m, endl(), error_m, float2string(), OpalData::getAuxiliaryOutputDirectory(), getEmittances(), OpalData::getInstance(), Cyclotron::getPHIinit(), gmsg, H_m, Hsc_m, initialize(), L2ErrorNorm(), m_m, N_m, niterations_m, nSectors_m, nSteps_m, nStepsPerSector_m, prinit_m, q_m, Attrib::Distribution::R, rinit_m, sigma_m, sigmas_m, truncOrder_m, updateInitialSigma(), updateSigma(), write_m, writeMatrix(), and writeOrbitOutput_m().

Here is the call graph for this function:

◆ reduce()

template<class matrix >
void SigmaGenerator::reduce ( matrix &  M)
private

Definition at line 319 of file SigmaGenerator.h.

◆ updateInitialSigma()

SigmaGenerator::matrix_t SigmaGenerator::updateInitialSigma ( const matrix_t M,
const vector_t eigen,
sparse_matrix_t R,
sparse_matrix_t invR 
)
private

Computes the new initial sigma matrix.

Parameters
Mis the 6x6 one turn transfer matrix
Ris the transformation matrix
invRis the inverse transformation matrix

Definition at line 619 of file SigmaGenerator.cpp.

References beta_m, Util::combineFilePath(), E_m, emittance_m, fabs(), float2string(), gamma_m, OpalData::getAuxiliaryOutputDirectory(), OpalData::getInstance(), prod(), Attrib::Distribution::R, sign(), sqrt(), write_m, and writeMatrix().

Referenced by match().

Here is the call graph for this function:

◆ updateSigma()

void SigmaGenerator::updateSigma ( const std::vector< matrix_t > &  Mscs,
const std::vector< matrix_t > &  Mcycs 
)
private

Computes new sigma matrices (one for each angle)

Mscs is a vector of all space charge maps Mcycs is a vector of all cyclotron maps

Definition at line 763 of file SigmaGenerator.cpp.

References Util::combineFilePath(), E_m, float2string(), OpalData::getAuxiliaryOutputDirectory(), OpalData::getInstance(), niterations_m, nSteps_m, prod(), sigmas_m, write_m, and writeMatrix().

Referenced by match().

Here is the call graph for this function:

◆ writeMatrix()

void SigmaGenerator::writeMatrix ( std::ofstream &  os,
const matrix_t m 
)
private

Definition at line 859 of file SigmaGenerator.cpp.

References endl(), and write_m.

Referenced by initialize(), match(), updateInitialSigma(), and updateSigma().

Here is the call graph for this function:

◆ writeOrbitOutput_m()

void SigmaGenerator::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 
)
private

Called within SigmaGenerator::match().

Parameters
r_turnis the radius [m]
peois the momentum
h_turnis the inverse bending radius
fidx_turnis the field index
ds_turnis the path length element

Definition at line 822 of file SigmaGenerator.cpp.

References Util::combineFilePath(), E_m, endl(), float2string(), OpalData::getAuxiliaryOutputDirectory(), OpalData::getInstance(), and write_m.

Referenced by match().

Here is the call graph for this function:

Member Data Documentation

◆ beta_m

double SigmaGenerator::beta_m
private

Velocity (c/v), \( \left[\beta\right] = 1 \).

Definition at line 165 of file SigmaGenerator.h.

Referenced by getEmittances(), initialize(), SigmaGenerator(), and updateInitialSigma().

◆ converged_m

bool SigmaGenerator::converged_m
private

Is true if converged, false otherwise.

Definition at line 173 of file SigmaGenerator.h.

Referenced by getIterations(), and match().

◆ delta_m

Series SigmaGenerator::delta_m
private

Definition at line 210 of file SigmaGenerator.h.

Referenced by SigmaGenerator().

◆ E_m

double SigmaGenerator::E_m
private

Stores the user-define energy, \( \left[E\right] = MeV \).

Definition at line 155 of file SigmaGenerator.h.

Referenced by initialize(), match(), updateInitialSigma(), updateSigma(), and writeOrbitOutput_m().

◆ Emax_m

double SigmaGenerator::Emax_m
private

Maximum energy reached in cyclotron, \( \left[E_{max}\right] = MeV \).

Definition at line 177 of file SigmaGenerator.h.

◆ Emin_m

double SigmaGenerator::Emin_m
private

Minimum energy needed in cyclotron, \( \left[E_{min}\right] = MeV \).

Definition at line 175 of file SigmaGenerator.h.

Referenced by SigmaGenerator().

◆ emittance_m

std::array<double,3> SigmaGenerator::emittance_m
private

Stores the desired emittances, \( \left[\varepsilon_{x}\right] = \left[\varepsilon_{y}\right] = \left[\varepsilon_{z}\right] = m \ rad \)

Definition at line 151 of file SigmaGenerator.h.

Referenced by getEmittances(), initialize(), SigmaGenerator(), and updateInitialSigma().

◆ error_m

double SigmaGenerator::error_m
private

Error of computation.

Definition at line 189 of file SigmaGenerator.h.

Referenced by getError(), and match().

◆ gamma2_m

double SigmaGenerator::gamma2_m
private

Relativistic factor squared.

Definition at line 161 of file SigmaGenerator.h.

Referenced by initialize(), and SigmaGenerator().

◆ gamma_m

double SigmaGenerator::gamma_m
private

Relativistic factor (which can be computed out ot the kinetic energy and rest mass (potential energy)), \( \left[\gamma\right] = 1 \)

Definition at line 159 of file SigmaGenerator.h.

Referenced by getEmittances(), initialize(), SigmaGenerator(), and updateInitialSigma().

◆ H_m

Hamiltonian SigmaGenerator::H_m
private

Stores the Hamiltonian of the cyclotron.

Definition at line 204 of file SigmaGenerator.h.

Referenced by match(), and SigmaGenerator().

◆ Hsc_m

SpaceCharge SigmaGenerator::Hsc_m
private

Stores the Hamiltonian for the space charge.

Definition at line 207 of file SigmaGenerator.h.

Referenced by match(), and SigmaGenerator().

◆ I_m

double SigmaGenerator::I_m
private

Stores the value of the current, \( \left[I\right] = A \).

Definition at line 147 of file SigmaGenerator.h.

Referenced by initialize(), and SigmaGenerator().

◆ l_m

Series SigmaGenerator::l_m
private

Definition at line 210 of file SigmaGenerator.h.

Referenced by SigmaGenerator().

◆ m_m

double SigmaGenerator::m_m
private

Is the mass of the particles, \( \left[m\right] = \frac{MeV}{c^{2}} \).

Definition at line 167 of file SigmaGenerator.h.

Referenced by initialize(), match(), and SigmaGenerator().

◆ N_m

unsigned int SigmaGenerator::N_m
private

Number of integration steps for closed orbit computation.

Definition at line 179 of file SigmaGenerator.h.

Referenced by match().

◆ nh_m

double SigmaGenerator::nh_m
private

Harmonic number, \( \left[N_{h}\right] = 1 \).

Definition at line 163 of file SigmaGenerator.h.

Referenced by initialize(), and SigmaGenerator().

◆ niterations_m

unsigned int SigmaGenerator::niterations_m
private

Is the number of iterations needed for convergence.

Definition at line 171 of file SigmaGenerator.h.

Referenced by getIterations(), match(), and updateSigma().

◆ nSectors_m

unsigned int SigmaGenerator::nSectors_m
private

Number of (symmetric) sectors the field is averaged over.

Definition at line 181 of file SigmaGenerator.h.

Referenced by match().

◆ nSteps_m

unsigned int SigmaGenerator::nSteps_m
private

Number of integration steps in total.

Definition at line 186 of file SigmaGenerator.h.

Referenced by initialize(), match(), and updateSigma().

◆ nStepsPerSector_m

unsigned int SigmaGenerator::nStepsPerSector_m
private

Number of integration steps per sector (--> also: number of maps per sector)

Definition at line 183 of file SigmaGenerator.h.

Referenced by match().

◆ prinit_m

double SigmaGenerator::prinit_m
private

Definition at line 213 of file SigmaGenerator.h.

Referenced by getInjectionMomentum(), and match().

◆ px_m

Series SigmaGenerator::px_m
private

Definition at line 210 of file SigmaGenerator.h.

Referenced by SigmaGenerator().

◆ pz_m

Series SigmaGenerator::pz_m
private

Definition at line 210 of file SigmaGenerator.h.

Referenced by SigmaGenerator().

◆ q_m

double SigmaGenerator::q_m
private

Is the particle charge [e].

Definition at line 169 of file SigmaGenerator.h.

Referenced by initialize(), match(), and SigmaGenerator().

◆ rdm_m

RealDiracMatrix SigmaGenerator::rdm_m
private

RealDiracMatrix-class member used for decoupling

Definition at line 273 of file SigmaGenerator.h.

Referenced by decouple().

◆ rinit_m

double SigmaGenerator::rinit_m
private

Definition at line 212 of file SigmaGenerator.h.

Referenced by getInjectionRadius(), and match().

◆ sigma_m

matrix_t SigmaGenerator::sigma_m
private

Stores the stationary distribution (sigma matrix)

Definition at line 198 of file SigmaGenerator.h.

Referenced by getSigma(), and match().

◆ sigmas_m

std::vector<matrix_t> SigmaGenerator::sigmas_m
private

Vector storing the second moment matrix for each angle.

Definition at line 201 of file SigmaGenerator.h.

Referenced by initialize(), match(), and updateSigma().

◆ truncOrder_m

unsigned int SigmaGenerator::truncOrder_m
private

Truncation order of the power series for the Hamiltonian (cyclotron and space charge)

Definition at line 192 of file SigmaGenerator.h.

Referenced by match(), and SigmaGenerator().

◆ wo_m

double SigmaGenerator::wo_m
private

Is the orbital frequency, \( \left[\omega_{o}\right] = \frac{1}{s} \).

Definition at line 153 of file SigmaGenerator.h.

Referenced by initialize(), and SigmaGenerator().

◆ write_m

bool SigmaGenerator::write_m
private

Decides for writing output (default: true)

Definition at line 195 of file SigmaGenerator.h.

Referenced by initialize(), match(), updateInitialSigma(), updateSigma(), writeMatrix(), and writeOrbitOutput_m().

◆ x_m

Series SigmaGenerator::x_m
private

All variables x, px, z, pz, l, delta.

Definition at line 210 of file SigmaGenerator.h.

Referenced by SigmaGenerator().

◆ z_m

Series SigmaGenerator::z_m
private

Definition at line 210 of file SigmaGenerator.h.

Referenced by SigmaGenerator().


The documentation for this class was generated from the following files: