pynkowski.data.healpix

This submodule contains the class for scalar fields in HEALPix format and an easy interface for $P^2$ maps in this format.

  1'''This submodule contains the class for scalar fields in HEALPix format and an easy interface for $P^2$ maps in this format.'''
  2import numpy as np
  3import healpy as hp
  4from .base_da import DataField
  5from .utils_da import get_theta, healpix_derivatives, healpix_second_derivatives
  6
  7def _hotspot(healpix_map):
  8    """Find the local maxima of the input HEALPix map.
  9
 10    Returns
 11    -------
 12    pixels : np.array
 13        Indices of the pixels which are local maxima.
 14
 15    values : np.array
 16        Values of input map which are local maxima.
 17
 18    """
 19    nside = hp.npix2nside(healpix_map.shape[0])
 20    neigh = hp.get_all_neighbours(nside, np.arange(12*nside**2))
 21    
 22    extT = np.concatenate([healpix_map, [np.min(healpix_map)-1.]])
 23    neigh_matrix = extT[neigh]
 24
 25    mask = np.all(healpix_map > neigh_matrix, axis=0)
 26    pixels = np.argwhere(mask).flatten()
 27    values = healpix_map[pixels]
 28    return(pixels, values)
 29
 30
 31class Healpix(DataField):
 32    """Class for spherical scalar fields in HEALPix format.
 33
 34    Parameters
 35    ----------
 36    field : np.array
 37        Values of the scalar field in HEALPix format in RING scheme.
 38        
 39    normalise : bool, optional
 40        If `True`, the map is normalised to have unit variance. Note: the mean of the map is not forced to be 0.
 41        Default: `True`.
 42
 43    mask : np.array or None, optional
 44        Mask where the field is considered. It is a bool array of the same shape that `field`.
 45        Default: all data is included.
 46        
 47    Attributes
 48    ----------
 49    field : np.array
 50        Data of the field as a HEALPix map in RING scheme.
 51        
 52    nside : int
 53        Parameter `nside` of the map. The number of pixels is `12*nside**2`.
 54    
 55    dim : int
 56        Dimension of the space where the field is defined. In this case, the space is the sphere and this is 2.
 57        
 58    name : str
 59        Name of the field. In this case, "HEALPix map"
 60    
 61    first_der : np.array or None
 62        First **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(dim, field.shape)`.
 63    
 64    second_der : np.array or None
 65        Second **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(dim*(dim+1)/2, field.shape)`.
 66        The order of the derivatives is diagonal first, e.g. in `dim=3`: `11`, `22`, `33`, `12`, `13`, `23`.
 67        
 68    mask : np.array
 69        Mask where the field is considered. It is a bool array of the same shape that `field`.
 70        
 71    """   
 72    def __init__(self, field, normalise=True, mask=None):
 73        super().__init__(field, 2, name="HEALPix map", mask=mask)
 74        self.nside = hp.get_nside(self.field)
 75        if hp.get_nside(self.mask) != self.nside:
 76            raise ValueError('The map and the mask have different nside')
 77        
 78        if normalise:
 79            # self.field -= self.field.mean()
 80            σ2 = self.get_variance()
 81            self.field /= np.sqrt(σ2)
 82            
 83    def get_variance(self):
 84        """Compute the variance of the Healpix scalar map within the mask. 
 85
 86        Returns
 87        -------
 88        var : float
 89            The variance of the Healpix map within the mask.
 90
 91        """    
 92        return (np.var(self.field[self.mask]))
 93    
 94    def get_first_der(self, lmax=None):
 95        """Compute the covariant first derivatives of the input Healpix scalar map. 
 96        It stores:
 97        
 98        - first covariant derivative wrt θ in self.first_der[0]
 99        - first covariant derivative wrt ϕ in self.first_der[1]
100        
101        Parameters
102        ----------
103        lmax : int or None, optional
104            Maximum multipole used in the computation of the derivatives.
105            Default: 3*nside - 1
106        """    
107        self.first_der = np.array(healpix_derivatives(self.field, gradient=True, lmax=lmax))  # order: θ, ϕ
108            
109    def get_second_der(self, lmax=None):
110        """Compute the covariant second derivatives of the input Healpix scalar map. 
111        It stores:
112        
113        - second covariant derivative wrt θθ in self.second_der[0]
114        - second covariant derivative wrt ϕϕ in self.second_der[1]
115        - second covariant derivative wrt θϕ in self.second_der[2]
116
117        Parameters
118        ----------
119        lmax : int or None, optional
120            Maximum multipole used in the computation of the derivatives.
121            Default: 3*nside - 1
122        """    
123        if self.first_der is None:
124            self.get_first_der()
125        theta = get_theta(self.nside)
126        
127        second_partial = healpix_second_derivatives(self.first_der[0], np.cos(theta) * self.first_der[1], lmax=lmax)  #order θθ, ϕϕ, θϕ
128        
129        self.second_der = np.array([second_partial[0],
130                                    second_partial[1]/np.cos(theta)**2. + self.first_der[0]*np.tan(theta),
131                                    (second_partial[2]/np.cos(theta) - self.first_der[1] * np.tan(theta))])  #order θθ, ϕϕ, θϕ
132
133    def maxima_list(self):
134        """Compute the values of the local maxima of the HEALPix map.
135
136        Returns
137        -------
138        values : np.array
139            Values of the map which are local maxima.
140
141        """
142        pixels, values = _hotspot(self.field)
143        return values[self.mask[pixels]]
144
145    def minima_list(self):
146        """Compute the values of the local minima of the HEALPix map.
147
148        Returns
149        -------
150        values : np.array
151            Values of the map which are local minima.
152
153        """
154        pixels, values = _hotspot(-self.field)
155        return -values[self.mask[pixels]]
156
157
158class HealpixP2(Healpix):
159    """Class for spherical scalar fields in HEALPix format.
160
161    Parameters
162    ----------
163    Q : np.array
164        Values of the Q component in HEALPix format in RING scheme.
165        
166    U : np.array
167        Values of the U component in HEALPix format in RING scheme.
168        
169    normalise : bool, optional
170        If `True`, Q and U are normalised to unit variance. Note: the normalisation is computed as the average of both variances, the mean of both maps are not forced to be 0.
171    
172    mask : np.array or None, optional
173        Mask where the field if considered. It is a bool array of the same shape that `field`.
174        Default: all data is included.
175        
176    Attributes
177    ----------
178    Q : np.array
179        Data of the $Q$ component as a HEALPix map in RING scheme.
180        
181    U : np.array
182        Data of the $U$ component as a HEALPix map in RING scheme.
183        
184    field : np.array
185        Data of the $P^2$ field as a HEALPix map in RING scheme.
186        
187    nside : int
188        Parameter `nside` of the map. The number of pixels is `12*nside**2`.
189    
190    dim : int
191        Dimension of the space where the field is defined. In this case, the space is the sphere and this is 2.
192        
193    name : str
194        Name of the field. In this case, "HEALPix map"
195    
196    first_der : np.array or None
197        First **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(dim, field.shape)`.
198    
199    second_der : np.array or None
200        Second **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(dim*(dim+1)/2, field.shape)`.
201        The order of the derivatives is diagonal first, e.g. in `dim=3`: `11`, `22`, `33`, `12`, `13`, `23`.
202        
203    mask : np.array
204        Mask where the field if considered. It is a bool array of the same shape that `field`.
205        
206    """   
207    def __init__(self, Q, U, normalise=True, mask=None):
208        self.Q = Q
209        self.U = U
210        self.mask = mask
211        if normalise:
212            # self.Q -= self.Q.mean()
213            # self.U -= self.U.mean()
214            σ2 = self.get_variance()
215            self.Q /= np.sqrt(σ2)
216            self.U /= np.sqrt(σ2)
217        
218        super().__init__(self.Q**2. + self.U**2., normalise=False, mask=mask)
219            
220    def get_variance(self):
221        """Compute the variance of the input Healpix scalar map within the input mask. 
222
223        Returns
224        -------
225        var : float
226            The variance of the input Healpix map within the input mask.
227
228        """    
229        return (np.var(self.Q[self.mask]) + np.var(self.U[self.mask]))/2.
230    
231
232__all__ = ["Healpix", "HealpixP2"]
233__docformat__ = "numpy"
class Healpix(pynkowski.data.base_da.DataField):
 32class Healpix(DataField):
 33    """Class for spherical scalar fields in HEALPix format.
 34
 35    Parameters
 36    ----------
 37    field : np.array
 38        Values of the scalar field in HEALPix format in RING scheme.
 39        
 40    normalise : bool, optional
 41        If `True`, the map is normalised to have unit variance. Note: the mean of the map is not forced to be 0.
 42        Default: `True`.
 43
 44    mask : np.array or None, optional
 45        Mask where the field is considered. It is a bool array of the same shape that `field`.
 46        Default: all data is included.
 47        
 48    Attributes
 49    ----------
 50    field : np.array
 51        Data of the field as a HEALPix map in RING scheme.
 52        
 53    nside : int
 54        Parameter `nside` of the map. The number of pixels is `12*nside**2`.
 55    
 56    dim : int
 57        Dimension of the space where the field is defined. In this case, the space is the sphere and this is 2.
 58        
 59    name : str
 60        Name of the field. In this case, "HEALPix map"
 61    
 62    first_der : np.array or None
 63        First **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(dim, field.shape)`.
 64    
 65    second_der : np.array or None
 66        Second **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(dim*(dim+1)/2, field.shape)`.
 67        The order of the derivatives is diagonal first, e.g. in `dim=3`: `11`, `22`, `33`, `12`, `13`, `23`.
 68        
 69    mask : np.array
 70        Mask where the field is considered. It is a bool array of the same shape that `field`.
 71        
 72    """   
 73    def __init__(self, field, normalise=True, mask=None):
 74        super().__init__(field, 2, name="HEALPix map", mask=mask)
 75        self.nside = hp.get_nside(self.field)
 76        if hp.get_nside(self.mask) != self.nside:
 77            raise ValueError('The map and the mask have different nside')
 78        
 79        if normalise:
 80            # self.field -= self.field.mean()
 81            σ2 = self.get_variance()
 82            self.field /= np.sqrt(σ2)
 83            
 84    def get_variance(self):
 85        """Compute the variance of the Healpix scalar map within the mask. 
 86
 87        Returns
 88        -------
 89        var : float
 90            The variance of the Healpix map within the mask.
 91
 92        """    
 93        return (np.var(self.field[self.mask]))
 94    
 95    def get_first_der(self, lmax=None):
 96        """Compute the covariant first derivatives of the input Healpix scalar map. 
 97        It stores:
 98        
 99        - first covariant derivative wrt θ in self.first_der[0]
100        - first covariant derivative wrt ϕ in self.first_der[1]
101        
102        Parameters
103        ----------
104        lmax : int or None, optional
105            Maximum multipole used in the computation of the derivatives.
106            Default: 3*nside - 1
107        """    
108        self.first_der = np.array(healpix_derivatives(self.field, gradient=True, lmax=lmax))  # order: θ, ϕ
109            
110    def get_second_der(self, lmax=None):
111        """Compute the covariant second derivatives of the input Healpix scalar map. 
112        It stores:
113        
114        - second covariant derivative wrt θθ in self.second_der[0]
115        - second covariant derivative wrt ϕϕ in self.second_der[1]
116        - second covariant derivative wrt θϕ in self.second_der[2]
117
118        Parameters
119        ----------
120        lmax : int or None, optional
121            Maximum multipole used in the computation of the derivatives.
122            Default: 3*nside - 1
123        """    
124        if self.first_der is None:
125            self.get_first_der()
126        theta = get_theta(self.nside)
127        
128        second_partial = healpix_second_derivatives(self.first_der[0], np.cos(theta) * self.first_der[1], lmax=lmax)  #order θθ, ϕϕ, θϕ
129        
130        self.second_der = np.array([second_partial[0],
131                                    second_partial[1]/np.cos(theta)**2. + self.first_der[0]*np.tan(theta),
132                                    (second_partial[2]/np.cos(theta) - self.first_der[1] * np.tan(theta))])  #order θθ, ϕϕ, θϕ
133
134    def maxima_list(self):
135        """Compute the values of the local maxima of the HEALPix map.
136
137        Returns
138        -------
139        values : np.array
140            Values of the map which are local maxima.
141
142        """
143        pixels, values = _hotspot(self.field)
144        return values[self.mask[pixels]]
145
146    def minima_list(self):
147        """Compute the values of the local minima of the HEALPix map.
148
149        Returns
150        -------
151        values : np.array
152            Values of the map which are local minima.
153
154        """
155        pixels, values = _hotspot(-self.field)
156        return -values[self.mask[pixels]]

