Module s2_sky_mod

module s2_sky_mod

        ! Uses
    use s2_types_mod, only: s2_sp, s2_spc, s2_dp, pi
    use s2_error_mod
    use s2_vect_mod
    use s2_pl_mod

        ! Types
    public type s2_sky

        ! Variables
    integer, public, parameter :: S2_SKY_RING = 1
    integer, public, parameter :: S2_SKY_NEST = 2
    integer, public, parameter :: S2_SKY_FUN_TYPE_SPHERE = 1
    integer, public, parameter :: S2_SKY_FUN_TYPE_PLANE = 2
    integer, public, parameter :: S2_SKY_DEFAULT_FTYPE = S2_SKY_FUN_TYPE_SPHERE
    integer, public, parameter :: S2_SKY_FILE_TYPE_MAP = 1
    integer, public, parameter :: S2_SKY_FILE_TYPE_SKY = 2
    integer, public, parameter :: S2_SKY_FILE_TYPE_ALM = 3

        ! Interfaces
    public interface s2_sky_init

        ! Subroutines and functions
    private function s2_sky_init_empty (nside, pix_scheme, lmax, mmax) result (sky)
    private function s2_sky_init_map (map, nside, pix_scheme, lmax, mmax) result (sky)
    private function s2_sky_init_alm (alm, lmax, mmax, nside, pix_scheme) result (sky)
    private function s2_sky_init_ab (x_ab, interp, nside, pix_scheme, lmax, mmax, sdw) result (sky)
    private function s2_sky_init_fun (fun, nside, pix_scheme, lmax, mmax, param, fun_type_in) result (sky)
    private function s2_sky_init_file (filename, file_type, extension) result (sky)
    private function s2_sky_init_file_map (filename, extension_in) result (sky)
    private function s2_sky_init_file_alm (filename) result (sky)
    private function s2_sky_init_file_sky (filename) result (sky)
    private function s2_sky_init_copy (orig) result (copy)
    public subroutine s2_sky_free (sky)
    public subroutine s2_sky_remove_map (sky)
    public subroutine s2_sky_compute_alm (sky, lmax, mmax, message)
    public subroutine s2_sky_compute_alm_iter (sky, iter_order, lmax, mmax, message)
    public subroutine s2_sky_compute_map (sky, nside, message)
    public subroutine s2_sky_map_convert (sky, pix_scheme)
    public subroutine s2_sky_conv (sky, beam)
    public subroutine s2_sky_offset (sky, offset)
    public subroutine s2_sky_scale (sky, scale)
    public function s2_sky_add (a, b, subtract) result (c)
    public function s2_sky_product (a, b, divide) result (c)
    public subroutine s2_sky_thres (sky, thres_lower, thres_upper)
    public subroutine s2_sky_thres_abs (sky, thres)
    public function s2_sky_error_twonorm (a, b) result (er)
    public function s2_sky_rms (a) result (rms)
    public subroutine s2_sky_dilate (sky, a, b, norm_preserve)
    public subroutine s2_sky_rotate (sky, alpha, beta, gamma)
    public function s2_sky_power_map (sky) result (power)
    public function s2_sky_power_alm (sky) result (power)
    public function s2_sky_azimuthal_bl (sky, cutoff_prop_in) result (mmax_min)
    public function s2_sky_admiss (sky) result (admiss)
    public subroutine s2_sky_admiss_dil (sky, dilation1, dilation2, admiss, norm_preserve_in)
    public subroutine s2_sky_extract_ab_fsht (sky, xtp, L)
    public subroutine s2_sky_extract_ab (sky, x_ab)
    private function s2_sky_interp_ab (x_ab, alpha, beta, interp, sdw) result (interp_val)
    public subroutine s2_sky_downsample (sky, nside_down, mask)
    private subroutine s2_sky_down_nsideby2 (sky, mask_in)
    public subroutine s2_sky_upsample (sky, nside_up)
    private subroutine s2_sky_up_nsideby2 (sky)
    private subroutine s2_sky_valid_sizes (sky)
    public subroutine s2_sky_draw_dot (sky, alpha, beta, large)
    public subroutine s2_sky_write_map_file (sky, filename, comment)
    public subroutine s2_sky_write_alm_file (sky, filename, comment)
    public subroutine s2_sky_io_fits_write (filename, sky, comment)
    private subroutine s2_sky_io_fits_read (filename, sky)
    private subroutine s2_sky_io_fits_error_check (status, halt)
    private subroutine s2_sky_io_fits_exists (filename, status, exists)
    private subroutine s2_sky_io_fits_del (filename, status)
    public subroutine s2_sky_set_lmax (sky, lmax, mmax_in)
    public subroutine s2_sky_set_nside (sky, nside)
    public function s2_sky_get_init (sky) result (init)
    public function s2_sky_get_nside (sky) result (nside)
    public function s2_sky_get_npix (sky) result (npix)
    public function s2_sky_get_lmax (sky) result (lmax)
    public function s2_sky_get_mmax (sky) result (mmax)
    public function s2_sky_get_pix_scheme (sky) result (pix_scheme)
    public subroutine s2_sky_get_map (sky, map)
    public function s2_sky_get_map_pix (sky, ipix) result (map_pix)
    public subroutine s2_sky_get_alm (sky, alm)
    public function s2_sky_get_cl (sky) result (cl)
    public subroutine s2_sky_get_cm (sky, cm)
    public function s2_sky_get_map_status (sky) result (map_status)
    public function s2_sky_get_alm_status (sky) result (alm_status)
    public function s2_sky_get_n_param (sky) result (n_param)
    public subroutine s2_sky_get_param (sky, param)

