S2LET library : documentation for IDL routines and interfaces

This page was created by the IDL library routine This page was created by the IDL library routine This page was created by the IDL library routine This page was created by the IDL library routine This page was created by the IDL library routine mk_html_help. For more information on this routine, refer to the IDL Online Help Navigator or type:

     ? mk_html_help

at the IDL command line prompt.

Last modified: Tue Sep 22 12:42:56 2015.


List of Routines


Routine Descriptions

S2LET_AXISYM_HPX_WAV_ANALYSIS

[Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_axisym_hpx_wav_analysis

 PURPOSE:
   Compute wavelet transform of a Healpix map

 CALLING SEQUENCE:
   f_wav = s2let_axisym_hpx_wav_analysis(f, B, L, J_min)

 INPUTS:
   f     - input healpix nside (number of pixels is 12*nside^2)
   B     - Wavelet parameter
   L     - Band-limit to be used for the spherical harmonic transforms
   J_min - First wavelet scale to be used
   wavtype - Wavelet type (1: scale-discretised, 2:needlets, 3: cubic splines)

 OUTPUTS:
   f_wav - Struc containing the wavelets

 EXAMPLE:

   see s2let_hpx_demo

   read_fits_map, file, f
   f_wav = s2let_axisym_hpx_wav_analysis(f, B, L, J_min)
   f_rec = s2let_axisym_hpx_wav_synthesis(f_wav)
   J_max = s2let_j_max(L, b)
   mollview, f_rec, title='Reconstructed map'
   mollview, f_wav.scal, title='Scaling map'
   for j=0, J_max-J_min do begin
       mollview, f_wav.(j), title='Wavelet map '
            +strtrim(j+1,2)+' on '+strtrim(J_max-J_min+1,2)
   endfor
   
 COMMENTS:
   Nside of the input map is automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_axisym_hpx_wav_analysis.pro)


S2LET_AXISYM_HPX_WAV_SYNTHESIS

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_axisym_hpx_wav_synthesis

 PURPOSE:
   Reconstruct a healpix map from its wavelet transform

 CALLING SEQUENCE:
   f = s2let_axisym_hpx_wav_synthesis(f_wav)

 INPUTS:
   f_wav - Struc containing the wavelets maps

 OUTPUTS:
   f - The reconstructed Healpix map

 EXAMPLE:

   see s2let_hpx_demo

   read_fits_map, file, f
   f_wav = s2let_axisym_hpx_wav_analysis(f, B, L, J_min)
   f_rec = s2let_axisym_hpx_wav_synthesis(f_wav)
   J_max = s2let_j_max(L, b)
   mollview, f_rec, title='Reconstructed map'
   mollview, f_wav.scal, title='Scaling map'
   for j=0, J_max-J_min do begin
       mollview, f_wav.(j), title='Wavelet map '
            +strtrim(j+1,2)+' on '+strtrim(J_max-J_min+1,2)
   endfor
   
 COMMENTS:
   B, L and J_min of the input wavelets are automatically detected 

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_axisym_hpx_wav_synthesis.pro)


S2LET_AXISYM_MW_WAV_ANALYSIS

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_axisym_mw_wav_analysis

 PURPOSE:
   Compute full resolution exact wavelet transform of a complex MW map

 CALLING SEQUENCE:
   f_wav = s2let_axisym_mw_wav_analysis(f, B, J_min)

 INPUTS:
   f     - input MW map (number of pixels is L*(2*L-1))
   B     - Wavelet parameter
   J_min - First wavelet scale to be used
   wavtype - Wavelet type (1: scale-discretised, 2:needlets, 3: cubic splines)

 OUTPUTS:
   f_wav - Struc containing the wavelets

 EXAMPLE:

   see s2let_demo and s2let_test

   f = s2let_read_mw_real_map(file)
   if multires eq 0 then f_wav = s2let_axisym_mw_wav_analysis_real(f, B, J_min) 
      else f_wav = s2let_axisym_mw_wav_analysis_multires_real(f, B, J_min)
   if multires eq 0 then f_rec = s2let_axisym_mw_wav_synthesis_real(f_wav) 
      else f_rec = s2let_axisym_mw_wav_synthesis_multires_real(f_wav)
   L = s2let_get_mw_bandlimit(f)
   J_max = s2let_j_max(L, B)
   s2let_plot_mollweide, f_rec, title='Band-limited map'
   s2let_plot_mollweide, f_wav.scal, title='Scaling map'
   for j=0, J_max-J_min do begin
      s2let_plot_mollweide, f_wav.(j), title='Wavelet map '
      +strtrim(j+1,2)+' on '+strtrim(J_max-J_min+1,2)
   endfor
   
 COMMENTS:
   The resolution/bandlimit L of the input map is automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_axisym_mw_wav_analysis.pro)


