ROOT logo
////////////////////////////////////////////////////////////////////////////////
//                                                                            //
// AliFemtoYlm - the class to calculate varous components of spherical        //
//  harmonics                                                                 //
//                                                                            //
// Authors: Adam Kisiel kisiel@mps.ohio-state.edu                             //
//                                                                            //
////////////////////////////////////////////////////////////////////////////////

#include "AliFemtoYlm.h"
#include <TMath.h>
#include <iostream>

double *AliFemtoYlm::fgPrefactors = 0x0;
int    *AliFemtoYlm::fgPrefshift = 0x0;
int    *AliFemtoYlm::fgPlmshift = 0x0;

AliFemtoYlm::AliFemtoYlm() {
  InitializeYlms();
}

AliFemtoYlm::~AliFemtoYlm() {
  free(fgPrefactors);
  free(fgPrefshift);
  free(fgPlmshift);
}


AliFemtoYlm::AliFemtoYlm(const AliFemtoYlm& aYlm){
  fgPrefshift = aYlm.fgPrefshift;
  InitializeYlms();
}

AliFemtoYlm& AliFemtoYlm::operator=(const AliFemtoYlm& aYlm){
  if (this == &aYlm)
    return *this;

  InitializeYlms();

  return *this;
}

std::complex<double> AliFemtoYlm::Ceiphi(double phi){
  return std::complex<double>(cos(phi),sin(phi));
}

double AliFemtoYlm::Legendre(int ell, int em, double ctheta){
  // Calculate a single Legendre value
  // *** Warning - NOT optimal - calculated all Plms up to L !!!
  double lbuf[36];
  AliFemtoYlm::LegendreUpToYlm(ell, ctheta, lbuf);

  return lbuf[fgPlmshift[ell]-abs(em)];
}

std::complex<double> AliFemtoYlm::Ylm(int ell,int m,double theta,double phi){
  // Calculate Ylm spherical input
  double ctheta;
  std::complex<double> answer;
  std::complex<double> ci(0.0,1.0);
  ctheta=cos(theta);
  answer = (fgPrefactors[fgPrefshift[ell]+m]*Legendre(ell,m,ctheta))*Ceiphi(m*phi);

  return answer;
}

std::complex<double> AliFemtoYlm::Ylm(int ell, int m, double x, double y, double z){
  // Calculate Ylm cartesian input
  std::complex<double> answer; 
  double ctheta,phi;
  double r = sqrt(x*x+y*y+z*z);
  if ( r < 1e-10 || fabs(z) < 1e-10 ) ctheta = 0.0;
  else ctheta=z/r;
  phi=atan2(y,x);
  answer = (fgPrefactors[fgPrefshift[ell]+m]*Legendre(ell,m,ctheta))*Ceiphi(m*phi);

  return answer;	
}

double AliFemtoYlm::ReYlm(int ell, int m, double theta, double phi){
	return real(AliFemtoYlm::Ylm(ell,m,theta,phi));
}

double AliFemtoYlm::ImYlm(int ell, int m, double theta, double phi){
	return imag(AliFemtoYlm::Ylm(ell,m,theta,phi));
}

double AliFemtoYlm::ReYlm(int ell, int m, double x,double y,double z){
	return real(AliFemtoYlm::Ylm(ell,m,x,y,z));
}

double AliFemtoYlm::ImYlm(int ell, int m, double x,double y,double z){
	return imag(AliFemtoYlm::Ylm(ell,m,x,y,z));
}