end module s2_sky_mod
Provides functionality to support and manipulate a function defined on the sky(/sphere). Support is provided to representing the sky in both real (map) space and in harmonic (alm) space, and also to convert between the two representations. One may also scale, dilate, rotate, downsample and add functions defined on the sky. The sky pixelisation used is currently based on HEALPix.

Author: J. D. McEwen (mcewen@mrao.cam.ac.uk)

Version: 0.1 August 2004


Description of Types

s2_sky

public type s2_sky
    private
    logical :: init = .false.
    integer :: nside = 0
    integer :: npix = 0
    integer :: lmax = 0
    integer :: mmax = 0
    integer :: pix_scheme = S2_SKY_RING
    real (kind=s2_sp), allocatable, dimension (:) :: map
    complex (kind=s2_spc), allocatable, dimension (:,:) :: alm
    logical :: map_status = .false.
    logical :: alm_status = .false.
    real (kind=s2_sp), allocatable, dimension (:) :: param
    integer :: n_param = 0
end type s2_sky

Description of Variables

S2_SKY_RING

integer, public, parameter :: S2_SKY_RING = 1
Healpix pixelisation ring type.

S2_SKY_NEST

integer, public, parameter :: S2_SKY_NEST = 2
Healpix pixelisation nest type.

S2_SKY_FUN_TYPE_SPHERE

integer, public, parameter :: S2_SKY_FUN_TYPE_SPHERE = 1
Specifier to indicate function passed to construct sky is defined on the sphere.

S2_SKY_FUN_TYPE_PLANE

integer, public, parameter :: S2_SKY_FUN_TYPE_PLANE = 2
Specifier to indicate function passed to construct sky is defined on the plane (and will be numerically projected onto the sphere).

S2_SKY_DEFAULT_FTYPE

integer, public, parameter :: S2_SKY_DEFAULT_FTYPE = S2_SKY_FUN_TYPE_SPHERE
Default function type if fun_type flag is not specified.

S2_SKY_FILE_TYPE_MAP

integer, public, parameter :: S2_SKY_FILE_TYPE_MAP = 1
Map fits file type.

S2_SKY_FILE_TYPE_SKY

integer, public, parameter :: S2_SKY_FILE_TYPE_SKY = 2
Full s2_sky file type.

S2_SKY_FILE_TYPE_ALM

integer, public, parameter :: S2_SKY_FILE_TYPE_ALM = 3
Alm fits file type.

Description of Interfaces

s2_sky_init

public interface s2_sky_init
    module procedure s2_sky_init_map
    module procedure s2_sky_init_alm
    module procedure s2_sky_init_ab
    module procedure s2_sky_init_fun
    module procedure s2_sky_init_file
    module procedure s2_sky_init_copy
end interface s2_sky_init

Description of Subroutines and Functions

s2_sky_init_empty

private function s2_sky_init_empty (nside, pix_scheme, lmax, mmax) result (sky)
    integer, optional, intent(in) :: nside
    integer, optional, intent(in) :: pix_scheme
    integer, optional, intent(in) :: lmax
    integer, optional, intent(in) :: mmax
    type (s2_sky) :: sky
    ! Calls: s2_error, s2_sky_valid_sizes