S2LET_AXISYM_MW_WAV_ANALYSIS_MULTIRES

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_axisym_mw_wav_analysis_multires

 PURPOSE:
   Compute multiresolution exact wavelet transform of a complex MW map

 CALLING SEQUENCE:
   f_wav = s2let_axisym_mw_wav_analysis_multires(f, B, J_min)

 INPUTS:
   f     - input MW map (number of pixels is L*(2*L-1))
   B     - Wavelet parameter
   J_min - First wavelet scale to be used
   wavtype - Wavelet type (1: scale-discretised, 2:needlets, 3: cubic splines)

 OUTPUTS:
   f_wav - Struc containing the wavelets

 EXAMPLE:

   see s2let_demo and s2let_test

   f = s2let_read_mw_real_map(file)
   if multires eq 0 then f_wav = s2let_axisym_mw_wav_analysis_real(f, B, J_min) 
      else f_wav = s2let_axisym_mw_wav_analysis_multires_real(f, B, J_min)
   if multires eq 0 then f_rec = s2let_axisym_mw_wav_synthesis_real(f_wav) 
      else f_rec = s2let_axisym_mw_wav_synthesis_multires_real(f_wav)
   L = s2let_get_mw_bandlimit(f)
   J_max = s2let_j_max(L, B)
   s2let_plot_mollweide, f_rec, title='Band-limited map'
   s2let_plot_mollweide, f_wav.scal, title='Scaling map'
   for j=0, J_max-J_min do begin
      s2let_plot_mollweide, f_wav.(j), title='Wavelet map '
      +strtrim(j+1,2)+' on '+strtrim(J_max-J_min+1,2)
   endfor
   
 COMMENTS:
   The resolution/bandlimit L of the input map is automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_axisym_mw_wav_analysis_multires.pro)


S2LET_AXISYM_MW_WAV_ANALYSIS_MULTIRES_REAL

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_axisym_mw_wav_analysis_multires_real

 PURPOSE:
   Compute multiresolution exact wavelet transform of a real MW map

 CALLING SEQUENCE:
   f_wav = s2let_axisym_mw_wav_analysis_multires_real(f, B, J_min)

 INPUTS:
   f     - input MW map (number of pixels is L*(2*L-1))
   B     - Wavelet parameter
   J_min - First wavelet scale to be used
   wavtype - Wavelet type (1: scale-discretised, 2:needlets, 3: cubic splines)

 OUTPUTS:
   f_wav - Struc containing the wavelets

 EXAMPLE:

   see s2let_demo and s2let_test

   f = s2let_read_mw_real_map(file)
   if multires eq 0 then f_wav = s2let_axisym_mw_wav_analysis_real(f, B, J_min) 
      else f_wav = s2let_axisym_mw_wav_analysis_multires_real(f, B, J_min)
   if multires eq 0 then f_rec = s2let_axisym_mw_wav_synthesis_real(f_wav) 
      else f_rec = s2let_axisym_mw_wav_synthesis_multires_real(f_wav)
   L = s2let_get_mw_bandlimit(f)
   J_max = s2let_j_max(L, B)
   s2let_plot_mollweide, f_rec, title='Band-limited map'
   s2let_plot_mollweide, f_wav.scal, title='Scaling map'
   for j=0, J_max-J_min do begin
      s2let_plot_mollweide, f_wav.(j), title='Wavelet map '
      +strtrim(j+1,2)+' on '+strtrim(J_max-J_min+1,2)
   endfor
   
 COMMENTS:
   The resolution/bandlimit L of the input map is automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_axisym_mw_wav_analysis_multires_real.pro)