Class for spherical scalar fields in HEALPix format.

Parameters
  • field (np.array): Values of the scalar field in HEALPix format in RING scheme.
  • normalise (bool, optional): If True, the map is normalised to have unit variance. Note: the mean of the map is not forced to be 0. Default: True.
  • mask (np.array or None, optional): Mask where the field is considered. It is a bool array of the same shape that field. Default: all data is included.
Attributes
  • field (np.array): Data of the field as a HEALPix map in RING scheme.
  • nside (int): Parameter nside of the map. The number of pixels is 12*nside**2.
  • dim (int): Dimension of the space where the field is defined. In this case, the space is the sphere and this is 2.
  • name (str): Name of the field. In this case, "HEALPix map"
  • first_der (np.array or None): First covariant derivatives of the field in an orthonormal basis of the space. Same structure as field, and shape (dim, field.shape).
  • second_der (np.array or None): Second covariant derivatives of the field in an orthonormal basis of the space. Same structure as field, and shape (dim*(dim+1)/2, field.shape). The order of the derivatives is diagonal first, e.g. in dim=3: 11, 22, 33, 12, 13, 23.
  • mask (np.array): Mask where the field is considered. It is a bool array of the same shape that field.
Healpix(field, normalise=True, mask=None)
73    def __init__(self, field, normalise=True, mask=None):
74        super().__init__(field, 2, name="HEALPix map", mask=mask)
75        self.nside = hp.get_nside(self.field)
76        if hp.get_nside(self.mask) != self.nside:
77            raise ValueError('The map and the mask have different nside')
78        
79        if normalise:
80            # self.field -= self.field.mean()
81            σ2 = self.get_variance()
82            self.field /= np.sqrt(σ2)
def get_variance(self):
84    def get_variance(self):
85        """Compute the variance of the Healpix scalar map within the mask. 
86
87        Returns
88        -------
89        var : float
90            The variance of the Healpix map within the mask.
91
92        """    
93        return (np.var(self.field[self.mask]))

