pynkowski.data.scalar

  1import numpy as np
  2import healpy as hp
  3
  4try:
  5    from tqdm.notebook import tqdm
  6except:
  7    tqdm = lambda x: x
  8    print('tqdm not loaded')
  9    
 10from .utils import get_theta, derivatives, second_derivatives
 11 
 12 
 13 
 14
 15class Scalar():
 16    """Class to compute Minkowski functionals (MFs) and extrema of Healpix scalar maps. It computes and stores spatial first and
 17    second derivatives of Healpix scalar maps.
 18
 19    Parameters
 20    ----------
 21    Smap : np.array
 22        The input Healpix map where all the statistics will be computed.
 23    
 24    normalise : bool, optional
 25        If True, divide input Smap by its standard deviation. Default: True.
 26    
 27    mask : np.array, optional
 28        An input Healpix mask. All the statiscal quantities will be computed only within this mask. Default: None (all sky is considered).
 29    
 30    Attributes
 31    ----------
 32    Smap : np.array
 33        The Healpix map where all the statistics are computed.
 34        
 35    nside : int
 36        The Nside parameter of Smap
 37    
 38    mask : np.array
 39        A Healpix mask whithin which all the statiscal quantities are computed.
 40        
 41    grad_phi : np.array
 42        The phi covariant first derivative of Smap (only if first or second MFs are computed).
 43
 44    grad_theta : np.array
 45        The theta covariant first derivative of Smap (only if first or second MFs are computed).
 46
 47    der_phi : np.array
 48        The phi partial first derivative of Smap (only if first or second MFs are computed). 
 49
 50    der_theta_phi : np.array
 51        The second covariant derivative wrt theta and phi of Smap (only if second MF is computed). 
 52
 53    der_phi_phi : np.array
 54        The second covariant derivative wrt phi of Smap (only if second MF is computed). 
 55
 56    der_theta_theta : np.array
 57        The second covariant derivative wrt theta of Smap (only if second MF is computed). 
 58
 59    Notes
 60    -----
 61    The derivatives are always computed full-sky regardless of input mask.
 62        
 63    """    
 64    def __init__(self, Smap, normalise=True, mask=None):
 65        """Initialise the class to compute Minkowski functionals (MFs) and extrema of Healpix scalar maps. 
 66
 67        Parameters
 68        ----------
 69        Smap : np.array
 70            The input Healpix map where all the statistics will be computed.
 71
 72        normalise : bool, optional
 73            If True, divide input Smap by its standard deviation. Default: True.
 74
 75        mask : np.array, optional
 76            An input Healpix mask. All the statiscal quantities will be computed only within this mask. Default: None (all sky is considered).
 77
 78
 79        """    
 80        self.Smap = Smap.copy()
 81        self.nside = hp.get_nside(Smap)
 82        if mask is None:
 83            self.mask = np.ones(Smap.shape, dtype='bool')
 84        else:
 85            self.mask = mask
 86            if hp.get_nside(mask) != self.nside:
 87                raise ValueError('The map and the mask have different nside')
 88        if normalise:
 89            σ = self.get_variance()
 90            self.Smap /= np.sqrt(σ)
 91        self.grad_phi = None
 92        self.der_phi_phi = None
 93    
 94    
 95    def __repr__(self):
 96        return(f'map = {self.Smap}')
 97    
 98    def get_variance(self):
 99        """compute the variance of the input Healpix scalar map within the input mask. 
100
101        Returns
102        -------
103        var : float
104            The variance of the input Healpix map within the input mask.
105
106        """    
107        return (np.var(self.Smap[self.mask]))
108
109    def set_pix(self, pixs):
110        """return the values of the input Healpix scalar map in pixels pixs. 
111
112        Parameters
113        ----------
114        pixs : np.array
115            The indices of the input map pixels whose values of the map are returned.
116
117        Returns
118        -------
119        values : np.array
120            The values of the input Healpix scalar map in pixels pixs.
121
122        """    
123        return self.Smap[pixs]
124
125    def get_gradient(self):
126        """Compute the covariant and partial first derivatives of the input Healpix scalar map. 
127        It stores:
128        
129        - first covariant derivative wrt theta in self.grad_theta
130        - first partial derivative wrt phi in self.der_phi
131        - first covariant derivative wrt phi in self.grad_phi
132
133        """    
134        S_grad = derivatives(self.Smap, gradient=True)
135        theta = get_theta(self.nside)
136        self.grad_theta = Scalar(S_grad[0], normalise=False)
137        self.der_phi = Scalar(np.cos(theta) * S_grad[1], normalise=False)
138        self.grad_phi = Scalar(S_grad[1], normalise=False)   
139        
140        
141    def get_hessian(self):
142        """compute the covariant second derivatives of the input Healpix scalar map. 
143        It stores:
144        
145        - second covariant derivative wrt theta in self.der_theta_theta
146        - second covariant derivative wrt phi in self.der_phi_phi
147        - second covariant derivative wrt theta and phi in self.der_theta_phi
148
149        """    
150        if self.grad_phi == None:
151            self.get_gradient()
152        theta = get_theta(self.nside)
153        
154        S_der_der = second_derivatives(self.grad_theta.Smap, self.der_phi.Smap)
155        
156        self.der_theta_theta = Scalar(S_der_der[0], normalise=False)
157        self.der_phi_phi = Scalar(S_der_der[1]/np.cos(theta)**2. + self.grad_theta.Smap*np.tan(theta), normalise=False)
158        self.der_theta_phi = Scalar((S_der_der[2]/np.cos(theta) - self.grad_phi.Smap * np.tan(theta)) , normalise=False)
159        
160
161    def get_κ(self, pixs):
162        """Compute the geodesic curvature multiplied by the modulus of the gradient in pixels pixs. If not already computed, it computes 
163        the first and second covariant derivatives of the input map.
164
165        Parameters
166        ----------
167        pixs : np.array
168            The indices of the input map pixels where geodesic curvature is computed.
169
170        Returns
171        -------
172        k : np.array
173            The geodesic curvature in pixels `pixs`.
174
175        """    
176        num = 2.*self.grad_theta.set_pix(pixs)*self.grad_phi.set_pix(pixs)*self.der_theta_phi.set_pix(pixs) - self.grad_phi.set_pix(pixs)**2. * self.der_theta_theta.set_pix(pixs) - self.grad_theta.set_pix(pixs)**2. * self.der_phi_phi.set_pix(pixs)
177        den = self.grad_theta.set_pix(pixs)**2. + self.grad_phi.set_pix(pixs)**2.
178        return num / den
179        
180
181    def V0_pixel(self, u):
182        """Determine where input Healpix scalar map Smap is greater than threshold u. 
183
184        Parameters
185        ----------
186        u : float
187            The threshold considered for the computation of first Minkowski functional V0.
188
189        Returns
190        -------
191        v0map : np.array
192            a bool array with the same shape as the input map, with False where input map values are lower than threshold u.
193
194        """    
195        return self.Smap>u
196    
197    def V0_iter(self, u):
198        """Compute the normalised first Minkowski functional v0 at the threshold u within the given mask. 
199
200        Parameters
201        ----------
202        u : float
203            The threshold considered for the computation of v0.
204
205        Returns
206        -------
207        v0 : np.array
208            First normalised Minkowski functional evaluated at threshold u within the given mask.
209
210        """    
211        return np.mean(self.V0_pixel(u)[self.mask])
212    
213    def V0(self, us):
214        """Compute the normalised first Minkowski functional v0 at the different thresholds us within the given mask.
215
216        Parameters
217        ----------
218        us : np.array
219            The thresholds considered for the computation of v0.
220
221        Returns
222        -------
223        v0s : np.array
224            First normalised Minkowski functional evaluated at thresholds us within the given mask.
225
226        """    
227        us = np.atleast_1d(us)
228        return np.array([self.V0_iter(u) for u in tqdm(us)])
229    
230    
231    def V1_pixel(self, u, du):
232        """Compute the modulus of the gradient where the values of the input map are between `u-du/2` and `u+du/2`. 
233
234        Parameters
235        ----------
236        u : float
237            The centered value of the bin considered for the computation.
238
239        du : float
240            The width of the bin considered for the computation.
241
242        Returns
243        -------
244        v1map : np.array
245            Modulus of the gradient where `u-du/2 < Smap < u+du/2`.
246
247        """    
248        theta = get_theta(self.nside)
249        
250        areas = np.zeros_like(theta)
251
252        mask = (u + du/2. > self.Smap) & (u - du/2. <= self.Smap) & ~(np.isclose(np.cos(theta),0, atol=1.e-2))
253        pixs = np.arange(12*self.nside**2)[mask]
254        
255        areas[pixs] = np.sqrt(self.grad_theta.set_pix(pixs)**2. + self.grad_phi.set_pix(pixs)**2.)
256        return areas/du / 4.
257    
258    def V1_iter(self, u, du):
259        """Compute the normalised second Minkowski functional v1 in the bin `u-du/2` and `u+du/2` within the given mask. 
260
261        Parameters
262        ----------
263        u : float
264            The centered value of the bin considered for the computation of v1.
265
266        du : float
267            The width of the bin considered for the computation of v1.
268
269        Returns
270        -------
271        v1 : float
272            Second normalised Minkowski functional evaluated in the bin u-du/2 and u+du/2 within the given mask.
273
274        """    
275        return np.mean(self.V1_pixel(u, du)[self.mask])
276    
277    def V1(self, us, edges=False):
278        """Compute the normalised second Minkowski functional v1 at the different thresholds us within the given mask.
279
280        Parameters
281        ----------
282        us : np.array
283            The thresholds considered for the computation of v2. See 'edges' for details.
284
285        edges : bool, optional
286            If False, us is considered as an array of uniformly distributed thresholds. 
287            If True, us is considered as a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform distributions of thresholds. 
288            In this last case, the thresholds are the central value of the given bins.
289            Default: False.
290
291        Returns
292        -------
293        v1s : np.array
294            Second normalised Minkowski functional evaluated at thresholds us within the given mask.
295
296        """
297        if self.grad_phi == None:
298            self.get_gradient()
299            
300        us = np.atleast_1d(us)
301        
302        if edges:
303            dus = (us[1:]-us[:-1])
304            us = (us[1:]+us[:-1])/2.
305        else:
306            if us.shape == (1,):
307                dus = np.array([0.1])
308            else:
309                dus = (us[1]-us[0])*np.ones(us.shape[0])           
310                if not (np.isclose(us[1:]-us[:-1], dus[0])).all():
311                    raise ValueError('The thresholds distributions is not uniform. Please set edges=True')
312            
313        return np.array([self.V1_iter(u, du) for (u, du) in zip(tqdm(us), dus)])
314
315    
316    def V2_pixel(self, u, du):
317        """Compute the geodesic curvature multiplied by the modulus of the gradient where the values of the input map are between `u-du/2` and `u+du/2`. 
318
319        Parameters
320        ----------
321        u : float
322            The centered value of the bin considered for the computation.
323
324        du : float
325            The width of the bin considered for the computation.
326
327        Returns
328        -------
329        v2map : np.array
330            Geodesic curvature multiplied by the modulus of the gradient where `u-du/2 < Smap < u+du/2`.
331
332        """           
333        theta = get_theta(self.nside)
334        
335        mask = (u + du/2. > self.Smap) & (u -du/2. <= self.Smap) & ~(np.isclose(np.cos(theta),0, atol=1.e-2))
336        pixs = np.arange(12*self.nside**2)[mask]
337        areas = np.zeros_like(theta)
338        areas[mask] = self.get_κ(pixs)
339        
340        return (areas/du / (2.*np.pi))
341    
342     
343    def V2_iter(self, u, du):
344        """Compute the normalised third Minkowski functional v2 in the bin `u-du/2` and `u+du/2` within the given mask. 
345
346        Parameters
347        ----------
348        u : float
349            The centered value of the bin considered for the computation of v2.
350
351        du : float
352            The width of the bin considered for the computation of v2.
353
354        Returns
355        -------
356        v2 : np.array
357            Third normalised Minkowski functional evaluated in the bin u-du/2 and u+du/2 within the given mask.
358
359        """    
360        return np.mean(self.V2_pixel(u, du)[self.mask])
361    
362    def V2(self, us, edges=False):
363        """Compute the normalised third Minkowski functional v2 at the different thresholds us within the given mask.
364
365        Parameters
366        ----------
367        us : np.array
368            The thresholds considered for the computation of v2. See 'edges' for details.
369
370        edges : bool, optional
371            If False, us is considered as an array of uniformly distributed thresholds. 
372            If True, us is considered as a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform distributions of thresholds. 
373            In this last case, the thresholds are the central value of the given bins.
374            Default: False.
375
376        Returns
377        -------
378        v2s : np.array
379            Third normalised Minkowski functional evaluated at thresholds us within the given mask.
380
381        """   
382        if self.grad_phi == None:
383            self.get_gradient()
384        if self.der_phi_phi == None:
385            self.get_hessian()
386            
387        us = np.atleast_1d(us)
388        
389        if edges:
390            dus = (us[1:]-us[:-1])
391            us = (us[1:]+us[:-1])/2.
392        else:
393            if us.shape == (1,):
394                dus = np.array([0.1])
395            else:
396                dus = (us[1]-us[0])*np.ones(us.shape[0])           
397                if not (np.isclose(us[1:]-us[:-1], dus[0])).all():
398                    raise ValueError('The thresholds distributions is not uniform. Please set edges=True')
399            
400        return np.array([self.V2_iter(u, du) for (u, du) in zip(tqdm(us), dus)])
401    
402
403    
404    def get_maxima(self):
405        """Find the local maxima of the input scalar map.
406
407        Returns
408        -------
409        pixels : np.array
410            Indices of the pixels which are local maxima.
411
412        values : np.array
413            Values of input map which are local maxima.
414
415        """    
416        neigh = hp.get_all_neighbours(self.nside, np.arange(12*self.nside**2))
417        
418        extT = np.concatenate([self.Smap, [np.min(self.Smap)-1.]])
419        neigh_matrix = extT[neigh]
420
421        mask = np.all(self.Smap > neigh_matrix, axis=0)
422        pixels = np.argwhere(mask).flatten()
423        values = self.Smap[pixels].flatten()
424
425        return(pixels, values)
426    
427    def get_minima(self):
428        """Find the local minima of the input scalar map.
429
430        Returns
431        -------
432        pixels : np.array
433            Indices of the pixels which are local minima
434
435        values : np.array
436            Values of input map which are local minima
437
438        """    
439        self.Smap = -self.Smap
440        pixels, values = self.get_maxima()
441        self.Smap = -self.Smap
442        
443        return(pixels, -values)
444        
445        
446        
447__all__ = ["Scalar"]
448__docformat__ = "numpy"
class Scalar:
 16class Scalar():
 17    """Class to compute Minkowski functionals (MFs) and extrema of Healpix scalar maps. It computes and stores spatial first and
 18    second derivatives of Healpix scalar maps.
 19
 20    Parameters
 21    ----------
 22    Smap : np.array
 23        The input Healpix map where all the statistics will be computed.
 24    
 25    normalise : bool, optional
 26        If True, divide input Smap by its standard deviation. Default: True.
 27    
 28    mask : np.array, optional
 29        An input Healpix mask. All the statiscal quantities will be computed only within this mask. Default: None (all sky is considered).
 30    
 31    Attributes
 32    ----------
 33    Smap : np.array
 34        The Healpix map where all the statistics are computed.
 35        
 36    nside : int
 37        The Nside parameter of Smap
 38    
 39    mask : np.array
 40        A Healpix mask whithin which all the statiscal quantities are computed.
 41        
 42    grad_phi : np.array
 43        The phi covariant first derivative of Smap (only if first or second MFs are computed).
 44
 45    grad_theta : np.array
 46        The theta covariant first derivative of Smap (only if first or second MFs are computed).
 47
 48    der_phi : np.array
 49        The phi partial first derivative of Smap (only if first or second MFs are computed). 
 50
 51    der_theta_phi : np.array
 52        The second covariant derivative wrt theta and phi of Smap (only if second MF is computed). 
 53
 54    der_phi_phi : np.array
 55        The second covariant derivative wrt phi of Smap (only if second MF is computed). 
 56
 57    der_theta_theta : np.array
 58        The second covariant derivative wrt theta of Smap (only if second MF is computed). 
 59
 60    Notes
 61    -----
 62    The derivatives are always computed full-sky regardless of input mask.
 63        
 64    """    
 65    def __init__(self, Smap, normalise=True, mask=None):
 66        """Initialise the class to compute Minkowski functionals (MFs) and extrema of Healpix scalar maps. 
 67
 68        Parameters
 69        ----------
 70        Smap : np.array
 71            The input Healpix map where all the statistics will be computed.
 72
 73        normalise : bool, optional
 74            If True, divide input Smap by its standard deviation. Default: True.
 75
 76        mask : np.array, optional
 77            An input Healpix mask. All the statiscal quantities will be computed only within this mask. Default: None (all sky is considered).
 78
 79
 80        """    
 81        self.Smap = Smap.copy()
 82        self.nside = hp.get_nside(Smap)
 83        if mask is None:
 84            self.mask = np.ones(Smap.shape, dtype='bool')
 85        else:
 86            self.mask = mask
 87            if hp.get_nside(mask) != self.nside:
 88                raise ValueError('The map and the mask have different nside')
 89        if normalise:
 90            σ = self.get_variance()
 91            self.Smap /= np.sqrt(σ)
 92        self.grad_phi = None
 93        self.der_phi_phi = None
 94    
 95    
 96    def __repr__(self):
 97        return(f'map = {self.Smap}')
 98    
 99    def get_variance(self):