S2LET_AXISYM_MW_WAV_ANALYSIS_REAL

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_axisym_mw_wav_analysis_real

 PURPOSE:
   Compute full resolution exact wavelet transform of a real MW map

 CALLING SEQUENCE:
   f_wav = s2let_axisym_mw_wav_analysis_real(f, B, J_min)

 INPUTS:
   f     - input MW map (number of pixels is L*(2*L-1))
   B     - Wavelet parameter
   J_min - First wavelet scale to be used
   wavtype - Wavelet type (1: scale-discretised, 2:needlets, 3: cubic splines)

 OUTPUTS:
   f_wav - Struc containing the wavelets

 EXAMPLE:

   see s2let_demo and s2let_test

   f = s2let_read_mw_real_map(file)
   if multires eq 0 then f_wav = s2let_axisym_mw_wav_analysis_real(f, B, J_min) 
      else f_wav = s2let_axisym_mw_wav_analysis_multires_real(f, B, J_min)
   if multires eq 0 then f_rec = s2let_axisym_mw_wav_synthesis_real(f_wav) 
      else f_rec = s2let_axisym_mw_wav_synthesis_multires_real(f_wav)
   L = s2let_get_mw_bandlimit(f)
   J_max = s2let_j_max(L, B)
   s2let_plot_mollweide, f_rec, title='Band-limited map'
   s2let_plot_mollweide, f_wav.scal, title='Scaling map'
   for j=0, J_max-J_min do begin
      s2let_plot_mollweide, f_wav.(j), title='Wavelet map '
      +strtrim(j+1,2)+' on '+strtrim(J_max-J_min+1,2)
   endfor
   
 COMMENTS:
   The resolution/bandlimit L of the input map is automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_axisym_mw_wav_analysis_real.pro)


S2LET_AXISYM_MW_WAV_SYNTHESIS

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_axisym_mw_wav_synthesis

 PURPOSE:
   Reconstruct a complex MW map from its full resolution exact wavelet transform

 CALLING SEQUENCE:
   f = s2let_axisym_mw_wav_synthesis(f_wav)

 INPUTS:
   f_wav - input exact wavelet transfrom 

 OUTPUTS:
   f     - Reconstructed MW map

 EXAMPLE:

   see s2let_demo and s2let_test

   f = s2let_read_mw_real_map(file)
   if multires eq 0 then f_wav = s2let_axisym_mw_wav_analysis_real(f, B, J_min) 
      else f_wav = s2let_axisym_mw_wav_analysis_multires_real(f, B, J_min)
   if multires eq 0 then f_rec = s2let_axisym_mw_wav_synthesis_real(f_wav) 
      else f_rec = s2let_axisym_mw_wav_synthesis_multires_real(f_wav)
   L = s2let_get_mw_bandlimit(f)
   J_max = s2let_j_max(L, B)
   s2let_plot_mollweide, f_rec, title='Band-limited map'
   s2let_plot_mollweide, f_wav.scal, title='Scaling map'
   for j=0, J_max-J_min do begin
      s2let_plot_mollweide, f_wav.(j), title='Wavelet map '
      +strtrim(j+1,2)+' on '+strtrim(J_max-J_min+1,2)
   endfor
   
 COMMENTS:
   The wavelet parameters B/J_min are automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_axisym_mw_wav_synthesis.pro)