Compute the variance of the Healpix scalar map within the mask.

Returns
  • var (float): The variance of the Healpix map within the mask.
def get_first_der(self, lmax=None):
 95    def get_first_der(self, lmax=None):
 96        """Compute the covariant first derivatives of the input Healpix scalar map. 
 97        It stores:
 98        
 99        - first covariant derivative wrt θ in self.first_der[0]
100        - first covariant derivative wrt ϕ in self.first_der[1]
101        
102        Parameters
103        ----------
104        lmax : int or None, optional
105            Maximum multipole used in the computation of the derivatives.
106            Default: 3*nside - 1
107        """    
108        self.first_der = np.array(healpix_derivatives(self.field, gradient=True, lmax=lmax))  # order: θ, ϕ

Compute the covariant first derivatives of the input Healpix scalar map. It stores:

  • first covariant derivative wrt θ in self.first_der[0]
  • first covariant derivative wrt ϕ in self.first_der[1]
Parameters
  • lmax (int or None, optional): Maximum multipole used in the computation of the derivatives. Default: 3*nside - 1
def get_second_der(self, lmax=None):
110    def get_second_der(self, lmax=None):
111        """Compute the covariant second derivatives of the input Healpix scalar map. 
112        It stores:
113        
114        - second covariant derivative wrt θθ in self.second_der[0]
115        - second covariant derivative wrt ϕϕ in self.second_der[1]
116        - second covariant derivative wrt θϕ in self.second_der[2]
117
118        Parameters
119        ----------
120        lmax : int or None, optional
121            Maximum multipole used in the computation of the derivatives.
122            Default: 3*nside - 1
123        """    
124        if self.first_der is None:
125            self.get_first_der()
126        theta = get_theta(self.nside)
127        
128        second_partial = healpix_second_derivatives(self.first_der[0], np.cos(theta) * self.first_der[1], lmax=lmax)  #order θθ, ϕϕ, θϕ
129        
130        self.second_der = np.array([second_partial[0],
131                                    second_partial[1]/np.cos(theta)**2. + self.first_der[0]*np.tan(theta),
132                                    (second_partial[2]/np.cos(theta) - self.first_der[1] * np.tan(theta))])  #order θθ, ϕϕ, θϕ