100        """compute the variance of the input Healpix scalar map within the input mask. 
101
102        Returns
103        -------
104        var : float
105            The variance of the input Healpix map within the input mask.
106
107        """    
108        return (np.var(self.Smap[self.mask]))
109
110    def set_pix(self, pixs):
111        """return the values of the input Healpix scalar map in pixels pixs. 
112
113        Parameters
114        ----------
115        pixs : np.array
116            The indices of the input map pixels whose values of the map are returned.
117
118        Returns
119        -------
120        values : np.array
121            The values of the input Healpix scalar map in pixels pixs.
122
123        """    
124        return self.Smap[pixs]
125
126    def get_gradient(self):
127        """Compute the covariant and partial first derivatives of the input Healpix scalar map. 
128        It stores:
129        
130        - first covariant derivative wrt theta in self.grad_theta
131        - first partial derivative wrt phi in self.der_phi
132        - first covariant derivative wrt phi in self.grad_phi
133
134        """    
135        S_grad = derivatives(self.Smap, gradient=True)
136        theta = get_theta(self.nside)
137        self.grad_theta = Scalar(S_grad[0], normalise=False)
138        self.der_phi = Scalar(np.cos(theta) * S_grad[1], normalise=False)
139        self.grad_phi = Scalar(S_grad[1], normalise=False)   
140        
141        
142    def get_hessian(self):
143        """compute the covariant second derivatives of the input Healpix scalar map. 
144        It stores:
145        
146        - second covariant derivative wrt theta in self.der_theta_theta
147        - second covariant derivative wrt phi in self.der_phi_phi
148        - second covariant derivative wrt theta and phi in self.der_theta_phi
149
150        """    
151        if self.grad_phi == None:
152            self.get_gradient()
153        theta = get_theta(self.nside)
154        
155        S_der_der = second_derivatives(self.grad_theta.Smap, self.der_phi.Smap)
156        
157        self.der_theta_theta = Scalar(S_der_der[0], normalise=False)
158        self.der_phi_phi = Scalar(S_der_der[1]/np.cos(theta)**2. + self.grad_theta.Smap*np.tan(theta), normalise=False)
159        self.der_theta_phi = Scalar((S_der_der[2]/np.cos(theta) - self.grad_phi.Smap * np.tan(theta)) , normalise=False)
160        
161
162    def get_κ(self, pixs):
163        """Compute the geodesic curvature multiplied by the modulus of the gradient in pixels pixs. If not already computed, it computes 
164        the first and second covariant derivatives of the input map.
165
166        Parameters
167        ----------
168        pixs : np.array
169            The indices of the input map pixels where geodesic curvature is computed.
170
171        Returns
172        -------
173        k : np.array
174            The geodesic curvature in pixels `pixs`.
175
176        """    
177        num = 2.*self.grad_theta.set_pix(pixs)*self.grad_phi.set_pix(pixs)*self.der_theta_phi.set_pix(pixs) - self.grad_phi.set_pix(pixs)**2. * self.der_theta_theta.set_pix(pixs) - self.grad_theta.set_pix(pixs)**2. * self.der_phi_phi.set_pix(pixs)
178        den = self.grad_theta.set_pix(pixs)**2. + self.grad_phi.set_pix(pixs)**2.
179        return num / den
180        
181
182    def V0_pixel(self, u):
183        """Determine where input Healpix scalar map Smap is greater than threshold u. 
184
185        Parameters
186        ----------
187        u : float
188            The threshold considered for the computation of first Minkowski functional V0.
189
190        Returns
191        -------
192        v0map : np.array
193            a bool array with the same shape as the input map, with False where input map values are lower than threshold u.
194
195        """    
196        return self.Smap>u
197    
198    def V0_iter(self, u):
199        """Compute the normalised first Minkowski functional v0 at the threshold u within the given mask. 
200
201        Parameters
202        ----------
203        u : float
204            The threshold considered for the computation of v0.
205
206        Returns
207        -------
208        v0 : np.array
209            First normalised Minkowski functional evaluated at threshold u within the given mask.
210
211        """    
212        return np.mean(self.V0_pixel(u)[self.mask])
213    
214    def V0(self, us):
215        """Compute the normalised first Minkowski functional v0 at the different thresholds us within the given mask.
216
217        Parameters
218        ----------
219        us : np.array
220            The thresholds considered for the computation of v0.
221
222        Returns
223        -------
224        v0s : np.array
225            First normalised Minkowski functional evaluated at thresholds us within the given mask.
226
227        """    
228        us = np.atleast_1d(us)
229        return np.array([self.V0_iter(u) for u in tqdm(us)])
230    
231    
232    def V1_pixel(self, u, du):
233        """Compute the modulus of the gradient where the values of the input map are between `u-du/2` and `u+du/2`. 
234
235        Parameters
236        ----------
237        u : float
238            The centered value of the bin considered for the computation.
239
240        du : float
241            The width of the bin considered for the computation.
242
243        Returns
244        -------
245        v1map : np.array
246            Modulus of the gradient where `u-du/2 < Smap < u+du/2`.
247
248        """    
249        theta = get_theta(self.nside)
250        
251        areas = np.zeros_like(theta)
252
253        mask = (u + du/2. > self.Smap) & (u - du/2. <= self.Smap) & ~(np.isclose(np.cos(theta),0, atol=1.e-2))
254        pixs = np.arange(12*self.nside**2)[mask]
255        
256        areas[pixs] = np.sqrt(self.grad_theta.set_pix(pixs)**2. + self.grad_phi.set_pix(pixs)**2.)
257        return areas/du / 4.
258    
259    def V1_iter(self, u, du):
260        """Compute the normalised second Minkowski functional v1 in the bin `u-du/2` and `u+du/2` within the given mask. 
261
262        Parameters
263        ----------
264        u : float
265            The centered value of the bin considered for the computation of v1.
266
267        du : float
268            The width of the bin considered for the computation of v1.
269
270        Returns
271        -------
272        v1 : float
273            Second normalised Minkowski functional evaluated in the bin u-du/2 and u+du/2 within the given mask.
274
275        """    
276        return np.mean(self.V1_pixel(u, du)[self.mask])
277    
278    def V1(self, us, edges=False):
279        """Compute the normalised second Minkowski functional v1 at the different thresholds us within the given mask.
280
281        Parameters
282        ----------
283        us : np.array
284            The thresholds considered for the computation of v2. See 'edges' for details.
285
286        edges : bool, optional
287            If False, us is considered as an array of uniformly distributed thresholds. 
288            If True, us is considered as a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform distributions of thresholds. 
289            In this last case, the thresholds are the central value of the given bins.
290            Default: False.
291
292        Returns
293        -------
294        v1s : np.array
295            Second normalised Minkowski functional evaluated at thresholds us within the given mask.
296
297        """
298        if self.grad_phi == None:
299            self.get_gradient()
300            
301        us = np.atleast_1d(us)
302        
303        if edges:
304            dus = (us[1:]-us[:-1])
305            us = (us[1:]+us[:-1])/2.
306        else:
307            if us.shape == (1,):
308                dus = np.array([0.1])
309            else:
310                dus = (us[1]-us[0])*np.ones(us.shape[0])           
311                if not (np.isclose(us[1:]-us[:-1], dus[0])).all():
312                    raise ValueError('The thresholds distributions is not uniform. Please set edges=True')
313            
314        return np.array([self.V1_iter(u, du) for (u, du) in zip(tqdm(us), dus)])
315
316    
317    def V2_pixel(self, u, du):
318        """Compute the geodesic curvature multiplied by the modulus of the gradient where the values of the input map are between `u-du/2` and `u+du/2`. 
319
320        Parameters
321        ----------
322        u : float
323            The centered value of the bin considered for the computation.
324
325        du : float
326            The width of the bin considered for the computation.
327
328        Returns
329        -------
330        v2map : np.array
331            Geodesic curvature multiplied by the modulus of the gradient where `u-du/2 < Smap < u+du/2`.
332
333        """           
334        theta = get_theta(self.nside)
335        
336        mask = (u + du/2. > self.Smap) & (u -du/2. <= self.Smap) & ~(np.isclose(np.cos(theta),0, atol=1.e-2))
337        pixs = np.arange(12*self.nside**2)[mask]
338        areas = np.zeros_like(theta)
339        areas[mask] = self.get_κ(pixs)
340        
341        return (areas/du / (2.*np.pi))
342    
343     
344    def V2_iter(self, u, du):
345        """Compute the normalised third Minkowski functional v2 in the bin `u-du/2` and `u+du/2` within the given mask. 
346
347        Parameters
348        ----------
349        u : float
350            The centered value of the bin considered for the computation of v2.
351
352        du : float
353            The width of the bin considered for the computation of v2.
354
355        Returns
356        -------
357        v2 : np.array
358            Third normalised Minkowski functional evaluated in the bin u-du/2 and u+du/2 within the given mask.
359
360        """    
361        return np.mean(self.V2_pixel(u, du)[self.mask])
362    
363    def V2(self, us, edges=False):
364        """Compute the normalised third Minkowski functional v2 at the different thresholds us within the given mask.
365
366        Parameters
367        ----------
368        us : np.array
369            The thresholds considered for the computation of v2. See 'edges' for details.
370
371        edges : bool, optional
372            If False, us is considered as an array of uniformly distributed thresholds. 
373            If True, us is considered as a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform distributions of thresholds. 
374            In this last case, the thresholds are the central value of the given bins.
375            Default: False.
376
377        Returns
378        -------
379        v2s : np.array
380            Third normalised Minkowski functional evaluated at thresholds us within the given mask.
381
382        """   
383        if self.grad_phi == None:
384            self.get_gradient()
385        if self.der_phi_phi == None:
386            self.get_hessian()
387            
388        us = np.atleast_1d(us)
389        
390        if edges:
391            dus = (us[1:]-us[:-1])
392            us = (us[1:]+us[:-1])/2.
393        else:
394            if us.shape == (1,):
395                dus = np.array([0.1])
396            else:
397                dus = (us[1]-us[0])*np.ones(us.shape[0])           
398                if not (np.isclose(us[1:]-us[:-1], dus[0])).all():
399                    raise ValueError('The thresholds distributions is not uniform. Please set edges=True')
400            
401        return np.array([self.V2_iter(u, du) for (u, du) in zip(tqdm(us), dus)])
402    
403
404    
405    def get_maxima(self):
406        """Find the local maxima of the input scalar map.
407
408        Returns
409        -------
410        pixels : np.array
411            Indices of the pixels which are local maxima.
412
413        values : np.array
414            Values of input map which are local maxima.
415
416        """    
417        neigh = hp.get_all_neighbours(self.nside, np.arange(12*self.nside**2))
418        
419        extT = np.concatenate([self.Smap, [np.min(self.Smap)-1.]])
420        neigh_matrix = extT[neigh]
421
422        mask = np.all(self.Smap > neigh_matrix, axis=0)
423        pixels = np.argwhere(mask).flatten()
424        values = self.Smap[pixels].flatten()
425
426        return(pixels, values)
427    
428    def get_minima(self):
429        """Find the local minima of the input scalar map.
430
431        Returns
432        -------
433        pixels : np.array
434            Indices of the pixels which are local minima
435
436        values : np.array
437            Values of input map which are local minima
438
439        """    
440        self.Smap = -self.Smap
441        pixels, values = self.get_maxima()
442        self.Smap = -self.Smap
443        
444        return(pixels, -values)