end function s2_sky_init_empty
Initialised an empty sky at a specified resolution. Neither the map or alm sky representations are allocated. They are allocated when required.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_init_map

private function s2_sky_init_map (map, nside, pix_scheme, lmax, mmax) result (sky)
    real (kind=s2_sp), intent(in), dimension (:) :: map
    integer, intent(in) :: nside
    integer, intent(in) :: pix_scheme
    integer, optional, intent(in) :: lmax
    integer, optional, intent(in) :: mmax
    type (s2_sky) :: sky
    ! Calls: s2_error
end function s2_sky_init_map
Initalise a sky from a specified map.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_init_alm

private function s2_sky_init_alm (alm, lmax, mmax, nside, pix_scheme) result (sky)
    complex (kind=s2_spc), intent(in), dimension (:,:) :: alm
    integer, intent(in) :: lmax
    integer, intent(in) :: mmax
    integer, optional, intent(in) :: nside
    integer, optional, intent(in) :: pix_scheme
    type (s2_sky) :: sky
    ! Calls: s2_error
end function s2_sky_init_alm
Initalise a sky from a specified alm.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_init_ab

private function s2_sky_init_ab (x_ab, interp, nside, pix_scheme, lmax, mmax, sdw) result (sky)
    real (kind=s2_sp), intent(in), dimension (:,:) :: x_ab
    logical, intent(in) :: interp
    integer, intent(in) :: nside
    integer, intent(in) :: pix_scheme
    integer, optional, intent(in) :: lmax
    integer, optional, intent(in) :: mmax
    logical, optional, intent(in) :: sdw
    type (s2_sky) :: sky
    ! Calls: pix2ang_nest, pix2ang_ring, s2_error
end function s2_sky_init_ab
Initialise a sky from the ecp (equi-sampled) alpha-beta array defined on the sphere.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 - November 2004

s2_sky_init_fun

private function s2_sky_init_fun (fun, nside, pix_scheme, lmax, mmax, param, fun_type_in) result (sky)
    interface fun
        function fun (theta, phi, param) result (val)
            real (kind=s2_sp), intent(in) :: theta
            real (kind=s2_sp), intent(in) :: phi
            real (kind=s2_sp), optional, intent(in), dimension (:) :: param
            real (kind=s2_sp) :: val
        end function fun
    end interface fun
    integer, intent(in) :: nside
    integer, intent(in) :: pix_scheme
    integer, optional, intent(in) :: lmax
    integer, optional, intent(in) :: mmax
    real (kind=s2_sp), optional, intent(in), dimension (:) :: param
    integer, optional, intent(in) :: fun_type_in
    type (s2_sky) :: sky
    ! Calls: pix2ang_nest, pix2ang_ring, s2_error
end function s2_sky_init_fun
Initialise a sky from a template function.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_init_file

private function s2_sky_init_file (filename, file_type, extension) result (sky)
    character (len=*), intent(in) :: filename
    integer, intent(in) :: file_type
    integer, optional, intent(in) :: extension
    type (s2_sky) :: sky
    ! Calls: s2_error
end function s2_sky_init_file
Wrapper to initialise a sky from either a map fits file or a full s2_sky fits file.

Variables:

Author: J. D. McEwen

Version: 0.1 April 2005

s2_sky_init_file_map

private function s2_sky_init_file_map (filename, extension_in) result (sky)
    character (len=*), intent(in) :: filename
    integer, optional, intent(in) :: extension_in
    type (s2_sky) :: sky
    ! Calls: input_map, s2_error
end function s2_sky_init_file_map
Initialise a sky from a map file.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_init_file_alm

private function s2_sky_init_file_alm (filename) result (sky)
    character (len=*), intent(in) :: filename
    type (s2_sky) :: sky
    ! Calls: fits2alms, s2_error
end function s2_sky_init_file_alm
Initialise a sky from a fits alm file.

Variables:

Author: J. D. McEwen

Version: 0.1 November 2005

s2_sky_init_file_sky

private function s2_sky_init_file_sky (filename) result (sky)
    character (len=*), intent(in) :: filename
    type (s2_sky) :: sky
    ! Calls: s2_error, s2_sky_io_fits_read
end function s2_sky_init_file_sky
Wrapper to initialise a sky data structure from a full s2_sky file. The sky structure is read and initialised by the routine s2_sky_io_fits_read.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 - April 2005

s2_sky_init_copy

private function s2_sky_init_copy (orig) result (copy)
    type (s2_sky), intent(in) :: orig
    type (s2_sky) :: copy
    ! Calls: s2_error