Compute the covariant second derivatives of the input Healpix scalar map. It stores:

  • second covariant derivative wrt θθ in self.second_der[0]
  • second covariant derivative wrt ϕϕ in self.second_der[1]
  • second covariant derivative wrt θϕ in self.second_der[2]
Parameters
  • lmax (int or None, optional): Maximum multipole used in the computation of the derivatives. Default: 3*nside - 1
def maxima_list(self):
134    def maxima_list(self):
135        """Compute the values of the local maxima of the HEALPix map.
136
137        Returns
138        -------
139        values : np.array
140            Values of the map which are local maxima.
141
142        """
143        pixels, values = _hotspot(self.field)
144        return values[self.mask[pixels]]

Compute the values of the local maxima of the HEALPix map.

Returns
  • values (np.array): Values of the map which are local maxima.
def minima_list(self):
146    def minima_list(self):
147        """Compute the values of the local minima of the HEALPix map.
148
149        Returns
150        -------
151        values : np.array
152            Values of the map which are local minima.
153
154        """
155        pixels, values = _hotspot(-self.field)
156        return -values[self.mask[pixels]]

Compute the values of the local minima of the HEALPix map.

Returns
  • values (np.array): Values of the map which are local minima.
class HealpixP2(Healpix):
159class HealpixP2(Healpix):
160    """Class for spherical scalar fields in HEALPix format.
161
162    Parameters
163    ----------
164    Q : np.array
165        Values of the Q component in HEALPix format in RING scheme.
166        
167    U : np.array
168        Values of the U component in HEALPix format in RING scheme.
169        
170    normalise : bool, optional
171        If `True`, Q and U are normalised to unit variance. Note: the normalisation is computed as the average of both variances, the mean of both maps are not forced to be 0.
172    
173    mask : np.array or None, optional
174        Mask where the field if considered. It is a bool array of the same shape that `field`.
175        Default: all data is included.
176        
177    Attributes
178    ----------
179    Q : np.array
180        Data of the $Q$ component as a HEALPix map in RING scheme.
181        
182    U : np.array
183        Data of the $U$ component as a HEALPix map in RING scheme.
184        
185    field : np.array
186        Data of the $P^2$ field as a HEALPix map in RING scheme.
187        
188    nside : int
189        Parameter `nside` of the map. The number of pixels is `12*nside**2`.
190    
191    dim : int
192        Dimension of the space where the field is defined. In this case, the space is the sphere and this is 2.
193        
194    name : str
195        Name of the field. In this case, "HEALPix map"
196    
197    first_der : np.array or None
198        First **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(dim, field.shape)`.
199    
200    second_der : np.array or None
201        Second **covariant** derivatives of the field in an orthonormal basis of the space. Same structure as `field`, and shape `(dim*(dim+1)/2, field.shape)`.
202        The order of the derivatives is diagonal first, e.g. in `dim=3`: `11`, `22`, `33`, `12`, `13`, `23`.
203        
204    mask : np.array
205        Mask where the field if considered. It is a bool array of the same shape that `field`.
206        
207    """   
208    def __init__(self, Q, U, normalise=True, mask=None):
209        self.Q = Q
210        self.U = U
211        self.mask = mask
212        if normalise:
213            # self.Q -= self.Q.mean()
214            # self.U -= self.U.mean()
215            σ2 = self.get_variance()
216            self.Q /= np.sqrt(σ2)
217            self.U /= np.sqrt(σ2)
218        
219        super().__init__(self.Q**2. + self.U**2., normalise=False, mask=mask)
220            
221    def get_variance(self):
222        """Compute the variance of the input Healpix scalar map within the input mask. 
223
224        Returns
225        -------
226        var : float
227            The variance of the input Healpix map within the input mask.
228
229        """    
230        return (np.var(self.Q[self.mask]) + np.var(self.U[self.mask]))/2.

