wigner.cc

Go to the documentation of this file.
00001 /*
00002  *  This file is part of libcxxsupport.
00003  *
00004  *  libcxxsupport is free software; you can redistribute it and/or modify
00005  *  it under the terms of the GNU General Public License as published by
00006  *  the Free Software Foundation; either version 2 of the License, or
00007  *  (at your option) any later version.
00008  *
00009  *  libcxxsupport is distributed in the hope that it will be useful,
00010  *  but WITHOUT ANY WARRANTY; without even the implied warranty of
00011  *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00012  *  GNU General Public License for more details.
00013  *
00014  *  You should have received a copy of the GNU General Public License
00015  *  along with libcxxsupport; if not, write to the Free Software
00016  *  Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00017  */
00018 
00019 /*
00020  *  libcxxsupport is being developed at the Max-Planck-Institut fuer Astrophysik
00021  *  and financially supported by the Deutsches Zentrum fuer Luft- und Raumfahrt
00022  *  (DLR).
00023  */
00024 
00025 /*! \file wigner.cc
00026  *  Several C++ classes for calculating Wigner matrices
00027  *
00028  *  Copyright (C) 2009-2015 Max-Planck-Society
00029  *  \author Martin Reinecke and others (see individual classes)
00030  */
00031 
00032 #include "wigner.h"
00033 #include "lsconstants.h"
00034 
00035 using namespace std;
00036 
00037 void wigner_d_halfpi_risbo_scalar::do_line0 (double *l1, int j)
00038   {
00039   double xj = pq/j;
00040   for (int i=n; i>=1; --i)
00041     l1[i] = xj*sqt[j]*(sqt[j-i]*l1[i] - sqt[i]*l1[i-1]);
00042   l1[0] = pq*l1[0];
00043   }
00044 void wigner_d_halfpi_risbo_scalar::do_line (const double *l1, double *l2,
00045   int j, int k)
00046   {
00047   double xj = pq/j;
00048   double t1 = xj*sqt[j-k];
00049   double t2 = xj*sqt[k];
00050   for (int i=n; i>=1; --i)
00051     l2[i] = t1 * (sqt[j-i]*l2[i] - sqt[i]*l2[i-1])
00052             +t2 * (sqt[j-i]*l1[i] + sqt[i]*l1[i-1]);
00053   l2[0] = sqt[j] * (t2*l1[0]+t1*l2[0]);
00054   }
00055 
00056 wigner_d_halfpi_risbo_scalar::wigner_d_halfpi_risbo_scalar(int lmax)
00057   : pq(.5*sqrt(2.)), sqt(2*lmax+1), d(lmax+2,lmax+2), n(-1)
00058   { for (tsize m=0; m<sqt.size(); ++m) sqt[m] = sqrt(double(m)); }
00059 
00060 const arr2<double> &wigner_d_halfpi_risbo_scalar::recurse ()
00061   {
00062   ++n;
00063   if (n==0)
00064     d[0][0] = 1;
00065   else if (n==1)
00066     {
00067     d[0][0] = .5; d[0][1] =-pq;
00068     d[1][0] = pq; d[1][1] = 0.;
00069     }
00070   else
00071     {
00072 //padding
00073     int flip = 1;
00074     for (int i=0; i<n; ++i)
00075       {
00076       d[i][n]=flip*d[i][n-2];
00077       d[n][i]=flip*d[n-2][i];
00078       flip=-flip;
00079       }
00080     d[n][n]=flip*d[n-2][n];
00081 
00082     do_line (d[n-1],d[n],2*n-1,n);
00083     for (int k=n; k>=2; --k)
00084       {
00085       do_line (d[k-2],d[k-1],2*n-1,k-1);
00086       do_line (d[k-1],d[k],2*n,k);
00087       }
00088     do_line0 (d[0],2*n-1);
00089     do_line (d[0],d[1],2*n,1);
00090     do_line0 (d[0],2*n);
00091     }
00092   return d;
00093   }
00094 
00095 void wigner_d_risbo_scalar::do_line0 (double *l1, int j)
00096   {
00097   double xj = 1./j;
00098   l1[j] = -p*l1[j-1];
00099   for (int i=j-1; i>=1; --i)
00100     l1[i] = xj*sqt[j]*(q*sqt[j-i]*l1[i] - p*sqt[i]*l1[i-1]);
00101   l1[0] = q*l1[0];
00102   }
00103 void wigner_d_risbo_scalar::do_line (const double *l1, double *l2, int j, int k)
00104   {
00105   double xj = 1./j;
00106   double t1 = xj*sqt[j-k]*q, t2 = xj*sqt[j-k]*p;
00107   double t3 = xj*sqt[k]*p,   t4 = xj*sqt[k]*q;
00108   l2[j] = sqt[j] * (t4*l1[j-1]-t2*l2[j-1]);
00109   for (int i=j-1; i>=1; --i)
00110     l2[i] = t1*sqt[j-i]*l2[i] - t2*sqt[i]*l2[i-1]
00111             +t3*sqt[j-i]*l1[i] + t4*sqt[i]*l1[i-1];
00112   l2[0] = sqt[j] * (t3*l1[0]+t1*l2[0]);
00113   }
00114 
00115 wigner_d_risbo_scalar::wigner_d_risbo_scalar(int lmax, double ang)
00116   : p(sin(ang/2)), q(cos(ang/2)), sqt(2*lmax+1),
00117     d(lmax+1,2*lmax+1), n(-1)
00118   { for (tsize m=0; m<sqt.size(); ++m) sqt[m] = sqrt(double(m)); }
00119 
00120 const arr2<double> &wigner_d_risbo_scalar::recurse ()
00121   {
00122   ++n;
00123   if (n==0)
00124     d[0][0] = 1;
00125   else if (n==1)
00126     {
00127     d[0][0] = q*q; d[0][1] = -p*q*sqt[2]; d[0][2] = p*p;
00128     d[1][0] = -d[0][1]; d[1][1] = q*q-p*p; d[1][2] = d[0][1];
00129     }
00130   else
00131     {
00132     // padding
00133     int sign = (n&1)? -1 : 1;
00134     for (int i=0; i<=2*n-2; ++i)
00135       {
00136       d[n][i] = sign*d[n-2][2*n-2-i];
00137       sign=-sign;
00138       }
00139     do_line (d[n-1],d[n],2*n-1,n);
00140     for (int k=n; k>=2; --k)
00141       {
00142       do_line (d[k-2],d[k-1],2*n-1,k-1);
00143       do_line (d[k-1],d[k],2*n,k);
00144       }
00145     do_line0 (d[0],2*n-1);
00146     do_line (d[0],d[1],2*n,1);
00147     do_line0 (d[0],2*n);
00148     }
00149   return d;
00150   }
00151 
00152 wigner_d_halfpi_risbo_openmp::wigner_d_halfpi_risbo_openmp(int lmax)
00153   : pq(.5*sqrt(2.)), sqt(2*lmax+1), d(lmax+2,lmax+2),
00154     dd(lmax+2,lmax+2), n(-1)
00155   { for (tsize m=0; m<sqt.size(); ++m) sqt[m] = sqrt(double(m)); }
00156 
00157 const arr2<double> &wigner_d_halfpi_risbo_openmp::recurse ()
00158   {
00159   ++n;
00160   if (n==0)
00161     d[0][0] = 1;
00162   else if (n==1)
00163     {
00164     d.fast_alloc(3,3);
00165     d[0][0] = .5; d[0][1] =-pq;
00166     d[1][0] = pq; d[1][1] = 0.;
00167     }
00168   else
00169     {
00170 //padding
00171     int flip = 1;
00172     for (int i=0; i<n; ++i)
00173       {
00174       d[i][n]=flip*d[i][n-2];
00175       d[n][i]=flip*d[n-2][i];
00176       flip=-flip;
00177       }
00178     d[n][n]=flip*d[n-2][n];
00179 
00180     for (int j=2*n-1; j<=2*n; ++j)
00181       {
00182       dd.fast_alloc(n+2,n+2);
00183       double tmpx1 = pq/j;
00184       dd[0][0] = pq*d[0][0];
00185       for (int i=1;i<=n; ++i)
00186         dd[0][i] = tmpx1*sqt[j]*(sqt[j-i]*d[0][i] - sqt[i]*d[0][i-1]);
00187 #pragma omp parallel
00188 {
00189       int k;
00190 #pragma omp for schedule(static)
00191       for (k=1; k<=n; ++k)
00192         {
00193         double stmp1=sqt[j-k]*tmpx1;
00194         double stmp2=sqt[k]*tmpx1;
00195         double save1 = stmp1*d[k][0], save2 = stmp2*d[k-1][0];
00196         dd[k][0] = sqt[j]*(save1+save2);
00197         for (int i=1; i<=n; ++i)
00198           {
00199           dd[k][i] = sqt[i]*(save2-save1);
00200           save1 = stmp1*d[k][i];
00201           save2 = stmp2*d[k-1][i];
00202           dd[k][i] += sqt[j-i]*(save1+save2);
00203           }
00204         }
00205 }
00206       dd.swap(d);
00207       }
00208     }
00209   return d;
00210   }
00211 
00212 wigner_d_risbo_openmp::wigner_d_risbo_openmp(int lmax, double ang)
00213   : p(sin(ang/2)), q(cos(ang/2)), sqt(2*lmax+1),
00214     d(lmax+1,2*lmax+1), dd(lmax+1,2*lmax+1), n(-1)
00215   { for (tsize m=0; m<sqt.size(); ++m) sqt[m] = sqrt(double(m)); }
00216 
00217 const arr2<double> &wigner_d_risbo_openmp::recurse ()
00218   {
00219   ++n;
00220   if (n==0)
00221     d[0][0] = 1;
00222   else if (n==1)
00223     {
00224     d[0][0] = q*q; d[0][1] = -p*q*sqt[2]; d[0][2] = p*p;
00225     d[1][0] = -d[0][1]; d[1][1] = q*q-p*p; d[1][2] = d[0][1];
00226     }
00227   else
00228     {
00229     // padding
00230     int sign = (n&1)? -1 : 1;
00231     for (int i=0; i<=2*n-2; ++i)
00232       {
00233       d[n][i] = sign*d[n-2][2*n-2-i];
00234       sign=-sign;
00235       }
00236     for (int j=2*n-1; j<=2*n; ++j)
00237       {
00238       double xj = 1./j;
00239       dd[0][0] = q*d[0][0];
00240       for (int i=1;i<j; ++i)
00241         dd[0][i] = xj*sqt[j]*(q*sqt[j-i]*d[0][i] - p*sqt[i]*d[0][i-1]);
00242       dd[0][j] = -p*d[0][j-1];
00243 #pragma omp parallel
00244 {
00245       int k;
00246 #pragma omp for schedule(static)
00247       for (k=1; k<=n; ++k)
00248         {
00249         double t1 = xj*sqt[j-k]*q, t2 = xj*sqt[j-k]*p;
00250         double t3 = xj*sqt[k  ]*p, t4 = xj*sqt[k  ]*q;
00251         dd[k][0] = xj*sqt[j]*(q*sqt[j-k]*d[k][0] + p*sqt[k]*d[k-1][0]);
00252         for (int i=1; i<j; ++i)
00253           dd[k][i] = t1*sqt[j-i]*d[k  ][i] - t2*sqt[i]*d[k  ][i-1]
00254                     + t3*sqt[j-i]*d[k-1][i] + t4*sqt[i]*d[k-1][i-1];
00255         dd[k][j] = -t2*sqt[j]*d[k][j-1] + t4*sqt[j]*d[k-1][j-1];
00256         }
00257 }
00258       dd.swap(d);
00259       }
00260     }
00261   return d;
00262   }
00263 
00264 
00265 wignergen_scalar::wignergen_scalar (int lmax_, const arr<double> &thetas,
00266   double epsilon)
00267   : eps(epsilon), lmax(lmax_),
00268     logsum(2*lmax+1), lc05(thetas.size()), ls05(thetas.size()),
00269     flm1(2*lmax+1), flm2(2*lmax+1),
00270     cf(maxscale+1-minscale), costh(thetas.size()), xl(lmax+1),
00271     thetaflip(thetas.size()),
00272     m1(-1234567890), m2(-1234567890), am1(-1234567890), am2(-1234567890),
00273     mlo(-1234567890), mhi(-1234567890),
00274     fx(lmax+2), result(lmax+1)
00275   {
00276   planck_assert(lmax>=0,"lmax too small");
00277   logsum[0] = 0.;
00278   for (tsize m=1; m<logsum.size(); ++m)
00279     logsum[m] = logsum[m-1]+log(static_cast<long double>(m));
00280   for (tsize lm=0; lm<flm1.size(); ++lm)
00281     {
00282     flm1[lm] = sqrt(1./(lm+1.));
00283     flm2[lm] = sqrt(lm/(lm+1.));
00284     }
00285   for (tsize i=0; i<cf.size(); ++i)
00286     cf[i] = ldexp(1.,(int(i)+minscale)*large_exponent2);
00287 
00288   fsmall = ldexp(1.,-large_exponent2);
00289   fbig = ldexp(1.,large_exponent2);
00290 
00291   for (tsize i=0; i<thetas.size(); ++i)
00292     {
00293     double theta=fmodulo(thetas[i],twopi);
00294     if (theta>pi) theta-=twopi;
00295     thetaflip[i]=(theta<0);
00296     theta=abs(theta); // now theta is in (0; pi)
00297     // tiny adjustments to make sure cos and sin (theta/2) are positive
00298     if (theta==0.) theta=1e-16;
00299     if (abs_approx(theta,pi,1e-15)) theta=pi-1e-15;
00300     costh[i]=cos(theta);
00301     lc05[i]=log(cos(0.5L*theta));
00302     ls05[i]=log(sin(0.5L*theta));
00303     }
00304   xl[0]=0;
00305   for (tsize l=1; l<xl.size(); ++l) xl[l]=1./l;
00306 
00307   for (tsize l=0; l<fx.size(); ++l)
00308     fx[l][0]=fx[l][1]=fx[l][2]=0.;
00309   }
00310 
00311 void wignergen_scalar::prepare (int m1_, int m2_)
00312   {
00313   if ((m1_==m1) && (m2_==m2)) return;
00314 
00315   int mlo_=abs(m1_), mhi_=abs(m2_);
00316   if (mhi_<mlo_) swap(mhi_,mlo_);
00317   bool ms_similar = ((mhi==mhi_) && (mlo==mlo_));
00318   bool flip_m_sign = ((m1*m2)!=(m1_*m2_));
00319 
00320   m1=m1_; m2=m2_;
00321   mlo=am1=abs(m1); mhi=am2=abs(m2);
00322   if (mhi<mlo) swap(mhi,mlo);
00323 
00324   if (ms_similar)
00325     {
00326     if (flip_m_sign)
00327       for (int l=mhi; l<lmax; ++l)
00328         fx[l+1][1]=-fx[l+1][1];
00329     }
00330   else
00331     {
00332     for (int l=mhi; l<lmax; ++l)
00333       {
00334       double t = flm1[l+m1]*flm1[l-m1]*flm1[l+m2]*flm1[l-m2];
00335       double lt = 2*l+1;
00336       double l1 = l+1;
00337       fx[l+1][0]=l1*lt*t;
00338       fx[l+1][1]=m1*m2*xl[l]*xl[l+1];
00339       t = flm2[l+m1]*flm2[l-m1]*flm2[l+m2]*flm2[l-m2];
00340       fx[l+1][2]=t*l1*xl[l];
00341       }
00342     }
00343 
00344   prefactor = 0.5L*(logsum[2*mhi]-logsum[mhi+mlo]-logsum[mhi-mlo]);
00345 
00346   preMinus = false;
00347   if (mhi==am1)
00348     {
00349     cosPow = mhi-m2; sinPow = mhi+m2;
00350     if (m1>=0)
00351       { swap(cosPow, sinPow); preMinus=((mhi-m2)&1); }
00352     }
00353   else
00354     {
00355     cosPow = mhi+m1; sinPow = mhi-m1;
00356     if (m2<0)
00357       { swap(cosPow, sinPow); preMinus=((mhi+m1)&1); }
00358     }
00359   }
00360 
00361 const arr<double> &wignergen_scalar::calc (int nth, int &firstl)
00362   {
00363   calc(nth, firstl, result);
00364   return result;
00365   }
00366 
00367 void wignergen_scalar::calc (int nth, int &firstl, arr<double> &resx) const
00368   {
00369   int l=mhi;
00370   const dbl3 *fy = &fx[0];
00371   const double cth = costh[nth];
00372   double *res = &resx[0];
00373   long double logval = prefactor + lc05[nth]*cosPow + ls05[nth]*sinPow;
00374   logval *= inv_ln2;
00375   int scale = int (logval/large_exponent2)-minscale;
00376   double rec1 = 0.;
00377   double rec2 = double(exp(ln2*(logval-(scale+minscale)*large_exponent2)));
00378   if (preMinus ^ (thetaflip[nth] && ((am1+am2)&1))) rec2 = -rec2;
00379 
00380   while(scale<0) // iterate until we reach the realm of IEEE numbers
00381     {
00382     if (++l>lmax) break;
00383     rec1 = (cth - fy[l][1])*fy[l][0]*rec2 - fy[l][2]*rec1;
00384     if (++l>lmax) break;
00385     rec2 = (cth - fy[l][1])*fy[l][0]*rec1 - fy[l][2]*rec2;
00386 
00387     while (abs(rec2)>fbig)
00388       {
00389       rec1 *= fsmall;
00390       rec2 *= fsmall;
00391       ++scale;
00392       }
00393     }
00394 
00395   if (scale<0) { firstl=lmax+1; return; }
00396   rec1 *= cf[scale];
00397   rec2 *= cf[scale];
00398 
00399   for (;l<lmax-1;l+=2) // iterate until we cross the eps threshold
00400     {
00401     if (abs(rec2)>eps) break;
00402     rec1 = (cth - fy[l+1][1])*fy[l+1][0]*rec2 - fy[l+1][2]*rec1;
00403     if (abs(rec1)>eps) { swap(rec1,rec2); ++l; break; }
00404     rec2 = (cth - fy[l+2][1])*fy[l+2][0]*rec1 - fy[l+2][2]*rec2;
00405     }
00406   if ((abs(rec2)<=eps) && (++l<=lmax))
00407     {
00408     rec1 = (cth - fy[l][1])*fy[l][0]*rec2 - fy[l][2]*rec1;
00409     swap (rec1,rec2);
00410     }
00411 
00412   firstl = l;
00413   if (l>lmax) return;
00414 
00415   res[l]=rec2;
00416 
00417   for (;l<lmax-1;l+=2)
00418     {
00419     res[l+1] = rec1 = (cth - fy[l+1][1])*fy[l+1][0]*rec2 - fy[l+1][2]*rec1;
00420     res[l+2] = rec2 = (cth - fy[l+2][1])*fy[l+2][0]*rec1 - fy[l+2][2]*rec2;
00421     }
00422   while (true)
00423     {
00424     if (++l>lmax) break;
00425     res[l] = rec1 = (cth - fy[l][1])*fy[l][0]*rec2 - fy[l][2]*rec1;
00426     if (++l>lmax) break;
00427     res[l] = rec2 = (cth - fy[l][1])*fy[l][0]*rec1 - fy[l][2]*rec2;
00428     }
00429   }
00430 
00431 #ifdef __SSE2__
00432 
00433 #define RENORMALIZE \
00434   do \
00435     { \
00436     double rec1a, rec1b, rec2a, rec2b, cfa, cfb; \
00437     rec1.writeTo(rec1a,rec1b); rec2.writeTo(rec2a,rec2b); \
00438     corfac.writeTo(cfa,cfb); \
00439     while (abs(rec2a)>fbig) \
00440       { \
00441       rec1a*=fsmall; rec2a*=fsmall; ++scale1; \
00442       cfa = (scale1<0) ? 0. : cf[scale1]; \
00443       } \
00444     while (abs(rec2b)>fbig) \
00445       { \
00446       rec1b*=fsmall; rec2b*=fsmall; ++scale2; \
00447       cfb = (scale2<0) ? 0. : cf[scale2]; \
00448       } \
00449     rec1.readFrom(rec1a,rec1b); rec2.readFrom(rec2a,rec2b); \
00450     corfac.readFrom(cfa,cfb); \
00451     } \
00452   while(0)
00453 
00454 #define GETPRE(prea,preb,lv) \
00455   prea=(cth-fy[lv][1])*fy[lv][0]; \
00456   preb=fy[lv][2];
00457 
00458 #define NEXTSTEP(prea,preb,prec,pred,reca,recb,lv) \
00459   { \
00460   prec = fy[lv][1]; \
00461   preb *= reca; \
00462   prea *= recb; \
00463   V2df t0 (fy[lv][0]); \
00464   prec = cth-prec; \
00465   pred = fy[lv][2]; \
00466   reca = prea-preb; \
00467   prec *= t0; \
00468   }
00469 
00470 const arr_align<V2df,16> &wignergen::calc (int nth1, int nth2, int &firstl)
00471   {
00472   calc(nth1, nth2, firstl, result2);
00473   return result2;
00474   }
00475 
00476 void wignergen::calc (int nth1, int nth2, int &firstl,
00477   arr_align<V2df,16> &resx) const
00478   {
00479   int l=mhi;
00480   const dbl3 *fy = &fx[0];
00481   const V2df cth(costh[nth1],costh[nth2]);
00482   V2df *res = &resx[0];
00483   long double logval1 = prefactor + lc05[nth1]*cosPow + ls05[nth1]*sinPow,
00484               logval2 = prefactor + lc05[nth2]*cosPow + ls05[nth2]*sinPow;
00485   logval1 *= inv_ln2;
00486   logval2 *= inv_ln2;
00487   int scale1 = int (logval1/large_exponent2)-minscale,
00488       scale2 = int (logval2/large_exponent2)-minscale;
00489   V2df rec1(0.);
00490   double tr1 = double(exp(ln2*(logval1-(scale1+minscale)*large_exponent2))),
00491          tr2 = double(exp(ln2*(logval2-(scale2+minscale)*large_exponent2)));
00492   if (preMinus ^ (thetaflip[nth1] && ((am1+am2)&1))) tr1 = -tr1;
00493   if (preMinus ^ (thetaflip[nth2] && ((am1+am2)&1))) tr2 = -tr2;
00494   V2df rec2(tr1,tr2);
00495   V2df corfac ((scale1<0) ? 0. : cf[scale1], (scale2<0) ? 0. : cf[scale2]);
00496 
00497   V2df eps2(eps);
00498   V2df fbig2(fbig);
00499 
00500   V2df pre0,pre1,pre2,pre3;
00501 
00502   GETPRE(pre0,pre1,l+1)
00503   if ((scale1<0) && (scale2<0))
00504     {
00505     while (true)
00506       {
00507       if (++l>lmax) break;
00508       NEXTSTEP(pre0,pre1,pre2,pre3,rec1,rec2,l+1)
00509       if (++l>lmax) break;
00510       NEXTSTEP(pre2,pre3,pre0,pre1,rec2,rec1,l+1)
00511       if (any(abs(rec2).gt(fbig2)))
00512         {
00513         RENORMALIZE;
00514         if ((scale1>=0) || (scale2>=0)) break;
00515         }
00516       }
00517     }
00518 
00519   if (l<=lmax)
00520     {
00521     GETPRE(pre0,pre1,l+1)
00522     while (true)
00523       {
00524       V2df t1;
00525       res[l]=t1=rec2*corfac;
00526       if (any(abs(t1).gt(eps2)))
00527         break;
00528 
00529       if (++l>lmax) break;
00530       NEXTSTEP(pre0,pre1,pre2,pre3,rec1,rec2,l+1)
00531 
00532       res[l]=t1=rec1*corfac;
00533       if (any(abs(t1).gt(eps2)))
00534         { swap(rec1,rec2); break; }
00535 
00536       if (++l>lmax) break;
00537       NEXTSTEP(pre2,pre3,pre0,pre1,rec2,rec1,l+1)
00538 
00539       if (any(abs(rec2).gt(fbig2)))
00540         RENORMALIZE;
00541       }
00542     }
00543   firstl=l;
00544   if (l>lmax) return;
00545 
00546   GETPRE(pre0,pre1,l+1)
00547   while (true)
00548     {
00549     V2df t1;
00550     res[l]=t1=rec2*corfac;
00551     if (all(abs(t1).ge(eps2)))
00552       break;
00553 
00554     if (++l>lmax) break;
00555     NEXTSTEP(pre0,pre1,pre2,pre3,rec1,rec2,l+1)
00556 
00557     res[l]=t1=rec1*corfac;
00558     if (all(abs(t1).ge(eps2)))
00559       { swap(rec1,rec2); break; }
00560 
00561     if (++l>lmax) break;
00562     NEXTSTEP(pre2,pre3,pre0,pre1,rec2,rec1,l+1)
00563 
00564     if (any(abs(rec2).gt(fbig2)))
00565       RENORMALIZE;
00566     }
00567 
00568   if (l>lmax) return;
00569   rec1*=corfac;
00570   rec2*=corfac;
00571 
00572   GETPRE(pre0,pre1,l+1)
00573   for (;l<lmax-1;l+=2)
00574     {
00575     res[l] = rec2;
00576     NEXTSTEP(pre0,pre1,pre2,pre3,rec1,rec2,l+2)
00577     res[l+1] = rec1;
00578     NEXTSTEP(pre2,pre3,pre0,pre1,rec2,rec1,l+3)
00579     }
00580 
00581   res[l] = rec2;
00582   if (++l<=lmax)
00583     {
00584     NEXTSTEP(pre0,pre1,pre2,pre3,rec1,rec2,l+1)
00585     res[l] = rec1;
00586     }
00587   }
00588 
00589 #endif /* __SSE2__ */
00590 
00591 wigner_estimator::wigner_estimator (int lmax_, double epsPow_)
00592   : lmax(lmax_), xlmax(1./lmax_), epsPow(epsPow_) {}
00593 
00594 void wigner_estimator::prepare_m (int m1_, int m2_)
00595   {
00596   m1=abs(m1_); m2=abs(m2_);
00597   mbig=max(m1,m2);
00598   double cos1=m1*xlmax, cos2=m2*xlmax;
00599   double s1s2=sqrt((1.-cos1*cos1)*(1.-cos2*cos2));
00600   cosm1m2=cos1*cos2+s1s2;
00601   }
00602 
00603 bool wigner_estimator::canSkip (double theta) const
00604   {
00605   if (mbig==lmax) return false; // don't have a good criterion for this case
00606   double delta = m1*m1 + m2*m2 - abs(2.*m1*m2*cos(theta));
00607   double sth = sin(theta);
00608   if (abs_approx(sth,0.,1e-7)) return (delta>1.); // close to a pole
00609   return (((sqrt(delta)-epsPow)*cosm1m2/abs(sth)) > lmax);
00610   }

Generated on Thu Oct 8 14:48:51 2015 for LevelS C++ support library