Class to compute Minkowski functionals (MFs) and extrema of Healpix scalar maps. It computes and stores spatial first and second derivatives of Healpix scalar maps.

Parameters
  • Smap (np.array): The input Healpix map where all the statistics will be computed.
  • normalise (bool, optional): If True, divide input Smap by its standard deviation. Default: True.
  • mask (np.array, optional): An input Healpix mask. All the statiscal quantities will be computed only within this mask. Default: None (all sky is considered).
Attributes
  • Smap (np.array): The Healpix map where all the statistics are computed.
  • nside (int): The Nside parameter of Smap
  • mask (np.array): A Healpix mask whithin which all the statiscal quantities are computed.
  • grad_phi (np.array): The phi covariant first derivative of Smap (only if first or second MFs are computed).
  • grad_theta (np.array): The theta covariant first derivative of Smap (only if first or second MFs are computed).
  • der_phi (np.array): The phi partial first derivative of Smap (only if first or second MFs are computed).
  • der_theta_phi (np.array): The second covariant derivative wrt theta and phi of Smap (only if second MF is computed).
  • der_phi_phi (np.array): The second covariant derivative wrt phi of Smap (only if second MF is computed).
  • der_theta_theta (np.array): The second covariant derivative wrt theta of Smap (only if second MF is computed).