S2LET_AXISYM_MW_WAV_SYNTHESIS_MULTIRES

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_axisym_mw_wav_synthesis_multires

 PURPOSE:
   Reconstruct a complex MW map from its multiresolution exact wavelet transform

 CALLING SEQUENCE:
   f = s2let_axisym_mw_wav_synthesis_multires(f_wav)

 INPUTS:
   f_wav - input exact wavelet transfrom 

 OUTPUTS:
   f     - Reconstructed MW map

 EXAMPLE:

   see s2let_demo and s2let_test

   f = s2let_read_mw_real_map(file)
   if multires eq 0 then f_wav = s2let_axisym_mw_wav_analysis_real(f, B, J_min) 
      else f_wav = s2let_axisym_mw_wav_analysis_multires_real(f, B, J_min)
   if multires eq 0 then f_rec = s2let_axisym_mw_wav_synthesis_real(f_wav) 
      else f_rec = s2let_axisym_mw_wav_synthesis_multires_real(f_wav)
   L = s2let_get_mw_bandlimit(f)
   J_max = s2let_j_max(L, B)
   s2let_plot_mollweide, f_rec, title='Band-limited map'
   s2let_plot_mollweide, f_wav.scal, title='Scaling map'
   for j=0, J_max-J_min do begin
      s2let_plot_mollweide, f_wav.(j), title='Wavelet map '
      +strtrim(j+1,2)+' on '+strtrim(J_max-J_min+1,2)
   endfor
   
 COMMENTS:
   The wavelet parameters B/J_min are automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_axisym_mw_wav_synthesis_multires.pro)


S2LET_AXISYM_MW_WAV_SYNTHESIS_MULTIRES_REAL

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_axisym_mw_wav_synthesis_multires_real

 PURPOSE:
   Reconstruct a complex MW map from its full resolution exact wavelet transform

 CALLING SEQUENCE:
   f = s2let_axisym_mw_wav_synthesis_multires_real(f_wav)

 INPUTS:
   f_wav - input exact wavelet transfrom 

 OUTPUTS:
   f     - Reconstructed MW map

 EXAMPLE:

   see s2let_demo and s2let_test

   f = s2let_read_mw_real_map(file)
   if multires eq 0 then f_wav = s2let_axisym_mw_wav_analysis_real(f, B, J_min) 
      else f_wav = s2let_axisym_mw_wav_analysis_multires_real(f, B, J_min)
   if multires eq 0 then f_rec = s2let_axisym_mw_wav_synthesis_real(f_wav) 
      else f_rec = s2let_axisym_mw_wav_synthesis_multires_real(f_wav)
   L = s2let_get_mw_bandlimit(f)
   J_max = s2let_j_max(L, B)
   s2let_plot_mollweide, f_rec, title='Band-limited map'
   s2let_plot_mollweide, f_wav.scal, title='Scaling map'
   for j=0, J_max-J_min do begin
      s2let_plot_mollweide, f_wav.(j), title='Wavelet map '
      +strtrim(j+1,2)+' on '+strtrim(J_max-J_min+1,2)
   endfor
   
 COMMENTS:
   The wavelet parameters B/J_min are automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_axisym_mw_wav_synthesis_multires_real.pro)


S2LET_AXISYM_MW_WAV_SYNTHESIS_REAL

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_axisym_mw_wav_synthesis_real

 PURPOSE:
   Reconstruct a real MW map from its full resolution exact wavelet transform

 CALLING SEQUENCE:
   f = s2let_axisym_mw_wav_synthesis_real(f_wav)

 INPUTS:
   f_wav - input exact wavelet transfrom 

 OUTPUTS:
   f     - Reconstructed MW map

 EXAMPLE:

   see s2let_demo and s2let_test

   f = s2let_read_mw_real_map(file)
   if multires eq 0 then f_wav = s2let_axisym_mw_wav_analysis_real(f, B, J_min) 
      else f_wav = s2let_axisym_mw_wav_analysis_multires_real(f, B, J_min)
   if multires eq 0 then f_rec = s2let_axisym_mw_wav_synthesis_real(f_wav) 
      else f_rec = s2let_axisym_mw_wav_synthesis_multires_real(f_wav)
   L = s2let_get_mw_bandlimit(f)
   J_max = s2let_j_max(L, B)
   s2let_plot_mollweide, f_rec, title='Band-limited map'
   s2let_plot_mollweide, f_wav.scal, title='Scaling map'
   for j=0, J_max-J_min do begin
      s2let_plot_mollweide, f_wav.(j), title='Wavelet map '
      +strtrim(j+1,2)+' on '+strtrim(J_max-J_min+1,2)
   endfor
   
 COMMENTS:
   The wavelet parameters B/J_min are automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_axisym_mw_wav_synthesis_real.pro)