end function s2_sky_init_copy
Initialise a sky as a copy of another sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_free

public subroutine s2_sky_free (sky)
    type (s2_sky), intent(inout) :: sky
    ! Calls: s2_error
end subroutine s2_sky_free
Free all data associated with an initialised sky and reset all other attributes.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_remove_map

public subroutine s2_sky_remove_map (sky)
    type (s2_sky), intent(inout) :: sky
    ! Calls: s2_error
end subroutine s2_sky_remove_map
Free all data associated with a map and reset all map related attributes.

Variables:

Author: J. D. McEwen

Version: 0.1 November 2007

s2_sky_compute_alm

public subroutine s2_sky_compute_alm (sky, lmax, mmax, message)
    type (s2_sky), intent(inout) :: sky
    integer, optional, intent(in) :: lmax
    integer, optional, intent(in) :: mmax
    logical, optional, intent(in) :: message
    ! Calls: map2alm, s2_error, s2_sky_map_convert, s2_sky_set_lmax, s2_sky_valid_sizes
end subroutine s2_sky_compute_alm
Compute the alms for a sky from the map representation. Errors occur if the initialised sky does not have a map defined.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_compute_alm_iter

public subroutine s2_sky_compute_alm_iter (sky, iter_order, lmax, mmax, message)
    type (s2_sky), intent(inout) :: sky
    integer, intent(in) :: iter_order
    integer, optional, intent(in) :: lmax
    integer, optional, intent(in) :: mmax
    logical, optional, intent(in) :: message
    ! Calls: alm2map, map2alm, s2_error, s2_sky_map_convert, s2_sky_set_lmax, s2_sky_valid_sizes
end subroutine s2_sky_compute_alm_iter
Compute the alms for a sky from the map representation using iteration. Errors occur if the initialised sky does not have a map defined.

Variables:

Author: J. D. McEwen

Version: 0.1 November 2007

s2_sky_compute_map

public subroutine s2_sky_compute_map (sky, nside, message)
    type (s2_sky), intent(inout) :: sky
    integer, optional, intent(in) :: nside
    logical, optional, intent(in) :: message
    ! Calls: alm2map, s2_error, s2_sky_set_nside, s2_sky_valid_sizes
end subroutine s2_sky_compute_map
Compute the map for a sky from the alm representation. Errors occur if the initialised sky does not have an alm defined.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_map_convert

public subroutine s2_sky_map_convert (sky, pix_scheme)
    type (s2_sky), intent(inout) :: sky
    integer, intent(in) :: pix_scheme
    ! Calls: convert_nest2ring, convert_ring2nest, s2_error
end subroutine s2_sky_map_convert
Convert the map pixelisation scheme to that specified. If the map is already in the required pixelisation scheme then do nothing.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_conv

public subroutine s2_sky_conv (sky, beam)
    type (s2_sky), intent(inout) :: sky
    type (s2_pl), intent(in) :: beam
    ! Calls: s2_error, s2_pl_conv, s2_sky_compute_map
end subroutine s2_sky_conv
Convolve a sky with a beam.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 May 2005

s2_sky_offset

public subroutine s2_sky_offset (sky, offset)
    type (s2_sky), intent(inout) :: sky
    real (kind=s2_sp), intent(in) :: offset
    ! Calls: s2_error, s2_sky_compute_alm, s2_sky_compute_map
end subroutine s2_sky_offset
Add an offset to the sky.

Variables:

Author: J. D. McEwen

Version: 0.1 June 2006

s2_sky_scale

public subroutine s2_sky_scale (sky, scale)
    type (s2_sky), intent(inout) :: sky
    real (kind=s2_sp), intent(in) :: scale
    ! Calls: s2_error
end subroutine s2_sky_scale
Scale both the sky map and alm.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_add

public function s2_sky_add (a, b, subtract) result (c)
    type (s2_sky), intent(inout) :: a
    type (s2_sky), intent(inout) :: b
    logical, optional, intent(in) :: subtract
    type (s2_sky) :: c
    ! Calls: s2_error, s2_sky_compute_map, s2_sky_map_convert
end function s2_sky_add
Add the maps of two skys.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_product

public function s2_sky_product (a, b, divide) result (c)
    type (s2_sky), intent(inout) :: a
    type (s2_sky), intent(inout) :: b
    logical, optional, intent(in) :: divide
    type (s2_sky) :: c
    ! Calls: s2_error, s2_sky_compute_map, s2_sky_map_convert