Notes

The derivatives are always computed full-sky regardless of input mask.

Scalar(Smap, normalise=True, mask=None)
65    def __init__(self, Smap, normalise=True, mask=None):
66        """Initialise the class to compute Minkowski functionals (MFs) and extrema of Healpix scalar maps. 
67
68        Parameters
69        ----------
70        Smap : np.array
71            The input Healpix map where all the statistics will be computed.
72
73        normalise : bool, optional
74            If True, divide input Smap by its standard deviation. Default: True.
75
76        mask : np.array, optional
77            An input Healpix mask. All the statiscal quantities will be computed only within this mask. Default: None (all sky is considered).
78
79
80        """    
81        self.Smap = Smap.copy()
82        self.nside = hp.get_nside(Smap)
83        if mask is None:
84            self.mask = np.ones(Smap.shape, dtype='bool')
85        else:
86            self.mask = mask
87            if hp.get_nside(mask) != self.nside:
88                raise ValueError('The map and the mask have different nside')
89        if normalise:
90            σ = self.get_variance()
91            self.Smap /= np.sqrt(σ)
92        self.grad_phi = None
93        self.der_phi_phi = None

Initialise the class to compute Minkowski functionals (MFs) and extrema of Healpix scalar maps.

Parameters
  • Smap (np.array): The input Healpix map where all the statistics will be computed.
  • normalise (bool, optional): If True, divide input Smap by its standard deviation. Default: True.
  • mask (np.array, optional): An input Healpix mask. All the statiscal quantities will be computed only within this mask. Default: None (all sky is considered).