S2LET_DEMO1

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_demo1

 PURPOSE:
   Demo 1 : test the exactness of all MW spherical harmonics and wavelet
   transforms for a tomography map of the earth and plot the wavelet maps

 OPTIONAL KEYWORDS:
   B - The wavelet parameter for the test
   J_min - The first wavelet scale to be used for the transform
   multires - multiresolution flag (1 if yes, else 0)
   wavtype - Wavelet type (1: scale-discretised, 2:needlets, 3: cubic splines)
   DEFAULT VALUES: B=3, J_min=2, multires=0, wavtype=1

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_demo1.pro)


S2LET_DEMO2

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_demo2

 PURPOSE:
   Demo : read Healpix map, convert it to MW map, then test the
   exactness of all MW spherical harmonics and wavelet
   transforms and plot the wavelet maps

 OPTIONAL KEYWORDS:
   L - The bandlimit for the spherical harmonic transforms
   B - The wavelet parameter for the test
   J_min - The first wavelet scale to be used for the transform
   multires - multiresolution flag (1 if yes, else 0)
   wavtype - Wavelet type (1: scale-discretised, 2:needlets, 3: cubic splines)
   DEFAULT VALUES: L=128, B=3, J_min=2, multires=0, wavtype=1

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_demo2.pro)


S2LET_DYLIB_EXISTS

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_dylib_exists

 PURPOSE:
   Check if the s2let dynamic library exists somewhere in the path

 CALLING SEQUENCE:
   status = s2let_dylib_exists()

 OUTPUT
   1 if found, 0 otherwise

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_dylib_exists.pro)


S2LET_GET_DYLIB

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_get_dylib

 PURPOSE:
   Get the location/name of the s2let dynamic library

 CALLING SEQUENCE:
   location = s2let_get_dylib()

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_get_dylib.pro)


S2LET_GET_MW_BANDLIMIT

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_get_mw_bandlimit

 PURPOSE:
   Detect the resolution/bandlimit L of an input MW map

 CALLING SEQUENCE:
   L = s2let_get_mw_bandlimit(f)

 INPUT:
   f - array/map of size L*(2*L-1)

 OUTPUT
   L - the resolution/bandlimit

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_get_mw_bandlimit.pro)


S2LET_GET_WAV_BANDLIMIT

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_get_wav_bandlimit

 PURPOSE:
   Compute the bandlimit of the j-th wavelet constructed with parameter B

 CALLING SEQUENCE:
   bl = s2let_get_wav_bandlimit(B, j)

 OUTPUT
   The bandlimit of the j-th wavelet constructed with parameter B

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_get_wav_bandlimit.pro)


S2LET_HEALPIX2MW

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_healpix2mw

 PURPOSE:
   Convert a real Healpix map to an MW map 
   through the spherical harmonic transform 

 CALLING SEQUENCE:
   mwmap = s2let_healpix2mw(hpxmap, lmax=lmax)

 INPUT:
   hpxmap - The Healpix map (npix = 12*nside*nside)
   lmax - The bandlimit of the SHA transform, thus the
       resolution of the MW map

 OUTPUT
   mwmap - The MW map (npix = lmax*(2*lmax-1))

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_healpix2mw.pro)