end function s2_sky_product
Multiply two sky maps together.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2005

s2_sky_thres

public subroutine s2_sky_thres (sky, thres_lower, thres_upper)
    type (s2_sky), intent(inout) :: sky
    real (kind=s2_sp), intent(in) :: thres_lower
    real (kind=s2_sp), intent(in) :: thres_upper
    ! Calls: s2_error, s2_sky_compute_alm, s2_sky_compute_map
end subroutine s2_sky_thres
Threshold all values outside range to specified limits.

Variables:

Author: J. D. McEwen

Version: 0.1 April 2007

s2_sky_thres_abs

public subroutine s2_sky_thres_abs (sky, thres)
    type (s2_sky), intent(inout) :: sky
    real (kind=s2_sp), intent(in) :: thres
    ! Calls: s2_error, s2_sky_compute_alm, s2_sky_compute_map
end subroutine s2_sky_thres_abs
Threshold all values below (in absolute value) specified threshold.

Variables:

Author: J. D. McEwen

Version: 0.1 April 2007

s2_sky_error_twonorm

public function s2_sky_error_twonorm (a, b) result (er)
    type (s2_sky), intent(inout) :: a
    type (s2_sky), intent(inout) :: b
    real (kind=s2_sp) :: er
    ! Calls: s2_error, s2_sky_compute_map, s2_sky_map_convert
end function s2_sky_error_twonorm
Compute two-norm error (i.e. mean square error) or two sky maps.

Variables:

Author: J. D. McEwen

Version: 0.1 April 2007

s2_sky_rms

public function s2_sky_rms (a) result (rms)
    type (s2_sky), intent(inout) :: a
    real (kind=s2_sp) :: rms
    ! Calls: s2_error, s2_sky_compute_map
end function s2_sky_rms
Compute root-mean squared value of sky (same as error_twonorm with second sky set to zero).

Variables:

Author: J. D. McEwen

Version: 0.1 April 2007

s2_sky_dilate

public subroutine s2_sky_dilate (sky, a, b, norm_preserve)
    type (s2_sky), intent(inout) :: sky
    real (kind=s2_sp), intent(in) :: a
    real (kind=s2_sp), intent(in) :: b
    logical, intent(in) :: norm_preserve
    ! Calls: ang2pix_nest, ang2pix_ring, pix2ang_nest, pix2ang_ring, s2_error
end subroutine s2_sky_dilate
Dilate the map of a sky. The dilation may be performed in two orthogonal directions. To perform the usual isotropic dilation set a=b.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_rotate

public subroutine s2_sky_rotate (sky, alpha, beta, gamma)
    type (s2_sky), intent(inout) :: sky
    real (kind=s2_sp) :: alpha
    real (kind=s2_sp) :: beta
    real (kind=s2_sp) :: gamma
    ! Calls: ang2pix_nest, ang2pix_ring, pix2ang_nest, pix2ang_ring, s2_error, s2_vect_free, s2_vect_rotate
end subroutine s2_sky_rotate
Rotate the map of a sky. The rotation is perform in real space. The rotation is defined by the rotation routine in the s2_vect_mod class.

Variables

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_power_map

public function s2_sky_power_map (sky) result (power)
    type (s2_sky), intent(in) :: sky
    real (kind=s2_sp) :: power
    ! Calls: s2_error
end function s2_sky_power_map
Compute power of function on sphere in map space.

Variables:

Author: J. D. McEwen

Version: 0.1 May 2006

s2_sky_power_alm

public function s2_sky_power_alm (sky) result (power)
    type (s2_sky), intent(in) :: sky
    real (kind=s2_sp) :: power
    ! Calls: s2_error, s2_pl_free
end function s2_sky_power_alm
Compute power of function on sphere in alm space.

Variables:

Author: J. D. McEwen

Version: 0.1 May 2006

s2_sky_azimuthal_bl

public function s2_sky_azimuthal_bl (sky, cutoff_prop_in) result (mmax_min)
    type (s2_sky), intent(in) :: sky
    real (kind=s2_sp), optional, intent(in) :: cutoff_prop_in
    integer :: mmax_min
    ! Calls: s2_error, s2_pl_free, s2_pl_get_spec, s2_sky_get_cm
end function s2_sky_azimuthal_bl
Find azimuthal band limit of sky. Finds lowest m' value such that cutoff_prop*100 percent of the cm power is contained in the alms with m index below m'.

Notes:

Variables