def get_variance(self)
 99    def get_variance(self):
100        """compute the variance of the input Healpix scalar map within the input mask. 
101
102        Returns
103        -------
104        var : float
105            The variance of the input Healpix map within the input mask.
106
107        """    
108        return (np.var(self.Smap[self.mask]))

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.
def set_pix(self, pixs)
110    def set_pix(self, pixs):
111        """return the values of the input Healpix scalar map in pixels pixs. 
112
113        Parameters
114        ----------
115        pixs : np.array
116            The indices of the input map pixels whose values of the map are returned.
117
118        Returns
119        -------
120        values : np.array
121            The values of the input Healpix scalar map in pixels pixs.
122
123        """    
124        return self.Smap[pixs]

return the values of the input Healpix scalar map in pixels pixs.

Parameters
  • pixs (np.array): The indices of the input map pixels whose values of the map are returned.
Returns
  • values (np.array): The values of the input Healpix scalar map in pixels pixs.
def get_gradient(self)
126    def get_gradient(self):
127        """Compute the covariant and partial first derivatives of the input Healpix scalar map. 
128        It stores:
129        
130        - first covariant derivative wrt theta in self.grad_theta
131        - first partial derivative wrt phi in self.der_phi
132        - first covariant derivative wrt phi in self.grad_phi
133
134        """    
135        S_grad = derivatives(self.Smap, gradient=True)
136        theta = get_theta(self.nside)
137        self.grad_theta = Scalar(S_grad[0], normalise=False)
138        self.der_phi = Scalar(np.cos(theta) * S_grad[1], normalise=False)
139        self.grad_phi = Scalar(S_grad[1], normalise=False)   

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

  • first covariant derivative wrt theta in self.grad_theta
  • first partial derivative wrt phi in self.der_phi
  • first covariant derivative wrt phi in self.grad_phi