S2LET_HPX_ALM2MAP_REAL

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_hpx_alm2map_real

 PURPOSE:
   Reconstruct a real Healpix map from its spherical harmonic
   transform (nb: not an exact transform due to healpix)

 CALLING SEQUENCE:
   f = s2let_hpx_alm2map_real(alm, nside)

 INPUT:
   alm - a complex array containing the spherical harmonic transform
         alm_{el em} can be accessed at alm(ind) with ind = el*el + el + em
   nside - resolution for the output map

 OUTPUT
   f - The Healpix map (npix = 12*nside*nside)

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_hpx_alm2map_real.pro)


S2LET_HPX_DEMO

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_hpx_demo

 PURPOSE:
   Demo : test all HPX spherical harmonics and wavelet transforms
   (not exact due to healpix) for a simulated CMB map and plot the wavelet maps

 OPTIONAL KEYWORDS:
   L - The bandlimit for the spherical harmonic transforms
   B - The wavelet parameter for the test
   J_min - The first wavelet scale to be used for the transform
   wavtype - Wavelet type (1: scale-discretised, 2:needlets, 3: cubic splines)
   DEFAULT VALUES: L=192, B=7, J_min=2, wavtype=1

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_hpx_demo.pro)


S2LET_J_MAX

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_j_max

 PURPOSE:
   Compute the maximum wavelet to be used, given a bandlimit L and a
   wavelet parameter B

 CALLING SEQUENCE:
   J_max = s2let_j_max(L, B)

 OUTPUT
   Compute the maximum wavelet to be used, given a bandlimit L and a
   wavelet parameter B;

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_j_max.pro)


S2LET_MW_ALM2MAP

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_mw_alm2map

 PURPOSE:
   Reconstruct a complex MW map from its exact spherical harmonic transform

 CALLING SEQUENCE:
   f = s2let_mw_alm2map(alm)

 INPUT:
   alm - a complex array containing the spherical harmonic transform
         alm_{el em} can be accessed at alm(ind) with ind = el*el + el + em

 OUTPUT
   f - The MW map (npix = L*(2*L-1), L is detected)

 COMMENT:
   The resolution/bandlimit L is automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_mw_alm2map.pro)


S2LET_MW_ALM2MAP_REAL

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_mw_alm2map_real

 PURPOSE:
   Reconstruct a real MW map from its exact spherical harmonic transform

 CALLING SEQUENCE:
   f = s2let_mw_alm2map_real(alm)

 INPUT:
   alm - a complex array containing the spherical harmonic transform
         alm_{el em} can be accessed at alm(ind) with ind = el*el + el + em

 OUTPUT
   f - The MW map (npix = L*(2*L-1), L is detected)

 COMMENT:
   The resolution/bandlimit L is automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_mw_alm2map_real.pro)


S2LET_MW_MAP2ALM_REAL[1]

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_mw_map2alm_real

 PURPOSE:
   Compute the spherical harmonic transform of a real Healpix map
   (nb: not an exact transform due to healpix)

 CALLING SEQUENCE:
   alm = s2let_hpx_map2alm_real(f, lmax)

 INPUTS
   f - The input healpix map (npix = 12*nside*nside, nside is detected)

 OUTPUT:
   alm - a complex array containing the spherical harmonic transform
         alm_{el em} can be accessed at alm(ind) with ind = el*el + el + em

 COMMENT:
   The resolution nside is automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_hpx_map2alm_real.pro)


S2LET_MW_MAP2ALM_REAL[2]

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_mw_map2alm_real

 PURPOSE:
   Compute exact spherical harmonic transform of a complex MW map

 CALLING SEQUENCE:
   alm = s2let_mw_map2alm(f)

 INPUTS
   f - The input MW map (npix = L*(2*L-1), L is detected)

 OUTPUT:
   alm - a complex array containing the spherical harmonic transform
         alm_{el em} can be accessed at alm(ind) with ind = el*el + el + em

 COMMENT:
   The resolution/bandlimit L is automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_mw_map2alm.pro)