Author: J. D. McEwen

Version: 0.1 May 2005

s2_sky_admiss

public function s2_sky_admiss (sky) result (admiss)
    type (s2_sky), intent(in) :: sky
    real (kind=s2_sp) :: admiss
    ! Calls: pix2ang_nest, pix2ang_ring, s2_error
end function s2_sky_admiss
Calculte the normalised numerical admissibility of a wavelet defined on the sphere.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_admiss_dil

public subroutine s2_sky_admiss_dil (sky, dilation1, dilation2, admiss, norm_preserve_in)
    type (s2_sky), intent(in) :: sky
    real (kind=s2_sp), intent(in), dimension (:) :: dilation1
    real (kind=s2_sp), intent(in), dimension (:) :: dilation2
    real (kind=s2_sp), intent(out), dimension (:) :: admiss
    logical, optional, intent(in) :: norm_preserve_in
    ! Calls: s2_error, s2_sky_dilate, s2_sky_free
end subroutine s2_sky_admiss_dil
Compute admissibility for function defined on sphere for a range of dilations.

Variables:

Author: J. D. McEwen

Version: 0.1 - November 2004

s2_sky_extract_ab_fsht

public subroutine s2_sky_extract_ab_fsht (sky, xtp, L)
    type (s2_sky), intent(in) :: sky
    real (kind=s2_sp), intent(out), dimension (0:L,0:2*L) :: xtp
    integer, intent(in) :: L
    ! Calls: ang2pix_nest, ang2pix_ring, s2_error
end subroutine s2_sky_extract_ab_fsht
Extra an ecp (equispaced) sampled theta-phi array over the sphere for the grid used for the Fast Spherical Harmonic Transform.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 - April 2008

s2_sky_extract_ab

public subroutine s2_sky_extract_ab (sky, x_ab)
    type (s2_sky), intent(in) :: sky
    real (kind=s2_sp), intent(out), dimension (:,:) :: x_ab
    ! Calls: ang2pix_nest, ang2pix_ring, s2_error
end subroutine s2_sky_extract_ab
Extra an ecp (equispaced) sampled alpha-beta array over the sphere corresponding to the given sky.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 - November 2004

s2_sky_interp_ab

private function s2_sky_interp_ab (x_ab, alpha, beta, interp, sdw) result (interp_val)
    real (kind=s2_sp), intent(in), dimension (:,:) :: x_ab
    real (kind=s2_dp), intent(inout) :: alpha
    real (kind=s2_dp), intent(inout) :: beta
    logical, intent(in) :: interp
    logical, optional, intent(in) :: sdw
    real (kind=s2_sp) :: interp_val
end function s2_sky_interp_ab
Interpolate a pixelised ab array for any continuous value in the appropriate angle ranges.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 - November 2004

s2_sky_downsample

public subroutine s2_sky_downsample (sky, nside_down, mask)
    type (s2_sky), intent(inout) :: sky
    integer, intent(in) :: nside_down
    logical, optional, intent(in) :: mask
    ! Calls: s2_error, s2_sky_down_nsideby2
end subroutine s2_sky_downsample
Downsample a sky map nside to the specified (lower!) nside.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_down_nsideby2

private subroutine s2_sky_down_nsideby2 (sky, mask_in)
    type (s2_sky), intent(inout) :: sky
    logical, optional, intent(in) :: mask_in
    ! Calls: s2_error, s2_sky_map_convert
end subroutine s2_sky_down_nsideby2
Downsample a sky map nside by a factor of two.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_upsample

public subroutine s2_sky_upsample (sky, nside_up)
    type (s2_sky), intent(inout) :: sky
    integer, intent(in) :: nside_up
    ! Calls: s2_error, s2_sky_up_nsideby2
end subroutine s2_sky_upsample
Up-sample a sky map nside to the specified (greater!) nside. Map is up-sampled by simply setting all sub pixels in a nested region to the value of the parent pixel).

Variables:

Author: J. D. McEwen

Version: 0.1 August 2005

s2_sky_up_nsideby2

private subroutine s2_sky_up_nsideby2 (sky)
    type (s2_sky), intent(inout) :: sky
    ! Calls: s2_error, s2_sky_map_convert
end subroutine s2_sky_up_nsideby2
Up-sample a sky map nside by a factor of two (simply set all sub pixels in a nested region to the value of the parent pixel).

Variables:

Author: J. D. McEwen

Version: 0.1 August 2005

s2_sky_valid_sizes