def get_hessian(self)
142    def get_hessian(self):
143        """compute the covariant second derivatives of the input Healpix scalar map. 
144        It stores:
145        
146        - second covariant derivative wrt theta in self.der_theta_theta
147        - second covariant derivative wrt phi in self.der_phi_phi
148        - second covariant derivative wrt theta and phi in self.der_theta_phi
149
150        """    
151        if self.grad_phi == None:
152            self.get_gradient()
153        theta = get_theta(self.nside)
154        
155        S_der_der = second_derivatives(self.grad_theta.Smap, self.der_phi.Smap)
156        
157        self.der_theta_theta = Scalar(S_der_der[0], normalise=False)
158        self.der_phi_phi = Scalar(S_der_der[1]/np.cos(theta)**2. + self.grad_theta.Smap*np.tan(theta), normalise=False)
159        self.der_theta_phi = Scalar((S_der_der[2]/np.cos(theta) - self.grad_phi.Smap * np.tan(theta)) , normalise=False)

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

  • second covariant derivative wrt theta in self.der_theta_theta
  • second covariant derivative wrt phi in self.der_phi_phi
  • second covariant derivative wrt theta and phi in self.der_theta_phi
def get_κ(self, pixs)
162    def get_κ(self, pixs):
163        """Compute the geodesic curvature multiplied by the modulus of the gradient in pixels pixs. If not already computed, it computes 
164        the first and second covariant derivatives of the input map.
165
166        Parameters
167        ----------
168        pixs : np.array
169            The indices of the input map pixels where geodesic curvature is computed.
170
171        Returns
172        -------
173        k : np.array
174            The geodesic curvature in pixels `pixs`.
175
176        """    
177        num = 2.*self.grad_theta.set_pix(pixs)*self.grad_phi.set_pix(pixs)*self.der_theta_phi.set_pix(pixs) - self.grad_phi.set_pix(pixs)**2. * self.der_theta_theta.set_pix(pixs) - self.grad_theta.set_pix(pixs)**2. * self.der_phi_phi.set_pix(pixs)
178        den = self.grad_theta.set_pix(pixs)**2. + self.grad_phi.set_pix(pixs)**2.
179        return num / den

Compute the geodesic curvature multiplied by the modulus of the gradient in pixels pixs. If not already computed, it computes the first and second covariant derivatives of the input map.

Parameters
  • pixs (np.array): The indices of the input map pixels where geodesic curvature is computed.