S2LET_MW_MAP2ALM_REAL[3]

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_mw_map2alm_real

 PURPOSE:
   Compute exact spherical harmonic transform of a real MW map

 CALLING SEQUENCE:
   alm = s2let_mw_map2alm_real(f)

 INPUTS
   f - The input MW map (npix = L*(2*L-1), L is detected)

 OUTPUT:
   alm - a complex array containing the spherical harmonic transform
         alm_{el em} can be accessed at alm(ind) with ind = el*el + el + em

 COMMENT:
   The resolution/bandlimit L is automatically detected

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_mw_map2alm_real.pro)


S2LET_MW_PIXEL_EDGES

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_mw_pixel_edges

 PURPOSE:
   Get (theta, phi) coordinates for the corners of the i-th pixel in
   the MW sampling

 CALLING SEQUENCE:
   arr = s2let_mw_pixel_edges(L, i)

 INPUTS
   L - resolution/bandlimit of the MW map
   i - pixel index

 OUTPUT:
   arr - an array containing four numbers: theta1, theta2, phi1 and phi2
       which are the locations of the corners of the pixel

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_mw_pixel_edges.pro)


S2LET_MW_PLOT_MOLLWEIDE

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_mw_plot_mollweide

 PURPOSE:
   Plot a real MW map (from data of FITS file) under Mollweide projection

 CALLING SEQUENCE:
   s2let_mw_plot_mollweide, mapfile, title='My map'

 INPUTS
   maporfile - filename for the FITS or map read with s2let_read_mw_real_map

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_mw_plot_mollweide.pro)


S2LET_MW_SAMPLING

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_mw_sampling

 PURPOSE:
   Compute the coordinates of the nodes of the MW sampling for a
   given bandlimit/resolution

 CALLING SEQUENCE:
   s2let_mw_sampling, L, thetas, phis

 INPUT:
   L - resolution/bandlimit

 OUTPUT
   thetas - The theta coordinates of the nodes
   phis - The phi coordinates of the nodes

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_mw_sampling.pro)


S2LET_MW_WRITE_REAL_MAP

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_mw_write_real_map

 PURPOSE:
   Write a real MW map to a FITS file

 CALLING SEQUENCE:
   s2let_mw_write_real_map, map, file

 INPUTS:
   map - input MW map (npix=L*(2*L-1), L is detected)
   file - filename for the FITS

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_mw_write_real_map.pro)


S2LET_READ_MW_REAL_MAP

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_read_mw_real_map

 PURPOSE:
   Read a real MW map from a FITS file

 CALLING SEQUENCE:
   f = s2let_mw_read_real_map(file)

 INPUTS
   file - filename for the FITS

 OUTPUT:
   map - MW map (npix=L*(2*L-1), L is detected)

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_mw_read_real_map.pro)


S2LET_TEST

[Previous Routine] [Next Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_test

 PURPOSE:
   Test the exactness of all MW spherical harmonics and wavelet
   transform for real and complex signals on the sphere

 OPTIONAL KEYWORDS:
   L - The bandlimit/resolution for the test
   B - The wavelet parameter for the test
   J_min - The first wavelet scale to be used for the transform
   multires - multiresolution flag (1 if yes, else 0)
   wavtype - Wavelet type (1: scale-discretised, 2:needlets, 3: cubic splines)
   DEFAULT VALUES: L=64, B=3, J_min=2, multires=0, wavtype=1

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_test.pro)


S2LET_VALID_WAV_PARAMETERS

[Previous Routine] [List of Routines]
 S2LET package - Copyright (C) 2012 
 Boris Leistedt & Jason McEwen

 NAME:
   s2let_valid_wav_parameters

 PURPOSE:
   Check and test a set of wavelet parameters

 CALLING SEQUENCE:
   status = s2let_valid_wav_parameters(B, L, J_min)

 INPUTS:
   B     - Wavelet parameter
   L     - Band-limit to be used for the spherical harmonic transforms
   J_min - First wavelet scale to be used
   wavtype - Wavelet type (1: scale-discretised, 2:needlets, 3: cubic splines)

(See /Users/bl/Dropbox/Wavelets/s2let/src/main/idl/s2let_valid_wav_parameters.pro)