Main Page | Modules | Class Hierarchy | Alphabetical List | Class List | File List | Class Members | File Members | Related Pages

ylmgen.h

00001 /*
00002  *  This file is part of Healpix_cxx.
00003  *
00004  *  Healpix_cxx 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  *  Healpix_cxx 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 Healpix_cxx; if not, write to the Free Software
00016  *  Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00017  *
00018  *  For more information about HEALPix, see http://healpix.jpl.nasa.gov
00019  */
00020 
00021 /*
00022  *  Healpix_cxx is being developed at the Max-Planck-Institut fuer Astrophysik
00023  *  and financially supported by the Deutsches Zentrum fuer Luft- und Raumfahrt
00024  *  (DLR).
00025  */
00026 
00027 /*
00028  *  Code for efficient calculation of Y_lm(phi=0,theta)
00029  *
00030  *  Copyright (C) 2005 Max-Planck-Society
00031  *  Author: Martin Reinecke
00032  */
00033 
00034 #ifndef HEALPIX_YLMGEN_H
00035 #define HEALPIX_YLMGEN_H
00036 
00037 #include <cmath>
00038 #include "arr.h"
00039 #include "constants.h"
00040 
00041 /*! Class for efficient calculation of Y_lm(theta,phi=0) */
00042 class Ylmgen
00043   {
00044   private:
00045     double fsmall, fbig, eps, cth_crit;
00046     int lmax, mmax, m_last, m_crit;
00047     arr<double> cf;
00048     arr<double[2]> recfac;
00049     arr<double> mfac;
00050 
00051     enum { large_exponent2 = 90, minscale=-4 };
00052 
00053     void recalc_recfac (int m)
00054       {
00055       using namespace std;
00056 
00057       if (m_last==m) return;
00058 
00059       int m2 = m*m;
00060       double f_old=1;
00061       for (int l=m; l<recfac.size(); ++l)
00062         {
00063         int l2 = (l+1)*(l+1);
00064         recfac[l][0] = sqrt(double(4*l2 - 1) / (l2-m2));
00065         recfac[l][1] = recfac[l][0]/f_old;
00066         f_old = recfac[l][0];
00067         }
00068 
00069       m_last=m;
00070       }
00071 
00072   public:
00073     /*! Creates a generator which will calculate Y_lm(theta,phi=0)
00074         up to \a l=l_max and \a m=m_max. It may regard Y_lm whose absolute
00075         magnitude is smaller than \a epsilon as zero. */
00076     Ylmgen (int l_max, int m_max, double epsilon=1e-30)
00077       : eps(epsilon), cth_crit(2.), lmax(l_max), mmax(m_max), m_last(-1),
00078         m_crit(mmax+1), cf(-minscale+11), recfac(lmax+1), mfac(mmax+1)
00079       {
00080       using namespace std;
00081 
00082       fsmall = ldexp(1.,-large_exponent2);
00083       fbig   = ldexp(1., large_exponent2);
00084       for (int m=0; m<cf.size(); ++m)
00085         cf[m] = ldexp(1.,(m+minscale)*large_exponent2);
00086 
00087       mfac[0] = 1;
00088       for (int m=1; m<mfac.size(); ++m)
00089         mfac[m] = mfac[m-1]*sqrt((2*m+1.)/(2*m));
00090       for (int m=0; m<mfac.size(); ++m)
00091         mfac[m] = inv_ln2*log(inv_sqrt4pi*mfac[m]);
00092       }
00093 
00094     /*! For a colatitude given by \a cth and \a sth (representing cos(theta)
00095         and sin(theta)) and a multipole moment \a m, calculate the
00096         Y_lm(theta,phi=0) for \a m<=l<=lmax and return in it \a result[l].
00097         On exit, \a firstl is the \a l index of the first Y_lm with an
00098         absolute magnitude larger than \a epsilon. If \a firstl>lmax, all
00099         absolute values are smaller than \a epsilon.
00100         \a result[l] is undefined for all \a l<firstl. */
00101     void get_Ylm (double cth, double sth, int m, arr<double> &result,
00102       int &firstl)
00103       {
00104       using namespace std;
00105 
00106       planck_assert (m<=mmax, "get_Ylm: m larger than mmax");
00107 
00108       if (((m>=m_crit)&&(abs(cth)>=cth_crit)) || ((m>0)&&(sth==0)))
00109         { firstl=lmax+1; return; }
00110 
00111       recalc_recfac(m);
00112       result.alloc(lmax+1);
00113 
00114       double logval = mfac[m];
00115       if (m>0) logval += m*inv_ln2*log(sth);
00116       int scale = int (logval/large_exponent2)-minscale;
00117       double corfac = (scale<0) ? 0. : cf[scale];
00118       double lam_1 = 0;
00119       double lam_2 = exp(ln2*(logval-(scale+minscale)*large_exponent2));
00120       if (m&1) lam_2 = -lam_2;
00121 
00122       int l=m;
00123       while (true)
00124         {
00125         if (abs(lam_2*corfac)>eps) break;
00126         if (++l>lmax) break;
00127         double lam_0 = cth*lam_2*recfac[l-1][0] - lam_1*recfac[l-1][1];
00128         if (abs(lam_0*corfac)>eps) { lam_1=lam_2; lam_2=lam_0; break; }
00129         if (++l>lmax) break;
00130         lam_1 = cth*lam_0*recfac[l-1][0] - lam_2*recfac[l-1][1];
00131         if (abs(lam_1*corfac)>eps) { lam_2=lam_1; lam_1=lam_0; break; }
00132         if (++l>lmax) break;
00133         lam_2 = cth*lam_1*recfac[l-1][0] - lam_0*recfac[l-1][1];
00134 
00135         while (abs(lam_2)>fbig)
00136           {
00137           lam_1 *= fsmall;
00138           lam_2 *= fsmall;
00139           ++scale;
00140           corfac = (scale<0) ? 0. : cf[scale];
00141           }
00142         }
00143 
00144       firstl=l;
00145       if (l>lmax)
00146         { m_crit=m; cth_crit=abs(cth); return; }
00147 
00148       lam_1*=corfac;
00149       lam_2*=corfac;
00150 
00151       for (;l<lmax-2;l+=3)
00152         {
00153         result[l]=lam_2;
00154         double lam_0 = cth*lam_2*recfac[l][0] - lam_1*recfac[l][1];
00155         result[l+1] = lam_0;
00156         lam_1 = cth*lam_0*recfac[l+1][0] - lam_2*recfac[l+1][1];
00157         result[l+2] = lam_1;
00158         lam_2 = cth*lam_1*recfac[l+2][0] - lam_0*recfac[l+2][1];
00159         }
00160       while (true)
00161         {
00162         result[l]=lam_2;
00163         if (++l>lmax) break;
00164         double lam_0 = cth*lam_2*recfac[l-1][0] - lam_1*recfac[l-1][1];
00165         result[l] = lam_0;
00166         if (++l>lmax) break;
00167         lam_1 = cth*lam_0*recfac[l-1][0] - lam_2*recfac[l-1][1];
00168         result[l] = lam_1;
00169         if (++l>lmax) break;
00170         lam_2 = cth*lam_1*recfac[l-1][0] - lam_0*recfac[l-1][1];
00171         }
00172       }
00173 
00174   };
00175 
00176 #endif

Generated on Fri Jul 8 09:37:15 2005 for Healpix C++