Returns
  • k (np.array): The geodesic curvature in pixels pixs.
def V0_pixel(self, u)
182    def V0_pixel(self, u):
183        """Determine where input Healpix scalar map Smap is greater than threshold u. 
184
185        Parameters
186        ----------
187        u : float
188            The threshold considered for the computation of first Minkowski functional V0.
189
190        Returns
191        -------
192        v0map : np.array
193            a bool array with the same shape as the input map, with False where input map values are lower than threshold u.
194
195        """    
196        return self.Smap>u

Determine where input Healpix scalar map Smap is greater than threshold u.

Parameters
  • u (float): The threshold considered for the computation of first Minkowski functional V0.
Returns
  • v0map (np.array): a bool array with the same shape as the input map, with False where input map values are lower than threshold u.
def V0_iter(self, u)
198    def V0_iter(self, u):
199        """Compute the normalised first Minkowski functional v0 at the threshold u within the given mask. 
200
201        Parameters
202        ----------
203        u : float
204            The threshold considered for the computation of v0.
205
206        Returns
207        -------
208        v0 : np.array
209            First normalised Minkowski functional evaluated at threshold u within the given mask.
210
211        """    
212        return np.mean(self.V0_pixel(u)[self.mask])

Compute the normalised first Minkowski functional v0 at the threshold u within the given mask.

Parameters
  • u (float): The threshold considered for the computation of v0.
Returns
  • v0 (np.array): First normalised Minkowski functional evaluated at threshold u within the given mask.
def V0(self, us)
214    def V0(self, us):
215        """Compute the normalised first Minkowski functional v0 at the different thresholds us within the given mask.
216
217        Parameters
218        ----------
219        us : np.array
220            The thresholds considered for the computation of v0.
221
222        Returns
223        -------
224        v0s : np.array
225            First normalised Minkowski functional evaluated at thresholds us within the given mask.
226
227        """    
228        us = np.atleast_1d(us)
229        return np.array([self.V0_iter(u) for u in tqdm(us)])

Compute the normalised first Minkowski functional v0 at the different thresholds us within the given mask.

Parameters
  • us (np.array): The thresholds considered for the computation of v0.
Returns
  • v0s (np.array): First normalised Minkowski functional evaluated at thresholds us within the given mask.
def V1_pixel(self, u, du)
232    def V1_pixel(self, u, du):
233        """Compute the modulus of the gradient where the values of the input map are between `u-du/2` and `u+du/2`. 
234
235        Parameters
236        ----------
237        u : float
238            The centered value of the bin considered for the computation.
239
240        du : float
241            The width of the bin considered for the computation.
242
243        Returns
244        -------
245        v1map : np.array
246            Modulus of the gradient where `u-du/2 < Smap < u+du/2`.
247
248        """    
249        theta = get_theta(self.nside)
250        
251        areas = np.zeros_like(theta)
252
253        mask = (u + du/2. > self.Smap) & (u - du/2. <= self.Smap) & ~(np.isclose(np.cos(theta),0, atol=1.e-2))
254        pixs = np.arange(12*self.nside**2)[mask]
255        
256        areas[pixs] = np.sqrt(self.grad_theta.set_pix(pixs)**2. + self.grad_phi.set_pix(pixs)**2.)
257        return areas/du / 4.

Compute the modulus of the gradient where the values of the input map are between u-du/2 and u+du/2.

Parameters
  • u (float): The centered value of the bin considered for the computation.
  • du (float): The width of the bin considered for the computation.
Returns
  • v1map (np.array): Modulus of the gradient where u-du/2 < Smap < u+du/2.
def V1_iter(self, u, du)
259    def V1_iter(self, u, du):
260        """Compute the normalised second Minkowski functional v1 in the bin `u-du/2` and `u+du/2` within the given mask. 
261
262        Parameters
263        ----------
264        u : float
265            The centered value of the bin considered for the computation of v1.
266
267        du : float
268            The width of the bin considered for the computation of v1.
269
270        Returns
271        -------
272        v1 : float
273            Second normalised Minkowski functional evaluated in the bin u-du/2 and u+du/2 within the given mask.
274
275        """    
276        return np.mean(self.V1_pixel(u, du)[self.mask])

Compute the normalised second Minkowski functional v1 in the bin u-du/2 and u+du/2 within the given mask.

Parameters
  • u (float): The centered value of the bin considered for the computation of v1.
  • du (float): The width of the bin considered for the computation of v1.
Returns
  • v1 (float): Second normalised Minkowski functional evaluated in the bin u-du/2 and u+du/2 within the given mask.
def V1(self, us, edges=False)
278    def V1(self, us, edges=False):
279        """Compute the normalised second Minkowski functional v1 at the different thresholds us within the given mask.
280
281        Parameters
282        ----------
283        us : np.array
284            The thresholds considered for the computation of v2. See 'edges' for details.
285
286        edges : bool, optional
287            If False, us is considered as an array of uniformly distributed thresholds. 
288            If True, us is considered as a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform distributions of thresholds. 
289            In this last case, the thresholds are the central value of the given bins.
290            Default: False.
291
292        Returns
293        -------
294        v1s : np.array
295            Second normalised Minkowski functional evaluated at thresholds us within the given mask.
296
297        """
298        if self.grad_phi == None:
299            self.get_gradient()
300            
301        us = np.atleast_1d(us)
302        
303        if edges:
304            dus = (us[1:]-us[:-1])
305            us = (us[1:]+us[:-1])/2.
306        else:
307            if us.shape == (1,):
308                dus = np.array([0.1])
309            else:
310                dus = (us[1]-us[0])*np.ones(us.shape[0])           
311                if not (np.isclose(us[1:]-us[:-1], dus[0])).all():
312                    raise ValueError('The thresholds distributions is not uniform. Please set edges=True')
313            
314        return np.array([self.V1_iter(u, du) for (u, du) in zip(tqdm(us), dus)])

Compute the normalised second Minkowski functional v1 at the different thresholds us within the given mask.

