This section contains a description of routines which were changed since TRACY's translation from PASCAL to C. For a description of the routines which were left untouched refer to the Tracy-2 manual by J. Bengtsson [2].
gemit in physlib.h
Void gemit();
Calculates radiation integrals (emittances), damping times, fractional tunes and generalized sigma matrices ala CHAO. Requires DA mode and lattice with cavity to run. The global 3D array globval.ElemMat contains the 6x6 sigma matrices for all positions. See also printsigma.
gcmat in lsoc.h
Void gcmat(long bpm, long bpmdis[mnp], long corr, long corrdis[mnp], long plane,
double (*A)[mnp], double (*U)[mnp], double (*V)[mnp], double *w,
long wdis[mnp], double (*Ai)[mnp])
Determine the corrector-bpm correlation matrix A:
-----------
\/beta beta
i j
A = ------------- cos(nu pi - 2 pi|nu - nu |)
ij 2 sin(pi nu) i j
and perform Singular Value Decomposition (SVD) of A for plane plane with A^-1 = Ai = V*w*U^T. Please refer to gcmatprint for the description of the remaining parameters.
lsoc in physlib.h
Static Void lsoc(long niter, long bpm, long bpmdis[mnp], long corr, long corrdis[mnp],
double maxkick, long maxbits, long plane,
(*U)[mnp], (*V)[mnp], double *w);
Perform an orbit correction based on the SVD algorithm. Iterate niter times for plane plane. Please refer to gcmatprint for the description of the remaining parameters. See
Numerical Recipes chapter 2.6 for an introduction to the SVD algorithm.
GirSet in bcosys.h
Void GirSet(double gfrac, double jfrac, double efrac, double bfrac);
The rms girder alignment error defined in input file "bcosys.dat" ist multiplied with bfrac. The rms girder joints error ist multiplied with jfrac and the rms element error ist multiplied with efrac. If bfrac is != 0
GirSet simulates a simple beam-based alignment (BBA) procedure of qudrupoles with respect to adjacent monitors (method=1) or sextupoles with respect to adjacent monitors (method=2). It simply correlates monitors-quadrupole (sextupole) combinations assuming a certain rms error of the aligment procedure. bfrac=0.1 would correspond to a BBA error which is 10 times smaller than the rms element error. The BBA settings are dumped to the file "bbaset.dat". The variable method is a local variable to GirSet. Angle errors of elements are now included for girders and elements. The rms angle is defined in "bcosys.dat" in micro degrees !-)
This section contains a description of routines which were added since TRACY's translation from PASCAL to C. For a description of the routines which were left untouched refer to the Tracy-2 manual by J. Bengtsson [2].
getbumprec4a in t2bump.c
Void getbumprec4a (long ncorr, long *corr, long plane,
long corr1, long corr2, long corr3, long corr4, bumprec *bump);
calculate asymmetric closed orbit bump for 4 correctors corr1-4 and store kick ratios in bump
getbumprec4s in t2bump.c
Void getbumprec4s (long ncorr, long *corr, long plane,
long corr1, long corr2, long corr3, long corr4, bumprec *bump);
calculate symmetric closed orbit bump for 4 correctors corr1-4 and store kick ratios in bump
SetUpBump4a in t2bump.c
Void SetUpBump4a (long ncorr, long *corr, double dnumin, long plane);
initialize ncorr asymmetric closed orbit bumps utilizing 4 adjacent coils around the machine with a minimal phase difference larger than dnumin for the plane plane.corr contains element numbers of the correction coils (See inibump in physlib.h for details).
SetUpBump4s in t2bump.c
Void SetUpBump4s (long ncorr, long *corr, double dnumin, long plane);
initialize ncorr symmetric closed orbit bumps utilizing 4 adjacent coils around the machine with a minimal phase difference larger than dnumin for the plane plane.corr contains element numbers of the correction coils (See inibump in physlib.h for details).
InitBUMP4a in t2bump.c
Void InitBUMP4a (long *ncorr, long *hcorr, long *vcorr, double dnuhmin, double dnuvmin);
initialize ncorr[0] horizontal and ncorr[1] vertical asymmetric closed orbit bumps utilizing 4 adjacent coils around the machine with a minimal phase difference larger than dnuhmin and dnuvmin for the horizontal and the
vertical plane respectively. hcorr and vcorr contain element numbers
of the correction coils (See inibump in physlib.h for details).
InitBUMP4s in t2bump.c
Void InitBUMP4s (long *ncorr, long *hcorr, long *vcorr, double dnuhmin, double dnuvmin);
initialize ncorr[0] horizontal and ncorr[1] vertical symmetric closed orbit bumps utilizing 4 adjacent coils around the machine with a minimal phase difference larger than dnuhmin and dnuvmin for the horizontal and the
vertical plane respectively. hcorr and vcorr contain element numbers
of the correction coils (See inibump in physlib.h for details).
EigenVal in fit.c
Void EigenVal (double (*Ai)[mnp], long n, double *wr, double *wi);
calculate real wr and imaginary wi part of the eigenvalues for the
matrix Ai. Uses routines from the Numerical Recipes library.
QuadFit in fit.c
Void QuadFit (double *X, double *Y, long ndata, double *A, double Covar[][3],
double *Chisq);
perform a least square quadratic fit to (X,Y) using ndata samples. A contains the coefficients of the fit and Chisq the 3x3 covariance matrix. Uses routines from the Numerical Recipes library.
NormEigenVec in eigenv.c
Void NormEigenVec (vector *Vr, vector *Vi, double *wr, double *wi, vector *t6a);
sort and normalize complex eigenvectors (Vr,Vi) with eigenvalues (wr,wi) and store the result in t6a. Used for the calculation of
generalized sigma matrices ala CHAO in emit.
ElemIndex in t2lat.c
long ElemIndex(Char *name);
return element family index. Note: in the PASCAL version the element family index could be comfortably accessed using the element name. This is no longer possible because we gave up on the interpretive PASCAL-S.
Mpole_GetdT in t2elem.c
double Mpole_GetdT(long Fnum1, long Knum1);
return total roll angle of the element Cell[ElemFam[Fnum1 - 1].KidList[Knum1 - 1]].Elem which is a sum of a design ,a systematic error and a random error part.
Mpole_DefdTpar in t2elem.c
double Mpole_DefdTpar(long Fnum1, long Knum1, double PdTpar);
Set design roll angle to PdTpar degrees.
Mpole_DefdTsys in t2elem.c
double Mpole_DefdTsys(long Fnum1, long Knum1, double PdTsys);
Set systematic roll angle error to PdTsys degrees.
fft_double in fourierd.c
void fft_double
( unsigned NumSamples, /* must be a power of 2 */
int InverseTransform, /* 0=forward FFT, 1=inverse FFT */
double *RealIn, /* array of input's real samples */
double *ImaginaryIn, /* array of input's imag samples */
double *RealOut, /* array of output's reals */
double *ImaginaryOut ); /* array of output's imaginaries */
performs a onedimensional fourier transform. Is an alternative to FFT.
Circumference in physlib.h
double Circumference();
returns the length of the ring.
printsigma() in physlib.h
Void printsigma();
dump beam ellipse sigmas and twists calculated from generalized sigma matrices ala CHAO at every element to the file "sigma.dat".
digitize in physlib.h
double digitize(double x, double maxkick, double maxsamp);
map x onto the integer interval (-maxsamp ... maxsamp) where maxsamp corresponds maxkick.
Dis_In in physlib.h
Void Dis_In(long *bpmdis, long *vcorrdis, long *hcorrdis, long *wvdis, long *whdis)
Read a number of flag arrays from a plain text file named "
dis.dat". The first line is always treated as a comment line. The first element on the following lines specifies the number of flags of type integer following on the same line. A space separated comment at the end of each line is allowed. bpmdis and v(h)corrdis are flag arrays marking bad bpm's and vertical (horizontal) correctors (0=ok, 1=faulty). wv(h)dis defines which SVD w values belonging to the vertical (horizontal) motion are explicitly set to zero (0=keep value, 1=set value to zero) in order to reduce the influence of bpm noise.
gcmatprint in lsoc.h
Void gcmatprint(long bpm, long bpmdis[mnp], long corr, long corrdis[mnp], long plane,
double (*A)[mnp], double (*U)[mnp], double (*V)[mnp], double *w,
long wdis[mnp], double (*Ai)[mnp])
dump result of Singular Value Decomposition (SVD) of the corrector-bpm correlation matrix A for plane plane performed by dsvdc with A^-1 = Ai = V*w*U^T to the file "svd.dat". bpm and corr are element family indices (see ElemIndex) of bpm's and correctors. bpmdis and corrdis are flag arrays marking bad bpm's and correctors (0=ok, 1=faulty). wdis defines which w values are set to zero (0=keep value, 1=set value to zero) in order to reduce the influence of bpm noise.
minsq in cernlib.h
Void minsq(long mfun, long nvar, double *fm, double *xn, double *en, long iprint,
long niter, long cov, double xstep)
Interface to the cernlib routine MINSQ which minimizes a sum of squares of functions:
mfun = number of functions
nvar = number of variables
fm = one-dim array for f1,f2,...,fm
xn = one-dim array for x1,x2,...,xn
en = starting incremint of xs
iprint = print after iprint iterations 0=final print only
ninter = number of iterations
cov = =1 approx. error matrix printed, =0 otherwise
xstep = initial step size searching for minimum along the line
xstep <=0.5 recommended
A one-dim dummy array wf[NW] with NW=200 is declared within minsq.
Make sure that nvar+mfun(nvar+1)+3nvar(nvar+1)/2 does not exceed
NW=200. The user has to define the following routine:
Void fcn(long mfun, long nvar, double *fm, double *xn, long iflag);
iflag is =1 for first entry, =4 for normal entry and =3 for the final entry after the minimum has been found. A maximum of MFUN=10 functions of NVAR=10 variables is supported. For detailed information refer to the cernlib
writeups.
minvar in cernlib.h
Void minvar(double *x,double *y, double *r, double *eps, double *step, long maxf,
double a, double b, char *f)
Interface to the cernlib routine MINVAR which calculates to an attempted specified accuracy the abscissa of a local minimum of real-valued function f(x) lying in a given interval a,b together with the value of the function at the minimum:
x = estimate of the abscissa on exit estimate
y = on exit f(x)
r = on exit |x-x0|<r x0 x0=exact abscissa of a minimum of f
eps = accuracy |x-x0|/(1+|x|)<=eps
step = initial search step
a,b = search interval
f = function f(x) to be minimized
The user has to define a function f and parse a pointer to it to minvar:
float f(float *x, long *iflag)
iflag is =0 for first entry, and =1 for subsequent entries. For detailed information refer to the cernlib writeups.
SkewCorr in corr.h
Void SkewCorr(long niter, double *xn, double *fm, double *rn, double eps, double step,
long maxf, double a, double b)
Minimizes the emittance coupling based on minvar utilizing 3 families of skew quadrupoles in the long straight sections of the SLS (dsl11[p,m], dsl12[p,m], dsl13[p,m]). Starting skew values may be assigned before SkewCorr is invoked. The families are assumed to be orthogonal. In case of a small nonorthogonality the minimzation can be repeated niter times in order find the minimum:
niter = number of iterations
xn[MFUN] = estimate of the abscissa on exit estimate
fm[MFUN] = on exit f(x)
rn[MFUN] = on exit |x-x0|<r x0 x0=exact abscissa of a minimum of f
eps = accuracy |x-x0|/(1+|x|)<=eps
step = initial search step
a,b = search interval