void AliFemtoYlm::InitializeYlms()
{
  // Calculate prefactors for fast Ylm calculation

  double oneoversqrtpi = 1.0/TMath::Sqrt(TMath::Pi());

  if(fgPrefactors!=NULL)
    free(fgPrefactors);
  if(fgPrefshift!=NULL)
    free(fgPrefshift);
  if(fgPlmshift!=NULL)
    free(fgPlmshift);

  fgPrefactors = (double *) malloc(sizeof(double) * 36);
  fgPrefshift  = (int *) malloc(sizeof(int) * 6);
  fgPlmshift   = (int *) malloc(sizeof(int) * 6);

  // l=0 prefactors
  fgPrefactors[0] = 0.5*oneoversqrtpi;

  // l=1 prefactors
  fgPrefactors[1] = 0.5*sqrt(3.0/2.0)*oneoversqrtpi;
  fgPrefactors[2] = 0.5*sqrt(3.0)*oneoversqrtpi;
  fgPrefactors[3] = -fgPrefactors[1];

  // l=2 prefactors
  fgPrefactors[4] = 0.25*sqrt(15.0/2.0)*oneoversqrtpi;
  fgPrefactors[5] = 0.5*sqrt(15.0/2.0)*oneoversqrtpi;
  fgPrefactors[6] = 0.25*sqrt(5.0)*oneoversqrtpi;
  fgPrefactors[7] = -fgPrefactors[5];
  fgPrefactors[8] = fgPrefactors[4];

  // l=3 prefactors
  fgPrefactors[9]  = 0.125*sqrt(35.0)*oneoversqrtpi;
  fgPrefactors[10] = 0.25*sqrt(105.0/2.0)*oneoversqrtpi;
  fgPrefactors[11] = 0.125*sqrt(21.0)*oneoversqrtpi;
  fgPrefactors[12] = 0.25*sqrt(7.0)*oneoversqrtpi;
  fgPrefactors[13] = -fgPrefactors[11];
  fgPrefactors[14] = fgPrefactors[10];
  fgPrefactors[15] = -fgPrefactors[9];

  // l=4 prefactors
  fgPrefactors[16] = 3.0/16.0*sqrt(35.0/2.0)*oneoversqrtpi;
  fgPrefactors[17] = 3.0/8.0*sqrt(35.0)*oneoversqrtpi;
  fgPrefactors[18] = 3.0/8.0*sqrt(5.0/2.0)*oneoversqrtpi;
  fgPrefactors[19] = 3.0/8.0*sqrt(5.0)*oneoversqrtpi;
  fgPrefactors[20] = 3.0/16.0*oneoversqrtpi;
  fgPrefactors[21] = -fgPrefactors[19];
  fgPrefactors[22] = fgPrefactors[18];
  fgPrefactors[23] = -fgPrefactors[17];
  fgPrefactors[24] = fgPrefactors[16];

  // l=5 prefactors
  fgPrefactors[25] = 3.0/32.0*sqrt(77.0)*oneoversqrtpi;
  fgPrefactors[26] = 3.0/16.0*sqrt(385.0/2.0)*oneoversqrtpi;
  fgPrefactors[27] = 1.0/32.0*sqrt(385.0)*oneoversqrtpi;
  fgPrefactors[28] = 1.0/8.0*sqrt(1155.0/2.0)*oneoversqrtpi;
  fgPrefactors[29] = 1.0/16.0*sqrt(165.0/2.0)*oneoversqrtpi;
  fgPrefactors[30] = 1.0/16.0*sqrt(11.0)*oneoversqrtpi;
  fgPrefactors[31] = -fgPrefactors[29];
  fgPrefactors[32] = fgPrefactors[28];
  fgPrefactors[33] = -fgPrefactors[27];
  fgPrefactors[34] = fgPrefactors[26];
  fgPrefactors[35] = -fgPrefactors[25];

  fgPrefshift[0] = 0;
  fgPrefshift[1] = 2;
  fgPrefshift[2] = 6;
  fgPrefshift[3] = 12;
  fgPrefshift[4] = 20;
  fgPrefshift[5] = 30;

  fgPlmshift[0] = 0;
  fgPlmshift[1] = 2;
  fgPlmshift[2] = 5;
  fgPlmshift[3] = 9;
  fgPlmshift[4] = 14;
  fgPlmshift[5] = 20;
}