Parameters
  • us (np.array): The thresholds considered for the computation of v2. See 'edges' for details.
  • edges (bool, optional): If False, us is considered as an array of uniformly distributed thresholds. If True, us is considered as a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform distributions of thresholds. In this last case, the thresholds are the central value of the given bins. Default: False.
Returns
  • v1s (np.array): Second normalised Minkowski functional evaluated at thresholds us within the given mask.
def V2_pixel(self, u, du)
317    def V2_pixel(self, u, du):
318        """Compute the geodesic curvature multiplied by the modulus of the gradient where the values of the input map are between `u-du/2` and `u+du/2`. 
319
320        Parameters
321        ----------
322        u : float
323            The centered value of the bin considered for the computation.
324
325        du : float
326            The width of the bin considered for the computation.
327
328        Returns
329        -------
330        v2map : np.array
331            Geodesic curvature multiplied by the modulus of the gradient where `u-du/2 < Smap < u+du/2`.
332
333        """           
334        theta = get_theta(self.nside)
335        
336        mask = (u + du/2. > self.Smap) & (u -du/2. <= self.Smap) & ~(np.isclose(np.cos(theta),0, atol=1.e-2))
337        pixs = np.arange(12*self.nside**2)[mask]
338        areas = np.zeros_like(theta)
339        areas[mask] = self.get_κ(pixs)
340        
341        return (areas/du / (2.*np.pi))

Compute the geodesic curvature multiplied by the modulus of the gradient where the values of the input map are between u-du/2 and u+du/2.

Parameters
  • u (float): The centered value of the bin considered for the computation.
  • du (float): The width of the bin considered for the computation.
Returns
  • v2map (np.array): Geodesic curvature multiplied by the modulus of the gradient where u-du/2 < Smap < u+du/2.
def V2_iter(self, u, du)
344    def V2_iter(self, u, du):
345        """Compute the normalised third Minkowski functional v2 in the bin `u-du/2` and `u+du/2` within the given mask. 
346
347        Parameters
348        ----------
349        u : float
350            The centered value of the bin considered for the computation of v2.
351
352        du : float
353            The width of the bin considered for the computation of v2.
354
355        Returns
356        -------
357        v2 : np.array
358            Third normalised Minkowski functional evaluated in the bin u-du/2 and u+du/2 within the given mask.
359
360        """    
361        return np.mean(self.V2_pixel(u, du)[self.mask])

Compute the normalised third Minkowski functional v2 in the bin u-du/2 and u+du/2 within the given mask.

Parameters
  • u (float): The centered value of the bin considered for the computation of v2.
  • du (float): The width of the bin considered for the computation of v2.
Returns
  • v2 (np.array): Third normalised Minkowski functional evaluated in the bin u-du/2 and u+du/2 within the given mask.
def V2(self, us, edges=False)
363    def V2(self, us, edges=False):
364        """Compute the normalised third Minkowski functional v2 at the different thresholds us within the given mask.
365
366        Parameters
367        ----------
368        us : np.array
369            The thresholds considered for the computation of v2. See 'edges' for details.
370
371        edges : bool, optional
372            If False, us is considered as an array of uniformly distributed thresholds. 
373            If True, us is considered as a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform distributions of thresholds. 
374            In this last case, the thresholds are the central value of the given bins.
375            Default: False.
376
377        Returns
378        -------
379        v2s : np.array
380            Third normalised Minkowski functional evaluated at thresholds us within the given mask.
381
382        """   
383        if self.grad_phi == None:
384            self.get_gradient()
385        if self.der_phi_phi == None:
386            self.get_hessian()
387            
388        us = np.atleast_1d(us)
389        
390        if edges:
391            dus = (us[1:]-us[:-1])
392            us = (us[1:]+us[:-1])/2.
393        else:
394            if us.shape == (1,):
395                dus = np.array([0.1])
396            else:
397                dus = (us[1]-us[0])*np.ones(us.shape[0])           
398                if not (np.isclose(us[1:]-us[:-1], dus[0])).all():
399                    raise ValueError('The thresholds distributions is not uniform. Please set edges=True')
400            
401        return np.array([self.V2_iter(u, du) for (u, du) in zip(tqdm(us), dus)])

Compute the normalised third Minkowski functional v2 at the different thresholds us within the given mask.

Parameters
  • us (np.array): The thresholds considered for the computation of v2. See 'edges' for details.
  • edges (bool, optional): If False, us is considered as an array of uniformly distributed thresholds. If True, us is considered as a monotonically increasing array of bin edges, including the rightmost edge, allowing for non-uniform distributions of thresholds. In this last case, the thresholds are the central value of the given bins. Default: False.
Returns
  • v2s (np.array): Third normalised Minkowski functional evaluated at thresholds us within the given mask.
def get_maxima(self)
405    def get_maxima(self):
406        """Find the local maxima of the input scalar map.
407
408        Returns
409        -------
410        pixels : np.array
411            Indices of the pixels which are local maxima.
412
413        values : np.array
414            Values of input map which are local maxima.
415
416        """    
417        neigh = hp.get_all_neighbours(self.nside, np.arange(12*self.nside**2))
418        
419        extT = np.concatenate([self.Smap, [np.min(self.Smap)-1.]])
420        neigh_matrix = extT[neigh]
421
422        mask = np.all(self.Smap > neigh_matrix, axis=0)
423        pixels = np.argwhere(mask).flatten()
424        values = self.Smap[pixels].flatten()
425
426        return(pixels, values)

Find the local maxima of the input scalar map.

Returns
  • pixels (np.array): Indices of the pixels which are local maxima.
  • values (np.array): Values of input map which are local maxima.
def get_minima(self)
428    def get_minima(self):
429        """Find the local minima of the input scalar map.
430
431        Returns
432        -------
433        pixels : np.array
434            Indices of the pixels which are local minima
435
436        values : np.array
437            Values of input map which are local minima
438
439        """    
440        self.Smap = -self.Smap
441        pixels, values = self.get_maxima()
442        self.Smap = -self.Smap
443        
444        return(pixels, -values)

Find the local minima of the input scalar map.

Returns
  • pixels (np.array): Indices of the pixels which are local minima
  • values (np.array): Values of input map which are local minima