private subroutine s2_sky_valid_sizes (sky)
    type (s2_sky), intent(in) :: sky
    ! Calls: s2_error
end subroutine s2_sky_valid_sizes
Check the nside, lmax and mmax resolutions are valid. Invalid if any are negative. Also warning given if lmax > 3*nside of mmax < lmax

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_draw_dot

public subroutine s2_sky_draw_dot (sky, alpha, beta, large)
    type (s2_sky), intent(inout) :: sky
    real (kind=s2_sp), intent(in), dimension (:) :: alpha
    real (kind=s2_sp), intent(in), dimension (:) :: beta
    logical, optional, intent(in) :: large
    ! Calls: ang2pix_nest, neighbours_nest, s2_error, s2_sky_map_convert
end subroutine s2_sky_draw_dot
Draw dots on sky map at specified positions.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_write_map_file

public subroutine s2_sky_write_map_file (sky, filename, comment)
    type (s2_sky), intent(in) :: sky
    character (len=*), intent(in) :: filename
    character (len=*), optional, intent(in) :: comment
    ! Calls: add_card, s2_error, write_bintab
end subroutine s2_sky_write_map_file
Write a sky map to a fits file.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_write_alm_file

public subroutine s2_sky_write_alm_file (sky, filename, comment)
    type (s2_sky), intent(in) :: sky
    character (len=*), intent(in) :: filename
    character (len=*), optional, intent(in) :: comment
    ! Calls: add_card, alms2fits, s2_error
end subroutine s2_sky_write_alm_file
Write sky alms to a fits file.

Variables:

Author: J. D. McEwen

Version: 0.1 November 2005

s2_sky_io_fits_write

public subroutine s2_sky_io_fits_write (filename, sky, comment)
    character (len=*), intent(in) :: filename
    type (s2_sky), intent(in) :: sky
    character (len=*), optional, intent(in) :: comment
    ! Calls: ftclos, ftfiou, ftgiou, ftibin, ftiimg, ftinit, ftp2de, ftpcle, ftpcom, ftpdat, ftphpr, ftpkyj, ftpkyl, ftpkys, s2_error, s2_sky_io_fits_error_check, s2_sky_io_fits_exists