void AliFemtoYlm::LegendreUpToYlm(int lmax, double ctheta, double *lbuf)
{
  // Calculate a set of legendre polynomials up to a given l
  // with spherical input
  double sins[6];
  double coss[6];
  sins[0] = 0.0;
  coss[0] = 1.0;
  sins[1] = sqrt(1-ctheta*ctheta);
  coss[1] = ctheta;
  for (int iter=2; iter<6; iter++) {
    sins[iter] = sins[iter-1]*sins[1];
    coss[iter] = coss[iter-1]*coss[1];
  }

  // Legendre polynomials l=0
  lbuf[0] = 1.0;

  // Legendre polynomials l=1
  if (lmax>0) {
    lbuf[1] = sins[1];
    lbuf[2] = coss[1];
  }

  // Legendre polynomials l=2
  if (lmax>1) {
    lbuf[3] = sins[2];
    lbuf[4] = sins[1]*coss[1];
    lbuf[5] = 3*coss[2]-1;
  }

  // Legendre polynomials l=3
  if (lmax>2) {
    lbuf[6] = sins[3];
    lbuf[7] = sins[2]*coss[1];
    lbuf[8] = (5*coss[2]-1)*sins[1];
    lbuf[9] = 5*coss[3]-3*coss[1];
  }

  // Legendre polynomials l=4
  if (lmax>3) {
    lbuf[10] = sins[4];
    lbuf[11] = sins[3]*coss[1];
    lbuf[12] = (7*coss[2]-1)*sins[2];
    lbuf[13] = (7*coss[3]-3*coss[1])*sins[1];
    lbuf[14] = 35*coss[4]-30*coss[2]+3;
  }

  // Legendre polynomials l=5
  if (lmax>4) {
    lbuf[15] = sins[5];
    lbuf[16] = sins[4]*coss[1];
    lbuf[17] = (9*coss[2]-1)*sins[3];
    lbuf[18] = (3*coss[3]-1*coss[1])*sins[2];
    lbuf[19] = (21*coss[4]-14*coss[2]+1)*sins[1];
    lbuf[20] = 63*coss[5]-70*coss[3]+15*coss[1];
  }
}

void AliFemtoYlm::YlmUpToL(int lmax, double x, double y, double z, std::complex<double> *ylms)
{
  // Calculate a set of Ylms up to a given l
  // with cartesian input
  double ctheta,phi;

  double r = sqrt(x*x+y*y+z*z);
  if ( r < 1e-10 || fabs(z) < 1e-10 ) ctheta = 0.0;
  else ctheta=z/r;
  phi=atan2(y,x);
  
  YlmUpToL(lmax, ctheta, phi, ylms);

}

