OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
MapAnalyser.cpp
Go to the documentation of this file.
1 // ------------------------------------------------------------------------
2 // Copyright: see Copyright.readme
3 // ------------------------------------------------------------------------
4 //
5 // Class: MapAnalyser
6 // Organizes the function for a linear map analysis from
7 // ThickTracker.
8 // Transfer map -> tunes, symplecticity and stability
9 // Sigma Matrix -> (not projected) beam emittance
10 //
11 // ------------------------------------------------------------------------
12 //
13 // $Author: ganz $
14 //
15 // ------------------------------------------------------------------------
16 
17 
18 #include "MapAnalyser.h"
19 
20 #include <fstream>
21 
22 
23 #include "Physics/Physics.h"
24 
25 #include <gsl/gsl_math.h>
26 #include <gsl/gsl_eigen.h>
27 #include <gsl/gsl_linalg.h>
28 
29 
31  : mapAnalysis_m(IpplTimings::getTimer("mapAnalysis"))
32  , bunchAnalysis_m(IpplTimings::getTimer("bunchAnalysis"))
33 { }
34 
35 //Analyzes a TransferMap for the tunes, symplecticity and stability
36 void MapAnalyser::linTAnalyze(const fMatrix_t& tMatrix){
37  std::ofstream tplot;
38  tplot.open ("tunePlot.txt",std::ios::app);
39  tplot << std::setprecision(16);
40  //================================
41 
43 
44 
45  fMatrix_t blocktMatrix, scalingMatrix, invScalingMatrix;
46 
47  //get the EigenValues/-Vectors
48  cfMatrix_t eigenValM, eigenVecM, invEigenVecM;
49  eigenDecomp_m(tMatrix, eigenValM, eigenVecM, invEigenVecM);
50 
51  //normalize eigen Vectors
52  for (int i=0; i <6 ; i++){
53  double temp=0;
54 
55  for (int j=0; j <6 ; j+=2){
56  temp += 2*(eigenVecM[j][i] * std::conj(eigenVecM[j+1][i])).imag();
57  }
58  temp=std::fabs(temp);
59 
60  if (temp > 1e-10){
61  //TODO is it neccessary to normalize the eigenValues?
62  //eigenValM[i][i] *= std::sqrt(temp);
63  for (int j=0; j <6 ; j++){
64  eigenVecM[j][i] /= std::sqrt(temp);
65  invEigenVecM[j][i] /= std::sqrt(temp);
66  }
67  }
68  }
69 
70  //TODO pack transposition in a neat function
71  cfMatrix_t eigenVecMT;
72  for (int i=0; i <6 ; i++){
73  for (int j=0; j <6 ; j++){
74  eigenVecMT[i][j]= eigenVecM[j][i];
75  }
76  }
77 
78  double prec =0.01;
79 
80  int idx=0;
81  std::array<int,6> pairs;
82 
83  fMatrix_t skewMatrix = createSkewMatrix_m();
84  cfMatrix_t cSkewMatrix = complexTypeMatrix_m(skewMatrix);
85 
86  fMatrix_t K = imagPartOfMatrix_m(eigenVecMT *cSkewMatrix* eigenVecM);
87 
88  for (int i=0; i<6;i++) pairs[i]=-1;
89 
90  //Printer
91  for (int i=0; i <6 ; i++){
92  for (int j=0; j <6 ; j++){
93  if (K[i][j]>1-prec && K[i][j]<1+prec){
94  pairs[idx]=i;
95  pairs[idx+1]=j;
96  idx+=2;
97 
98  if (idx==6) break;
99  }
100  }
101  }
102 
103 
104  //Fill empty elements in pairs
106 
107  //Compare eigenvalues
108  int negidx=6-1;
109  for (int i=0; i<6; i++) {
110  pairsIt = std::find (pairs.begin(), pairs.end(), i);
111  if (pairsIt != pairs.end()){
112  continue;
113  }
114 
115  for (int j=0; i!=j && j <6 ; j++){
116  double diff = std::abs( eigenValM[i][i] - eigenValM[j][j]);
117  if (diff< 0.001){
118  pairs[negidx]=i;
119  pairs[negidx-1]=j;
120  negidx-=2;
121  }
122  }
123  }
124 
125  //Fill the paris vector
126  for (int i=0; i < 6 && idx < 6; i++){
127  pairsIt = std::find (pairs.begin(), pairs.end(), i);
128 
129  if (pairsIt != pairs.end()){
130  continue;
131  } else{
132  pairs[idx]=i;
133  idx++;
134  }
135  }
136 
137  cfMatrix_t tempM = eigenVecM ;
138  //cfMatrix_t tempInvM = eigenVecM ;
139  cfMatrix_t tempValM = eigenValM ;
140 
141 
142  for (int i=0; i <6 ; i++){
143  eigenValM[i][i]=tempValM[pairs[i]][pairs[i]];
144  for (int j=0; j <6 ; j++){
145  eigenVecM[j][i]=tempM[j][pairs[i]];
146  //eigenVecMT[i][j]=eigenVecM[j][i];
147  //invEigenVecM[i][j]=tempInvM[pairs[i]][j];
148  }
149  }
150 
151  invEigenVecM= invertMatrix_m(eigenVecM);
152  cfMatrix_t cblocktMatrix = getBlockDiagonal_m(tMatrix, eigenVecM, invEigenVecM);
153 
154 
155 
156  for (int i=0; i<6; i++){
157 
158  }
159 
160  FVector<std::complex<double>, 3> betaTunes, betaTunes2;
161  FVector<double, 3> betaTunes3;
162 
163 
164  //rearrangeEigen(eigenValM, eigenVecM);
165 
166  for (int i = 0; i < 3; i++){
167  betaTunes[i]=std::log(eigenValM[i*2][i*2])/ (2*Physics::pi * std::complex<double>(0, 1));
168  betaTunes[i]= std::asin(cblocktMatrix[i*2][i*2+1])/(std::complex<double>(2*Physics::pi, 0));
169 
170  betaTunes2[i]= std::acos(cblocktMatrix[i*2][i*2])/(std::complex<double>(2*Physics::pi, 0));
171  betaTunes2[i]= std::acos(cblocktMatrix[i*2][i*2].real())/(std::complex<double>(2*Physics::pi, 0));
172 
173  betaTunes3[i]= std::acos(eigenValM[i*2][i*2].real()) / (2*Physics::pi);
174 
175  double lenEigenVec = 0;
176  for (int j = 0; j < 3; j++){
177  lenEigenVec += std::abs(eigenVecM[i][j]);
178 
179  }
180  lenEigenVec = std::sqrt(lenEigenVec);
181  tplot<<"1: "<<betaTunes[i] <<std::endl;
182  tplot<<"2: "<<betaTunes2[i]<< std::endl;
183  //tplot<<"3: "<<lenEigenVec<< std::endl;
184  tplot<<"3: ("<<betaTunes3[i]<< ",0)"<< std::endl;
185 
187  //TODO: do something with the tunes etc
188  }
189 
190 }
191 
192 
193 //Analyses the Sigma Matrix
196  fMatrix_t sigmaBlockM;
197 
198  fMatrix_t skewMatrix;
199  for (int i=0; i <6 ; i=i+2){
200  skewMatrix[0+i][1+i]=1.;
201  skewMatrix[1+i][0+i]=-1.;
202  }
203 
204  sigMatrix=sigMatrix*skewMatrix;
205 
206  cfMatrix_t eigenValM, eigenVecM, invEigenVecM;
207 
208  eigenDecomp_m(sigMatrix,eigenValM, eigenVecM, invEigenVecM);
209 
210  cfMatrix_t cblocksigM = getBlockDiagonal_m(sigMatrix, eigenVecM,invEigenVecM);
211 
212 
213  for (int i=0; i<6; i++){
214  for (int j =0; j<6; j++){
215  sigmaBlockM[i][j]=cblocksigM[i][j].real();
216  }
217 
218  }
219 
220  cfMatrix_t teigenValM, teigenVecM, tinvEigenVecM;
221  eigenDecomp_m(sigmaBlockM, teigenValM, teigenVecM, tinvEigenVecM);
222 
223  for (int i=0; i<6; i++){
224 
225  }
226 
227 
228  //TODO: Where to go with EigenValues?
230 }
231 
232 
233 //Returns the eigenDecomositon for the fMatrix M = eigenVec * eigenVal * invEigenVec
234 void MapAnalyser::eigenDecomp_m(const fMatrix_t& M, cfMatrix_t& eigenVal, cfMatrix_t& eigenVec, cfMatrix_t& invEigenVec){
235 
236  double data[36];
237  int idx, s;
238  for (int i = 0; i < 6; i++) {
239  for (int j = 0; j < 6; j++) {
240  idx = i * 6 + j;
241  data[idx] = M[i][j];
242  }
243  }
244 
245 
246  gsl_matrix_view m = gsl_matrix_view_array(data, 6, 6);
247  gsl_vector_complex *eval = gsl_vector_complex_alloc(6);
248  gsl_matrix_complex *evec = gsl_matrix_complex_alloc(6, 6);
249  gsl_matrix_complex *eveci = gsl_matrix_complex_alloc(6, 6);
250  gsl_permutation * p = gsl_permutation_alloc(6);
251  gsl_eigen_nonsymmv_workspace * w = gsl_eigen_nonsymmv_alloc(6);
252 
253  //get Eigenvalues and Eigenvectors
254  gsl_eigen_nonsymmv(&m.matrix, eval, evec, w);
255  gsl_eigen_nonsymmv_free(w);
256  gsl_eigen_nonsymmv_sort(eval, evec, GSL_EIGEN_SORT_ABS_DESC);
257 
258  for (int i = 0; i < 6; ++i) {
259  eigenVal[i][i] = std::complex<double>(
260  GSL_REAL(gsl_vector_complex_get(eval, i)),
261  GSL_IMAG(gsl_vector_complex_get(eval, i)));
262  for (int j = 0; j < 6; ++j) {
263  eigenVec[i][j] = std::complex<double>(
264  GSL_REAL(gsl_matrix_complex_get(evec, i, j)),
265  GSL_IMAG(gsl_matrix_complex_get(evec, i, j)));
266  }
267  }
268 
269  //invert Eigenvectormatrix
270  gsl_linalg_complex_LU_decomp(evec, p, &s);
271  gsl_linalg_complex_LU_invert(evec, p, eveci);
272 
273  //Create invEigenVecMatrix
274  for (int i = 0; i < 6; ++i) {
275  for (int j = 0; j < 6; ++j) {
276  invEigenVec[i][j] = std::complex<double>(
277  GSL_REAL(gsl_matrix_complex_get(eveci, i, j)),
278  GSL_IMAG(gsl_matrix_complex_get(eveci, i, j)));
279  }
280  }
281 
282 
283  //free space
284  gsl_vector_complex_free(eval);
285  gsl_matrix_complex_free(evec);
286  gsl_matrix_complex_free(eveci);
287 
288 
289 }
290 
291 
292 //Transforms the Matirx to a block diagonal rotation Matrix
294  cfMatrix_t& eigenVecM, cfMatrix_t& invEigenVecM){
295 
296  cfMatrix_t cM, qMatrix, invqMatrix, nMatrix, invnMatrix, rMatrix;
297 
298  for (int i=0; i<6; i++){
299  for (int j =0; j<6; j++){
300  cM[i][j]=std::complex<double>(M[i][j],0);
301  }
302  }
303 
304  for (int i=0; i <6 ; i=i+2){
305  qMatrix[0+i][0+i]=std::complex<double>(1.,0);
306  qMatrix[0+i][1+i]=std::complex<double>(0,1.);
307  qMatrix[1+i][0+i]=std::complex<double>(1.,0);
308  qMatrix[1+i][1+i]=std::complex<double>(0,-1);
309 
310  invqMatrix[0+i][0+i]=std::complex<double>(1.,0);
311  invqMatrix[0+i][1+i]=std::complex<double>(1.,0);
312  invqMatrix[1+i][0+i]=std::complex<double>(0.,-1.);
313  invqMatrix[1+i][1+i]=std::complex<double>(0,1.);
314  }
315  qMatrix/=std::sqrt(2.);
316  invqMatrix/=std::sqrt(2);
317 
318  nMatrix=eigenVecM*qMatrix;
319  invnMatrix= invqMatrix* invEigenVecM;
320 
321 
322  rMatrix= invnMatrix * cM * nMatrix;
323 
324  return rMatrix;
325 
326 }
327 
328 
330  cfMatrix_t N1, cinvN, cR, ctM, N2;
331  fMatrix_t R1, S, sigmaS;
332 
333  for (int i=0; i<6; i++){
334  for (int j =0; j<6; j++){
335  ctM[i][j]=std::complex<double>(tM[i][j],0);
336  }
337  }
338 
339  S = createSkewMatrix_m();
340  sigmaS = Sigma*S;
341 
342  setNMatrix_m(sigmaS, N2, cinvN);
343 
344  std::array<double, 3> phi;
345 
346  for (int i = 0; i < 3; i++){
347  phi[i] = std::atan(oldN[2*i+1][i].real()/oldN[2*i][2*i].real());
348  }
349 
350 
351  R1=createRotMatrix_m(phi);
352 
353  for (int i=0; i<6; i++){
354  for (int j =0; j<6; j++){
355  cR[i][j]=std::complex<double>(R1[i][j],0);
356  N1[i][j]=oldN[i][j].real();
357  }
358  }
359 }
360 
361 
363 
364  cfMatrix_t eigenValM, eigenVecM, invEigenVecM, eigenVecMT;
365 
366 
367  eigenDecomp_m(M, eigenValM, eigenVecM, invEigenVecM);
368 
369  cfMatrix_t cM, qMatrix, invqMatrix, nMatrix, invnMatrix, rMatrix;
370 
371  //std::ofstream tmap;
372  //tmap.open ("TransferMap.txt",std::ios::app);
373  //tmap << std::setprecision(16);
374 
375 
376  for (int i=0; i<6; i++){
377  for (int j =0; j<6; j++){
378  cM[i][j]=std::complex<double>(M[i][j],0);
379  }
380  }
381 
382  for (int i=0; i <6 ; i=i+2){
383  qMatrix[0+i][0+i]=std::complex<double>(1.,0);
384  qMatrix[0+i][1+i]=std::complex<double>(0,1.);
385  qMatrix[1+i][0+i]=std::complex<double>(1.,0);
386  qMatrix[1+i][1+i]=std::complex<double>(0,-1);
387 
388  invqMatrix[0+i][0+i]=std::complex<double>(1.,0);
389  invqMatrix[0+i][1+i]=std::complex<double>(1.,0);
390  invqMatrix[1+i][0+i]=std::complex<double>(0.,-1.);
391  invqMatrix[1+i][1+i]=std::complex<double>(0,1.);
392  }
393  qMatrix/=std::sqrt(2.);
394  invqMatrix/=std::sqrt(2);
395 
396 
397  N=eigenVecM*qMatrix;
398  invN= invqMatrix* invEigenVecM;
399 }
400 
401 
403  fMatrix_t R;
404 
405  for (int i = 0; i < 3; i++){
406  R[2*i][2*i] = std::cos(phi[1]);
407  R[2*i+1][2*i+1] = R[2*i][2*i];
408  R[2*i][2*i+1] = std::sin(phi[1]);
409  R[2*i+1][2*i] = -R[2*i][2*i+1];
410  }
411  return R;
412 }
413 
414 
416  fMatrix_t S;
417 
418  for (int i = 0; i < 3; i++){
419  S[2*i][2*i+1] = 1;
420  S[2*i+1][2*i] = -1;
421  }
422  return S;
423 }
424 
425 
427  fMatrix_t M;
428  for (int i=0; i<6; i++){
429  for (int j =0; j<6; j++){
430  M[i][j]=cM[i][j].real();
431  }
432  }
433  return M;
434 }
435 
437  fMatrix_t M;
438  for (int i=0; i<6; i++){
439  for (int j =0; j<6; j++){
440  M[i][j]=cM[i][j].imag();
441  }
442  }
443  return M;
444 }
445 
447  cfMatrix_t cM;
448  for (int i=0; i<6; i++){
449  for (int j =0; j<6; j++){
450  cM[i][j]=std::complex<double>(M[i][j],0);
451  }
452  }
453  return cM;
454 }
455 
456 
458 
459  gsl_set_error_handler_off();
460  //gsl_vector_complex *m = gsl_vector_complex_alloc(6);
461  gsl_matrix_complex *m = gsl_matrix_complex_alloc(6, 6);
462  gsl_matrix_complex *invm = gsl_matrix_complex_alloc(6, 6);
463  gsl_permutation * p = gsl_permutation_alloc(6);
464  gsl_complex temp;
465  int s;
466 
467 
468  //Create invEigenVecMatrix
469  for (int i = 0; i < 6; ++i) {
470  for (int j = 0; j < 6; ++j) {
471  GSL_SET_COMPLEX(&temp,std::real(M[i][j]),std::imag(M[i][j]));
472  gsl_matrix_complex_set( m, i, j, temp);
473  }
474  }
475 
476  //invert Eigenvectormatrix
477  int eigenDecompStatus = gsl_linalg_complex_LU_decomp(m, p, &s);
478  if (eigenDecompStatus != 0){
479  std::cout<< "This actually works" << std::endl;
480  //gsl_set_error_handler (NULL);
481 
482  }
483 
484  int invertStatus = gsl_linalg_complex_LU_invert(m, p, invm);
485 
486  if ( invertStatus ) {
487  std::cout << "Error" << std::endl;
488  std::exit(1);
489  }
490 
491 
492  if (invertStatus != 0){
493  std::cout<< "This actually works" << std::endl;
494  //gsl_set_error_handler (NULL);
495 
496  }
497 
498  cfMatrix_t invM;
499  //Create invEigenVecMatrix
500  for (int i = 0; i < 6; ++i) {
501  for (int j = 0; j < 6; ++j) {
502  invM[i][j] = std::complex<double>(
503  GSL_REAL(gsl_matrix_complex_get(invm, i, j)),
504  GSL_IMAG(gsl_matrix_complex_get(invm, i, j)));
505  }
506  }
507 
508  //free space
509  gsl_matrix_complex_free(m);
510  gsl_matrix_complex_free(invm);
511  gsl_permutation_free(p);
512 
513 
514  return invM;
515 
516 }
517 
518 
519 //TODO Work in progress
521 /*#ifdef PHIL_WRITE
522  std::ofstream out;
523  out.open("OUT.txt", std::ios::app);
524 #endif
525  double precision = 1e-9;
526  int pair = 0;
527  std::complex<double> mult;
528  std::pair<int, int> pairs[3];
529  cfMatrix_t eigVa, eigVe, temp;
530 
531  out << "+++++++++++++++++++++++++++++++" << std::endl;
532  //out << EigenVec << std::endl;
533  for (int i = 0; pair < 3 && i <2*DIM; i++){
534  for (int j = 0; pair < 3 && j < i; j++){
535  mult=0;
536  if (true){//std::abs(EigenVal[i][i].real() - EigenVal[j][j].real()) < precision
537  //&& std::abs(EigenVal[i][i].imag() + EigenVal[j][j].imag() ) < precision ){
538  out << eigenVal[i][i] << " " << eigenVal[j][j] << std::endl;
539  out << "/////////////Calc///////////" << std::endl;
540  for (int k = 0; k <DIM; k++){
541 
542  out << -1. * EigenVec[2*k+1][i] << EigenVec[2*k][j] << std::endl;
543  out << EigenVec[2*k][i] << EigenVec[2*k+1][j] << std::endl;
544 
545 
546  mult += -1. * EigenVec[2*k+1][i] * EigenVec[2*k][j];
547  mult += EigenVec[2*k][i] * EigenVec[2*k+1][j];
548  out <<"in between mult: " << mult << std::endl;
549  }
550 
551  out << mult << std::endl;
552  out << "Criteria " << std::abs(mult.imag()-1.) << std::endl;
553  if (std::abs(mult.imag()-1.) < precision){
554  pairs[pair]= std::pair<int,int>(i,j);
555  out << i << " " << j << std::endl;
556  pair++;
557  }
558 
559  }
560 
561  }
562 
563  }*/
564 }
IpplTimings::TimerRef bunchAnalysis_m
Definition: MapAnalyser.h:116
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:808
constexpr double e
The value of .
Definition: Physics.h:40
IpplTimings::TimerRef mapAnalysis_m
Definition: MapAnalyser.h:115
void setNMatrix_m(fMatrix_t &M, cfMatrix_t &N, cfMatrix_t &invN)
A templated representation for vectors.
Definition: PartBunchBase.h:26
PETE_TUTree< FnFabs, typename T::PETE_Expr_t > fabs(const PETE_Expr< T > &l)
Definition: PETE.h:815
Tps< T > sin(const Tps< T > &x)
Sine.
Definition: TpsMath.h:111
fMatrix_t realPartOfMatrix_m(cfMatrix_t cM)
fMatrix_t imagPartOfMatrix_m(cfMatrix_t cM)
Tps< T > log(const Tps< T > &x)
Natural logarithm.
Definition: TpsMath.h:182
FLieGenerator< T, N > real(const FLieGenerator< std::complex< T >, N > &)
Take real part of a complex generator.
fMatrix_t createSkewMatrix_m()
void rearrangeEigen_m(cfMatrix_t &EigenVal, cfMatrix_t &EigenVec)
cfMatrix_t complexTypeMatrix_m(fMatrix_t M)
void linSigAnalyze(fMatrix_t &sigMatrix)
constexpr double pi
The value of .
Definition: Physics.h:31
m_complex conj(const m_complex &c)
Definition: MVector.h:105
void printPhaseShift_m(fMatrix_t &Sigma, fMatrix_t tM, cfMatrix_t &oldN)
PETE_TUTree< FnArcTan, typename T::PETE_Expr_t > atan(const PETE_Expr< T > &l)
Definition: PETE.h:810
static void startTimer(TimerRef t)
Definition: IpplTimings.h:187
void eigenDecomp_m(const fMatrix_t &M, cfMatrix_t &eigenVal, cfMatrix_t &eigenVec, cfMatrix_t &invEigenVec)
fMatrix_t createRotMatrix_m(std::array< double, 3 > phi)
const T * find(const T table[], const std::string &name)
Look up name.
Definition: TFind.h:34
Tps< T > sqrt(const Tps< T > &x)
Square root.
Definition: TpsMath.h:91
cfMatrix_t invertMatrix_m(const cfMatrix_t &M)
Tps< T > cos(const Tps< T > &x)
Cosine.
Definition: TpsMath.h:129
cfMatrix_t getBlockDiagonal_m(const fMatrix_t &M, cfMatrix_t &eigenVecM, cfMatrix_t &invEigenVecM)
PETE_TUTree< FnArcSin, typename T::PETE_Expr_t > asin(const PETE_Expr< T > &l)
Definition: PETE.h:809
void linTAnalyze(const fMatrix_t &tMatrix)
Definition: MapAnalyser.cpp:36
std::string::iterator iterator
Definition: MSLang.h:16
static void stopTimer(TimerRef t)
Definition: IpplTimings.h:192
#define K
Definition: integrate.cpp:118
FLieGenerator< T, N > imag(const FLieGenerator< std::complex< T >, N > &)
Take imaginary part of a complex generator.
Inform & endl(Inform &inf)
Definition: Inform.cpp:42