end subroutine s2_sky_io_fits_write
Write a s2_sky object to a fits file. All s2_sky data is writtend to the file (i.e. both the map and alms are written if present, as well as any sky function parameters and all other s2_sky variables.

Variables:

Author: J. D. McEwen

Version: 0.1 - April 2005

s2_sky_io_fits_read

private subroutine s2_sky_io_fits_read (filename, sky)
    character (len=*), intent(in) :: filename
    type (s2_sky), intent(out) :: sky
    ! Calls: ftclos, ftfiou, ftg2de, ftgcve, ftgiou, ftgkyj, ftgkyl, ftmahd, ftopen, ftthdu, s2_error, s2_sky_io_fits_error_check, s2_sky_io_fits_exists
end subroutine s2_sky_io_fits_read
Reads a full s2_sky fits file and allocates a new sky structure with the data read.

Variables:

Author: J. D. McEwen

Version: 0.1 - April 2005

s2_sky_io_fits_error_check

private subroutine s2_sky_io_fits_error_check (status, halt)
    integer, intent(inout) :: status
    logical, intent(in) :: halt
    ! Calls: ftgerr, ftgmsg
end subroutine s2_sky_io_fits_error_check
Check if a fits error has occured and print error message. Halt program execution if halt flag is set.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 - November 2004

s2_sky_io_fits_exists

private subroutine s2_sky_io_fits_exists (filename, status, exists)
    character (len=*), intent(in) :: filename
    integer, intent(inout) :: status
    logical, intent(out) :: exists
    ! Calls: ftclos, ftfiou, ftgiou, ftopen, s2_sky_io_fits_error_check
end subroutine s2_sky_io_fits_exists
Check if a fits file exists.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 - November 2004

s2_sky_io_fits_del

private subroutine s2_sky_io_fits_del (filename, status)
    character (len=*), intent(in) :: filename
    integer, intent(inout) :: status
    ! Calls: ftcmsg, ftdelt, ftfiou, ftgiou, ftopen
end subroutine s2_sky_io_fits_del
Delete a fits file.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 - November 2004

s2_sky_set_lmax

public subroutine s2_sky_set_lmax (sky, lmax, mmax_in)
    type (s2_sky), intent(inout) :: sky
    integer, intent(in) :: lmax
    integer, optional, intent(in) :: mmax_in
    ! Calls: s2_sky_valid_sizes
end subroutine s2_sky_set_lmax
Set a sky lmax. If alm is calculated it is removed.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_set_nside

public subroutine s2_sky_set_nside (sky, nside)
    type (s2_sky), intent(inout) :: sky
    integer, intent(in) :: nside
    ! Calls: s2_sky_valid_sizes
end subroutine s2_sky_set_nside
Set a sky nside. If map is calculated it is removed.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_init

public function s2_sky_get_init (sky) result (init)
    type (s2_sky), intent(in) :: sky
    logical :: init
end function s2_sky_get_init
Get init variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_nside

public function s2_sky_get_nside (sky) result (nside)
    type (s2_sky), intent(in) :: sky
    integer :: nside
    ! Calls: s2_error
end function s2_sky_get_nside
Get nside variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_npix

public function s2_sky_get_npix (sky) result (npix)
    type (s2_sky), intent(in) :: sky
    integer :: npix
    ! Calls: s2_error
end function s2_sky_get_npix
Get nside variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_lmax

public function s2_sky_get_lmax (sky) result (lmax)
    type (s2_sky), intent(in) :: sky
    integer :: lmax
    ! Calls: s2_error
end function s2_sky_get_lmax
Get lmax variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_mmax

public function s2_sky_get_mmax (sky) result (mmax)
    type (s2_sky), intent(in) :: sky
    integer :: mmax
    ! Calls: s2_error
end function s2_sky_get_mmax
Get mmax variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_pix_scheme

public function s2_sky_get_pix_scheme (sky) result (pix_scheme)
    type (s2_sky), intent(in) :: sky
    integer :: pix_scheme
    ! Calls: s2_error
end function s2_sky_get_pix_scheme
Get pix_scheme variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_map

public subroutine s2_sky_get_map (sky, map)
    type (s2_sky), intent(in) :: sky
    real (kind=s2_sp), intent(out), dimension (:) :: map
    ! Calls: s2_error
end subroutine s2_sky_get_map
Get map variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_map_pix

public function s2_sky_get_map_pix (sky, ipix) result (map_pix)
    type (s2_sky), intent(in) :: sky
    integer, intent(in) :: ipix
    real (kind=s2_sp) :: map_pix
    ! Calls: s2_error
end function s2_sky_get_map_pix
Get map pixel from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_alm

public subroutine s2_sky_get_alm (sky, alm)
    type (s2_sky), intent(in) :: sky
    complex (kind=s2_spc), intent(out), dimension (:,:) :: alm
    ! Calls: s2_error
end subroutine s2_sky_get_alm
Get alm variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_cl

public function s2_sky_get_cl (sky) result (cl)
    type (s2_sky), intent(in) :: sky
    type (s2_pl) :: cl
    ! Calls: s2_error
end function s2_sky_get_cl
Compute cl for the given sky.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 January 2005

s2_sky_get_cm

public subroutine s2_sky_get_cm (sky, cm)
    type (s2_sky), intent(in) :: sky
    real (kind=s2_sp), intent(out), dimension (0:) :: cm
    ! Calls: s2_error
end subroutine s2_sky_get_cm
Compute cm for the given sky.

Notes:

Variables:

Author: J. D. McEwen

Version: 0.1 April 2005

s2_sky_get_map_status

public function s2_sky_get_map_status (sky) result (map_status)
    type (s2_sky), intent(in) :: sky
    logical :: map_status
    ! Calls: s2_error
end function s2_sky_get_map_status
Get map_status variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_alm_status

public function s2_sky_get_alm_status (sky) result (alm_status)
    type (s2_sky), intent(in) :: sky
    logical :: alm_status
    ! Calls: s2_error
end function s2_sky_get_alm_status
Get map_status variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 August 2004

s2_sky_get_n_param

public function s2_sky_get_n_param (sky) result (n_param)
    type (s2_sky), intent(in) :: sky
    integer :: n_param
    ! Calls: s2_error
end function s2_sky_get_n_param
Get n_param variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 November 2004

s2_sky_get_param

public subroutine s2_sky_get_param (sky, param)
    type (s2_sky), intent(in) :: sky
    real (kind=s2_sp), intent(out), dimension (:) :: param
    ! Calls: s2_error
end subroutine s2_sky_get_param
Get param variable from the passed sky.

Variables:

Author: J. D. McEwen

Version: 0.1 November 2004