Class for spherical scalar fields in HEALPix format.

Parameters
  • Q (np.array): Values of the Q component in HEALPix format in RING scheme.
  • U (np.array): Values of the U component in HEALPix format in RING scheme.
  • normalise (bool, optional): If True, Q and U are normalised to unit variance. Note: the normalisation is computed as the average of both variances, the mean of both maps are not forced to be 0.
  • mask (np.array or None, optional): Mask where the field if considered. It is a bool array of the same shape that field. Default: all data is included.
Attributes
  • Q (np.array): Data of the $Q$ component as a HEALPix map in RING scheme.
  • U (np.array): Data of the $U$ component as a HEALPix map in RING scheme.
  • field (np.array): Data of the $P^2$ field as a HEALPix map in RING scheme.
  • nside (int): Parameter nside of the map. The number of pixels is 12*nside**2.
  • dim (int): Dimension of the space where the field is defined. In this case, the space is the sphere and this is 2.
  • name (str): Name of the field. In this case, "HEALPix map"
  • first_der (np.array or None): First covariant derivatives of the field in an orthonormal basis of the space. Same structure as field, and shape (dim, field.shape).
  • second_der (np.array or None): Second covariant derivatives of the field in an orthonormal basis of the space. Same structure as field, and shape (dim*(dim+1)/2, field.shape). The order of the derivatives is diagonal first, e.g. in dim=3: 11, 22, 33, 12, 13, 23.
  • mask (np.array): Mask where the field if considered. It is a bool array of the same shape that field.
HealpixP2(Q, U, normalise=True, mask=None)
208    def __init__(self, Q, U, normalise=True, mask=None):
209        self.Q = Q
210        self.U = U
211        self.mask = mask
212        if normalise:
213            # self.Q -= self.Q.mean()
214            # self.U -= self.U.mean()
215            σ2 = self.get_variance()
216            self.Q /= np.sqrt(σ2)
217            self.U /= np.sqrt(σ2)
218        
219        super().__init__(self.Q**2. + self.U**2., normalise=False, mask=mask)
def get_variance(self):
221    def get_variance(self):
222        """Compute the variance of the input Healpix scalar map within the input mask. 
223
224        Returns
225        -------
226        var : float
227            The variance of the input Healpix map within the input mask.
228
229        """    
230        return (np.var(self.Q[self.mask]) + np.var(self.U[self.mask]))/2.

Compute the variance of the input Healpix scalar map within the input mask.

Returns
  • var (float): The variance of the input Healpix map within the input mask.