void AliFemtoYlm::YlmUpToL(int lmax, double ctheta, double phi, std::complex<double> *ylms)
{
  // Calculate a set of Ylms up to a given l
  // with spherical input
  int lcur = 0;  
  double lpol;
  
  double coss[6];
  double sins[6];

  double lbuf[36];
  LegendreUpToYlm(lmax, ctheta, lbuf);

  for (int iter=1; iter<=lmax; iter++) {
    coss[iter-1] = cos(iter*phi);
    sins[iter-1] = sin(iter*phi);
  }
  ylms[lcur++] = fgPrefactors[0]*lbuf[0] * std::complex<double>(1,0);
  
  for (int il = 1; il<=lmax; il++) {
    // First im = 0
    ylms[lcur+il] = fgPrefactors[fgPrefshift[il]]*lbuf[fgPlmshift[il]]*std::complex<double>(1.0,0.0);
    // Im != 0
    for (int im=1; im<=il; im++) {
      lpol = lbuf[fgPlmshift[il]-im];
      ylms[lcur+il-im] = fgPrefactors[fgPrefshift[il]-im]*lpol*std::complex<double>(coss[im-1],-sins[im-1]);
      ylms[lcur+il+im] = fgPrefactors[fgPrefshift[il]+im]*lpol*std::complex<double>(coss[im-1],sins[im-1]);
    }
    lcur += 2*il + 1;
  }
}
 AliFemtoYlm.cxx:1
 AliFemtoYlm.cxx:2
 AliFemtoYlm.cxx:3
 AliFemtoYlm.cxx:4
 AliFemtoYlm.cxx:5
 AliFemtoYlm.cxx:6
 AliFemtoYlm.cxx:7
 AliFemtoYlm.cxx:8
 AliFemtoYlm.cxx:9
 AliFemtoYlm.cxx:10
 AliFemtoYlm.cxx:11
 AliFemtoYlm.cxx:12
 AliFemtoYlm.cxx:13
 AliFemtoYlm.cxx:14
 AliFemtoYlm.cxx:15
 AliFemtoYlm.cxx:16
 AliFemtoYlm.cxx:17
 AliFemtoYlm.cxx:18
 AliFemtoYlm.cxx:19
 AliFemtoYlm.cxx:20
 AliFemtoYlm.cxx:21
 AliFemtoYlm.cxx:22
 AliFemtoYlm.cxx:23
 AliFemtoYlm.cxx:24
 AliFemtoYlm.cxx:25
 AliFemtoYlm.cxx:26
 AliFemtoYlm.cxx:27
 AliFemtoYlm.cxx:28
 AliFemtoYlm.cxx:29
 AliFemtoYlm.cxx:30
 AliFemtoYlm.cxx:31
 AliFemtoYlm.cxx:32
 AliFemtoYlm.cxx:33
 AliFemtoYlm.cxx:34
 AliFemtoYlm.cxx:35
 AliFemtoYlm.cxx:36
 AliFemtoYlm.cxx:37
 AliFemtoYlm.cxx:38
 AliFemtoYlm.cxx:39
 AliFemtoYlm.cxx:40
 AliFemtoYlm.cxx:41
 AliFemtoYlm.cxx:42
 AliFemtoYlm.cxx:43
 AliFemtoYlm.cxx:44
 AliFemtoYlm.cxx:45
 AliFemtoYlm.cxx:46
 AliFemtoYlm.cxx:47
 AliFemtoYlm.cxx:48
 AliFemtoYlm.cxx:49
 AliFemtoYlm.cxx:50
 AliFemtoYlm.cxx:51
 AliFemtoYlm.cxx:52
 AliFemtoYlm.cxx:53
 AliFemtoYlm.cxx:54
 AliFemtoYlm.cxx:55
 AliFemtoYlm.cxx:56
 AliFemtoYlm.cxx:57
 AliFemtoYlm.cxx:58
 AliFemtoYlm.cxx:59
 AliFemtoYlm.cxx:60
 AliFemtoYlm.cxx:61
 AliFemtoYlm.cxx:62
 AliFemtoYlm.cxx:63
 AliFemtoYlm.cxx:64
 AliFemtoYlm.cxx:65
 AliFemtoYlm.cxx:66
 AliFemtoYlm.cxx:67
 AliFemtoYlm.cxx:68
 AliFemtoYlm.cxx:69
 AliFemtoYlm.cxx:70
 AliFemtoYlm.cxx:71
 AliFemtoYlm.cxx:72
 AliFemtoYlm.cxx:73
 AliFemtoYlm.cxx:74
 AliFemtoYlm.cxx:75
 AliFemtoYlm.cxx:76
 AliFemtoYlm.cxx:77
 AliFemtoYlm.cxx:78
 AliFemtoYlm.cxx:79
 AliFemtoYlm.cxx:80
 AliFemtoYlm.cxx:81
 AliFemtoYlm.cxx:82
 AliFemtoYlm.cxx:83
 AliFemtoYlm.cxx:84
 AliFemtoYlm.cxx:85
 AliFemtoYlm.cxx:86
 AliFemtoYlm.cxx:87
 AliFemtoYlm.cxx:88
 AliFemtoYlm.cxx:89
 AliFemtoYlm.cxx:90
 AliFemtoYlm.cxx:91
 AliFemtoYlm.cxx:92
 AliFemtoYlm.cxx:93
 AliFemtoYlm.cxx:94
 AliFemtoYlm.cxx:95
 AliFemtoYlm.cxx:96
 AliFemtoYlm.cxx:97
 AliFemtoYlm.cxx:98
 AliFemtoYlm.cxx:99
 AliFemtoYlm.cxx:100
 AliFemtoYlm.cxx:101
 AliFemtoYlm.cxx:102
 AliFemtoYlm.cxx:103
 AliFemtoYlm.cxx:104
 AliFemtoYlm.cxx:105
 AliFemtoYlm.cxx:106
 AliFemtoYlm.cxx:107
 AliFemtoYlm.cxx:108
 AliFemtoYlm.cxx:109
 AliFemtoYlm.cxx:110
 AliFemtoYlm.cxx:111
 AliFemtoYlm.cxx:112
 AliFemtoYlm.cxx:113
 AliFemtoYlm.cxx:114
 AliFemtoYlm.cxx:115
 AliFemtoYlm.cxx:116
 AliFemtoYlm.cxx:117
 AliFemtoYlm.cxx:118
 AliFemtoYlm.cxx:119
 AliFemtoYlm.cxx:120
 AliFemtoYlm.cxx:121
 AliFemtoYlm.cxx:122
 AliFemtoYlm.cxx:123
 AliFemtoYlm.cxx:124
 AliFemtoYlm.cxx:125
 AliFemtoYlm.cxx:126
 AliFemtoYlm.cxx:127
 AliFemtoYlm.cxx:128
 AliFemtoYlm.cxx:129
 AliFemtoYlm.cxx:130
 AliFemtoYlm.cxx:131
 AliFemtoYlm.cxx:132
 AliFemtoYlm.cxx:133
 AliFemtoYlm.cxx:134
 AliFemtoYlm.cxx:135
 AliFemtoYlm.cxx:136
 AliFemtoYlm.cxx:137
 AliFemtoYlm.cxx:138
 AliFemtoYlm.cxx:139
 AliFemtoYlm.cxx:140
 AliFemtoYlm.cxx:141
 AliFemtoYlm.cxx:142
 AliFemtoYlm.cxx:143
 AliFemtoYlm.cxx:144
 AliFemtoYlm.cxx:145
 AliFemtoYlm.cxx:146
 AliFemtoYlm.cxx:147
 AliFemtoYlm.cxx:148
 AliFemtoYlm.cxx:149
 AliFemtoYlm.cxx:150
 AliFemtoYlm.cxx:151
 AliFemtoYlm.cxx:152
 AliFemtoYlm.cxx:153
 AliFemtoYlm.cxx:154
 AliFemtoYlm.cxx:155
 AliFemtoYlm.cxx:156
 AliFemtoYlm.cxx:157
 AliFemtoYlm.cxx:158
 AliFemtoYlm.cxx:159
 AliFemtoYlm.cxx:160
 AliFemtoYlm.cxx:161
 AliFemtoYlm.cxx:162
 AliFemtoYlm.cxx:163
 AliFemtoYlm.cxx:164
 AliFemtoYlm.cxx:165
 AliFemtoYlm.cxx:166
 AliFemtoYlm.cxx:167
 AliFemtoYlm.cxx:168
 AliFemtoYlm.cxx:169
 AliFemtoYlm.cxx:170
 AliFemtoYlm.cxx:171
 AliFemtoYlm.cxx:172
 AliFemtoYlm.cxx:173
 AliFemtoYlm.cxx:174
 AliFemtoYlm.cxx:175
 AliFemtoYlm.cxx:176
 AliFemtoYlm.cxx:177
 AliFemtoYlm.cxx:178
 AliFemtoYlm.cxx:179
 AliFemtoYlm.cxx:180
 AliFemtoYlm.cxx:181
 AliFemtoYlm.cxx:182
 AliFemtoYlm.cxx:183
 AliFemtoYlm.cxx:184
 AliFemtoYlm.cxx:185
 AliFemtoYlm.cxx:186
 AliFemtoYlm.cxx:187
 AliFemtoYlm.cxx:188
 AliFemtoYlm.cxx:189
 AliFemtoYlm.cxx:190
 AliFemtoYlm.cxx:191
 AliFemtoYlm.cxx:192
 AliFemtoYlm.cxx:193
 AliFemtoYlm.cxx:194
 AliFemtoYlm.cxx:195
 AliFemtoYlm.cxx:196
 AliFemtoYlm.cxx:197
 AliFemtoYlm.cxx:198
 AliFemtoYlm.cxx:199
 AliFemtoYlm.cxx:200
 AliFemtoYlm.cxx:201
 AliFemtoYlm.cxx:202
 AliFemtoYlm.cxx:203
 AliFemtoYlm.cxx:204
 AliFemtoYlm.cxx:205
 AliFemtoYlm.cxx:206
 AliFemtoYlm.cxx:207
 AliFemtoYlm.cxx:208
 AliFemtoYlm.cxx:209
 AliFemtoYlm.cxx:210
 AliFemtoYlm.cxx:211
 AliFemtoYlm.cxx:212
 AliFemtoYlm.cxx:213
 AliFemtoYlm.cxx:214
 AliFemtoYlm.cxx:215
 AliFemtoYlm.cxx:216
 AliFemtoYlm.cxx:217
 AliFemtoYlm.cxx:218
 AliFemtoYlm.cxx:219
 AliFemtoYlm.cxx:220
 AliFemtoYlm.cxx:221
 AliFemtoYlm.cxx:222
 AliFemtoYlm.cxx:223
 AliFemtoYlm.cxx:224
 AliFemtoYlm.cxx:225
 AliFemtoYlm.cxx:226
 AliFemtoYlm.cxx:227
 AliFemtoYlm.cxx:228
 AliFemtoYlm.cxx:229
 AliFemtoYlm.cxx:230
 AliFemtoYlm.cxx:231
 AliFemtoYlm.cxx:232
 AliFemtoYlm.cxx:233
 AliFemtoYlm.cxx:234
 AliFemtoYlm.cxx:235
 AliFemtoYlm.cxx:236
 AliFemtoYlm.cxx:237
 AliFemtoYlm.cxx:238
 AliFemtoYlm.cxx:239
 AliFemtoYlm.cxx:240
 AliFemtoYlm.cxx:241
 AliFemtoYlm.cxx:242
 AliFemtoYlm.cxx:243
 AliFemtoYlm.cxx:244
 AliFemtoYlm.cxx:245
 AliFemtoYlm.cxx:246
 AliFemtoYlm.cxx:247
 AliFemtoYlm.cxx:248
 AliFemtoYlm.cxx:249
 AliFemtoYlm.cxx:250
 AliFemtoYlm.cxx:251
 AliFemtoYlm.cxx:252
 AliFemtoYlm.cxx:253
 AliFemtoYlm.cxx:254
 AliFemtoYlm.cxx:255
 AliFemtoYlm.cxx:256
 AliFemtoYlm.cxx:257
 AliFemtoYlm.cxx:258
 AliFemtoYlm.cxx:259
 AliFemtoYlm.cxx:260
 AliFemtoYlm.cxx:261
 AliFemtoYlm.cxx:262
 AliFemtoYlm.cxx:263
 AliFemtoYlm.cxx:264
 AliFemtoYlm.cxx:265
 AliFemtoYlm.cxx:266
 AliFemtoYlm.cxx:267
 AliFemtoYlm.cxx:268
 AliFemtoYlm.cxx:269
 AliFemtoYlm.cxx:270
 AliFemtoYlm.cxx:271
 AliFemtoYlm.cxx:272
 AliFemtoYlm.cxx:273
 AliFemtoYlm.cxx:274
 AliFemtoYlm.cxx:275
 AliFemtoYlm.cxx:276
 AliFemtoYlm.cxx:277
 AliFemtoYlm.cxx:278
 AliFemtoYlm.cxx:279
 AliFemtoYlm.cxx:280